Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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// [1] A.F.Bielajew, NIMB, 111 (1996) 195-208
65// [2] I.Kawrakow, A.F.Bielajew, NIMB 134(1998) 325-336
66// 28.04.2017 M. Novak: New representation of the angular distribution data with
67// significantly reduced data size.
68// 23.08.2017 M. Novak: Added funtionality to handle Mott-correction to the
69// base GS angular distributions and some other factors (screening
70// parameter, first and second moments) when Mott-correction is
71// activated in the GS-MSC model.
72// 26.10.2025 M. Novak: added the related technical note as the proper reference.
73//
74// References:
75// M. Novak: https://arxiv.org/abs/2410.13361
76// -----------------------------------------------------------------------------
77
79
80
82#include "Randomize.hh"
83#include "G4Log.hh"
84#include "G4Exp.hh"
85
86#include "G4GSMottCorrection.hh"
87#include "G4MaterialTable.hh"
88#include "G4Material.hh"
91#include "G4EmParameters.hh"
92
93#include "G4String.hh"
94
95#include <fstream>
96#include <cstdlib>
97#include <cmath>
98
99#include <iostream>
100#include <iomanip>
101
102// perecomputed GS angular distributions, based on the Screened-Rutherford DCS
103// are the same for e- and e+ so make sure we load them only onece
104G4bool G4GoudsmitSaundersonTable::gIsInitialised = false;
105//
106std::vector<G4GoudsmitSaundersonTable::GSMSCAngularDtr*> G4GoudsmitSaundersonTable::gGSMSCAngularDistributions1;
107std::vector<G4GoudsmitSaundersonTable::GSMSCAngularDtr*> G4GoudsmitSaundersonTable::gGSMSCAngularDistributions2;
108//
109std::vector<double> G4GoudsmitSaundersonTable::gMoliereBc;
110std::vector<double> G4GoudsmitSaundersonTable::gMoliereXc2;
111
112
114 fIsElectron = iselectron;
115 // set initial values: final values will be set in the Initialize method
116 fLogLambda0 = 0.; // will be set properly at init.
117 fLogDeltaLambda = 0.; // will be set properly at init.
118 fInvLogDeltaLambda = 0.; // will be set properly at init.
119 fInvDeltaQ1 = 0.; // will be set properly at init.
120 fDeltaQ2 = 0.; // will be set properly at init.
121 fInvDeltaQ2 = 0.; // will be set properly at init.
122 //
123 fLowEnergyLimit = 0.1*CLHEP::keV; // will be set properly at init.
124 fHighEnergyLimit = 100.0*CLHEP::MeV; // will be set properly at init.
125 //
126 fIsMottCorrection = false; // will be set properly at init.
127 fIsPWACorrection = false; // will be set properly at init.
128 fMottCorrection = nullptr;
129 //
130 fNumSPCEbinPerDec = 3;
131}
132
134 for (std::size_t i=0; i<gGSMSCAngularDistributions1.size(); ++i) {
135 if (gGSMSCAngularDistributions1[i]) {
136 delete [] gGSMSCAngularDistributions1[i]->fUValues;
137 delete [] gGSMSCAngularDistributions1[i]->fParamA;
138 delete [] gGSMSCAngularDistributions1[i]->fParamB;
139 delete gGSMSCAngularDistributions1[i];
140 }
141 }
142 gGSMSCAngularDistributions1.clear();
143 for (std::size_t i=0; i<gGSMSCAngularDistributions2.size(); ++i) {
144 if (gGSMSCAngularDistributions2[i]) {
145 delete [] gGSMSCAngularDistributions2[i]->fUValues;
146 delete [] gGSMSCAngularDistributions2[i]->fParamA;
147 delete [] gGSMSCAngularDistributions2[i]->fParamB;
148 delete gGSMSCAngularDistributions2[i];
149 }
150 }
151 gGSMSCAngularDistributions2.clear();
152 if (fMottCorrection) {
153 delete fMottCorrection;
154 fMottCorrection = nullptr;
155 }
156 // clear scp correction data
157 for (std::size_t imc=0; imc<fSCPCPerMatCuts.size(); ++imc) {
158 if (fSCPCPerMatCuts[imc]) {
159 fSCPCPerMatCuts[imc]->fVSCPC.clear();
160 delete fSCPCPerMatCuts[imc];
161 }
162 }
163 fSCPCPerMatCuts.clear();
164 //
165 gIsInitialised = false;
166}
167
168void G4GoudsmitSaundersonTable::Initialise(G4double lownergylimit, G4double highenergylimit) {
169 fLowEnergyLimit = lownergylimit;
170 fHighEnergyLimit = highenergylimit;
171 G4double lLambdaMin = G4Log(gLAMBMIN);
172 G4double lLambdaMax = G4Log(gLAMBMAX);
173 fLogLambda0 = lLambdaMin;
174 fLogDeltaLambda = (lLambdaMax-lLambdaMin)/(gLAMBNUM-1.);
175 fInvLogDeltaLambda = 1./fLogDeltaLambda;
176 fInvDeltaQ1 = 1./((gQMAX1-gQMIN1)/(gQNUM1-1.));
177 fDeltaQ2 = (gQMAX2-gQMIN2)/(gQNUM2-1.);
178 fInvDeltaQ2 = 1./fDeltaQ2;
179 // load precomputed angular distributions and set up several values used during the sampling
180 // these are particle independet => they go to static container: load them only onece
181 if (!gIsInitialised) {
182 // load pre-computed GS angular distributions (computed based on Screened-Rutherford DCS)
183 LoadMSCData();
184 gIsInitialised = true;
185 }
186 InitMoliereMSCParams();
187 // Mott-correction: particle(e- or e+) dependet so init them
188 if (fIsMottCorrection) {
189 if (!fMottCorrection) {
190 fMottCorrection = new G4GSMottCorrection(fIsElectron);
191 }
192 fMottCorrection->Initialise();
193 }
194 // init scattering power correction data; used only together with Mott-correction
195 // (Moliere's parameters must be initialised before)
196 if (fMottCorrection) {
198 }
199}
200
201
202// samplig multiple scattering angles cos(theta) and sin(thata)
203// - including no-scattering, single, "few" scattering cases as well
204// - Mott-correction will be included if it was requested by the user (i.e. if fIsMottCorrection=true)
205// lambdaval : s/lambda_el
206// qval : s/lambda_el G1
207// scra : screening parameter
208// cost : will be the smapled cos(theta)
209// sint : will be the smapled sin(theta)
210// lekin : logarithm of the current kinetic energy
211// beta2 : the corresponding beta square
212// matindx : index of the current material
213// returns true if it was msc
215 G4double &sint, G4double lekin, G4double beta2, G4int matindx,
216 GSMSCAngularDtr **gsDtr, G4int &mcekini, G4int &mcdelti,
217 G4double &transfPar, G4bool isfirst) {
218 G4double rand0 = G4UniformRand();
219 G4double expn = G4Exp(-lambdaval);
220 //
221 // no scattering case
222 if (rand0<expn) {
223 cost = 1.0;
224 sint = 0.0;
225 return false;
226 }
227 //
228 // single scattering case : sample from the single scattering PDF
229 // - Mott-correction will be included if it was requested by the user (i.e. if fIsMottCorrection=true)
230 if (rand0<(1.+lambdaval)*expn) {
231 // cost is sampled in SingleScattering()
232 cost = SingleScattering(lambdaval, scra, lekin, beta2, matindx);
233 // add protections
234 if (cost<-1.0) cost = -1.0;
235 if (cost>1.0) cost = 1.0;
236 // compute sin(theta) from the sampled cos(theta)
237 G4double dum0 = 1.-cost;
238 sint = std::sqrt(dum0*(2.0-dum0));
239 return false;
240 }
241 //
242 // handle this case:
243 // -lambdaval < 1 i.e. mean #elastic events along the step is < 1 but
244 // the currently sampled case is not 0 or 1 scattering. [Our minimal
245 // lambdaval (that we have precomputed, transformed angular distributions
246 // stored in a form of equally probabe intervalls together with rational
247 // interp. parameters) is 1.]
248 // -probability of having n elastic events follows Poisson stat. with
249 // lambdaval parameter.
250 // -the max. probability (when lambdaval=1) of having more than one
251 // elastic events is 0.2642411 and the prob of having 2,3,..,n elastic
252 // events decays rapidly with n. So set a max n to 10.
253 // -sampling of this cases is done in a one-by-one single elastic event way
254 // where the current #elastic event is sampled from the Poisson distr.
255 if (lambdaval<1.0) {
256 G4double prob, cumprob;
257 prob = cumprob = expn;
258 G4double curcost,cursint;
259 // init cos(theta) and sin(theta) to the zero scattering values
260 cost = 1.0;
261 sint = 0.0;
262 for (G4int iel=1; iel<10; ++iel) {
263 // prob of having iel scattering from Poisson
264 prob *= lambdaval/(G4double)iel;
265 cumprob += prob;
266 //
267 //sample cos(theta) from the singe scattering pdf:
268 // - Mott-correction will be included if it was requested by the user (i.e. if fIsMottCorrection=true)
269 curcost = SingleScattering(lambdaval, scra, lekin, beta2, matindx);
270 G4double dum0 = 1.-curcost;
271 cursint = dum0*(2.0-dum0); // sin^2(theta)
272 //
273 // if we got current deflection that is not too small
274 // then update cos(theta) sin(theta)
275 if (cursint>1.0e-20) {
276 cursint = std::sqrt(cursint);
277 G4double curphi = CLHEP::twopi*G4UniformRand();
278 cost = cost*curcost-sint*cursint*std::cos(curphi);
279 sint = std::sqrt(std::max(0.0, (1.0-cost)*(1.0+cost)));
280 }
281 //
282 // check if we have done enough scattering i.e. sampling from the Poisson
283 if (rand0<cumprob) {
284 return false;
285 }
286 }
287 // if reached the max iter i.e. 10
288 return false;
289 }
290 //
291 // multiple scattering case with lambdavalue >= 1:
292 // - use the precomputed and transformed Goudsmit-Saunderson angular
293 // distributions to sample cos(theta)
294 // - Mott-correction will be included if it was requested by the user (i.e. if fIsMottCorrection=true)
295 cost = SampleCosTheta(lambdaval, qval, scra, lekin, beta2, matindx, gsDtr, mcekini, mcdelti, transfPar, isfirst);
296 // add protections
297 if (cost<-1.0) cost = -1.0;
298 if (cost> 1.0) cost = 1.0;
299 // compute cos(theta) and sin(theta) from the sampled 1-cos(theta)
300 G4double dum0 = 1.0-cost;
301 sint = std::sqrt(dum0*(2.0-dum0));
302 // return true if it was msc
303 return true;
304}
305
306
308 G4double lekin, G4double beta2, G4int matindx,
309 GSMSCAngularDtr **gsDtr, G4int &mcekini,G4int &mcdelti,
310 G4double &transfPar, G4bool isfirst) {
311 G4double cost = 1.;
312 // determine the base GS angular distribution if it is the first call (when sub-step sampling is used)
313 if (isfirst) {
314 *gsDtr = GetGSAngularDtr(scra, lambdaval, qval, transfPar);
315 }
316 // sample cost from the GS angular distribution (computed based on Screened-Rutherford DCS)
317 cost = SampleGSSRCosTheta(*gsDtr, transfPar);
318 // Mott-correction if it was requested by the user
319 if (fIsMottCorrection && *gsDtr) { // no Mott-correction in case of izotropic theta
320 static const G4int nlooplim = 1000;
321 G4int nloop = 0 ; // rejection loop counter
322// G4int ekindx = -1; // evaluate only in the first call
323// G4int deltindx = -1 ; // evaluate only in the first call
324 G4double val = fMottCorrection->GetMottRejectionValue(lekin, beta2, qval, cost, matindx, mcekini, mcdelti);
325 while (G4UniformRand()>val && ++nloop<nlooplim) {
326 // sampling cos(theta)
327 cost = SampleGSSRCosTheta(*gsDtr, transfPar);
328 val = fMottCorrection->GetMottRejectionValue(lekin, beta2, qval, cost, matindx, mcekini, mcdelti);
329 };
330 }
331 return cost;
332}
333
334
335// returns with cost sampled from the GS angular distribution computed based on Screened-Rutherford DCS
337 // check if isotropic theta (i.e. cost is uniform on [-1:1])
338 if (!gsDtr) {
339 return 1.-2.0*G4UniformRand();
340 }
341 //
342 // sampling form the selected distribution
343 G4double ndatm1 = gsDtr->fNumData-1.;
344 G4double delta = 1.0/ndatm1;
345 // determine lower cumulative bin inidex
346 G4double rndm = G4UniformRand();
347 G4int indxl = rndm*ndatm1;
348 G4double aval = rndm-indxl*delta;
349 G4double dum0 = delta*aval;
350
351 G4double dum1 = (1.0+gsDtr->fParamA[indxl]+gsDtr->fParamB[indxl])*dum0;
352 G4double dum2 = delta*delta + gsDtr->fParamA[indxl]*dum0 + gsDtr->fParamB[indxl]*aval*aval;
353 G4double sample = gsDtr->fUValues[indxl] + dum1/dum2 *(gsDtr->fUValues[indxl+1]-gsDtr->fUValues[indxl]);
354 // transform back u to cos(theta) :
355 // this is the sampled cos(theta) = (2.0*para*sample)/(1.0-sample+para)
356 return 1.-(2.0*transfpar*sample)/(1.0-sample+transfpar);
357}
358
359
360// determine the GS angular distribution we need to sample from: will set other things as well ...
362 G4double &lambdaval, G4double &qval, G4double &transfpar) {
363 GSMSCAngularDtr *dtr = nullptr;
364 G4bool first = false;
365 // isotropic cost above gQMAX2 (i.e. dtr stays nullptr)
366 if (qval<gQMAX2) {
367 G4int lamIndx = -1; // lambda value index
368 G4int qIndx = -1; // lambda value index
369 // init to second grid Q values
370 G4int numQVal = gQNUM2;
371 G4double minQVal = gQMIN2;
372 G4double invDelQ = fInvDeltaQ2;
373 G4double pIndxH = 0.; // probability of taking higher index
374 // check if first or second grid needs to be used
375 if (qval<gQMIN2) { // first grid
376 first = true;
377 // protect against qval<gQMIN1
378 if (qval<gQMIN1) {
379 qval = gQMIN1;
380 qIndx = 0;
381 //pIndxH = 0.;
382 }
383 // set to first grid Q values
384 numQVal = gQNUM1;
385 minQVal = gQMIN1;
386 invDelQ = fInvDeltaQ1;
387 }
388 // make sure that lambda = s/lambda_el is in [gLAMBMIN,gLAMBMAX)
389 // lambda<gLAMBMIN=1 is already handeled before so lambda>= gLAMBMIN for sure
390 if (lambdaval>=gLAMBMAX) {
391 lambdaval = gLAMBMAX-1.e-8;
392 lamIndx = gLAMBNUM-1;
393 }
394 G4double lLambda = G4Log(lambdaval);
395 //
396 // determine lower lambda (=s/lambda_el) index: linear interp. on log(lambda) scale
397 if (lamIndx<0) {
398 pIndxH = (lLambda-fLogLambda0)*fInvLogDeltaLambda;
399 lamIndx = (G4int)(pIndxH); // lower index of the lambda bin
400 pIndxH = pIndxH-lamIndx; // probability of taking the higher index distribution
401 if (G4UniformRand()<pIndxH) {
402 ++lamIndx;
403 }
404 }
405 //
406 // determine lower Q (=s/lambda_el G1) index: linear interp. on Q
407 if (qIndx<0) {
408 pIndxH = (qval-minQVal)*invDelQ;
409 qIndx = (G4int)(pIndxH); // lower index of the Q bin
410 pIndxH = pIndxH-qIndx;
411 if (G4UniformRand()<pIndxH) {
412 ++qIndx;
413 }
414 }
415 // set indx
416 G4int indx = lamIndx*numQVal+qIndx;
417 if (first) {
418 dtr = gGSMSCAngularDistributions1[indx];
419 } else {
420 dtr = gGSMSCAngularDistributions2[indx];
421 }
422 // dtr might be nullptr that indicates isotropic cot distribution because:
423 // - if the selected lamIndx, qIndx correspond to L(=s/lambda_el) and Q(=s/lambda_el G1) such that G1(=Q/L) > 1
424 // G1 should always be < 1 and if G1 is ~1 -> the dtr is isotropic (this can only happen in case of the 2. grid)
425 //
426 // compute the transformation parameter
427 if (lambdaval>10.0) {
428 transfpar = 0.5*(-2.77164+lLambda*( 2.94874-lLambda*(0.1535754-lLambda*0.00552888) ));
429 } else {
430 transfpar = 0.5*(1.347+lLambda*(0.209364-lLambda*(0.45525-lLambda*(0.50142-lLambda*0.081234))));
431 }
432 transfpar *= (lambdaval+4.0)*scra;
433 }
434 // return with the selected GS angular distribution that we need to sample cost from (if nullptr => isotropic cost)
435 return dtr;
436}
437
438
440 gGSMSCAngularDistributions1.resize(gLAMBNUM*gQNUM1,nullptr);
441 const G4String str1 = G4EmParameters::Instance()->GetDirLEDATA() + "/msc_GS/GSGrid_1/gsDistr_";
442 for (G4int il=0; il<gLAMBNUM; ++il) {
443 G4String fname = str1 + std::to_string(il);
444 std::ifstream infile(fname,std::ios::in);
445 if (!infile.is_open()) {
446 G4String msgc = "Cannot open file: " + fname;
447 G4Exception("G4GoudsmitSaundersonTable::LoadMSCData()","em0006",
448 FatalException, msgc.c_str());
449 return;
450 }
451 for (G4int iq=0; iq<gQNUM1; ++iq) {
452 auto gsd = new GSMSCAngularDtr();
453 infile >> gsd->fNumData;
454 gsd->fUValues = new G4double[gsd->fNumData]();
455 gsd->fParamA = new G4double[gsd->fNumData]();
456 gsd->fParamB = new G4double[gsd->fNumData]();
457 G4double ddummy;
458 infile >> ddummy; infile >> ddummy;
459 for (G4int i=0; i<gsd->fNumData; ++i) {
460 infile >> gsd->fUValues[i];
461 infile >> gsd->fParamA[i];
462 infile >> gsd->fParamB[i];
463 }
464 gGSMSCAngularDistributions1[il*gQNUM1+iq] = gsd;
465 }
466 infile.close();
467 }
468 //
469 // second grid
470 gGSMSCAngularDistributions2.resize(gLAMBNUM*gQNUM2,nullptr);
471 const G4String str2 = G4EmParameters::Instance()->GetDirLEDATA() + "/msc_GS/GSGrid_2/gsDistr_";
472 for (G4int il=0; il<gLAMBNUM; ++il) {
473 G4String fname = str2 + std::to_string(il);
474 std::ifstream infile(fname,std::ios::in);
475 if (!infile.is_open()) {
476 G4String msgc = "Cannot open file: " + fname;
477 G4Exception("G4GoudsmitSaundersonTable::LoadMSCData()","em0006",
478 FatalException, msgc.c_str());
479 return;
480 }
481 for (G4int iq=0; iq<gQNUM2; ++iq) {
482 G4int numData;
483 infile >> numData;
484 if (numData>1) {
485 auto gsd = new GSMSCAngularDtr();
486 gsd->fNumData = numData;
487 gsd->fUValues = new G4double[gsd->fNumData]();
488 gsd->fParamA = new G4double[gsd->fNumData]();
489 gsd->fParamB = new G4double[gsd->fNumData]();
490 double ddummy;
491 infile >> ddummy; infile >> ddummy;
492 for (G4int i=0; i<gsd->fNumData; ++i) {
493 infile >> gsd->fUValues[i];
494 infile >> gsd->fParamA[i];
495 infile >> gsd->fParamB[i];
496 }
497 gGSMSCAngularDistributions2[il*gQNUM2+iq] = gsd;
498 } else {
499 gGSMSCAngularDistributions2[il*gQNUM2+iq] = nullptr;
500 }
501 }
502 infile.close();
503 }
504}
505
506// samples cost in single scattering based on Screened-Rutherford DCS
507// (with Mott-correction if it was requested)
509 G4double lekin, G4double beta2,
510 G4int matindx) {
511 G4double rand1 = G4UniformRand();
512 // sample cost from the Screened-Rutherford DCS
513 G4double cost = 1.-2.0*scra*rand1/(1.0-rand1+scra);
514 // Mott-correction if it was requested by the user
515 if (fIsMottCorrection) {
516 static const G4int nlooplim = 1000; // rejection loop limit
517 G4int nloop = 0 ; // loop counter
518 G4int ekindx = -1 ; // evaluate only in the first call
519 G4int deltindx = 0 ; // single scattering case
520 G4double q1 = 0.; // not used when deltindx = 0;
521 // computing Mott rejection function value
522 G4double val = fMottCorrection->GetMottRejectionValue(lekin, beta2, q1, cost,
523 matindx, ekindx, deltindx);
524 while (G4UniformRand()>val && ++nloop<nlooplim) {
525 // sampling cos(theta) from the Screened-Rutherford DCS
526 rand1 = G4UniformRand();
527 cost = 1.-2.0*scra*rand1/(1.0-rand1+scra);
528 // computing Mott rejection function value
529 val = fMottCorrection->GetMottRejectionValue(lekin, beta2, q1, cost, matindx,
530 ekindx, deltindx);
531 };
532 }
533 return cost;
534}
535
536
538 G4int matindx, G4double &mcToScr,
539 G4double &mcToQ1, G4double &mcToG2PerG1) {
540 if (fIsMottCorrection) {
541 fMottCorrection->GetMottCorrectionFactors(logekin, beta2, matindx, mcToScr, mcToQ1, mcToG2PerG1);
542 }
543}
544
545
546// compute material dependent Moliere MSC parameters at initialisation
547void G4GoudsmitSaundersonTable::InitMoliereMSCParams() {
548 const G4double const1 = 7821.6; // [cm2/g]
549 const G4double const2 = 0.1569; // [cm2 MeV2 / g]
550 const G4double finstrc2 = 5.325135453E-5; // fine-structure const. square
551
552 G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
553 // get number of materials in the table
554 std::size_t numMaterials = theMaterialTable->size();
555 // make sure that we have long enough vectors
556 if(gMoliereBc.size()<numMaterials) {
557 gMoliereBc.resize(numMaterials);
558 gMoliereXc2.resize(numMaterials);
559 }
560 G4double xi = 1.0;
561 G4int maxZ = 200;
562 if (fIsMottCorrection || fIsPWACorrection) {
563 // xi = 1.0; <= always set to 1 from now on
565 }
566 //
567 for (std::size_t imat=0; imat<numMaterials; ++imat) {
568 const G4Material* theMaterial = (*theMaterialTable)[imat];
569 const G4ElementVector* theElemVect = theMaterial->GetElementVector();
570 const G4int numelems = (G4int)theMaterial->GetNumberOfElements();
571 //
572 const G4double* theNbAtomsPerVolVect = theMaterial->GetVecNbOfAtomsPerVolume();
573 G4double theTotNbAtomsPerVol = theMaterial->GetTotNbOfAtomsPerVolume();
574 //
575 G4double zs = 0.0;
576 G4double zx = 0.0;
577 G4double ze = 0.0;
578 G4double sa = 0.0;
579 //
580 for(G4int ielem = 0; ielem < numelems; ielem++) {
581 G4double zet = (*theElemVect)[ielem]->GetZ();
582 if (zet>maxZ) {
583 zet = (G4double)maxZ;
584 }
585 G4double iwa = (*theElemVect)[ielem]->GetN();
586 G4double ipz = theNbAtomsPerVolVect[ielem]/theTotNbAtomsPerVol;
587 G4double dum = ipz*zet*(zet+xi);
588 zs += dum;
589 ze += dum*(-2.0/3.0)*G4Log(zet);
590 zx += dum*G4Log(1.0+3.34*finstrc2*zet*zet);
591 sa += ipz*iwa;
592 }
593 G4double density = theMaterial->GetDensity()*CLHEP::cm3/CLHEP::g; // [g/cm3]
594 //
595 gMoliereBc[theMaterial->GetIndex()] = const1*density*zs/sa*G4Exp(ze/zs)/G4Exp(zx/zs); //[1/cm]
596 gMoliereXc2[theMaterial->GetIndex()] = const2*density*zs/sa; // [MeV2/cm]
597 // change to Geant4 internal units of 1/length and energ2/length
598 gMoliereBc[theMaterial->GetIndex()] *= 1.0/CLHEP::cm;
599 gMoliereXc2[theMaterial->GetIndex()] *= CLHEP::MeV*CLHEP::MeV/CLHEP::cm;
600 }
601}
602
603
604// this method is temporary, will be removed/replaced with a more effictien solution after 10.3.ref09
606 G4int imc = matcut->GetIndex();
607 G4double corFactor = 1.0;
608 if (!(fSCPCPerMatCuts[imc]->fIsUse) || ekin<=fSCPCPerMatCuts[imc]->fPrCut) {
609 return corFactor;
610 }
611 // get the scattering power correction factor
612 G4double lekin = G4Log(ekin);
613 G4double remaining = (lekin-fSCPCPerMatCuts[imc]->fLEmin)*fSCPCPerMatCuts[imc]->fILDel;
614 std::size_t lindx = (std::size_t)remaining;
615 remaining -= lindx;
616 std::size_t imax = fSCPCPerMatCuts[imc]->fVSCPC.size()-1;
617 if (lindx>=imax) {
618 corFactor = fSCPCPerMatCuts[imc]->fVSCPC[imax];
619 } else {
620 corFactor = fSCPCPerMatCuts[imc]->fVSCPC[lindx] + remaining*(fSCPCPerMatCuts[imc]->fVSCPC[lindx+1]-fSCPCPerMatCuts[imc]->fVSCPC[lindx]);
621 }
622 return corFactor;
623}
624
625
627 // get the material-cuts table
629 std::size_t numMatCuts = thePCTable->GetTableSize();
630 // clear container if any
631 for (std::size_t imc=0; imc<fSCPCPerMatCuts.size(); ++imc) {
632 if (fSCPCPerMatCuts[imc]) {
633 fSCPCPerMatCuts[imc]->fVSCPC.clear();
634 delete fSCPCPerMatCuts[imc];
635 fSCPCPerMatCuts[imc] = nullptr;
636 }
637 }
638 //
639 // set size of the container and create the corresponding data structures
640 fSCPCPerMatCuts.resize(numMatCuts,nullptr);
641 // loop over the material-cuts and create scattering power correction data structure for each
642 for (G4int imc=0; imc<(G4int)numMatCuts; ++imc) {
643 const G4MaterialCutsCouple *matCut = thePCTable->GetMaterialCutsCouple(imc);
644 // caclulate sm := \sum_i n_i Z_i and sz := \sum_i ni Z_i(Z_i + \xi_0) with \xi_0=1
645 const G4Material* theMaterial = matCut->GetMaterial();
646 const G4ElementVector* theElemVect = theMaterial->GetElementVector();
647 const G4int numelems = (G4int)theMaterial->GetNumberOfElements();
648 const G4double* theNbAtomsPerVolVect = theMaterial->GetVecNbOfAtomsPerVolume();
649 G4double theTotNbAtomsPerVol = theMaterial->GetTotNbOfAtomsPerVolume();
650 G4double zs = 0.0;
651 G4double sm = 0.0;
652 for(G4int ielem = 0; ielem < numelems; ielem++) {
653 G4double zet = (*theElemVect)[ielem]->GetZ();
654 G4double ipz = theNbAtomsPerVolVect[ielem]/theTotNbAtomsPerVol;
655 zs = zs + ipz*zet*(zet + 1.0); // \xi_0 = 1
656 sm = sm + ipz*zet;
657 }
658 // ready to calculate the scattering power correction
659 // get e- production cut in the current material-cuts in energy
660 G4double limit;
661 G4double ecut;
662 if (fIsElectron) {
663 ecut = (*(thePCTable->GetEnergyCutsVector(idxG4ElectronCut)))[matCut->GetIndex()];
664 limit = 2.*ecut;
665 } else {
666 ecut = (*(thePCTable->GetEnergyCutsVector(idxG4PositronCut)))[matCut->GetIndex()];
667 limit = ecut;
668 }
669 G4double min = std::max(limit,fLowEnergyLimit);
670 G4double max = fHighEnergyLimit;
671 if (min>=max) {
672 fSCPCPerMatCuts[imc] = new SCPCorrection();
673 fSCPCPerMatCuts[imc]->fIsUse = false;
674 fSCPCPerMatCuts[imc]->fPrCut = min;
675 continue;
676 }
677 G4int numEbins = fNumSPCEbinPerDec*G4lrint(std::log10(max/min));
678 numEbins = std::max(numEbins,3);
679 G4double lmin = G4Log(min);
680 G4double ldel = G4Log(max/min)/(numEbins-1.0);
681 fSCPCPerMatCuts[imc] = new SCPCorrection();
682 fSCPCPerMatCuts[imc]->fVSCPC.resize(numEbins,1.0);
683 fSCPCPerMatCuts[imc]->fIsUse = true;
684 fSCPCPerMatCuts[imc]->fPrCut = min;
685 fSCPCPerMatCuts[imc]->fLEmin = lmin;
686 fSCPCPerMatCuts[imc]->fILDel = 1./ldel;
687 for (G4int ie=0; ie<numEbins; ++ie) {
688 G4double ekin = G4Exp(lmin+ie*ldel);
689 G4double scpCorr = 1.0;
690 // compute correction factor: I.Kawrakow NIMB 114(1996)307-326 (Eqs(32-37))
691 if (ie>0) {
692 G4double tau = ekin/CLHEP::electron_mass_c2;
693 G4double tauCut = ecut/CLHEP::electron_mass_c2;
694 // Moliere's screening parameter
695 G4int matindx = (G4int)matCut->GetMaterial()->GetIndex();
696 G4double A = GetMoliereXc2(matindx)/(4.0*tau*(tau+2.)*GetMoliereBc(matindx));
697 G4double gr = (1.+2.*A)*G4Log(1.+1./A)-2.;
698 G4double dum0 = (tau+2.)/(tau+1.);
699 G4double dum1 = tau+1.;
700 G4double gm = G4Log(0.5*tau/tauCut) + (1.+dum0*dum0)*G4Log(2.*(tau-tauCut+2.)/(tau+4.))
701 - 0.25*(tau+2.)*( tau+2.+2.*(2.*tau+1.)/(dum1*dum1))*
702 G4Log((tau+4.)*(tau-tauCut)/tau/(tau-tauCut+2.))
703 + 0.5*(tau-2*tauCut)*(tau+2.)*(1./(tau-tauCut)-1./(dum1*dum1));
704 if (gm<gr) {
705 gm = gm/gr;
706 } else {
707 gm = 1.;
708 }
709 scpCorr = 1.0 - gm*sm/zs;
710 }
711 fSCPCPerMatCuts[imc]->fVSCPC[ie] = scpCorr;
712 }
713 }
714}
std::vector< const G4Element * > G4ElementVector
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition G4Exp.hh:132
G4double G4Log(G4double x)
Definition G4Log.hh:169
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 G4EmParameters * Instance()
const G4String & GetDirLEDATA() const
G4double SampleGSSRCosTheta(const GSMSCAngularDtr *gsDrt, G4double transfpar)
G4double SampleCosTheta(G4double lambdaval, G4double qval, G4double scra, G4double lekin, G4double beta2, G4int matindx, GSMSCAngularDtr **gsDtr, G4int &mcekini, G4int &mcdelti, G4double &transfPar, G4bool isfirst)
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)
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)
GSMSCAngularDtr * GetGSAngularDtr(G4double scra, G4double &lambdaval, G4double &qval, G4double &transfpar)
void Initialise(G4double lownergylimit, G4double highenergylimit)
G4double GetMoliereXc2(G4int matindx)
const G4Material * GetMaterial() const
G4double GetDensity() const
const G4ElementVector * GetElementVector() const
G4double GetTotNbOfAtomsPerVolume() const
const G4double * GetVecNbOfAtomsPerVolume() const
std::size_t GetIndex() const
static G4MaterialTable * GetMaterialTable()
std::size_t GetNumberOfElements() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
const std::vector< G4double > * GetEnergyCutsVector(std::size_t pcIdx) const
static G4ProductionCutsTable * GetProductionCutsTable()
std::string to_string(G4FermiAtomicMass mass)
int G4lrint(double ad)
Definition templates.hh:134