Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PairProductionRelModel.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 file
30//
31//
32// File name: G4PairProductionRelModel
33//
34// Author: Andreas Schaelicke
35//
36// Creation date: 02.04.2009
37//
38// Modifications:
39// 20.03.17 Change LPMconstant such that it gives suppression variable 's'
40// that consistent to Migdal's one; fix a small bug in 'logTS1'
41// computation; suppression is consistent now with the one in the
42// brem. model (F.Hariri)
43// 28-05-18 New version with improved screening function approximation, improved
44// LPM function approximation, efficiency, documentation and cleanup.
45// Corrected call to selecting target atom in the final state sampling.
46// (M. Novak)
47//
48// Class Description:
49//
50// Main References:
51// J.W.Motz et.al., Rev. Mod. Phys. 41 (1969) 581.
52// S.Klein, Rev. Mod. Phys. 71 (1999) 1501.
53// T.Stanev et.al., Phys. Rev. D25 (1982) 1291.
54// M.L.Ter-Mikaelian, High-energy Electromagnetic Processes in Condensed Media,
55// Wiley, 1972.
56//
57// -------------------------------------------------------------------
58
61#include "G4SystemOfUnits.hh"
62#include "G4Gamma.hh"
63#include "G4Electron.hh"
64#include "G4Positron.hh"
66#include "G4LossTableManager.hh"
67#include "G4ModifiedTsai.hh"
68#include "G4Exp.hh"
69#include "G4Pow.hh"
70#include "G4AutoLock.hh"
71
73
74// LPM constant: \alpha(mc^2)^2/(4\pi*\hbar c)
76 CLHEP::fine_structure_const*CLHEP::electron_mass_c2*CLHEP::electron_mass_c2
77 /(4.*CLHEP::pi*CLHEP::hbarc);
78
79// abscissas and weights of an 8 point Gauss-Legendre quadrature
80// for numerical integration on [0,1]
82 1.98550718e-02, 1.01666761e-01, 2.37233795e-01, 4.08282679e-01,
83 5.91717321e-01, 7.62766205e-01, 8.98333239e-01, 9.80144928e-01
84};
86 5.06142681e-02, 1.11190517e-01, 1.56853323e-01, 1.81341892e-01,
87 1.81341892e-01, 1.56853323e-01, 1.11190517e-01, 5.06142681e-02
88};
89
90// elastic and inelatic radiation logarithms for light elements (where the
91// Thomas-Fermi model doesn't work): computed by using Dirac-Fock model of atom.
93 0.0, 5.3104, 4.7935, 4.7402, 4.7112, 4.6694, 4.6134, 4.5520
94};
96 0.0, 5.9173, 5.6125, 5.5377, 5.4728, 5.4174, 5.3688, 5.3236
97};
98
99// constant cross section factor
101 4.*CLHEP::fine_structure_const*CLHEP::classic_electr_radius
102 *CLHEP::classic_electr_radius;
103
104// gamma energy limit above which LPM suppression will be applied (if the
105// fIsUseLPMCorrection flag is true)
107
108// special data structure per element i.e. per Z
109std::vector<G4PairProductionRelModel::ElementData*> G4PairProductionRelModel::gElementData;
110
111// LPM supression functions evaluated at initialisation time
112G4PairProductionRelModel::LPMFuncs G4PairProductionRelModel::gLPMFuncs;
113
114namespace
115{
116 G4Mutex thePairProdRelMutex = G4MUTEX_INITIALIZER;
117}
118
119// CTR
121 const G4String& nam)
122 : G4VEmModel(nam), fIsUseLPMCorrection(true), fIsUseCompleteScreening(false),
123 fLPMEnergy(0.), fG4Calc(G4Pow::GetInstance()), fTheGamma(G4Gamma::Gamma()),
124 fTheElectron(G4Electron::Electron()), fThePositron(G4Positron::Positron()),
125 fParticleChange(nullptr)
126{
127 // gamma energy below which the parametrized atomic x-section is used (30 GeV)
128 fParametrizedXSectionThreshold = 30.0*CLHEP::GeV;
129 // gamma energy below the Coulomb correction is turned off (50 MeV)
130 fCoulombCorrectionThreshold = 50.0*CLHEP::MeV;
131 // set angular generator used in the final state kinematics computation
133}
134
135// DTR
137{
138 if (isFirstInstance) {
139 // clear ElementData container
140 for (auto const & ptr : gElementData) { delete ptr; }
141 gElementData.clear();
142 // clear LPMFunctions (if any)
144 gLPMFuncs.fLPMFuncG.clear();
145 gLPMFuncs.fLPMFuncPhi.clear();
146 gLPMFuncs.fIsInitialized = false;
147 }
148 }
149}
150
152 const G4DataVector& cuts)
153{
155
156 if (isFirstInstance || gElementData.empty()) {
157 // init element data and LPM funcs
158 G4AutoLock l(&thePairProdRelMutex);
159 if (gElementData.empty()) {
160 isFirstInstance = true;
161 gElementData.resize(gMaxZet+1, nullptr);
162 }
163 // static data should be initialised only in the one instance
164 InitialiseElementData();
166 InitLPMFunctions();
167 }
168 l.unlock();
169 }
170 // element selectors should be initialised in the master thread
171 if (IsMaster()) {
173 }
174}
175
181
183 G4double Z)
184{
185 G4double xSection = 0.0;
186 // check if LPM suppression needs to be used
187 const G4bool isLPM = (fIsUseLPMCorrection && gammaEnergy>gEgLPMActivation);
188 // determine the kinematical limits (taken into account the correction due to
189 // the way in which the Coulomb correction is applied i.e. avoid negative DCS)
190 const G4int iz = std::min(gMaxZet, G4lrint(Z));
191 const G4double eps0 = CLHEP::electron_mass_c2/gammaEnergy;
192 // Coulomb correction is always included in the DCS even below 50 MeV (note:
193 // that this DCS is only used to get the integrated x-section)
194 const G4double dmax = gElementData[iz]->fDeltaMaxHigh;
195 const G4double dmin = 4.*eps0*gElementData[iz]->fDeltaFactor;
196 const G4double eps1 = 0.5 - 0.5*std::sqrt(1.-dmin/dmax);
197 const G4double epsMin = std::max(eps0, eps1);
198 const G4double epsMax = 0.5; // DCS is symmetric around eps=0.5
199 // let Et be the total energy transferred to the e- or to the e+
200 // the [Et-min, Et-max] interval will be divided into i=1,2,..,n subintervals
201 // with width of dInterv = (Et-max - Et-min)/n and numerical integration will
202 // be done in each sub-inteval using the xi = (Et - Et_i-min)/dInterv variable
203 // that is in [0,1]. The 8-point GL q. is used for the integration on [0,1].
204 const G4int numSub = 2;
205 const G4double dInterv= (epsMax - epsMin)*gammaEnergy/G4double(numSub);
206 G4double minEti = epsMin*gammaEnergy; // Et-min i.e. Et_0-min
207 for (G4int i = 0; i < numSub; ++i) {
208 for (G4int ngl = 0; ngl < 8; ++ngl) {
209 const G4double Et = (minEti + gXGL[ngl]*dInterv);
210 const G4double xs = isLPM ? ComputeRelDXSectionPerAtom(Et, gammaEnergy, Z)
211 : ComputeDXSectionPerAtom(Et, gammaEnergy, Z);
212 xSection += gWGL[ngl]*xs;
213 }
214 // update minimum Et of the sub-inteval
215 minEti += dInterv;
216 }
217 // apply corrections of variable transformation and half interval integration
218 xSection = std::max(2.*xSection*dInterv, 0.);
219 return xSection;
220}
221
222// DCS WITHOUT LPM SUPPRESSION
223// Computes DCS value for a given target element (Z), initial gamma energy (Eg),
224// total energy transferred to one of the e-/e+ pair(Et) WITHOUT LPM suppression
225// The constant factor 4 \alpha r_0^2 Z (Z +\eta(Z)) is not included here and
226// the returned value will be differential in total energy transfer instead of
227// the eps=Et/Eg. The computed part of the DCS
228// NORMAL CASE: DEFAULT STTING (i.e. fIsUseCompleteScreening = FALSE)
229// ds/deps(Et,Eg,Z) = ds/deps(eps,Z) = (eps^2+(1-eps)^2)*[phi1(d)/4-ln(Z)/3-fc]
230// + 2*eps(1-eps)*[phi2(d)/4-ln(Z)/3-fc]/3 where the universal (in the TF model)
231// screening variable d=d(eps)=136Z^(-1/3)eps0/[eps*(1-eps)] with eps0=mc^2/Eg.
232// COMPLETE SCREENING (when d(eps) approx-equal-to 0) : NEED TO BE SET BY USER
233// ds/deps(Et,Eg,Z) = ds/deps(eps,Z) = (eps^2+(1-eps)^2+eps*(1-eps)/3)*[Lel-fc]
234// -eps(1-eps)/9 where Lel=phi1(0)/4-ln(Z)/3 is the elastic(coherent) radiation
235// logarithm, fc is the Coulomb correction and the relation phi2(0)/4-ln(Z)/3 =
236// phi1(0)/4-1/6-ln(Z)/3 = Lel-1/6 (due to phi2(0)=phi1(0)-2/3) was used.
238 G4double gammaEnergy,
239 G4double Z)
240{
241 G4double xSection = 0.;
242 const G4int iz = std::min(gMaxZet, G4lrint(Z));
243 const G4double eps = pEnergy/gammaEnergy;
244 const G4double epsm = 1.-eps;
245 const G4double dum = eps*epsm;
247 // complete screening:
248 const G4double Lel = gElementData[iz]->fLradEl;
249 const G4double fc = gElementData[iz]->fCoulomb;
250 xSection = (eps*eps + epsm*epsm + 2.*dum/3.)*(Lel-fc) - dum/9.;
251 } else {
252 // normal case:
253 const G4double eps0 = CLHEP::electron_mass_c2/gammaEnergy;
254 const G4double fc = gElementData[iz]->fCoulomb;
255 const G4double lnZ13 = gElementData[iz]->fLogZ13;
256 const G4double delta = gElementData[iz]->fDeltaFactor*eps0/dum;
257 G4double phi1, phi2;
258 ComputePhi12(delta, phi1, phi2);
259 xSection = (eps*eps + epsm*epsm)*(0.25*phi1-lnZ13-fc)
260 + 2.*dum*(0.25*phi2-lnZ13-fc)/3.;
261 }
262 // non-const. part of the DCS differential in total energy transfer not in eps
263 // ds/dEt=ds/deps deps/dEt with deps/dEt=1/Eg
264 return std::max(xSection, 0.0)/gammaEnergy;
265}
266
267// DCS WITH POSSIBLE LPM SUPPRESSION
268// Computes DCS value for a given target element (Z), initial gamma energy (Eg),
269// total energy transferred to one of the e-/e+ pair(Et) WITH LPM suppression.
270// For a given Z, the LPM suppression will depend on the material through the
271// LMP-Energy. This will determine the suppression variable s and the LPM sup-
272// pression functions xi(s), fi(s) and G(s).
273// The constant factor 4 \alpha r_0^2 Z (Z +\eta(Z)) is not included here and
274// the returned value will be differential in total energy transfer instead of
275// the eps=Et/Eg. The computed part of the DCS
276// NORMAL CASE: DEFAULT STTING (i.e. fIsUseCompleteScreening = FALSE)
277// ds/deps(Et,Eg,Z)=ds/deps(eps,Z) = xi(s)*{ (eps^2+(1-eps)^2)*[2fi(s)/3+G(s)/3]
278// *[phi1(d)/4-ln(Z)/3-fc] + 2*eps(1-eps)*G(s)*[phi2(d)/4-ln(Z)/3-fc]/3 } where
279// the universal (in the TF model) screening variable d=d(eps)=136Z^(-1/3)eps0
280// /[eps*(1-eps)] with eps0=mc^2/Eg.
281// COMPLETE SCREENING (when d(eps) approx-equal-to 0) : NEED TO BE SET BY USER
282// ds/deps(Et,Eg,Z) = ds/deps(eps,Z) = xi(s)*{ [Lel-fc]*[ (eps^2+(1-eps)^2+eps
283// *(1-eps)/3)*2fi(s)/3 + G(s)/3] - eps(1-eps)*G(s)/9 }
284// Note, that when the LPM suppression is absent i.e. xi(s)=fi(s)=G(s)=1, both
285// the normal and the complete screening DCS give back the NO-LMP case above.
287 G4double gammaEnergy,
288 G4double Z)
289{
290 G4double xSection = 0.;
291 const G4int iz = std::min(gMaxZet, G4lrint(Z));
292 const G4double eps = pEnergy/gammaEnergy;
293 const G4double epsm = 1.-eps;
294 const G4double dum = eps*epsm;
295 // evaluate LPM suppression functions
296 G4double fXiS, fGS, fPhiS;
297 ComputeLPMfunctions(fXiS, fGS, fPhiS, eps, gammaEnergy, iz);
299 // complete screening:
300 const G4double Lel = gElementData[iz]->fLradEl;
301 const G4double fc = gElementData[iz]->fCoulomb;
302 xSection = (Lel-fc)*((eps*eps+epsm*epsm)*2.*fPhiS + fGS)/3. - dum*fGS/9.;
303 } else {
304 // normal case:
305 const G4double eps0 = CLHEP::electron_mass_c2/gammaEnergy;
306 const G4double fc = gElementData[iz]->fCoulomb;
307 const G4double lnZ13 = gElementData[iz]->fLogZ13;
308 const G4double delta = gElementData[iz]->fDeltaFactor*eps0/dum;
309 G4double phi1, phi2;
310 ComputePhi12(delta, phi1, phi2);
311 xSection = (eps*eps + epsm*epsm)*(2.*fPhiS+fGS)*(0.25*phi1-lnZ13-fc)/3.
312 + 2.*dum*fGS*(0.25*phi2-lnZ13-fc)/3.;
313 }
314 // non-const. part of the DCS differential in total energy transfer not in eps
315 // ds/dEt=ds/deps deps/dEt with deps/dEt=1/Eg
316 return std::max(fXiS*xSection, 0.0)/gammaEnergy;
317}
318
321 G4double gammaEnergy, G4double Z, G4double, G4double, G4double)
322{
323 G4double crossSection = 0.0 ;
324 // check kinematical limit
325 if ( gammaEnergy <= 2.0*electron_mass_c2 ) { return crossSection; }
326 // compute the atomic cross section either by using x-section parametrization
327 // or by numerically integrationg the DCS (with or without LPM)
328 if ( gammaEnergy < fParametrizedXSectionThreshold) {
329 // using the parametrized cross sections (max up to 80 GeV)
330 crossSection = ComputeParametrizedXSectionPerAtom(gammaEnergy, Z);
331 } else {
332 // by numerical integration of the DCS:
333 // Computes the cross section with or without LPM suppression depending on
334 // settings (by default with if the gamma energy is above a given threshold)
335 // and using or not using complete sreening approximation (by default not).
336 // Only the dependent part is computed in the numerical integration of the DCS
337 // i.e. the result must be multiplied here with 4 \alpha r_0^2 Z(Z+\eta(Z))
338 crossSection = ComputeXSectionPerAtom(gammaEnergy, Z);
339 // apply the constant factors:
340 // - eta(Z) is a correction to account interaction in the field of e-
341 // - gXSecFactor = 4 \alpha r_0^2
342 const G4int iz = std::min(gMaxZet, G4lrint(Z));
343 const G4double eta = gElementData[iz]->fEtaValue;
344 crossSection *= gXSecFactor*Z*(Z+eta);
345 }
346 // final protection
347 return std::max(crossSection, 0.);
348}
349
355
356void
357G4PairProductionRelModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
358 const G4MaterialCutsCouple* couple,
359 const G4DynamicParticle* aDynamicGamma,
360 G4double,
361 G4double)
362// The secondaries e+e- energies are sampled using the Bethe - Heitler
363// cross sections with Coulomb correction.
364// A modified version of the random number techniques of Butcher & Messel
365// is used (Nuc Phys 20(1960),15).
366//
367// GEANT4 internal units.
368//
369// Note 1 : Effects due to the breakdown of the Born approximation at
370// low energy are ignored.
371// Note 2 : The differential cross section implicitly takes account of
372// pair creation in both nuclear and atomic electron fields.
373// However triplet prodution is not generated.
374{
375 const G4Material* mat = couple->GetMaterial();
376 const G4double gammaEnergy = aDynamicGamma->GetKineticEnergy();
377 const G4double eps0 = CLHEP::electron_mass_c2/gammaEnergy ;
378 //
379 // check kinematical limit: gamma energy(Eg) must be at least 2 e- rest mass
380 // (but the model should be used at higher energies above 100 MeV)
381 if (eps0 > 0.5) { return; }
382 //
383 // select target atom of the material
384 const G4Element* anElement = SelectTargetAtom(couple, fTheGamma, gammaEnergy,
385 aDynamicGamma->GetLogKineticEnergy());
386 CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
387 //
388 // 'eps' is the total energy transferred to one of the e-/e+ pair in initial
389 // gamma energy units Eg. Since the corresponding DCS is symmetric on eps=0.5,
390 // the kinematical limits for eps0=mc^2/Eg <= eps <= 0.5
391 // 1. 'eps' is sampled uniformly on the [eps0, 0.5] inteval if Eg<Egsmall
392 // 2. otherwise, on the [eps_min, 0.5] interval according to the DCS (case 2.)
393 G4double eps;
394 // case 1.
395 static const G4double Egsmall = 2.*CLHEP::MeV;
396 if (gammaEnergy < Egsmall) {
397 eps = eps0 + (0.5-eps0)*rndmEngine->flat();
398 } else {
399 // case 2.
400 // get the Coulomb factor for the target element (Z) and gamma energy (Eg)
401 // F(Z) = 8*ln(Z)/3 if Eg <= 50 [MeV] => no Coulomb correction
402 // F(Z) = 8*ln(Z)/3 + 8*fc(Z) if Eg > 50 [MeV] => fc(Z) is the Coulomb cor.
403 //
404 // The screening variable 'delta(eps)' = 136*Z^{-1/3}*eps0/[eps(1-eps)]
405 // Due to the Coulomb correction, the DCS can go below zero even at
406 // kinematicaly allowed eps > eps0 values. In order to exclude this eps
407 // range with negative DCS, the minimum eps value will be set to eps_min =
408 // max[eps0, epsp] with epsp is the solution of SF(delta(epsp)) - F(Z)/2 = 0
409 // with SF being the screening function (SF1=SF2 at high value of delta).
410 // The solution is epsp = 0.5 - 0.5*sqrt[ 1 - 4*136*Z^{-1/3}eps0/deltap]
411 // with deltap = Exp[(42.038-F(Z))/8.29]-0.958. So the limits are:
412 // - when eps=eps_max = 0.5 => delta_min = 136*Z^{-1/3}*eps0/4
413 // - epsp = 0.5 - 0.5*sqrt[ 1 - delta_min/deltap]
414 // - and eps_min = max[eps0, epsp]
415 const G4int iZet = std::min(gMaxZet, anElement->GetZasInt());
416 const G4double deltaFactor = gElementData[iZet]->fDeltaFactor*eps0;
417 const G4double deltaMin = 4.*deltaFactor;
418 G4double deltaMax = gElementData[iZet]->fDeltaMaxLow;
419 G4double FZ = 8.*gElementData[iZet]->fLogZ13;
420 if ( gammaEnergy > fCoulombCorrectionThreshold ) { // Eg > 50 MeV ?
421 FZ += 8.*gElementData[iZet]->fCoulomb;
422 deltaMax = gElementData[iZet]->fDeltaMaxHigh;
423 }
424 // compute the limits of eps
425 const G4double epsp = 0.5 - 0.5*std::sqrt(1. - deltaMin/deltaMax) ;
426 const G4double epsMin = std::max(eps0,epsp);
427 const G4double epsRange = 0.5 - epsMin;
428 //
429 // sample the energy rate (eps) of the created electron (or positron)
431 ScreenFunction12(deltaMin, F10, F20);
432 F10 -= FZ;
433 F20 -= FZ;
434 const G4double NormF1 = std::max(F10 * epsRange * epsRange, 0.);
435 const G4double NormF2 = std::max(1.5 * F20 , 0.);
436 const G4double NormCond = NormF1/(NormF1 + NormF2);
437 // check if LPM correction is active
438 const G4bool isLPM = (fIsUseLPMCorrection && gammaEnergy>gEgLPMActivation);
440 // we will need 3 uniform random number for each trial of sampling
441 G4double rndmv[3];
442 G4double greject = 0.;
443 do {
444 rndmEngine->flatArray(3, rndmv);
445 if (NormCond > rndmv[0]) {
446 eps = 0.5 - epsRange * fG4Calc->A13(rndmv[1]);
447 const G4double delta = deltaFactor/(eps*(1.-eps));
448 if (isLPM) {
449 G4double lpmXiS, lpmGS, lpmPhiS, phi1, phi2;
450 ComputePhi12(delta, phi1, phi2);
451 ComputeLPMfunctions(lpmXiS, lpmGS, lpmPhiS, eps, gammaEnergy, iZet);
452 greject = lpmXiS*((2.*lpmPhiS+lpmGS)*phi1-lpmGS*phi2-lpmPhiS*FZ)/F10;
453 } else {
454 greject = (ScreenFunction1(delta)-FZ)/F10;
455 }
456 } else {
457 eps = epsMin + epsRange*rndmv[1];
458 const G4double delta = deltaFactor/(eps*(1.-eps));
459 if (isLPM) {
460 G4double lpmXiS, lpmGS, lpmPhiS, phi1, phi2;
461 ComputePhi12(delta, phi1, phi2);
462 ComputeLPMfunctions(lpmXiS, lpmGS, lpmPhiS, eps, gammaEnergy, iZet);
463 greject = lpmXiS*( (lpmPhiS+0.5*lpmGS)*phi1 + 0.5*lpmGS*phi2
464 -0.5*(lpmGS+lpmPhiS)*FZ )/F20;
465 } else {
466 greject = (ScreenFunction2(delta)-FZ)/F20;
467 }
468 }
469 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
470 } while (greject < rndmv[2]);
471 // end of eps sampling
472 }
473 //
474 // select charges randomly
475 G4double eTotEnergy, pTotEnergy;
476 if (rndmEngine->flat() > 0.5) {
477 eTotEnergy = (1.-eps)*gammaEnergy;
478 pTotEnergy = eps*gammaEnergy;
479 } else {
480 pTotEnergy = (1.-eps)*gammaEnergy;
481 eTotEnergy = eps*gammaEnergy;
482 }
483 //
484 // sample pair kinematics
485 //
486 const G4double eKinEnergy = std::max(0.,eTotEnergy - CLHEP::electron_mass_c2);
487 const G4double pKinEnergy = std::max(0.,pTotEnergy - CLHEP::electron_mass_c2);
488 //
489 G4ThreeVector eDirection, pDirection;
490 //
492 eKinEnergy, pKinEnergy, eDirection, pDirection);
493 // create G4DynamicParticle object for the particle1
494 auto aParticle1 = new G4DynamicParticle(fTheElectron,eDirection,eKinEnergy);
495
496 // create G4DynamicParticle object for the particle2
497 auto aParticle2 = new G4DynamicParticle(fThePositron,pDirection,pKinEnergy);
498 // Fill output vector
499 fvect->push_back(aParticle1);
500 fvect->push_back(aParticle2);
501 // kill incident photon
504}
505
506// should be called only by the master and at initialisation
507void G4PairProductionRelModel::InitialiseElementData()
508{
509 // create for all elements that are in the detector
510 auto elemTable = G4Element::GetElementTable();
511 for (auto const & elem : *elemTable) {
512 const G4int iz = std::min(gMaxZet, elem->GetZasInt());
513 if (nullptr == gElementData[iz]) { // create it if doesn't exist yet
514 const G4double logZ13 = elem->GetIonisation()->GetlogZ3();
515 const G4double Z13 = elem->GetIonisation()->GetZ3();
516 const G4double fc = elem->GetfCoulomb();
517 const G4double FZLow = 8.*logZ13;
518 const G4double FZHigh = 8.*(logZ13 + fc);
519 G4double Fel;
520 G4double Finel;
521 if (iz<5) { // use data from Dirac-Fock atomic model
522 Fel = gFelLowZet[iz];
523 Finel = gFinelLowZet[iz];
524 } else { // use the results of the Thomas-Fermi-Moliere model
525 Fel = G4Log(184.) - logZ13;
526 Finel = G4Log(1194.) - 2.*logZ13;
527 }
528 auto elD = new ElementData();
529 elD->fLogZ13 = logZ13;
530 elD->fCoulomb = fc;
531 elD->fLradEl = Fel;
532 elD->fDeltaFactor = 136./Z13;
533 elD->fDeltaMaxLow = G4Exp((42.038 - FZLow)/8.29) - 0.958;
534 elD->fDeltaMaxHigh = G4Exp((42.038 - FZHigh)/8.29) - 0.958;
535 elD->fEtaValue = Finel/(Fel-fc);
536 elD->fLPMVarS1Cond = std::sqrt(2.)*Z13*Z13/(184.*184.);
537 elD->fLPMILVarS1Cond = 1./G4Log(elD->fLPMVarS1Cond);
538 gElementData[iz] = elD;
539 }
540 }
541}
542
543// s goes up to 2 with ds = 0.01 be default
544void G4PairProductionRelModel::InitLPMFunctions() {
545 if (!gLPMFuncs.fIsInitialized) {
546 const G4int num = gLPMFuncs.fSLimit*gLPMFuncs.fISDelta+1;
547 gLPMFuncs.fLPMFuncG.resize(num);
548 gLPMFuncs.fLPMFuncPhi.resize(num);
549 for (G4int i=0; i<num; ++i) {
550 const G4double sval = i/gLPMFuncs.fISDelta;
551 ComputeLPMGsPhis(gLPMFuncs.fLPMFuncG[i],gLPMFuncs.fLPMFuncPhi[i],sval);
552 }
553 gLPMFuncs.fIsInitialized = true;
554 }
555}
556
557// used only at initialisation time
558void G4PairProductionRelModel::ComputeLPMGsPhis(G4double &funcGS, G4double &funcPhiS, const G4double varShat) {
559 if (varShat < 0.01) {
560 funcPhiS = 6.0*varShat*(1.0-CLHEP::pi*varShat);
561 funcGS = 12.0*varShat-2.0*funcPhiS;
562 } else {
563 const G4double varShat2 = varShat*varShat;
564 const G4double varShat3 = varShat*varShat2;
565 const G4double varShat4 = varShat2*varShat2;
566 if (varShat < 0.415827397755) { // Stanev ap.: for \psi(s) and compute G(s)
567 funcPhiS = 1.0-G4Exp( -6.0*varShat*(1.0+varShat*(3.0-CLHEP::pi))
568 + varShat3/(0.623+0.796*varShat+0.658*varShat2));
569 // 1-\exp \left\{-4s-\frac{8s^2}{1+3.936s+4.97s^2-0.05s^3+7.5s^4} \right\}
570 const G4double funcPsiS = 1.0-G4Exp( -4.0*varShat - 8.0*varShat2/(1.0
571 + 3.936*varShat+4.97*varShat2-0.05*varShat3+7.5*varShat4));
572 // G(s) = 3 \psi(s) - 2 \phi(s)
573 funcGS = 3.0*funcPsiS - 2.0*funcPhiS;
574 } else if (varShat < 1.55) {
575 funcPhiS = 1.0-G4Exp( -6.0*varShat*(1.0+varShat*(3.0-CLHEP::pi))
576 + varShat3/(0.623+0.796*varShat+0.658*varShat2));
577 const G4double dum0 = -0.16072300849123999+3.7550300067531581*varShat
578 -1.7981383069010097 *varShat2
579 +0.67282686077812381*varShat3
580 -0.1207722909879257 *varShat4;
581 funcGS = std::tanh(dum0);
582 } else {
583 funcPhiS = 1.0-0.01190476/varShat4;
584 if (varShat < 1.9156) {
585 const G4double dum0 = -0.16072300849123999+3.7550300067531581*varShat
586 -1.7981383069010097 *varShat2
587 +0.67282686077812381*varShat3
588 -0.1207722909879257 *varShat4;
589 funcGS = std::tanh(dum0);
590 } else {
591 funcGS = 1.0-0.0230655/varShat4;
592 }
593 }
594 }
595}
596
597// used at run-time to get some pre-computed LPM function values
598void G4PairProductionRelModel::GetLPMFunctions(G4double &lpmGs,
599 G4double &lpmPhis,
600 const G4double sval) {
601 if (sval < gLPMFuncs.fSLimit) {
602 G4double val = sval*gLPMFuncs.fISDelta;
603 const G4int ilow = (G4int)val;
604 val -= ilow;
605 lpmGs = (gLPMFuncs.fLPMFuncG[ilow+1]-gLPMFuncs.fLPMFuncG[ilow])*val
606 + gLPMFuncs.fLPMFuncG[ilow];
607 lpmPhis = (gLPMFuncs.fLPMFuncPhi[ilow+1]-gLPMFuncs.fLPMFuncPhi[ilow])*val
608 + gLPMFuncs.fLPMFuncPhi[ilow];
609 } else {
610 G4double ss = sval*sval;
611 ss *= ss;
612 lpmPhis = 1.0-0.01190476/ss;
613 lpmGs = 1.0-0.0230655/ss;
614 }
615}
616
617void G4PairProductionRelModel::ComputeLPMfunctions(G4double &funcXiS,
618 G4double &funcGS, G4double &funcPhiS, const G4double eps,
619 const G4double egamma, const G4int izet)
620{
621 // 1. y = E_+/E_{\gamma} with E_+ being the total energy transfered
622 // to one of the e-/e+ pair
623 // s' = \sqrt{ \frac{1}{8} \frac{1}{y(1-y)} \frac{E^{KL}_{LPM}}{E_{\gamma}} }
624 const G4double varSprime = std::sqrt(0.125*fLPMEnergy/(eps*egamma*(1.0-eps)));
625 const G4double condition = gElementData[izet]->fLPMVarS1Cond;
626 funcXiS = 2.0;
627 if (varSprime > 1.0) {
628 funcXiS = 1.0;
629 } else if (varSprime > condition) {
630 const G4double dum = gElementData[izet]->fLPMILVarS1Cond;
631 const G4double funcHSprime = G4Log(varSprime)*dum;
632 funcXiS = 1.0 + funcHSprime
633 - 0.08*(1.0-funcHSprime)*funcHSprime*(2.0-funcHSprime)*dum;
634 }
635 // 2. s=\frac{s'}{\sqrt{\xi(s')}}
636 const G4double varShat = varSprime / std::sqrt(funcXiS);
637 GetLPMFunctions(funcGS, funcPhiS, varShat);
638 // MAKE SURE SUPPRESSION IS SMALLER THAN 1: due to Migdal's approximation on xi
639 if (funcXiS * funcPhiS > 1. || varShat > 0.57) {
640 funcXiS = 1. / funcPhiS;
641 }
642}
643
644// Calculates the microscopic cross section in GEANT4 internal units. Same as in
645// G4BetheHeitlerModel and should be used below 80 GeV since it start to deverge
646// from the cross section data above 80-90 GeV:
647// Parametrized formula (L. Urban) is used to estimate the atomic cross sections
648// given numerically in the table of [Hubbell, J. H., Heinz Albert Gimm, and I.
649// Overbo: "Pair, Triplet, and Total Atomic Cross Sections (and Mass Attenuation
650// Coefficients) for 1 MeV‐100 GeV Photons in Elements Z= 1 to 100." Journal of
651// physical and chemical reference data 9.4 (1980): 1023-1148.]
652//
653// The formula gives a good approximation of the data from 1.5 MeV to 100 GeV.
654// below 1.5 MeV: sigma=sigma(1.5MeV)*(GammaEnergy-2electronmass)
655// *(GammaEnergy-2electronmass)
658 G4double Z)
659{
660 G4double xSection = 0.0 ;
661 // short versions
662 static const G4double kMC2 = CLHEP::electron_mass_c2;
663 // zero cross section below the kinematical limit: Eg<2mc^2
664 if (Z < 0.9 || gammaE <= 2.0*kMC2) { return xSection; }
665 //
666 static const G4double gammaEnergyLimit = 1.5*CLHEP::MeV;
667 // set coefficients a, b c
668 static const G4double a0 = 8.7842e+2*CLHEP::microbarn;
669 static const G4double a1 = -1.9625e+3*CLHEP::microbarn;
670 static const G4double a2 = 1.2949e+3*CLHEP::microbarn;
671 static const G4double a3 = -2.0028e+2*CLHEP::microbarn;
672 static const G4double a4 = 1.2575e+1*CLHEP::microbarn;
673 static const G4double a5 = -2.8333e-1*CLHEP::microbarn;
674
675 static const G4double b0 = -1.0342e+1*CLHEP::microbarn;
676 static const G4double b1 = 1.7692e+1*CLHEP::microbarn;
677 static const G4double b2 = -8.2381 *CLHEP::microbarn;
678 static const G4double b3 = 1.3063 *CLHEP::microbarn;
679 static const G4double b4 = -9.0815e-2*CLHEP::microbarn;
680 static const G4double b5 = 2.3586e-3*CLHEP::microbarn;
681
682 static const G4double c0 = -4.5263e+2*CLHEP::microbarn;
683 static const G4double c1 = 1.1161e+3*CLHEP::microbarn;
684 static const G4double c2 = -8.6749e+2*CLHEP::microbarn;
685 static const G4double c3 = 2.1773e+2*CLHEP::microbarn;
686 static const G4double c4 = -2.0467e+1*CLHEP::microbarn;
687 static const G4double c5 = 6.5372e-1*CLHEP::microbarn;
688 // check low energy limit of the approximation (1.5 MeV)
689 G4double gammaEnergyOrg = gammaE;
690 if (gammaE < gammaEnergyLimit) { gammaE = gammaEnergyLimit; }
691 // compute gamma energy variables
692 const G4double x = G4Log(gammaE/kMC2);
693 const G4double x2 = x *x;
694 const G4double x3 = x2*x;
695 const G4double x4 = x3*x;
696 const G4double x5 = x4*x;
697 //
698 const G4double F1 = a0 + a1*x + a2*x2 + a3*x3 + a4*x4 + a5*x5;
699 const G4double F2 = b0 + b1*x + b2*x2 + b3*x3 + b4*x4 + b5*x5;
700 const G4double F3 = c0 + c1*x + c2*x2 + c3*x3 + c4*x4 + c5*x5;
701 // compute the approximated cross section
702 xSection = (Z + 1.)*(F1*Z + F2*Z*Z + F3);
703 // check if we are below the limit of the approximation and apply correction
704 if (gammaEnergyOrg < gammaEnergyLimit) {
705 const G4double dum = (gammaEnergyOrg-2.*kMC2)/(gammaEnergyLimit-2.*kMC2);
706 xSection *= dum*dum;
707 }
708 return xSection;
709}
710
711
#define F10
#define F20
G4double condition(const G4ErrorSymMatrix &m)
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition G4Exp.hh:180
const G4double a0
G4double G4Log(G4double x)
Definition G4Log.hh:227
#define elem(i, j)
#define G4MUTEX_INITIALIZER
std::mutex G4Mutex
@ fStopAndKill
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
virtual double flat()=0
virtual void flatArray(const int size, double *vect)=0
G4double GetLogKineticEnergy() const
G4double GetKineticEnergy() const
static G4ElementTable * GetElementTable()
Definition G4Element.cc:389
G4int GetZasInt() const
Definition G4Element.hh:120
const G4Material * GetMaterial() const
G4double GetRadlen() const
G4double ComputeRelDXSectionPerAtom(G4double eplusEnergy, G4double gammaEnergy, G4double Z)
static const G4double gFelLowZet[8]
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4double ComputeParametrizedXSectionPerAtom(G4double gammaEnergy, G4double Z)
static const G4double gEgLPMActivation
void ScreenFunction12(const G4double delta, G4double &f1, G4double &f2)
static const G4double gFinelLowZet[8]
static std::vector< ElementData * > gElementData
G4PairProductionRelModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="BetheHeitlerLPM")
static const G4double gXGL[8]
G4double ScreenFunction1(const G4double delta)
void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double) override
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cut=0., G4double emax=DBL_MAX) override
G4ParticleDefinition * fThePositron
void ComputePhi12(const G4double delta, G4double &phi1, G4double &phi2)
G4double ComputeXSectionPerAtom(G4double gammaEnergy, G4double Z)
G4double ScreenFunction2(const G4double delta)
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4double ComputeDXSectionPerAtom(G4double eplusEnergy, G4double gammaEnergy, G4double Z)
void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
static const G4double gWGL[8]
G4ParticleChangeForGamma * fParticleChange
G4ParticleDefinition * fTheElectron
void SetProposedKineticEnergy(G4double proposedKinEnergy)
Definition G4Pow.hh:49
G4double A13(G4double A) const
Definition G4Pow.cc:116
virtual void SamplePairDirections(const G4DynamicParticle *dp, G4double elecKinEnergy, G4double posiKinEnergy, G4ThreeVector &dirElectron, G4ThreeVector &dirPositron, G4int Z=0, const G4Material *mat=nullptr)
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
G4VEmAngularDistribution * GetAngularDistribution()
G4ParticleChangeForGamma * GetParticleChangeForGamma()
std::vector< G4EmElementSelector * > * GetElementSelectors()
G4bool IsMaster() const
void SetAngularDistribution(G4VEmAngularDistribution *)
const G4Element * SelectTargetAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
void ProposeTrackStatus(G4TrackStatus status)
int G4lrint(double ad)
Definition templates.hh:134