Geant4 11.4.0
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 "G4LPMFunction.hh"
69#include "G4Exp.hh"
70#include "G4Pow.hh"
71#include "G4AutoLock.hh"
72
74
75// LPM constant: \alpha(mc^2)^2/(4\pi*\hbar c)
77 CLHEP::fine_structure_const*CLHEP::electron_mass_c2*CLHEP::electron_mass_c2
78 /(4.*CLHEP::pi*CLHEP::hbarc);
79
80// abscissas and weights of an 8 point Gauss-Legendre quadrature
81// for numerical integration on [0,1]
83 1.98550718e-02, 1.01666761e-01, 2.37233795e-01, 4.08282679e-01,
84 5.91717321e-01, 7.62766205e-01, 8.98333239e-01, 9.80144928e-01
85};
87 5.06142681e-02, 1.11190517e-01, 1.56853323e-01, 1.81341892e-01,
88 1.81341892e-01, 1.56853323e-01, 1.11190517e-01, 5.06142681e-02
89};
90
91// elastic and inelatic radiation logarithms for light elements (where the
92// Thomas-Fermi model doesn't work): computed by using Dirac-Fock model of atom.
94 0.0, 5.3104, 4.7935, 4.7402, 4.7112, 4.6694, 4.6134, 4.5520
95};
97 0.0, 5.9173, 5.6125, 5.5377, 5.4728, 5.4174, 5.3688, 5.3236
98};
99
100// constant cross section factor
102 4.*CLHEP::fine_structure_const*CLHEP::classic_electr_radius
103 *CLHEP::classic_electr_radius;
104
105// gamma energy limit above which LPM suppression will be applied (if the
106// fIsUseLPMCorrection flag is true)
108
109// special data structure per element i.e. per Z
110std::vector<G4PairProductionRelModel::ElementData*> G4PairProductionRelModel::gElementData;
111
112namespace
113{
114 G4Mutex thePairProdRelMutex = G4MUTEX_INITIALIZER;
115}
116
117// CTR
119 const G4String& nam)
121 fLPMEnergy(0.), fG4Calc(G4Pow::GetInstance()), fTheGamma(G4Gamma::Gamma()),
122 fTheElectron(G4Electron::Electron()), fThePositron(G4Positron::Positron()),
123 fParticleChange(nullptr)
124{
125 // gamma energy below which the parametrized atomic x-section is used (30 GeV)
126 fParametrizedXSectionThreshold = 30.0*CLHEP::GeV;
127 // gamma energy below the Coulomb correction is turned off (50 MeV)
128 fCoulombCorrectionThreshold = 50.0*CLHEP::MeV;
129 // set angular generator used in the final state kinematics computation
131}
132
133// DTR
135{
136 if (isFirstInstance) {
137 // clear ElementData container
138 for (auto const & ptr : gElementData) { delete ptr; }
139 gElementData.clear();
140 }
141}
142
144 const G4DataVector& cuts)
145{
147
148 if (isFirstInstance || gElementData.empty()) {
149 // init element data and LPM funcs
150 G4AutoLock l(&thePairProdRelMutex);
151 if (gElementData.empty()) {
152 isFirstInstance = true;
153 gElementData.resize(gMaxZet+1, nullptr);
154 }
155 // static data should be initialised only in the one instance
156 InitialiseElementData();
157 l.unlock();
158 }
159 // element selectors should be initialised in the master thread
160 if (IsMaster()) {
162 }
163}
164
170
172 G4double Z)
173{
174 G4double xSection = 0.0;
175 // check if LPM suppression needs to be used
176 const G4bool isLPM = (fIsUseLPMCorrection && gammaEnergy>gEgLPMActivation);
177 // determine the kinematical limits (taken into account the correction due to
178 // the way in which the Coulomb correction is applied i.e. avoid negative DCS)
179 const G4int iz = std::min(gMaxZet, G4lrint(Z));
180 const G4double eps0 = CLHEP::electron_mass_c2/gammaEnergy;
181 // Coulomb correction is always included in the DCS even below 50 MeV (note:
182 // that this DCS is only used to get the integrated x-section)
183 const G4double dmax = gElementData[iz]->fDeltaMaxHigh;
184 const G4double dmin = 4.*eps0*gElementData[iz]->fDeltaFactor;
185 const G4double eps1 = 0.5 - 0.5*std::sqrt(1.-dmin/dmax);
186 const G4double epsMin = std::max(eps0, eps1);
187 const G4double epsMax = 0.5; // DCS is symmetric around eps=0.5
188 // let Et be the total energy transferred to the e- or to the e+
189 // the [Et-min, Et-max] interval will be divided into i=1,2,..,n subintervals
190 // with width of dInterv = (Et-max - Et-min)/n and numerical integration will
191 // be done in each sub-inteval using the xi = (Et - Et_i-min)/dInterv variable
192 // that is in [0,1]. The 8-point GL q. is used for the integration on [0,1].
193 const G4int numSub = 2;
194 const G4double dInterv= (epsMax - epsMin)*gammaEnergy/G4double(numSub);
195 G4double minEti = epsMin*gammaEnergy; // Et-min i.e. Et_0-min
196 for (G4int i = 0; i < numSub; ++i) {
197 for (G4int ngl = 0; ngl < 8; ++ngl) {
198 const G4double Et = (minEti + gXGL[ngl]*dInterv);
199 const G4double xs = isLPM ? ComputeRelDXSectionPerAtom(Et, gammaEnergy, Z)
200 : ComputeDXSectionPerAtom(Et, gammaEnergy, Z);
201 xSection += gWGL[ngl]*xs;
202 }
203 // update minimum Et of the sub-inteval
204 minEti += dInterv;
205 }
206 // apply corrections of variable transformation and half interval integration
207 xSection = std::max(2.*xSection*dInterv, 0.);
208 return xSection;
209}
210
211// DCS WITHOUT LPM SUPPRESSION
212// Computes DCS value for a given target element (Z), initial gamma energy (Eg),
213// total energy transferred to one of the e-/e+ pair(Et) WITHOUT LPM suppression
214// The constant factor 4 \alpha r_0^2 Z (Z +\eta(Z)) is not included here and
215// the returned value will be differential in total energy transfer instead of
216// the eps=Et/Eg. The computed part of the DCS
217// NORMAL CASE: DEFAULT STTING (i.e. fIsUseCompleteScreening = FALSE)
218// ds/deps(Et,Eg,Z) = ds/deps(eps,Z) = (eps^2+(1-eps)^2)*[phi1(d)/4-ln(Z)/3-fc]
219// + 2*eps(1-eps)*[phi2(d)/4-ln(Z)/3-fc]/3 where the universal (in the TF model)
220// screening variable d=d(eps)=136Z^(-1/3)eps0/[eps*(1-eps)] with eps0=mc^2/Eg.
221// COMPLETE SCREENING (when d(eps) approx-equal-to 0) : NEED TO BE SET BY USER
222// ds/deps(Et,Eg,Z) = ds/deps(eps,Z) = (eps^2+(1-eps)^2+eps*(1-eps)/3)*[Lel-fc]
223// -eps(1-eps)/9 where Lel=phi1(0)/4-ln(Z)/3 is the elastic(coherent) radiation
224// logarithm, fc is the Coulomb correction and the relation phi2(0)/4-ln(Z)/3 =
225// phi1(0)/4-1/6-ln(Z)/3 = Lel-1/6 (due to phi2(0)=phi1(0)-2/3) was used.
227 G4double gammaEnergy,
228 G4double Z)
229{
230 G4double xSection = 0.;
231 const G4int iz = std::min(gMaxZet, G4lrint(Z));
232 const G4double eps = pEnergy/gammaEnergy;
233 const G4double epsm = 1.-eps;
234 const G4double dum = eps*epsm;
236 // complete screening:
237 const G4double Lel = gElementData[iz]->fLradEl;
238 const G4double fc = gElementData[iz]->fCoulomb;
239 xSection = (eps*eps + epsm*epsm + 2.*dum/3.)*(Lel-fc) - dum/9.;
240 } else {
241 // normal case:
242 const G4double eps0 = CLHEP::electron_mass_c2/gammaEnergy;
243 const G4double fc = gElementData[iz]->fCoulomb;
244 const G4double lnZ13 = gElementData[iz]->fLogZ13;
245 const G4double delta = gElementData[iz]->fDeltaFactor*eps0/dum;
246 G4double phi1, phi2;
247 ComputePhi12(delta, phi1, phi2);
248 xSection = (eps*eps + epsm*epsm)*(0.25*phi1-lnZ13-fc)
249 + 2.*dum*(0.25*phi2-lnZ13-fc)/3.;
250 }
251 // non-const. part of the DCS differential in total energy transfer not in eps
252 // ds/dEt=ds/deps deps/dEt with deps/dEt=1/Eg
253 return std::max(xSection, 0.0)/gammaEnergy;
254}
255
256// DCS WITH POSSIBLE LPM SUPPRESSION
257// Computes DCS value for a given target element (Z), initial gamma energy (Eg),
258// total energy transferred to one of the e-/e+ pair(Et) WITH LPM suppression.
259// For a given Z, the LPM suppression will depend on the material through the
260// LMP-Energy. This will determine the suppression variable s and the LPM sup-
261// pression functions xi(s), fi(s) and G(s).
262// The constant factor 4 \alpha r_0^2 Z (Z +\eta(Z)) is not included here and
263// the returned value will be differential in total energy transfer instead of
264// the eps=Et/Eg. The computed part of the DCS
265// NORMAL CASE: DEFAULT STTING (i.e. fIsUseCompleteScreening = FALSE)
266// ds/deps(Et,Eg,Z)=ds/deps(eps,Z) = xi(s)*{ (eps^2+(1-eps)^2)*[2fi(s)/3+G(s)/3]
267// *[phi1(d)/4-ln(Z)/3-fc] + 2*eps(1-eps)*G(s)*[phi2(d)/4-ln(Z)/3-fc]/3 } where
268// the universal (in the TF model) screening variable d=d(eps)=136Z^(-1/3)eps0
269// /[eps*(1-eps)] with eps0=mc^2/Eg.
270// COMPLETE SCREENING (when d(eps) approx-equal-to 0) : NEED TO BE SET BY USER
271// ds/deps(Et,Eg,Z) = ds/deps(eps,Z) = xi(s)*{ [Lel-fc]*[ (eps^2+(1-eps)^2+eps
272// *(1-eps)/3)*2fi(s)/3 + G(s)/3] - eps(1-eps)*G(s)/9 }
273// Note, that when the LPM suppression is absent i.e. xi(s)=fi(s)=G(s)=1, both
274// the normal and the complete screening DCS give back the NO-LMP case above.
276 G4double gammaEnergy,
277 G4double Z)
278{
279 G4double xSection = 0.;
280 const G4int iz = std::min(gMaxZet, G4lrint(Z));
281 const G4double eps = pEnergy/gammaEnergy;
282 const G4double epsm = 1.-eps;
283 const G4double dum = eps*epsm;
284 // evaluate LPM suppression functions
285 G4double fXiS, fGS, fPhiS;
286 ComputeLPMfunctions(fXiS, fGS, fPhiS, eps, gammaEnergy, iz);
288 // complete screening:
289 const G4double Lel = gElementData[iz]->fLradEl;
290 const G4double fc = gElementData[iz]->fCoulomb;
291 xSection = (Lel-fc)*((eps*eps+epsm*epsm)*2.*fPhiS + fGS)/3. - dum*fGS/9.;
292 } else {
293 // normal case:
294 const G4double eps0 = CLHEP::electron_mass_c2/gammaEnergy;
295 const G4double fc = gElementData[iz]->fCoulomb;
296 const G4double lnZ13 = gElementData[iz]->fLogZ13;
297 const G4double delta = gElementData[iz]->fDeltaFactor*eps0/dum;
298 G4double phi1, phi2;
299 ComputePhi12(delta, phi1, phi2);
300 xSection = (eps*eps + epsm*epsm)*(2.*fPhiS+fGS)*(0.25*phi1-lnZ13-fc)/3.
301 + 2.*dum*fGS*(0.25*phi2-lnZ13-fc)/3.;
302 }
303 // non-const. part of the DCS differential in total energy transfer not in eps
304 // ds/dEt=ds/deps deps/dEt with deps/dEt=1/Eg
305 return std::max(fXiS*xSection, 0.0)/gammaEnergy;
306}
307
310 G4double gammaEnergy, G4double Z, G4double, G4double, G4double)
311{
312 G4double crossSection = 0.0 ;
313 // check kinematical limit
314 if ( gammaEnergy <= 2.0*electron_mass_c2 ) { return crossSection; }
315 // compute the atomic cross section either by using x-section parametrization
316 // or by numerically integrationg the DCS (with or without LPM)
317 if ( gammaEnergy < fParametrizedXSectionThreshold) {
318 // using the parametrized cross sections (max up to 80 GeV)
319 crossSection = ComputeParametrizedXSectionPerAtom(gammaEnergy, Z);
320 } else {
321 // by numerical integration of the DCS:
322 // Computes the cross section with or without LPM suppression depending on
323 // settings (by default with if the gamma energy is above a given threshold)
324 // and using or not using complete sreening approximation (by default not).
325 // Only the dependent part is computed in the numerical integration of the DCS
326 // i.e. the result must be multiplied here with 4 \alpha r_0^2 Z(Z+\eta(Z))
327 crossSection = ComputeXSectionPerAtom(gammaEnergy, Z);
328 // apply the constant factors:
329 // - eta(Z) is a correction to account interaction in the field of e-
330 // - gXSecFactor = 4 \alpha r_0^2
331 const G4int iz = std::min(gMaxZet, G4lrint(Z));
332 const G4double eta = gElementData[iz]->fEtaValue;
333 crossSection *= gXSecFactor*Z*(Z+eta);
334 }
335 // final protection
336 return std::max(crossSection, 0.);
337}
338
344
345void
346G4PairProductionRelModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
347 const G4MaterialCutsCouple* couple,
348 const G4DynamicParticle* aDynamicGamma,
349 G4double,
350 G4double)
351// The secondaries e+e- energies are sampled using the Bethe - Heitler
352// cross sections with Coulomb correction.
353// A modified version of the random number techniques of Butcher & Messel
354// is used (Nuc Phys 20(1960),15).
355//
356// GEANT4 internal units.
357//
358// Note 1 : Effects due to the breakdown of the Born approximation at
359// low energy are ignored.
360// Note 2 : The differential cross section implicitly takes account of
361// pair creation in both nuclear and atomic electron fields.
362// However triplet prodution is not generated.
363{
364 const G4Material* mat = couple->GetMaterial();
365 const G4double gammaEnergy = aDynamicGamma->GetKineticEnergy();
366 const G4double eps0 = CLHEP::electron_mass_c2/gammaEnergy ;
367 //
368 // check kinematical limit: gamma energy(Eg) must be at least 2 e- rest mass
369 // (but the model should be used at higher energies above 100 MeV)
370 if (eps0 > 0.5) { return; }
371 //
372 // select target atom of the material
373 const G4Element* anElement = SelectTargetAtom(couple, fTheGamma, gammaEnergy,
374 aDynamicGamma->GetLogKineticEnergy());
375 CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
376 //
377 // 'eps' is the total energy transferred to one of the e-/e+ pair in initial
378 // gamma energy units Eg. Since the corresponding DCS is symmetric on eps=0.5,
379 // the kinematical limits for eps0=mc^2/Eg <= eps <= 0.5
380 // 1. 'eps' is sampled uniformly on the [eps0, 0.5] inteval if Eg<Egsmall
381 // 2. otherwise, on the [eps_min, 0.5] interval according to the DCS (case 2.)
382 G4double eps;
383 // case 1.
384 static const G4double Egsmall = 2.*CLHEP::MeV;
385 if (gammaEnergy < Egsmall) {
386 eps = eps0 + (0.5-eps0)*rndmEngine->flat();
387 } else {
388 // case 2.
389 // get the Coulomb factor for the target element (Z) and gamma energy (Eg)
390 // F(Z) = 8*ln(Z)/3 if Eg <= 50 [MeV] => no Coulomb correction
391 // F(Z) = 8*ln(Z)/3 + 8*fc(Z) if Eg > 50 [MeV] => fc(Z) is the Coulomb cor.
392 //
393 // The screening variable 'delta(eps)' = 136*Z^{-1/3}*eps0/[eps(1-eps)]
394 // Due to the Coulomb correction, the DCS can go below zero even at
395 // kinematicaly allowed eps > eps0 values. In order to exclude this eps
396 // range with negative DCS, the minimum eps value will be set to eps_min =
397 // max[eps0, epsp] with epsp is the solution of SF(delta(epsp)) - F(Z)/2 = 0
398 // with SF being the screening function (SF1=SF2 at high value of delta).
399 // The solution is epsp = 0.5 - 0.5*sqrt[ 1 - 4*136*Z^{-1/3}eps0/deltap]
400 // with deltap = Exp[(42.038-F(Z))/8.29]-0.958. So the limits are:
401 // - when eps=eps_max = 0.5 => delta_min = 136*Z^{-1/3}*eps0/4
402 // - epsp = 0.5 - 0.5*sqrt[ 1 - delta_min/deltap]
403 // - and eps_min = max[eps0, epsp]
404 const G4int iZet = std::min(gMaxZet, anElement->GetZasInt());
405 const G4double deltaFactor = gElementData[iZet]->fDeltaFactor*eps0;
406 const G4double deltaMin = 4.*deltaFactor;
407 G4double deltaMax = gElementData[iZet]->fDeltaMaxLow;
408 G4double FZ = 8.*gElementData[iZet]->fLogZ13;
409 if ( gammaEnergy > fCoulombCorrectionThreshold ) { // Eg > 50 MeV ?
410 FZ += 8.*gElementData[iZet]->fCoulomb;
411 deltaMax = gElementData[iZet]->fDeltaMaxHigh;
412 }
413 // compute the limits of eps
414 const G4double epsp = 0.5 - 0.5*std::sqrt(1. - deltaMin/deltaMax) ;
415 const G4double epsMin = std::max(eps0,epsp);
416 const G4double epsRange = 0.5 - epsMin;
417 //
418 // sample the energy rate (eps) of the created electron (or positron)
420 ScreenFunction12(deltaMin, F10, F20);
421 F10 -= FZ;
422 F20 -= FZ;
423 const G4double NormF1 = std::max(F10 * epsRange * epsRange, 0.);
424 const G4double NormF2 = std::max(1.5 * F20 , 0.);
425 const G4double NormCond = NormF1/(NormF1 + NormF2);
426 // check if LPM correction is active
427 const G4bool isLPM = (fIsUseLPMCorrection && gammaEnergy>gEgLPMActivation);
429 // we will need 3 uniform random number for each trial of sampling
430 G4double rndmv[3];
431 G4double greject = 0.;
432 do {
433 rndmEngine->flatArray(3, rndmv);
434 if (NormCond > rndmv[0]) {
435 eps = 0.5 - epsRange * fG4Calc->A13(rndmv[1]);
436 const G4double delta = deltaFactor/(eps*(1.-eps));
437 if (isLPM) {
438 G4double lpmXiS, lpmGS, lpmPhiS, phi1, phi2;
439 ComputePhi12(delta, phi1, phi2);
440 ComputeLPMfunctions(lpmXiS, lpmGS, lpmPhiS, eps, gammaEnergy, iZet);
441 greject = lpmXiS*((2.*lpmPhiS+lpmGS)*phi1-lpmGS*phi2-lpmPhiS*FZ)/F10;
442 } else {
443 greject = (ScreenFunction1(delta)-FZ)/F10;
444 }
445 } else {
446 eps = epsMin + epsRange*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*( (lpmPhiS+0.5*lpmGS)*phi1 + 0.5*lpmGS*phi2
453 -0.5*(lpmGS+lpmPhiS)*FZ )/F20;
454 } else {
455 greject = (ScreenFunction2(delta)-FZ)/F20;
456 }
457 }
458 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
459 } while (greject < rndmv[2]);
460 // end of eps sampling
461 }
462 //
463 // select charges randomly
464 G4double eTotEnergy, pTotEnergy;
465 if (rndmEngine->flat() > 0.5) {
466 eTotEnergy = (1.-eps)*gammaEnergy;
467 pTotEnergy = eps*gammaEnergy;
468 } else {
469 pTotEnergy = (1.-eps)*gammaEnergy;
470 eTotEnergy = eps*gammaEnergy;
471 }
472 //
473 // sample pair kinematics
474 //
475 const G4double eKinEnergy = std::max(0.,eTotEnergy - CLHEP::electron_mass_c2);
476 const G4double pKinEnergy = std::max(0.,pTotEnergy - CLHEP::electron_mass_c2);
477 //
478 G4ThreeVector eDirection, pDirection;
479 //
481 eKinEnergy, pKinEnergy, eDirection, pDirection);
482 // create G4DynamicParticle object for the particle1
483 auto aParticle1 = new G4DynamicParticle(fTheElectron,eDirection,eKinEnergy);
484
485 // create G4DynamicParticle object for the particle2
486 auto aParticle2 = new G4DynamicParticle(fThePositron,pDirection,pKinEnergy);
487 // Fill output vector
488 fvect->push_back(aParticle1);
489 fvect->push_back(aParticle2);
490 // kill incident photon
491 fParticleChange->SetProposedKineticEnergy(0.);
492 fParticleChange->ProposeTrackStatus(fStopAndKill);
493}
494
495// should be called only by the master and at initialisation
496void G4PairProductionRelModel::InitialiseElementData()
497{
498 // create for all elements that are in the detector
499 auto elemTable = G4Element::GetElementTable();
500 for (auto const & elem : *elemTable) {
501 const G4int iz = std::min(gMaxZet, elem->GetZasInt());
502 if (nullptr == gElementData[iz]) { // create it if doesn't exist yet
503 const G4double logZ13 = elem->GetIonisation()->GetlogZ3();
504 const G4double Z13 = elem->GetIonisation()->GetZ3();
505 const G4double fc = elem->GetfCoulomb();
506 const G4double FZLow = 8.*logZ13;
507 const G4double FZHigh = 8.*(logZ13 + fc);
508 G4double Fel;
509 G4double Finel;
510 if (iz<5) { // use data from Dirac-Fock atomic model
511 Fel = gFelLowZet[iz];
512 Finel = gFinelLowZet[iz];
513 } else { // use the results of the Thomas-Fermi-Moliere model
514 Fel = G4Log(184.) - logZ13;
515 Finel = G4Log(1194.) - 2.*logZ13;
516 }
517 auto elD = new ElementData();
518 elD->fLogZ13 = logZ13;
519 elD->fCoulomb = fc;
520 elD->fLradEl = Fel;
521 elD->fDeltaFactor = 136./Z13;
522 elD->fDeltaMaxLow = G4Exp((42.038 - FZLow)/8.29) - 0.958;
523 elD->fDeltaMaxHigh = G4Exp((42.038 - FZHigh)/8.29) - 0.958;
524 elD->fEtaValue = Finel/(Fel-fc);
525 elD->fLPMVarS1Cond = std::sqrt(2.)*Z13*Z13/(184.*184.);
526 elD->fLPMILVarS1Cond = 1./G4Log(elD->fLPMVarS1Cond);
527 gElementData[iz] = elD;
528 }
529 }
530}
531
532void G4PairProductionRelModel::ComputeLPMfunctions(G4double &funcXiS,
533 G4double &funcGS, G4double &funcPhiS, const G4double eps,
534 const G4double egamma, const G4int izet)
535{
536 // 1. y = E_+/E_{\gamma} with E_+ being the total energy transfered
537 // to one of the e-/e+ pair
538 // s' = \sqrt{ \frac{1}{8} \frac{1}{y(1-y)} \frac{E^{KL}_{LPM}}{E_{\gamma}} }
539 const G4double varSprime = std::sqrt(0.125*fLPMEnergy/(eps*egamma*(1.0-eps)));
540 const G4double condition = gElementData[izet]->fLPMVarS1Cond;
541 funcXiS = 2.0;
542 if (varSprime > 1.0) {
543 funcXiS = 1.0;
544 } else if (varSprime > condition) {
545 const G4double dum = gElementData[izet]->fLPMILVarS1Cond;
546 const G4double funcHSprime = G4Log(varSprime)*dum;
547 funcXiS = 1.0 + funcHSprime
548 - 0.08*(1.0-funcHSprime)*funcHSprime*(2.0-funcHSprime)*dum;
549 }
550 // 2. s=\frac{s'}{\sqrt{\xi(s')}}
551 const G4double varShat = varSprime / std::sqrt(funcXiS);
552 G4LPMFunction::GetLPMFunctions(funcGS, funcPhiS, varShat);
553 // MAKE SURE SUPPRESSION IS SMALLER THAN 1: due to Migdal's approximation on xi
554 if (funcXiS * funcPhiS > 1. || varShat > 0.57) {
555 funcXiS = 1. / funcPhiS;
556 }
557}
558
559// Calculates the microscopic cross section in GEANT4 internal units. Same as in
560// G4BetheHeitlerModel and should be used below 80 GeV since it start to deverge
561// from the cross section data above 80-90 GeV:
562// Parametrized formula (L. Urban) is used to estimate the atomic cross sections
563// given numerically in the table of [Hubbell, J. H., Heinz Albert Gimm, and I.
564// Overbo: "Pair, Triplet, and Total Atomic Cross Sections (and Mass Attenuation
565// Coefficients) for 1 MeV‐100 GeV Photons in Elements Z= 1 to 100." Journal of
566// physical and chemical reference data 9.4 (1980): 1023-1148.]
567//
568// The formula gives a good approximation of the data from 1.5 MeV to 100 GeV.
569// below 1.5 MeV: sigma=sigma(1.5MeV)*(GammaEnergy-2electronmass)
570// *(GammaEnergy-2electronmass)
573 G4double Z)
574{
575 G4double xSection = 0.0 ;
576 // short versions
577 static const G4double kMC2 = CLHEP::electron_mass_c2;
578 // zero cross section below the kinematical limit: Eg<2mc^2
579 if (Z < 0.9 || gammaE <= 2.0*kMC2) { return xSection; }
580 //
581 static const G4double gammaEnergyLimit = 1.5*CLHEP::MeV;
582 // set coefficients a, b c
583 static const G4double a0 = 8.7842e+2*CLHEP::microbarn;
584 static const G4double a1 = -1.9625e+3*CLHEP::microbarn;
585 static const G4double a2 = 1.2949e+3*CLHEP::microbarn;
586 static const G4double a3 = -2.0028e+2*CLHEP::microbarn;
587 static const G4double a4 = 1.2575e+1*CLHEP::microbarn;
588 static const G4double a5 = -2.8333e-1*CLHEP::microbarn;
589
590 static const G4double b0 = -1.0342e+1*CLHEP::microbarn;
591 static const G4double b1 = 1.7692e+1*CLHEP::microbarn;
592 static const G4double b2 = -8.2381 *CLHEP::microbarn;
593 static const G4double b3 = 1.3063 *CLHEP::microbarn;
594 static const G4double b4 = -9.0815e-2*CLHEP::microbarn;
595 static const G4double b5 = 2.3586e-3*CLHEP::microbarn;
596
597 static const G4double c0 = -4.5263e+2*CLHEP::microbarn;
598 static const G4double c1 = 1.1161e+3*CLHEP::microbarn;
599 static const G4double c2 = -8.6749e+2*CLHEP::microbarn;
600 static const G4double c3 = 2.1773e+2*CLHEP::microbarn;
601 static const G4double c4 = -2.0467e+1*CLHEP::microbarn;
602 static const G4double c5 = 6.5372e-1*CLHEP::microbarn;
603 // check low energy limit of the approximation (1.5 MeV)
604 G4double gammaEnergyOrg = gammaE;
605 if (gammaE < gammaEnergyLimit) { gammaE = gammaEnergyLimit; }
606 // compute gamma energy variables
607 const G4double x = G4Log(gammaE/kMC2);
608 const G4double x2 = x *x;
609 const G4double x3 = x2*x;
610 const G4double x4 = x3*x;
611 const G4double x5 = x4*x;
612 //
613 const G4double F1 = a0 + a1*x + a2*x2 + a3*x3 + a4*x4 + a5*x5;
614 const G4double F2 = b0 + b1*x + b2*x2 + b3*x3 + b4*x4 + b5*x5;
615 const G4double F3 = c0 + c1*x + c2*x2 + c3*x3 + c4*x4 + c5*x5;
616 // compute the approximated cross section
617 xSection = (Z + 1.)*(F1*Z + F2*Z*Z + F3);
618 // check if we are below the limit of the approximation and apply correction
619 if (gammaEnergyOrg < gammaEnergyLimit) {
620 const G4double dum = (gammaEnergyOrg-2.*kMC2)/(gammaEnergyLimit-2.*kMC2);
621 xSection *= dum*dum;
622 }
623 return xSection;
624}
G4TemplateAutoLock< G4Mutex > G4AutoLock
#define F10
#define F20
G4double condition(const G4ErrorSymMatrix &m)
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition G4Exp.hh:132
const G4double a0
G4double G4Log(G4double x)
Definition G4Log.hh:169
#define elem(i, j)
#define G4MUTEX_INITIALIZER
std::mutex G4Mutex
CLHEP::Hep3Vector G4ThreeVector
@ 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 const G4ElementTable * GetElementTable()
Definition G4Element.cc:401
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
Definition G4Pow.hh:49
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
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67
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 GetLPMFunctions(G4double &lpmFuncG, G4double &lpmFuncPhi, G4double sVar)
int G4lrint(double ad)
Definition templates.hh:134