Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 // $Id: G4PairProductionRelModel.cc 96934 2016-05-18 09:10:41Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4PairProductionRelModel
34 //
35 // Author: Andreas Schaelicke
36 //
37 // Creation date: 02.04.2009
38 //
39 // Modifications:
40 //
41 // Class Description:
42 //
43 // Main References:
44 // J.W.Motz et.al., Rev. Mod. Phys. 41 (1969) 581.
45 // S.Klein, Rev. Mod. Phys. 71 (1999) 1501.
46 // T.Stanev et.al., Phys. Rev. D25 (1982) 1291.
47 // M.L.Ter-Mikaelian, High-energy Electromagnetic Processes in Condensed Media,
48 // Wiley, 1972.
49 //
50 // -------------------------------------------------------------------
51 //
52 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
53 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
54 
56 #include "G4PhysicalConstants.hh"
57 #include "G4SystemOfUnits.hh"
58 #include "G4Gamma.hh"
59 #include "G4Electron.hh"
60 #include "G4Positron.hh"
61 #include "G4Log.hh"
62 
64 #include "G4LossTableManager.hh"
65 #include "G4Exp.hh"
66 
67 
68 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
69 
70 using namespace std;
71 
73 const G4double G4PairProductionRelModel::facFinel = G4Log(1194.); // 1440.
74 
75 const G4double G4PairProductionRelModel::preS1 = 1./(184.15*184.15);
77 
78 const G4double G4PairProductionRelModel::xgi[]={ 0.0199, 0.1017, 0.2372, 0.4083,
79  0.5917, 0.7628, 0.8983, 0.9801 };
80 const G4double G4PairProductionRelModel::wgi[]={ 0.0506, 0.1112, 0.1569, 0.1813,
81  0.1813, 0.1569, 0.1112, 0.0506 };
82 const G4double G4PairProductionRelModel::Fel_light[] = {0., 5.31 , 4.79 , 4.74 , 4.71};
83 const G4double G4PairProductionRelModel::Finel_light[] = {0., 6.144 , 5.621 , 5.805 , 5.924};
84 
88 
89 
91  const G4String& nam)
92  : G4VEmModel(nam),
94  fLPMflag(true),
95  lpmEnergy(0.),
96  use_completescreening(false)
97 {
98  fParticleChange = nullptr;
102 
104 
105  currentZ = z13 = z23 = lnZ = Fel = Finel = fCoulomb = phiLPM = gLPM = xiLPM = 0;
106 
107 }
108 
109 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
110 
112 {}
113 
114 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
115 
117  const G4DataVector& cuts)
118 {
120  if(IsMaster() && LowEnergyLimit() < HighEnergyLimit()) {
121  InitialiseElementSelectors(p, cuts);
122  }
123 }
124 
125 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
126 
128  G4VEmModel* masterModel)
129 {
130  if(LowEnergyLimit() < HighEnergyLimit()) {
132  }
133 }
134 
135 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
136 
138 {
139  G4double cross = 0.0;
140 
141  // number of intervals and integration step
142  G4double vcut = electron_mass_c2/totalEnergy ;
143 
144  // limits by the screening variable
145  G4double dmax = DeltaMax();
146  G4double dmin = min(DeltaMin(totalEnergy),dmax);
147  G4double vcut1 = 0.5 - 0.5*sqrt(1. - dmin/dmax) ;
148  vcut = max(vcut, vcut1);
149 
150 
151  G4double vmax = 0.5;
152  G4int n = 1; // needs optimisation
153 
154  G4double delta = (vmax - vcut)*totalEnergy/G4double(n);
155 
156  G4double e0 = vcut*totalEnergy;
157  G4double xs;
158 
159  // simple integration
160  for(G4int l=0; l<n; l++,e0 += delta) {
161  for(G4int i=0; i<8; i++) {
162 
163  G4double eg = (e0 + xgi[i]*delta);
164  if (fLPMflag && totalEnergy>100.*GeV)
165  xs = ComputeRelDXSectionPerAtom(eg,totalEnergy,Z);
166  else
167  xs = ComputeDXSectionPerAtom(eg,totalEnergy,Z);
168  cross += wgi[i]*xs;
169 
170  }
171  }
172 
173  cross *= delta*2.;
174 
175  return cross;
176 }
177 
178 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
179 
180 G4double
182  G4double totalEnergy,
183  G4double /*Z*/)
184 {
185  // most simple case - complete screening:
186 
187  // dsig/dE+ = 4 * alpha * Z**2 * r0**2 / k
188  // * [ (y**2 + (1-y**2) + 2/3*y*(1-y) ) * ( log (183 * Z**-1/3) + 1/9 * y*(1-y) ]
189  // y = E+/k
190  G4double yp=eplusEnergy/totalEnergy;
191  G4double ym=1.-yp;
192 
193  G4double cross = 0.;
195  cross = (yp*yp + ym*ym + 2./3.*ym*yp)*(Fel - fCoulomb) + yp*ym/9.;
196  else {
197  G4double delta = 0.25*DeltaMin(totalEnergy)/(yp*ym);
198  cross = (yp*yp + ym*ym)*(0.25*Phi1(delta) - lnZ/3. - fCoulomb)
199  + 2./3.*ym*yp*(0.25*Phi2(delta) - lnZ/3. - fCoulomb);
200  }
201  return cross/totalEnergy;
202 
203 }
204 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
205 
206 G4double
208  G4double totalEnergy,
209  G4double /*Z*/)
210 {
211  // most simple case - complete screening:
212 
213  // dsig/dE+ = 4 * alpha * Z**2 * r0**2 / k
214  // * [ (y**2 + (1-y**2) + 2/3*y*(1-y) ) * ( log (183 * Z**-1/3) + 1/9 * y*(1-y) ]
215  // y = E+/k
216  G4double yp=eplusEnergy/totalEnergy;
217  G4double ym=1.-yp;
218 
219  CalcLPMFunctions(totalEnergy,eplusEnergy); // gamma
220 
221  G4double cross = 0.;
223  cross = xiLPM*(2./3.*phiLPM*(yp*yp + ym*ym) + gLPM)*(Fel - fCoulomb);
224  else {
225  G4double delta = 0.25*DeltaMin(totalEnergy)/(yp*ym);
226  cross = (1./3.*gLPM + 2./3.*phiLPM)*(yp*yp + ym*ym)
227  *(0.25*Phi1(delta) - lnZ/3. - fCoulomb)
228  + 2./3.*gLPM*ym*yp*(0.25*Phi2(delta) - lnZ/3. - fCoulomb);
229  cross *= xiLPM;
230  }
231  return cross/totalEnergy;
232 
233 }
234 
235 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
236 
237 void
239 {
240  // *** calculate lpm variable s & sprime ***
241  // Klein eqs. (78) & (79)
242  G4double sprime = sqrt(0.125*k*lpmEnergy/(eplusEnergy*(k-eplusEnergy)));
243 
244  G4double s1 = preS1*z23;
245  G4double logS1 = 2./3.*lnZ-2.*facFel;
246  G4double logTS1 = logTwo+logS1;
247 
248  xiLPM = 2.;
249 
250  if (sprime>1)
251  xiLPM = 1.;
252  else if (sprime>sqrt(2.)*s1) {
253  G4double h = G4Log(sprime)/logTS1;
254  xiLPM = 1+h-0.08*(1-h)*(1-sqr(1-h))/logTS1;
255  }
256 
257  G4double s0 = sprime/sqrt(xiLPM);
258  // G4cout<<"k="<<k<<" y="<<eplusEnergy/k<<G4endl;
259  // G4cout<<"s0="<<s0<<G4endl;
260 
261  // *** calculate supression functions phi and G ***
262  // Klein eqs. (77)
263  G4double s2=s0*s0;
264  G4double s3=s0*s2;
265  G4double s4=s2*s2;
266 
267  if (s0<0.1) {
268  // high suppression limit
269  phiLPM = 6.*s0 - 18.84955592153876*s2 + 39.47841760435743*s3
270  - 57.69873135166053*s4;
271  gLPM = 37.69911184307752*s2 - 236.8705056261446*s3 + 807.7822389*s4;
272  }
273  else if (s0<1.9516) {
274  // intermediate suppression
275  // using eq.77 approxim. valid s0<2.
276  phiLPM = 1.-G4Exp(-6.*s0*(1.+(3.-pi)*s0)
277  +s3/(0.623+0.795*s0+0.658*s2));
278  if (s0<0.415827397755) {
279  // using eq.77 approxim. valid 0.07<s<2
280  G4double psiLPM = 1-G4Exp(-4*s0-8*s2/(1+3.936*s0+4.97*s2-0.05*s3+7.50*s4));
281  gLPM = 3*psiLPM-2*phiLPM;
282  }
283  else {
284  // using alternative parametrisiation
285  G4double pre = -0.16072300849123999 + s0*3.7550300067531581 + s2*-1.7981383069010097
286  + s3*0.67282686077812381 + s4*-0.1207722909879257;
287  gLPM = tanh(pre);
288  }
289  }
290  else {
291  // low suppression limit valid s>2.
292  phiLPM = 1. - 0.0119048/s4;
293  gLPM = 1. - 0.0230655/s4;
294  }
295 
296  // *** make sure suppression is smaller than 1 ***
297  // *** caused by Migdal approximation in xi ***
298  if (xiLPM*phiLPM>1. || s0>0.57) xiLPM=1./phiLPM;
299 }
300 
301 
302 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
303 
304 G4double
306  G4double gammaEnergy, G4double Z, G4double, G4double, G4double)
307 {
308  G4double crossSection = 0.0 ;
309  if ( gammaEnergy <= 2.0*electron_mass_c2 ) return crossSection;
310 
312  // choose calculator according to parameters and switches
313  // in the moment only one calculator:
314  crossSection=ComputeXSectionPerAtom(gammaEnergy,Z);
315 
316  G4double xi = Finel/(Fel - fCoulomb); // inelastic contribution
317  crossSection *= xsfactor*Z*(Z+xi);
318 
319  return crossSection;
320 }
321 
322 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
323 
324 void
325 G4PairProductionRelModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
326  const G4MaterialCutsCouple* couple,
327  const G4DynamicParticle* aDynamicGamma,
328  G4double,
329  G4double)
330 // The secondaries e+e- energies are sampled using the Bethe - Heitler
331 // cross sections with Coulomb correction.
332 // A modified version of the random number techniques of Butcher & Messel
333 // is used (Nuc Phys 20(1960),15).
334 //
335 // GEANT4 internal units.
336 //
337 // Note 1 : Effects due to the breakdown of the Born approximation at
338 // low energy are ignored.
339 // Note 2 : The differential cross section implicitly takes account of
340 // pair creation in both nuclear and atomic electron fields.
341 // However triplet prodution is not generated.
342 {
343  const G4Material* aMaterial = couple->GetMaterial();
344 
345  G4double GammaEnergy = aDynamicGamma->GetKineticEnergy();
346  G4ParticleMomentum GammaDirection = aDynamicGamma->GetMomentumDirection();
347 
348  G4double epsil ;
349  G4double epsil0 = electron_mass_c2/GammaEnergy ;
350  if(epsil0 > 1.0) { return; }
351 
352  // do it fast if GammaEnergy < 2. MeV
353  SetupForMaterial(theGamma, aMaterial, GammaEnergy);
354 
355  // select randomly one element constituing the material
356  const G4Element* anElement =
357  SelectRandomAtom(aMaterial, theGamma, GammaEnergy);
358 
359  CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
360 
361  if (GammaEnergy < Egsmall) {
362 
363  epsil = epsil0 + (0.5-epsil0)*rndmEngine->flat();
364 
365  } else {
366  // now comes the case with GammaEnergy >= 2. MeV
367 
368  // Extract Coulomb factor for this Element
369  G4double FZ = 8.*(anElement->GetIonisation()->GetlogZ3());
370  if (GammaEnergy > 50.*MeV) { FZ += 8.*(anElement->GetfCoulomb()); }
371 
372  // limits of the screening variable
373  G4double screenfac = 136.*epsil0/(anElement->GetIonisation()->GetZ3());
374  G4double screenmax = G4Exp ((42.24 - FZ)/8.368) - 0.952 ;
375  G4double screenmin = min(4.*screenfac,screenmax);
376 
377  // limits of the energy sampling
378  G4double epsil1 = 0.5 - 0.5*sqrt(1. - screenmin/screenmax) ;
379  G4double epsilmin = max(epsil0,epsil1) , epsilrange = 0.5 - epsilmin;
380 
381  //
382  // sample the energy rate of the created electron (or positron)
383  //
384  //G4double epsil, screenvar, greject ;
385  G4double screenvar, greject ;
386 
387  G4double F10 = ScreenFunction1(screenmin) - FZ;
388  G4double F20 = ScreenFunction2(screenmin) - FZ;
389  G4double NormF1 = max(F10*epsilrange*epsilrange,0.);
390  G4double NormF2 = max(1.5*F20,0.);
391 
392  do {
393  if ( NormF1/(NormF1+NormF2) > rndmEngine->flat() ) {
394  epsil = 0.5 - epsilrange*nist->GetZ13(rndmEngine->flat());
395  screenvar = screenfac/(epsil*(1-epsil));
396  if (fLPMflag && GammaEnergy>100.*GeV) {
397  CalcLPMFunctions(GammaEnergy,GammaEnergy*epsil);
398  greject = xiLPM*((gLPM+2.*phiLPM)*Phi1(screenvar) -
399  gLPM*Phi2(screenvar) - phiLPM*FZ)/F10;
400  }
401  else {
402  greject = (ScreenFunction1(screenvar) - FZ)/F10;
403  }
404 
405  } else {
406  epsil = epsilmin + epsilrange*rndmEngine->flat();
407  screenvar = screenfac/(epsil*(1-epsil));
408  if (fLPMflag && GammaEnergy>100.*GeV) {
409  CalcLPMFunctions(GammaEnergy,GammaEnergy*epsil);
410  greject = xiLPM*((0.5*gLPM+phiLPM)*Phi1(screenvar) +
411  0.5*gLPM*Phi2(screenvar) - 0.5*(gLPM+phiLPM)*FZ)/F20;
412  }
413  else {
414  greject = (ScreenFunction2(screenvar) - FZ)/F20;
415  }
416  }
417 
418  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
419  } while( greject < rndmEngine->flat());
420 
421  } // end of epsil sampling
422 
423  //
424  // fixe charges randomly
425  //
426 
427  G4double ElectTotEnergy, PositTotEnergy;
428  if (rndmEngine->flat() > 0.5) {
429 
430  ElectTotEnergy = (1.-epsil)*GammaEnergy;
431  PositTotEnergy = epsil*GammaEnergy;
432 
433  } else {
434 
435  PositTotEnergy = (1.-epsil)*GammaEnergy;
436  ElectTotEnergy = epsil*GammaEnergy;
437  }
438 
439  //
440  // scattered electron (positron) angles. ( Z - axis along the parent photon)
441  //
442  // universal distribution suggested by L. Urban
443  // (Geant3 manual (1993) Phys211),
444  // derived from Tsai distribution (Rev Mod Phys 49,421(1977))
445 
446  G4double u = -G4Log(rndmEngine->flat()*rndmEngine->flat());
447 
448  if (9. > 36.*rndmEngine->flat()) { u *= 1.6; }
449  else { u *= 0.53333; }
450 
451  G4double TetEl = u*electron_mass_c2/ElectTotEnergy;
452  G4double TetPo = u*electron_mass_c2/PositTotEnergy;
453  G4double Phi = twopi * rndmEngine->flat();
454  G4double dxEl= sin(TetEl)*cos(Phi),dyEl= sin(TetEl)*sin(Phi),dzEl=cos(TetEl);
455  G4double dxPo=-sin(TetPo)*cos(Phi),dyPo=-sin(TetPo)*sin(Phi),dzPo=cos(TetPo);
456 
457  //
458  // kinematic of the created pair
459  //
460  // the electron and positron are assumed to have a symetric
461  // angular distribution with respect to the Z axis along the parent photon.
462 
463  G4double ElectKineEnergy = max(0.,ElectTotEnergy - electron_mass_c2);
464 
465  G4ThreeVector ElectDirection (dxEl, dyEl, dzEl);
466  ElectDirection.rotateUz(GammaDirection);
467 
468  // create G4DynamicParticle object for the particle1
469  G4DynamicParticle* aParticle1= new G4DynamicParticle(
470  theElectron,ElectDirection,ElectKineEnergy);
471 
472  // the e+ is always created (even with Ekine=0) for further annihilation.
473 
474  G4double PositKineEnergy = max(0.,PositTotEnergy - electron_mass_c2);
475 
476  G4ThreeVector PositDirection (dxPo, dyPo, dzPo);
477  PositDirection.rotateUz(GammaDirection);
478 
479  // create G4DynamicParticle object for the particle2
480  G4DynamicParticle* aParticle2= new G4DynamicParticle(
481  thePositron,PositDirection,PositKineEnergy);
482 
483  // Fill output vector
484  fvect->push_back(aParticle1);
485  fvect->push_back(aParticle2);
486 
487  // kill incident photon
490 }
491 
492 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
493 
494 
496  const G4Material* mat, G4double)
497 {
499  // G4cout<<" lpmEnergy="<<lpmEnergy<<G4endl;
500 }
501 
502 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double) override
G4double Phi1(G4double delta) const
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:643
G4double ComputeXSectionPerAtom(G4double totalEnergy, G4double Z)
G4double GetKineticEnergy() const
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:146
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:636
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
const char * p
Definition: xmltok.h:285
G4double Phi2(G4double delta) const
G4PairProductionRelModel(const G4ParticleDefinition *p=0, const G4String &nam="BetheHeitlerLPM")
G4bool IsMaster() const
Definition: G4VEmModel.hh:725
G4double GetfCoulomb() const
Definition: G4Element.hh:191
static const G4double Finel_light[5]
virtual double flat()=0
#define F20
G4ParticleDefinition * thePositron
int G4int
Definition: G4Types.hh:78
static G4NistManager * Instance()
static constexpr double twopi
Definition: G4SIunits.hh:76
static const G4double xgi[8]
G4double ScreenFunction1(G4double ScreenVariable)
G4double GetZ13(G4double Z) const
static const G4double wgi[8]
const G4ThreeVector & GetMomentumDirection() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4double DeltaMin(G4double) const
void CalcLPMFunctions(G4double k, G4double eplus)
float electron_mass_c2
Definition: hepunit.py:274
G4double ScreenFunction2(G4double ScreenVariable)
#define F10
const G4int n
double flat()
Definition: G4AblaRandom.cc:47
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cut=0., G4double emax=DBL_MAX) override
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:809
G4double GetRadlen() const
Definition: G4Material.hh:220
G4double GetlogZ3() const
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
static G4Positron * Positron()
Definition: G4Positron.cc:94
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:817
G4ParticleChangeForGamma * fParticleChange
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4IonisParamElm * GetIonisation() const
Definition: G4Element.hh:199
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
static constexpr double GeV
Definition: G4SIunits.hh:217
static G4Electron * Electron()
Definition: G4Electron.cc:94
void SetProposedKineticEnergy(G4double proposedKinEnergy)
static constexpr double MeV
Definition: G4SIunits.hh:214
static constexpr double pi
Definition: G4SIunits.hh:75
G4double ComputeRelDXSectionPerAtom(G4double eplusEnergy, G4double totalEnergy, G4double Z)
T sqr(const T &x)
Definition: templates.hh:145
G4ParticleDefinition * theElectron
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
G4double GetZ3() const
static const G4double Fel_light[5]
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:544
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:132
const G4Material * GetMaterial() const
G4double ComputeDXSectionPerAtom(G4double eplusEnergy, G4double totalEnergy, G4double Z)