64const std::string G4GSMottCorrection::gElemSymbols[] = {
"H",
"He",
"Li",
"Be",
"B" ,
65 "C" ,
"N" ,
"O" ,
"F" ,
"Ne",
"Na",
"Mg",
"Al",
"Si",
"P" ,
"S",
"Cl",
"Ar",
"K" ,
"Ca",
"Sc",
66 "Ti",
"V" ,
"Cr",
"Mn",
"Fe",
"Co",
"Ni",
"Cu",
"Zn",
"Ga",
"Ge",
"As",
"Se",
"Br",
"Kr",
"Rb",
67 "Sr",
"Y" ,
"Zr",
"Nb",
"Mo",
"Tc",
"Ru",
"Rh",
"Pd",
"Ag",
"Cd",
"In",
"Sn",
"Sb",
"Te",
"I" ,
68 "Xe",
"Cs",
"Ba",
"La",
"Ce",
"Pr",
"Nd",
"Pm",
"Sm",
"Eu",
"Gd",
"Tb",
"Dy",
"Ho",
"Er",
"Tm",
69 "Yb",
"Lu",
"Hf",
"Ta",
"W" ,
"Re",
"Os",
"Ir",
"Pt",
"Au",
"Hg",
"Tl",
"Pb",
"Bi",
"Po",
"At",
70 "Rn",
"Fr",
"Ra",
"Ac",
"Th",
"Pa",
"U" ,
"Np",
"Pu",
"Am",
"Cm",
"Bk",
"Cf"};
74 fMaxEkin = CLHEP::electron_mass_c2*(1./std::sqrt(1.-gMaxBeta2)-1.);
75 fLogMinEkin =
G4Log(gMinEkin);
76 fInvLogDelEkin = (gNumEkin-gNumBeta2)/
G4Log(gMidEkin/gMinEkin);
77 G4double pt2 = gMidEkin*(gMidEkin+2.0*CLHEP::electron_mass_c2);
78 fMinBeta2 = pt2/(pt2+CLHEP::electron_mass_c2*CLHEP::electron_mass_c2);
79 fInvDelBeta2 = (gNumBeta2-1.)/(gMaxBeta2-fMinBeta2);
80 fInvDelDelta = (gNumDelta-1.)/gMaxDelta;
81 fInvDelAngle = gNumAngle-1.;
86 ClearMCDataPerElement();
87 ClearMCDataPerMaterial();
93 G4int ekinIndxLow = 0;
95 if (beta2>=gMaxBeta2) {
96 ekinIndxLow = gNumEkin - 1;
98 }
else if (beta2>=fMinBeta2) {
99 remRfaction = (beta2 - fMinBeta2) * fInvDelBeta2;
100 ekinIndxLow = (
G4int)remRfaction;
101 remRfaction -= ekinIndxLow;
102 ekinIndxLow += (gNumEkin - gNumBeta2);
103 }
else if (logekin>=fLogMinEkin) {
104 remRfaction = (logekin - fLogMinEkin) * fInvLogDelEkin;
105 ekinIndxLow = (
G4int)remRfaction;
106 remRfaction -= ekinIndxLow;
109 DataPerEkin *perEkinLow = fMCDataPerMaterial[matindx]->fDataPerEkin[ekinIndxLow];
110 mcToScr = perEkinLow->fMCScreening;
111 mcToQ1 = perEkinLow->fMCFirstMoment;
112 mcToG2PerG1 = perEkinLow->fMCSecondMoment;
113 if (remRfaction>0.) {
114 DataPerEkin *perEkinHigh = fMCDataPerMaterial[matindx]->fDataPerEkin[ekinIndxLow+1];
115 mcToScr += remRfaction*(perEkinHigh->fMCScreening - perEkinLow->fMCScreening);
116 mcToQ1 += remRfaction*(perEkinHigh->fMCFirstMoment - perEkinLow->fMCFirstMoment);
117 mcToG2PerG1 += remRfaction*(perEkinHigh->fMCSecondMoment - perEkinLow->fMCSecondMoment);
128 if (delta>=gMaxDelta) {
134 G4int ekinIndxLow = 0;
136 if (beta2>gMaxBeta2) {
137 ekinIndxLow = gNumEkin - 1;
139 }
else if (beta2>=fMinBeta2) {
140 probIndxHigh = (beta2 - fMinBeta2) * fInvDelBeta2;
141 ekinIndxLow = (
G4int)probIndxHigh;
142 probIndxHigh -= ekinIndxLow;
143 ekinIndxLow += (gNumEkin - gNumBeta2);
144 }
else if (logekin>fLogMinEkin) {
145 probIndxHigh = (logekin - fLogMinEkin) * fInvLogDelEkin;
146 ekinIndxLow = (
G4int)probIndxHigh;
147 probIndxHigh -= ekinIndxLow;
155 ekindx = ekinIndxLow;
161 G4double probIndxHigh = delta*fInvDelDelta;
163 probIndxHigh -= deltIndxLow;
169 deltindx = deltIndxLow;
173 DataPerDelta *perDelta = fMCDataPerMaterial[matindx]->fDataPerEkin[ekindx]->fDataPerDelta[deltindx];
176 G4double ang = std::sqrt(0.5*(1.-cost));
177 G4double remRfaction = ang*fInvDelAngle;
179 remRfaction -= angIndx;
180 if (angIndx<gNumAngle-2) {
181 val = remRfaction*(perDelta->fRejFuntion[angIndx+1]-perDelta->fRejFuntion[angIndx]) + perDelta->fRejFuntion[angIndx];
183 G4double dum = ang-1.+1./fInvDelAngle;
184 val = perDelta->fSA + dum*(perDelta->fSB + dum*(perDelta->fSC + dum*perDelta->fSD));
192 InitMCDataPerElement();
194 ClearMCDataPerMaterial();
196 InitMCDataPerMaterials();
200void G4GSMottCorrection::InitMCDataPerElement() {
202 if (fMCDataPerElement.size()<gMaxZet+1) {
203 fMCDataPerElement.resize(gMaxZet+1,
nullptr);
209 for (
G4int imc=0; imc<numMatCuts; ++imc) {
217 std::size_t numElems = elemVect->size();
218 for (std::size_t ielem=0; ielem<numElems; ++ielem) {
224 if (!fMCDataPerElement[izet]) {
225 LoadMCDataElement(
elem);
232void G4GSMottCorrection::InitMCDataPerMaterials() {
235 if (fMCDataPerMaterial.size()!=numMaterials) {
236 fMCDataPerMaterial.resize(numMaterials);
241 for (
G4int imc=0; imc<numMatCuts; ++imc) {
247 if (!fMCDataPerMaterial[mat->
GetIndex()]) {
248 InitMCDataMaterial(mat);
255void G4GSMottCorrection::LoadMCDataElement(
const G4Element *
elem) {
261 auto perElem =
new DataPerMaterial();
262 AllocateDataPerMaterial(perElem);
263 fMCDataPerElement[izet] = perElem;
268 path +=
"/msc_GS/MottCor/el/";
270 path +=
"/msc_GS/MottCor/pos/";
272 std::string fname = path+
"rej_"+gElemSymbols[izet-1];
273 std::istringstream infile(std::ios::in);
274 ReadCompressedFile(fname, infile);
276 for (
G4int iek=0; iek<gNumEkin; ++iek) {
277 DataPerEkin *perEkin = perElem->fDataPerEkin[iek];
279 infile >> perEkin->fMCScreening;
280 infile >> perEkin->fMCFirstMoment;
281 infile >> perEkin->fMCSecondMoment;
283 for (
G4int idel=0; idel<gNumDelta; ++idel) {
284 DataPerDelta *perDelta = perEkin->fDataPerDelta[idel];
286 for (
G4int iang=0; iang<gNumAngle; ++iang) {
287 infile >> perDelta->fRejFuntion[iang];
290 infile >> perDelta->fSA;
291 infile >> perDelta->fSB;
292 infile >> perDelta->fSC;
293 infile >> perDelta->fSD;
299void G4GSMottCorrection::ReadCompressedFile(std::string fname, std::istringstream &iss) {
300 std::string *dataString =
nullptr;
301 std::string compfilename(fname+
".z");
303 std::ifstream in(compfilename, std::ios::binary | std::ios::ate);
306 std::streamoff fileSize = in.tellg();
308 in.seekg(0,std::ios::beg);
310 Bytef *compdata =
new Bytef[fileSize];
312 in.read((
char*)compdata, fileSize);
315 uLongf complen = (uLongf)(fileSize*4);
316 Bytef *uncompdata =
new Bytef[complen];
317 while (
Z_OK!=
uncompress(uncompdata, &complen, compdata, fileSize)) {
321 uncompdata =
new Bytef[complen];
326 dataString =
new std::string((
char*)uncompdata, (
long)complen);
328 delete [] uncompdata;
330 std::string msg =
" Problem while trying to read " + compfilename +
" data file.\n";
336 iss.str(*dataString);
343void G4GSMottCorrection::InitMCDataMaterial(
const G4Material *mat) {
346 constexpr G4double finstrc2 = 5.325135453E-5;
348 G4double constFactor = CLHEP::electron_mass_c2*CLHEP::fine_structure_const/0.88534;
349 constFactor *= constFactor;
351 auto perMat =
new DataPerMaterial();
352 AllocateDataPerMaterial(perMat);
353 fMCDataPerMaterial[mat->
GetIndex()] = perMat;
369 for (
G4int ielem=0; ielem<numElems; ++ielem) {
370 G4double zet = (*elemVect)[ielem]->GetZ();
374 G4double iwa = (*elemVect)[ielem]->GetN();
375 G4double ipz = nbAtomsPerVolVect[ielem]/totNbAtomsPerVol;
378 ze += dum*(-2.0/3.0)*
G4Log(zet);
379 zx += dum*
G4Log(1.0+3.34*finstrc2*zet*zet);
384 moliereBc = const1*density*zs/sa*
G4Exp(ze/zs)/
G4Exp(zx/zs);
385 moliereXc2 = const2*density*zs/sa;
387 moliereBc *= 1.0/CLHEP::cm;
388 moliereXc2 *= CLHEP::MeV*CLHEP::MeV/CLHEP::cm;
391 for (
G4int iek=0; iek<gNumEkin; ++iek) {
394 G4double pt2 = ekin*(ekin+2.0*CLHEP::electron_mass_c2);
396 G4double b2 = fMinBeta2+(iek-(gNumEkin-gNumBeta2))/fInvDelBeta2;
397 ekin = CLHEP::electron_mass_c2*(1./std::sqrt(1.-b2)-1.);
398 pt2 = ekin*(ekin+2.0*CLHEP::electron_mass_c2);
401 for (
G4int ielem=0; ielem<numElems; ++ielem) {
409 G4double nZZPlus1 = nbAtomsPerVolVect[ielem]*zet*(zet+1.0)/totNbAtomsPerVol;
412 DataPerEkin *perElemPerEkin = fMCDataPerElement[izet]->fDataPerEkin[iek];
413 DataPerEkin *perMatPerEkin = perMat->fDataPerEkin[iek];
416 G4double mcScrCF = perElemPerEkin->fMCScreening;
421 perMatPerEkin->fMCScreening += nZZPlus1*
G4Log(Z23*mcScrCF);
424 mcScrCF *= constFactor*Z23/(4.*pt2);
434 perMatPerEkin->fMCFirstMoment += nZZPlus1*(
G4Log(1.+1./mcScrCF)-1./(1.+mcScrCF))*perElemPerEkin->fMCFirstMoment;
439 perMatPerEkin->fMCSecondMoment += nZZPlus1*perElemPerEkin->fMCSecondMoment;
443 for (
G4int idel=0; idel<gNumDelta; ++idel) {
444 DataPerDelta *perMatPerDelta = perMatPerEkin->fDataPerDelta[idel];
445 DataPerDelta *perElemPerDelta = perElemPerEkin->fDataPerDelta[idel];
447 for (
G4int iang=0; iang<gNumAngle; ++iang) {
448 perMatPerDelta->fRejFuntion[iang] += nZZPlus1*perElemPerDelta->fRejFuntion[iang];
451 perMatPerDelta->fSA += nZZPlus1*perElemPerDelta->fSA;
452 perMatPerDelta->fSB += nZZPlus1*perElemPerDelta->fSB;
453 perMatPerDelta->fSC += nZZPlus1*perElemPerDelta->fSC;
454 perMatPerDelta->fSD += nZZPlus1*perElemPerDelta->fSD;
458 if (ielem==numElems-1) {
463 perMatPerEkin->fMCScreening = constFactor*dumScr*moliereBc/moliereXc2;
467 G4double scrCorTed = constFactor*dumScr/(4.*pt2);
469 perMatPerEkin->fMCFirstMoment = perMatPerEkin->fMCFirstMoment/(zs*(dum0-1./(1.+scrCorTed)));
473 G4double G2PerG1 = 3.*(1.+scrCorTed)*((1.+2.*scrCorTed)*dum0-2.)/((1.+scrCorTed)*dum0-1.);
474 perMatPerEkin->fMCSecondMoment = perMatPerEkin->fMCSecondMoment/(zs*G2PerG1);
478 for (
G4int idel=0; idel<gNumDelta; ++idel) {
479 DataPerDelta *perMatPerDelta = perMatPerEkin->fDataPerDelta[idel];
482 for (
G4int iang=0; iang<gNumAngle; ++iang) {
483 if (perMatPerDelta->fRejFuntion[iang]>maxVal)
484 maxVal = perMatPerDelta->fRejFuntion[iang];
486 for (
G4int iang=0; iang<gNumAngle; ++iang) {
487 perMatPerDelta->fRejFuntion[iang] /=maxVal;
489 perMatPerDelta->fSA /= maxVal;
490 perMatPerDelta->fSB /= maxVal;
491 perMatPerDelta->fSC /= maxVal;
492 perMatPerDelta->fSD /= maxVal;
500void G4GSMottCorrection::AllocateDataPerMaterial(DataPerMaterial *data) {
501 data->fDataPerEkin =
new DataPerEkin*[gNumEkin]();
502 for (
G4int iek=0; iek<gNumEkin; ++iek) {
503 auto perEkin =
new DataPerEkin();
504 perEkin->fDataPerDelta =
new DataPerDelta*[gNumDelta]();
505 for (
G4int idel=0; idel<gNumDelta; ++idel) {
506 auto perDelta =
new DataPerDelta();
507 perDelta->fRejFuntion =
new double[gNumAngle]();
508 perEkin->fDataPerDelta[idel] = perDelta;
510 data->fDataPerEkin[iek] = perEkin;
514void G4GSMottCorrection::DeAllocateDataPerMaterial(DataPerMaterial *data) {
515 for (
G4int iek=0; iek<gNumEkin; ++iek) {
516 DataPerEkin *perEkin = data->fDataPerEkin[iek];
517 for (
G4int idel=0; idel<gNumDelta; ++idel) {
518 DataPerDelta *perDelta = perEkin->fDataPerDelta[idel];
519 delete [] perDelta->fRejFuntion;
522 delete [] perEkin->fDataPerDelta;
525 delete [] data->fDataPerEkin;
529void G4GSMottCorrection::ClearMCDataPerElement() {
530 for (std::size_t i=0; i<fMCDataPerElement.size(); ++i) {
531 if (fMCDataPerElement[i]) {
532 DeAllocateDataPerMaterial(fMCDataPerElement[i]);
533 delete fMCDataPerElement[i];
536 fMCDataPerElement.clear();
539void G4GSMottCorrection::ClearMCDataPerMaterial() {
540 for (std::size_t i=0; i<fMCDataPerMaterial.size(); ++i) {
541 if (fMCDataPerMaterial[i]) {
542 DeAllocateDataPerMaterial(fMCDataPerMaterial[i]);
543 delete fMCDataPerMaterial[i];
546 fMCDataPerMaterial.clear();
std::vector< const G4Element * > G4ElementVector
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
G4double G4Log(G4double x)
static G4EmParameters * Instance()
const G4String & GetDirLEDATA() const
G4GSMottCorrection(G4bool iselectron=true)
void GetMottCorrectionFactors(G4double logekin, G4double beta2, G4int matindx, G4double &mcToScr, G4double &mcToQ1, G4double &mcToG2PerG1)
G4double GetMottRejectionValue(G4double logekin, G4double G4beta2, G4double q1, G4double cost, G4int matindx, G4int &ekindx, G4int &deltindx)
const G4Material * GetMaterial() const
G4double GetDensity() const
const G4ElementVector * GetElementVector() const
G4double GetTotNbOfAtomsPerVolume() const
static std::size_t GetNumberOfMaterials()
const G4double * GetVecNbOfAtomsPerVolume() const
std::size_t GetIndex() const
std::size_t GetNumberOfElements() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
int ZEXPORT uncompress(Bytef *dest, uLongf *destLen, const Bytef *source, uLong sourceLen)