Geant4-11
G4GoudsmitSaundersonTable.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27// -----------------------------------------------------------------------------
28//
29// GEANT4 Class implementation file
30//
31// File name: G4GoudsmitSaundersonTable
32//
33// Author: Mihaly Novak / (Omrane Kadri)
34//
35// Creation date: 20.02.2009
36//
37// Class description:
38// Class to handle multiple scattering angular distributions precomputed by
39// using Kawrakow-Bielajew Goudsmit-Saunderson MSC model based on the screened
40// Rutherford DCS for elastic scattering of electrons/positrons [1,2]. This
41// class is used by G4GoudsmitSaundersonMscModel to sample the angular
42// deflection of electrons/positrons after travelling a given path.
43//
44// Modifications:
45// 04.03.2009 V.Ivanchenko cleanup and format according to Geant4 EM style
46// 26.08.2009 O.Kadri: avoiding unuseful calculations and optimizing the root
47// finding parameter error's within SampleTheta method
48// 08.02.2010 O.Kadri: reduce delared variables; reduce error of finding root
49// in secant method
50// 26.03.2010 O.Kadri: minimum of used arrays in computation within the dichotomie
51// finding method the error was the lowest value of uvalues
52// 12.05.2010 O.Kadri: changing of sqrt((b-a)*(b-a)) with fabs(b-a)
53// 18.05.2015 M. Novak This class has been completely replaced (only the original
54// class name was kept; class description was also inserted):
55// A new version of Kawrakow-Bielajew Goudsmit-Saunderson MSC model
56// based on the screened Rutherford DCS for elastic scattering of
57// electrons/positrons has been introduced[1,2]. The corresponding MSC
58// angular distributions over a 2D parameter grid have been recomputed
59// and the CDFs are now stored in a variable transformed (smooth) form
60// together with the corresponding rational interpolation parameters.
61// The new version is several times faster, more robust and accurate
62// compared to the earlier version (G4GoudsmitSaundersonMscModel class
63// that use these data has been also completely replaced)
64// 28.04.2017 M. Novak: New representation of the angular distribution data with
65// significantly reduced data size.
66// 23.08.2017 M. Novak: Added funtionality to handle Mott-correction to the
67// base GS angular distributions and some other factors (screening
68// parameter, first and second moments) when Mott-correction is
69// activated in the GS-MSC model.
70//
71// References:
72// [1] A.F.Bielajew, NIMB, 111 (1996) 195-208
73// [2] I.Kawrakow, A.F.Bielajew, NIMB 134(1998) 325-336
74//
75// -----------------------------------------------------------------------------
76
78
79
81#include "Randomize.hh"
82#include "G4Log.hh"
83#include "G4Exp.hh"
84
85#include "G4GSMottCorrection.hh"
86#include "G4MaterialTable.hh"
87#include "G4Material.hh"
90
91#include "G4String.hh"
92
93#include <fstream>
94#include <cstdlib>
95#include <cmath>
96
97#include <iostream>
98#include <iomanip>
99
100// perecomputed GS angular distributions, based on the Screened-Rutherford DCS
101// are the same for e- and e+ so make sure we load them only onece
103//
104std::vector<G4GoudsmitSaundersonTable::GSMSCAngularDtr*> G4GoudsmitSaundersonTable::gGSMSCAngularDistributions1;
105std::vector<G4GoudsmitSaundersonTable::GSMSCAngularDtr*> G4GoudsmitSaundersonTable::gGSMSCAngularDistributions2;
106//
107std::vector<double> G4GoudsmitSaundersonTable::gMoliereBc;
108std::vector<double> G4GoudsmitSaundersonTable::gMoliereXc2;
109
110
112 fIsElectron = iselectron;
113 // set initial values: final values will be set in the Initialize method
114 fLogLambda0 = 0.; // will be set properly at init.
115 fLogDeltaLambda = 0.; // will be set properly at init.
116 fInvLogDeltaLambda = 0.; // will be set properly at init.
117 fInvDeltaQ1 = 0.; // will be set properly at init.
118 fDeltaQ2 = 0.; // will be set properly at init.
119 fInvDeltaQ2 = 0.; // will be set properly at init.
120 //
121 fLowEnergyLimit = 0.1*CLHEP::keV; // will be set properly at init.
122 fHighEnergyLimit = 100.0*CLHEP::MeV; // will be set properly at init.
123 //
124 fIsMottCorrection = false; // will be set properly at init.
125 fIsPWACorrection = false; // will be set properly at init.
126 fMottCorrection = nullptr;
127 //
129}
130
132 for (size_t i=0; i<gGSMSCAngularDistributions1.size(); ++i) {
134 delete [] gGSMSCAngularDistributions1[i]->fUValues;
135 delete [] gGSMSCAngularDistributions1[i]->fParamA;
136 delete [] gGSMSCAngularDistributions1[i]->fParamB;
138 }
139 }
141 for (size_t i=0; i<gGSMSCAngularDistributions2.size(); ++i) {
143 delete [] gGSMSCAngularDistributions2[i]->fUValues;
144 delete [] gGSMSCAngularDistributions2[i]->fParamA;
145 delete [] gGSMSCAngularDistributions2[i]->fParamB;
147 }
148 }
150 if (fMottCorrection) {
151 delete fMottCorrection;
152 fMottCorrection = nullptr;
153 }
154 // clear scp correction data
155 for (size_t imc=0; imc<fSCPCPerMatCuts.size(); ++imc) {
156 if (fSCPCPerMatCuts[imc]) {
157 fSCPCPerMatCuts[imc]->fVSCPC.clear();
158 delete fSCPCPerMatCuts[imc];
159 }
160 }
161 fSCPCPerMatCuts.clear();
162 //
163 gIsInitialised = false;
164}
165
166void G4GoudsmitSaundersonTable::Initialise(G4double lownergylimit, G4double highenergylimit) {
167 fLowEnergyLimit = lownergylimit;
168 fHighEnergyLimit = highenergylimit;
169 G4double lLambdaMin = G4Log(gLAMBMIN);
170 G4double lLambdaMax = G4Log(gLAMBMAX);
171 fLogLambda0 = lLambdaMin;
172 fLogDeltaLambda = (lLambdaMax-lLambdaMin)/(gLAMBNUM-1.);
174 fInvDeltaQ1 = 1./((gQMAX1-gQMIN1)/(gQNUM1-1.));
175 fDeltaQ2 = (gQMAX2-gQMIN2)/(gQNUM2-1.);
177 // load precomputed angular distributions and set up several values used during the sampling
178 // these are particle independet => they go to static container: load them only onece
179 if (!gIsInitialised) {
180 // load pre-computed GS angular distributions (computed based on Screened-Rutherford DCS)
181 LoadMSCData();
182 gIsInitialised = true;
183 }
185 // Mott-correction: particle(e- or e+) dependet so init them
186 if (fIsMottCorrection) {
187 if (!fMottCorrection) {
188 fMottCorrection = new G4GSMottCorrection(fIsElectron);
189 }
190 fMottCorrection->Initialise();
191 }
192 // init scattering power correction data; used only together with Mott-correction
193 // (Moliere's parameters must be initialised before)
194 if (fMottCorrection) {
196 }
197}
198
199
200// samplig multiple scattering angles cos(theta) and sin(thata)
201// - including no-scattering, single, "few" scattering cases as well
202// - Mott-correction will be included if it was requested by the user (i.e. if fIsMottCorrection=true)
203// lambdaval : s/lambda_el
204// qval : s/lambda_el G1
205// scra : screening parameter
206// cost : will be the smapled cos(theta)
207// sint : will be the smapled sin(theta)
208// lekin : logarithm of the current kinetic energy
209// beta2 : the corresponding beta square
210// matindx : index of the current material
211// returns true if it was msc
213 G4double &sint, G4double lekin, G4double beta2, G4int matindx,
214 GSMSCAngularDtr **gsDtr, G4int &mcekini, G4int &mcdelti,
215 G4double &transfPar, G4bool isfirst) {
216 G4double rand0 = G4UniformRand();
217 G4double expn = G4Exp(-lambdaval);
218 //
219 // no scattering case
220 if (rand0<expn) {
221 cost = 1.0;
222 sint = 0.0;
223 return false;
224 }
225 //
226 // single scattering case : sample from the single scattering PDF
227 // - Mott-correction will be included if it was requested by the user (i.e. if fIsMottCorrection=true)
228 if (rand0<(1.+lambdaval)*expn) {
229 // cost is sampled in SingleScattering()
230 cost = SingleScattering(lambdaval, scra, lekin, beta2, matindx);
231 // add protections
232 if (cost<-1.0) cost = -1.0;
233 if (cost>1.0) cost = 1.0;
234 // compute sin(theta) from the sampled cos(theta)
235 G4double dum0 = 1.-cost;
236 sint = std::sqrt(dum0*(2.0-dum0));
237 return false;
238 }
239 //
240 // handle this case:
241 // -lambdaval < 1 i.e. mean #elastic events along the step is < 1 but
242 // the currently sampled case is not 0 or 1 scattering. [Our minimal
243 // lambdaval (that we have precomputed, transformed angular distributions
244 // stored in a form of equally probabe intervalls together with rational
245 // interp. parameters) is 1.]
246 // -probability of having n elastic events follows Poisson stat. with
247 // lambdaval parameter.
248 // -the max. probability (when lambdaval=1) of having more than one
249 // elastic events is 0.2642411 and the prob of having 2,3,..,n elastic
250 // events decays rapidly with n. So set a max n to 10.
251 // -sampling of this cases is done in a one-by-one single elastic event way
252 // where the current #elastic event is sampled from the Poisson distr.
253 if (lambdaval<1.0) {
254 G4double prob, cumprob;
255 prob = cumprob = expn;
256 G4double curcost,cursint;
257 // init cos(theta) and sin(theta) to the zero scattering values
258 cost = 1.0;
259 sint = 0.0;
260 for (G4int iel=1; iel<10; ++iel) {
261 // prob of having iel scattering from Poisson
262 prob *= lambdaval/(G4double)iel;
263 cumprob += prob;
264 //
265 //sample cos(theta) from the singe scattering pdf:
266 // - Mott-correction will be included if it was requested by the user (i.e. if fIsMottCorrection=true)
267 curcost = SingleScattering(lambdaval, scra, lekin, beta2, matindx);
268 G4double dum0 = 1.-curcost;
269 cursint = dum0*(2.0-dum0); // sin^2(theta)
270 //
271 // if we got current deflection that is not too small
272 // then update cos(theta) sin(theta)
273 if (cursint>1.0e-20) {
274 cursint = std::sqrt(cursint);
276 cost = cost*curcost-sint*cursint*std::cos(curphi);
277 sint = std::sqrt(std::max(0.0, (1.0-cost)*(1.0+cost)));
278 }
279 //
280 // check if we have done enough scattering i.e. sampling from the Poisson
281 if (rand0<cumprob) {
282 return false;
283 }
284 }
285 // if reached the max iter i.e. 10
286 return false;
287 }
288 //
289 // multiple scattering case with lambdavalue >= 1:
290 // - use the precomputed and transformed Goudsmit-Saunderson angular
291 // distributions to sample cos(theta)
292 // - Mott-correction will be included if it was requested by the user (i.e. if fIsMottCorrection=true)
293 cost = SampleCosTheta(lambdaval, qval, scra, lekin, beta2, matindx, gsDtr, mcekini, mcdelti, transfPar, isfirst);
294 // add protections
295 if (cost<-1.0) cost = -1.0;
296 if (cost> 1.0) cost = 1.0;
297 // compute cos(theta) and sin(theta) from the sampled 1-cos(theta)
298 G4double dum0 = 1.0-cost;
299 sint = std::sqrt(dum0*(2.0-dum0));
300 // return true if it was msc
301 return true;
302}
303
304
306 G4double lekin, G4double beta2, G4int matindx,
307 GSMSCAngularDtr **gsDtr, G4int &mcekini,G4int &mcdelti,
308 G4double &transfPar, G4bool isfirst) {
309 G4double cost = 1.;
310 // determine the base GS angular distribution if it is the first call (when sub-step sampling is used)
311 if (isfirst) {
312 *gsDtr = GetGSAngularDtr(scra, lambdaval, qval, transfPar);
313 }
314 // sample cost from the GS angular distribution (computed based on Screened-Rutherford DCS)
315 cost = SampleGSSRCosTheta(*gsDtr, transfPar);
316 // Mott-correction if it was requested by the user
317 if (fIsMottCorrection && *gsDtr) { // no Mott-correction in case of izotropic theta
318 static const G4int nlooplim = 1000;
319 G4int nloop = 0 ; // rejection loop counter
320// G4int ekindx = -1; // evaluate only in the first call
321// G4int deltindx = -1 ; // evaluate only in the first call
322 G4double val = fMottCorrection->GetMottRejectionValue(lekin, beta2, qval, cost, matindx, mcekini, mcdelti);
323 while (G4UniformRand()>val && ++nloop<nlooplim) {
324 // sampling cos(theta)
325 cost = SampleGSSRCosTheta(*gsDtr, transfPar);
326 val = fMottCorrection->GetMottRejectionValue(lekin, beta2, qval, cost, matindx, mcekini, mcdelti);
327 };
328 }
329 return cost;
330}
331
332
333// returns with cost sampled from the GS angular distribution computed based on Screened-Rutherford DCS
335 // check if isotropic theta (i.e. cost is uniform on [-1:1])
336 if (!gsDtr) {
337 return 1.-2.0*G4UniformRand();
338 }
339 //
340 // sampling form the selected distribution
341 G4double ndatm1 = gsDtr->fNumData-1.;
342 G4double delta = 1.0/ndatm1;
343 // determine lower cumulative bin inidex
344 G4double rndm = G4UniformRand();
345 G4int indxl = rndm*ndatm1;
346 G4double aval = rndm-indxl*delta;
347 G4double dum0 = delta*aval;
348
349 G4double dum1 = (1.0+gsDtr->fParamA[indxl]+gsDtr->fParamB[indxl])*dum0;
350 G4double dum2 = delta*delta + gsDtr->fParamA[indxl]*dum0 + gsDtr->fParamB[indxl]*aval*aval;
351 G4double sample = gsDtr->fUValues[indxl] + dum1/dum2 *(gsDtr->fUValues[indxl+1]-gsDtr->fUValues[indxl]);
352 // transform back u to cos(theta) :
353 // this is the sampled cos(theta) = (2.0*para*sample)/(1.0-sample+para)
354 return 1.-(2.0*transfpar*sample)/(1.0-sample+transfpar);
355}
356
357
358// determine the GS angular distribution we need to sample from: will set other things as well ...
360 G4double &lambdaval, G4double &qval, G4double &transfpar) {
361 GSMSCAngularDtr *dtr = nullptr;
362 G4bool first = false;
363 // isotropic cost above gQMAX2 (i.e. dtr stays nullptr)
364 if (qval<gQMAX2) {
365 G4int lamIndx = -1; // lambda value index
366 G4int qIndx = -1; // lambda value index
367 // init to second grid Q values
368 G4int numQVal = gQNUM2;
369 G4double minQVal = gQMIN2;
370 G4double invDelQ = fInvDeltaQ2;
371 G4double pIndxH = 0.; // probability of taking higher index
372 // check if first or second grid needs to be used
373 if (qval<gQMIN2) { // first grid
374 first = true;
375 // protect against qval<gQMIN1
376 if (qval<gQMIN1) {
377 qval = gQMIN1;
378 qIndx = 0;
379 //pIndxH = 0.;
380 }
381 // set to first grid Q values
382 numQVal = gQNUM1;
383 minQVal = gQMIN1;
384 invDelQ = fInvDeltaQ1;
385 }
386 // make sure that lambda = s/lambda_el is in [gLAMBMIN,gLAMBMAX)
387 // lambda<gLAMBMIN=1 is already handeled before so lambda>= gLAMBMIN for sure
388 if (lambdaval>=gLAMBMAX) {
389 lambdaval = gLAMBMAX-1.e-8;
390 lamIndx = gLAMBNUM-1;
391 }
392 G4double lLambda = G4Log(lambdaval);
393 //
394 // determine lower lambda (=s/lambda_el) index: linear interp. on log(lambda) scale
395 if (lamIndx<0) {
396 pIndxH = (lLambda-fLogLambda0)*fInvLogDeltaLambda;
397 lamIndx = (G4int)(pIndxH); // lower index of the lambda bin
398 pIndxH = pIndxH-lamIndx; // probability of taking the higher index distribution
399 if (G4UniformRand()<pIndxH) {
400 ++lamIndx;
401 }
402 }
403 //
404 // determine lower Q (=s/lambda_el G1) index: linear interp. on Q
405 if (qIndx<0) {
406 pIndxH = (qval-minQVal)*invDelQ;
407 qIndx = (G4int)(pIndxH); // lower index of the Q bin
408 pIndxH = pIndxH-qIndx;
409 if (G4UniformRand()<pIndxH) {
410 ++qIndx;
411 }
412 }
413 // set indx
414 G4int indx = lamIndx*numQVal+qIndx;
415 if (first) {
416 dtr = gGSMSCAngularDistributions1[indx];
417 } else {
418 dtr = gGSMSCAngularDistributions2[indx];
419 }
420 // dtr might be nullptr that indicates isotropic cot distribution because:
421 // - if the selected lamIndx, qIndx correspond to L(=s/lambda_el) and Q(=s/lambda_el G1) such that G1(=Q/L) > 1
422 // G1 should always be < 1 and if G1 is ~1 -> the dtr is isotropic (this can only happen in case of the 2. grid)
423 //
424 // compute the transformation parameter
425 if (lambdaval>10.0) {
426 transfpar = 0.5*(-2.77164+lLambda*( 2.94874-lLambda*(0.1535754-lLambda*0.00552888) ));
427 } else {
428 transfpar = 0.5*(1.347+lLambda*(0.209364-lLambda*(0.45525-lLambda*(0.50142-lLambda*0.081234))));
429 }
430 transfpar *= (lambdaval+4.0)*scra;
431 }
432 // return with the selected GS angular distribution that we need to sample cost from (if nullptr => isotropic cost)
433 return dtr;
434}
435
436
438 char* path = std::getenv("G4LEDATA");
439 if (!path) {
440 G4Exception("G4GoudsmitSaundersonTable::LoadMSCData()","em0006",
442 "Environment variable G4LEDATA not defined");
443 return;
444 }
445 //
447 const G4String str1 = G4String(path) + "/msc_GS/GSGrid_1/gsDistr_";
448 for (G4int il=0; il<gLAMBNUM; ++il) {
449 G4String fname = str1 + std::to_string(il);
450 std::ifstream infile(fname,std::ios::in);
451 if (!infile.is_open()) {
452 G4String msgc = "Cannot open file: " + fname;
453 G4Exception("G4GoudsmitSaundersonTable::LoadMSCData()","em0006",
454 FatalException, msgc.c_str());
455 return;
456 }
457 for (G4int iq=0; iq<gQNUM1; ++iq) {
458 GSMSCAngularDtr *gsd = new GSMSCAngularDtr();
459 infile >> gsd->fNumData;
460 gsd->fUValues = new G4double[gsd->fNumData]();
461 gsd->fParamA = new G4double[gsd->fNumData]();
462 gsd->fParamB = new G4double[gsd->fNumData]();
463 G4double ddummy;
464 infile >> ddummy; infile >> ddummy;
465 for (G4int i=0; i<gsd->fNumData; ++i) {
466 infile >> gsd->fUValues[i];
467 infile >> gsd->fParamA[i];
468 infile >> gsd->fParamB[i];
469 }
471 }
472 infile.close();
473 }
474 //
475 // second grid
477 const G4String str2 = G4String(path) + "/msc_GS/GSGrid_2/gsDistr_";
478 for (G4int il=0; il<gLAMBNUM; ++il) {
479 G4String fname = str2 + std::to_string(il);
480 std::ifstream infile(fname,std::ios::in);
481 if (!infile.is_open()) {
482 G4String msgc = "Cannot open file: " + fname;
483 G4Exception("G4GoudsmitSaundersonTable::LoadMSCData()","em0006",
484 FatalException, msgc.c_str());
485 return;
486 }
487 for (G4int iq=0; iq<gQNUM2; ++iq) {
488 G4int numData;
489 infile >> numData;
490 if (numData>1) {
491 GSMSCAngularDtr *gsd = new GSMSCAngularDtr();
492 gsd->fNumData = numData;
493 gsd->fUValues = new G4double[gsd->fNumData]();
494 gsd->fParamA = new G4double[gsd->fNumData]();
495 gsd->fParamB = new G4double[gsd->fNumData]();
496 double ddummy;
497 infile >> ddummy; infile >> ddummy;
498 for (G4int i=0; i<gsd->fNumData; ++i) {
499 infile >> gsd->fUValues[i];
500 infile >> gsd->fParamA[i];
501 infile >> gsd->fParamB[i];
502 }
504 } else {
505 gGSMSCAngularDistributions2[il*gQNUM2+iq] = nullptr;
506 }
507 }
508 infile.close();
509 }
510}
511
512// samples cost in single scattering based on Screened-Rutherford DCS
513// (with Mott-correction if it was requested)
515 G4double lekin, G4double beta2,
516 G4int matindx) {
517 G4double rand1 = G4UniformRand();
518 // sample cost from the Screened-Rutherford DCS
519 G4double cost = 1.-2.0*scra*rand1/(1.0-rand1+scra);
520 // Mott-correction if it was requested by the user
521 if (fIsMottCorrection) {
522 static const G4int nlooplim = 1000; // rejection loop limit
523 G4int nloop = 0 ; // loop counter
524 G4int ekindx = -1 ; // evaluate only in the first call
525 G4int deltindx = 0 ; // single scattering case
526 G4double q1 = 0.; // not used when deltindx = 0;
527 // computing Mott rejection function value
528 G4double val = fMottCorrection->GetMottRejectionValue(lekin, beta2, q1, cost,
529 matindx, ekindx, deltindx);
530 while (G4UniformRand()>val && ++nloop<nlooplim) {
531 // sampling cos(theta) from the Screened-Rutherford DCS
532 rand1 = G4UniformRand();
533 cost = 1.-2.0*scra*rand1/(1.0-rand1+scra);
534 // computing Mott rejection function value
535 val = fMottCorrection->GetMottRejectionValue(lekin, beta2, q1, cost, matindx,
536 ekindx, deltindx);
537 };
538 }
539 return cost;
540}
541
542
544 G4int matindx, G4double &mcToScr,
545 G4double &mcToQ1, G4double &mcToG2PerG1) {
546 if (fIsMottCorrection) {
547 fMottCorrection->GetMottCorrectionFactors(logekin, beta2, matindx, mcToScr, mcToQ1, mcToG2PerG1);
548 }
549}
550
551
552// compute material dependent Moliere MSC parameters at initialisation
554 const G4double const1 = 7821.6; // [cm2/g]
555 const G4double const2 = 0.1569; // [cm2 MeV2 / g]
556 const G4double finstrc2 = 5.325135453E-5; // fine-structure const. square
557
558 G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
559 // get number of materials in the table
560 size_t numMaterials = theMaterialTable->size();
561 // make sure that we have long enough vectors
562 if(gMoliereBc.size()<numMaterials) {
563 gMoliereBc.resize(numMaterials);
564 gMoliereXc2.resize(numMaterials);
565 }
566 G4double xi = 1.0;
567 G4int maxZ = 200;
569 // xi = 1.0; <= always set to 1 from now on
571 }
572 //
573 for (size_t imat=0; imat<numMaterials; ++imat) {
574 const G4Material* theMaterial = (*theMaterialTable)[imat];
575 const G4ElementVector* theElemVect = theMaterial->GetElementVector();
576 const G4int numelems = theMaterial->GetNumberOfElements();
577 //
578 const G4double* theNbAtomsPerVolVect = theMaterial->GetVecNbOfAtomsPerVolume();
579 G4double theTotNbAtomsPerVol = theMaterial->GetTotNbOfAtomsPerVolume();
580 //
581 G4double zs = 0.0;
582 G4double zx = 0.0;
583 G4double ze = 0.0;
584 G4double sa = 0.0;
585 //
586 for(G4int ielem = 0; ielem < numelems; ielem++) {
587 G4double zet = (*theElemVect)[ielem]->GetZ();
588 if (zet>maxZ) {
589 zet = (G4double)maxZ;
590 }
591 G4double iwa = (*theElemVect)[ielem]->GetN();
592 G4double ipz = theNbAtomsPerVolVect[ielem]/theTotNbAtomsPerVol;
593 G4double dum = ipz*zet*(zet+xi);
594 zs += dum;
595 ze += dum*(-2.0/3.0)*G4Log(zet);
596 zx += dum*G4Log(1.0+3.34*finstrc2*zet*zet);
597 sa += ipz*iwa;
598 }
599 G4double density = theMaterial->GetDensity()*CLHEP::cm3/CLHEP::g; // [g/cm3]
600 //
601 gMoliereBc[theMaterial->GetIndex()] = const1*density*zs/sa*G4Exp(ze/zs)/G4Exp(zx/zs); //[1/cm]
602 gMoliereXc2[theMaterial->GetIndex()] = const2*density*zs/sa; // [MeV2/cm]
603 // change to Geant4 internal units of 1/length and energ2/length
604 gMoliereBc[theMaterial->GetIndex()] *= 1.0/CLHEP::cm;
605 gMoliereXc2[theMaterial->GetIndex()] *= CLHEP::MeV*CLHEP::MeV/CLHEP::cm;
606 }
607}
608
609
610// this method is temporary, will be removed/replaced with a more effictien solution after 10.3.ref09
612 G4int imc = matcut->GetIndex();
613 G4double corFactor = 1.0;
614 if (!(fSCPCPerMatCuts[imc]->fIsUse) || ekin<=fSCPCPerMatCuts[imc]->fPrCut) {
615 return corFactor;
616 }
617 // get the scattering power correction factor
618 G4double lekin = G4Log(ekin);
619 G4double remaining = (lekin-fSCPCPerMatCuts[imc]->fLEmin)*fSCPCPerMatCuts[imc]->fILDel;
620 G4int lindx = (G4int)remaining;
621 remaining -= lindx;
622 G4int imax = fSCPCPerMatCuts[imc]->fVSCPC.size()-1;
623 if (lindx>=imax) {
624 corFactor = fSCPCPerMatCuts[imc]->fVSCPC[imax];
625 } else {
626 corFactor = fSCPCPerMatCuts[imc]->fVSCPC[lindx] + remaining*(fSCPCPerMatCuts[imc]->fVSCPC[lindx+1]-fSCPCPerMatCuts[imc]->fVSCPC[lindx]);
627 }
628 return corFactor;
629}
630
631
633 // get the material-cuts table
635 size_t numMatCuts = thePCTable->GetTableSize();
636 // clear container if any
637 for (size_t imc=0; imc<fSCPCPerMatCuts.size(); ++imc) {
638 if (fSCPCPerMatCuts[imc]) {
639 fSCPCPerMatCuts[imc]->fVSCPC.clear();
640 delete fSCPCPerMatCuts[imc];
641 fSCPCPerMatCuts[imc] = nullptr;
642 }
643 }
644 //
645 // set size of the container and create the corresponding data structures
646 fSCPCPerMatCuts.resize(numMatCuts,nullptr);
647 // loop over the material-cuts and create scattering power correction data structure for each
648 for (size_t imc=0; imc<numMatCuts; ++imc) {
649 const G4MaterialCutsCouple *matCut = thePCTable->GetMaterialCutsCouple(imc);
650 // get e- production cut in the current material-cuts in energy
651 G4double limit;
652 G4double ecut;
653 if (fIsElectron) {
654 ecut = (*(thePCTable->GetEnergyCutsVector(idxG4ElectronCut)))[matCut->GetIndex()];
655 limit = 2.*ecut;
656 } else {
657 ecut = (*(thePCTable->GetEnergyCutsVector(idxG4PositronCut)))[matCut->GetIndex()];
658 limit = ecut;
659 }
662 if (min>=max) {
663 fSCPCPerMatCuts[imc] = new SCPCorrection();
664 fSCPCPerMatCuts[imc]->fIsUse = false;
665 fSCPCPerMatCuts[imc]->fPrCut = min;
666 continue;
667 }
668 G4int numEbins = fNumSPCEbinPerDec*G4lrint(std::log10(max/min));
669 numEbins = std::max(numEbins,3);
670 G4double lmin = G4Log(min);
671 G4double ldel = G4Log(max/min)/(numEbins-1.0);
672 fSCPCPerMatCuts[imc] = new SCPCorrection();
673 fSCPCPerMatCuts[imc]->fVSCPC.resize(numEbins,1.0);
674 fSCPCPerMatCuts[imc]->fIsUse = true;
675 fSCPCPerMatCuts[imc]->fPrCut = min;
676 fSCPCPerMatCuts[imc]->fLEmin = lmin;
677 fSCPCPerMatCuts[imc]->fILDel = 1./ldel;
678 for (G4int ie=0; ie<numEbins; ++ie) {
679 G4double ekin = G4Exp(lmin+ie*ldel);
680 G4double scpCorr = 1.0;
681 // compute correction factor: I.Kawrakow NIMB 114(1996)307-326 (Eqs(32-37))
682 if (ie>0) {
684 G4double tauCut = ecut/CLHEP::electron_mass_c2;
685 // Moliere's screening parameter
686 G4int matindx = matCut->GetMaterial()->GetIndex();
687 G4double A = GetMoliereXc2(matindx)/(4.0*tau*(tau+2.)*GetMoliereBc(matindx));
688 G4double gr = (1.+2.*A)*G4Log(1.+1./A)-2.;
689 G4double dum0 = (tau+2.)/(tau+1.);
690 G4double dum1 = tau+1.;
691 G4double gm = G4Log(0.5*tau/tauCut) + (1.+dum0*dum0)*G4Log(2.*(tau-tauCut+2.)/(tau+4.))
692 - 0.25*(tau+2.)*( tau+2.+2.*(2.*tau+1.)/(dum1*dum1))*
693 G4Log((tau+4.)*(tau-tauCut)/tau/(tau-tauCut+2.))
694 + 0.5*(tau-2*tauCut)*(tau+2.)*(1./(tau-tauCut)-1./(dum1*dum1));
695 if (gm<gr) {
696 gm = gm/gr;
697 } else {
698 gm = 1.;
699 }
701 scpCorr = 1.-gm*z0/(z0*(z0+1.));
702 }
703 fSCPCPerMatCuts[imc]->fVSCPC[ie] = scpCorr;
704 }
705 }
706}
static const G4int imax
std::vector< const G4Element * > G4ElementVector
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
static const G4int maxZ
G4double G4Log(G4double x)
Definition: G4Log.hh:226
std::vector< G4Material * > G4MaterialTable
@ idxG4ElectronCut
@ idxG4PositronCut
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4double A[17]
#define G4UniformRand()
Definition: Randomize.hh:52
static G4int GetMaxZet()
G4double SampleGSSRCosTheta(const GSMSCAngularDtr *gsDrt, G4double transfpar)
static std::vector< GSMSCAngularDtr * > gGSMSCAngularDistributions2
G4double SampleCosTheta(G4double lambdaval, G4double qval, G4double scra, G4double lekin, G4double beta2, G4int matindx, GSMSCAngularDtr **gsDtr, G4int &mcekini, G4int &mcdelti, G4double &transfPar, G4bool isfirst)
static constexpr G4double gQMAX1
static std::vector< GSMSCAngularDtr * > gGSMSCAngularDistributions1
G4double ComputeScatteringPowerCorrection(const G4MaterialCutsCouple *matcut, G4double ekin)
void GetMottCorrectionFactors(G4double logekin, G4double beta2, G4int matindx, G4double &mcToScr, G4double &mcToQ1, G4double &mcToG2PerG1)
G4double SingleScattering(G4double lambdaval, G4double scra, G4double lekin, G4double beta2, G4int matindx)
static constexpr G4double gQMIN1
G4bool Sampling(G4double lambdaval, G4double qval, G4double scra, G4double &cost, G4double &sint, G4double lekin, G4double beta2, G4int matindx, GSMSCAngularDtr **gsDtr, G4int &mcekini, G4int &mcdelti, G4double &transfPar, G4bool isfirst)
G4GoudsmitSaundersonTable(G4bool iselectron)
GSMSCAngularDtr * GetGSAngularDtr(G4double scra, G4double &lambdaval, G4double &qval, G4double &transfpar)
static constexpr G4double gLAMBMAX
G4double GetMoliereBc(G4int matindx)
std::vector< SCPCorrection * > fSCPCPerMatCuts
void Initialise(G4double lownergylimit, G4double highenergylimit)
static constexpr G4double gQMIN2
static constexpr G4double gQMAX2
G4double GetMoliereXc2(G4int matindx)
static constexpr G4double gLAMBMIN
G4double GetZeffective() const
const G4Material * GetMaterial() const
G4double GetDensity() const
Definition: G4Material.hh:176
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:186
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:205
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:222
size_t GetNumberOfElements() const
Definition: G4Material.hh:182
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:202
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:672
size_t GetIndex() const
Definition: G4Material.hh:256
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
const std::vector< G4double > * GetEnergyCutsVector(std::size_t pcIdx) const
static G4ProductionCutsTable * GetProductionCutsTable()
static constexpr double electron_mass_c2
static constexpr double cm3
static constexpr double keV
static constexpr double twopi
Definition: SystemOfUnits.h:56
static constexpr double MeV
static constexpr double g
static constexpr double cm
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
string fname
Definition: test.py:308
int G4lrint(double ad)
Definition: templates.hh:134