Geant4.10
 All Data Structures 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 74581 2013-10-15 12:03:25Z 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 
85 static const G4double xsfactor =
87 static const G4double Egsmall=2.*MeV;
88 
89 
91  const G4String& nam)
92  : G4VEmModel(nam),
94  fLPMflag(true),
95  lpmEnergy(0.),
96  use_completescreening(false)
97 {
98  fParticleChange = 0;
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()) { InitialiseElementSelectors(p, cuts); }
121 }
122 
123 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
124 
126  G4VEmModel* masterModel)
127 {
129 }
130 
131 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
132 
134 {
135  G4double cross = 0.0;
136 
137  // number of intervals and integration step
138  G4double vcut = electron_mass_c2/totalEnergy ;
139 
140  // limits by the screening variable
141  G4double dmax = DeltaMax();
142  G4double dmin = min(DeltaMin(totalEnergy),dmax);
143  G4double vcut1 = 0.5 - 0.5*sqrt(1. - dmin/dmax) ;
144  vcut = max(vcut, vcut1);
145 
146 
147  G4double vmax = 0.5;
148  G4int n = 1; // needs optimisation
149 
150  G4double delta = (vmax - vcut)*totalEnergy/G4double(n);
151 
152  G4double e0 = vcut*totalEnergy;
153  G4double xs;
154 
155  // simple integration
156  for(G4int l=0; l<n; l++,e0 += delta) {
157  for(G4int i=0; i<8; i++) {
158 
159  G4double eg = (e0 + xgi[i]*delta);
160  if (fLPMflag && totalEnergy>100.*GeV)
161  xs = ComputeRelDXSectionPerAtom(eg,totalEnergy,Z);
162  else
163  xs = ComputeDXSectionPerAtom(eg,totalEnergy,Z);
164  cross += wgi[i]*xs;
165 
166  }
167  }
168 
169  cross *= delta*2.;
170 
171  return cross;
172 }
173 
174 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
175 
176 G4double
178  G4double totalEnergy,
179  G4double /*Z*/)
180 {
181  // most simple case - complete screening:
182 
183  // dsig/dE+ = 4 * alpha * Z**2 * r0**2 / k
184  // * [ (y**2 + (1-y**2) + 2/3*y*(1-y) ) * ( log (183 * Z**-1/3) + 1/9 * y*(1-y) ]
185  // y = E+/k
186  G4double yp=eplusEnergy/totalEnergy;
187  G4double ym=1.-yp;
188 
189  G4double cross = 0.;
191  cross = (yp*yp + ym*ym + 2./3.*ym*yp)*(Fel - fCoulomb) + 1./9.*yp*ym;
192  else {
193  G4double delta = 0.25*DeltaMin(totalEnergy)/(yp*ym);
194  cross = (yp*yp + ym*ym)*(0.25*Phi1(delta) - lnZ/3. - fCoulomb)
195  + 2./3.*ym*yp*(0.25*Phi2(delta) - lnZ/3. - fCoulomb);
196  }
197  return cross/totalEnergy;
198 
199 }
200 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
201 
202 G4double
204  G4double totalEnergy,
205  G4double /*Z*/)
206 {
207  // most simple case - complete screening:
208 
209  // dsig/dE+ = 4 * alpha * Z**2 * r0**2 / k
210  // * [ (y**2 + (1-y**2) + 2/3*y*(1-y) ) * ( log (183 * Z**-1/3) + 1/9 * y*(1-y) ]
211  // y = E+/k
212  G4double yp=eplusEnergy/totalEnergy;
213  G4double ym=1.-yp;
214 
215  CalcLPMFunctions(totalEnergy,eplusEnergy); // gamma
216 
217  G4double cross = 0.;
219  cross = xiLPM*(2./3.*phiLPM*(yp*yp + ym*ym) + gLPM)*(Fel - fCoulomb);
220  else {
221  G4double delta = 0.25*DeltaMin(totalEnergy)/(yp*ym);
222  cross = (1./3.*gLPM + 2./3.*phiLPM)*(yp*yp + ym*ym)
223  *(0.25*Phi1(delta) - lnZ/3. - fCoulomb)
224  + 2./3.*gLPM*ym*yp*(0.25*Phi2(delta) - lnZ/3. - fCoulomb);
225  cross *= xiLPM;
226  }
227  return cross/totalEnergy;
228 
229 }
230 
231 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
232 
233 void
235 {
236  // *** calculate lpm variable s & sprime ***
237  // Klein eqs. (78) & (79)
238  G4double sprime = sqrt(0.125*k*lpmEnergy/(eplusEnergy*(k-eplusEnergy)));
239 
240  G4double s1 = preS1*z23;
241  G4double logS1 = 2./3.*lnZ-2.*facFel;
242  G4double logTS1 = logTwo+logS1;
243 
244  xiLPM = 2.;
245 
246  if (sprime>1)
247  xiLPM = 1.;
248  else if (sprime>sqrt(2.)*s1) {
249  G4double h = G4Log(sprime)/logTS1;
250  xiLPM = 1+h-0.08*(1-h)*(1-sqr(1-h))/logTS1;
251  }
252 
253  G4double s0 = sprime/sqrt(xiLPM);
254  // G4cout<<"k="<<k<<" y="<<eplusEnergy/k<<G4endl;
255  // G4cout<<"s0="<<s0<<G4endl;
256 
257  // *** calculate supression functions phi and G ***
258  // Klein eqs. (77)
259  G4double s2=s0*s0;
260  G4double s3=s0*s2;
261  G4double s4=s2*s2;
262 
263  if (s0<0.1) {
264  // high suppression limit
265  phiLPM = 6.*s0 - 18.84955592153876*s2 + 39.47841760435743*s3
266  - 57.69873135166053*s4;
267  gLPM = 37.69911184307752*s2 - 236.8705056261446*s3 + 807.7822389*s4;
268  }
269  else if (s0<1.9516) {
270  // intermediate suppression
271  // using eq.77 approxim. valid s0<2.
272  phiLPM = 1.-G4Exp(-6.*s0*(1.+(3.-pi)*s0)
273  +s3/(0.623+0.795*s0+0.658*s2));
274  if (s0<0.415827397755) {
275  // using eq.77 approxim. valid 0.07<s<2
276  G4double psiLPM = 1-G4Exp(-4*s0-8*s2/(1+3.936*s0+4.97*s2-0.05*s3+7.50*s4));
277  gLPM = 3*psiLPM-2*phiLPM;
278  }
279  else {
280  // using alternative parametrisiation
281  G4double pre = -0.16072300849123999 + s0*3.7550300067531581 + s2*-1.7981383069010097
282  + s3*0.67282686077812381 + s4*-0.1207722909879257;
283  gLPM = tanh(pre);
284  }
285  }
286  else {
287  // low suppression limit valid s>2.
288  phiLPM = 1. - 0.0119048/s4;
289  gLPM = 1. - 0.0230655/s4;
290  }
291 
292  // *** make sure suppression is smaller than 1 ***
293  // *** caused by Migdal approximation in xi ***
294  if (xiLPM*phiLPM>1. || s0>0.57) xiLPM=1./phiLPM;
295 }
296 
297 
298 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
299 
300 G4double
302  G4double gammaEnergy, G4double Z,
304 {
305  G4double crossSection = 0.0 ;
306  // if ( Z < 0.9 ) return crossSection;
307  if ( gammaEnergy <= 2.0*electron_mass_c2 ) return crossSection;
308 
310  // choose calculator according to parameters and switches
311  // in the moment only one calculator:
312  crossSection=ComputeXSectionPerAtom(gammaEnergy,Z);
313 
314  G4double xi = Finel/(Fel - fCoulomb); // inelastic contribution
315  crossSection *= xsfactor*Z*(Z+xi);
316 
317  return crossSection;
318 }
319 
320 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
321 
322 void
323 G4PairProductionRelModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
324  const G4MaterialCutsCouple* couple,
325  const G4DynamicParticle* aDynamicGamma,
326  G4double,
327  G4double)
328 // The secondaries e+e- energies are sampled using the Bethe - Heitler
329 // cross sections with Coulomb correction.
330 // A modified version of the random number techniques of Butcher & Messel
331 // is used (Nuc Phys 20(1960),15).
332 //
333 // GEANT4 internal units.
334 //
335 // Note 1 : Effects due to the breakdown of the Born approximation at
336 // low energy are ignored.
337 // Note 2 : The differential cross section implicitly takes account of
338 // pair creation in both nuclear and atomic electron fields.
339 // However triplet prodution is not generated.
340 {
341  const G4Material* aMaterial = couple->GetMaterial();
342 
343  G4double GammaEnergy = aDynamicGamma->GetKineticEnergy();
344  G4ParticleMomentum GammaDirection = aDynamicGamma->GetMomentumDirection();
345 
346  G4double epsil ;
347  G4double epsil0 = electron_mass_c2/GammaEnergy ;
348  if(epsil0 > 1.0) { return; }
349 
350  // do it fast if GammaEnergy < 2. MeV
351  SetupForMaterial(theGamma, aMaterial, GammaEnergy);
352 
353  // select randomly one element constituing the material
354  const G4Element* anElement =
355  SelectRandomAtom(aMaterial, theGamma, GammaEnergy);
356 
357  if (GammaEnergy < Egsmall) {
358 
359  epsil = epsil0 + (0.5-epsil0)*G4UniformRand();
360 
361  } else {
362  // now comes the case with GammaEnergy >= 2. MeV
363 
364  // Extract Coulomb factor for this Element
365  G4double FZ = 8.*(anElement->GetIonisation()->GetlogZ3());
366  if (GammaEnergy > 50.*MeV) { FZ += 8.*(anElement->GetfCoulomb()); }
367 
368  // limits of the screening variable
369  G4double screenfac = 136.*epsil0/(anElement->GetIonisation()->GetZ3());
370  G4double screenmax = G4Exp ((42.24 - FZ)/8.368) - 0.952 ;
371  G4double screenmin = min(4.*screenfac,screenmax);
372 
373  // limits of the energy sampling
374  G4double epsil1 = 0.5 - 0.5*sqrt(1. - screenmin/screenmax) ;
375  G4double epsilmin = max(epsil0,epsil1) , epsilrange = 0.5 - epsilmin;
376 
377  //
378  // sample the energy rate of the created electron (or positron)
379  //
380  //G4double epsil, screenvar, greject ;
381  G4double screenvar, greject ;
382 
383  G4double F10 = ScreenFunction1(screenmin) - FZ;
384  G4double F20 = ScreenFunction2(screenmin) - FZ;
385  G4double NormF1 = max(F10*epsilrange*epsilrange,0.);
386  G4double NormF2 = max(1.5*F20,0.);
387 
388  do {
389  if ( NormF1/(NormF1+NormF2) > G4UniformRand() ) {
390  epsil = 0.5 - epsilrange*nist->GetZ13(G4UniformRand());
391  screenvar = screenfac/(epsil*(1-epsil));
392  if (fLPMflag && GammaEnergy>100.*GeV) {
393  CalcLPMFunctions(GammaEnergy,GammaEnergy*epsil);
394  greject = xiLPM*((gLPM+2.*phiLPM)*Phi1(screenvar) -
395  gLPM*Phi2(screenvar) - phiLPM*FZ)/F10;
396  }
397  else {
398  greject = (ScreenFunction1(screenvar) - FZ)/F10;
399  }
400 
401  } else {
402  epsil = epsilmin + epsilrange*G4UniformRand();
403  screenvar = screenfac/(epsil*(1-epsil));
404  if (fLPMflag && GammaEnergy>100.*GeV) {
405  CalcLPMFunctions(GammaEnergy,GammaEnergy*epsil);
406  greject = xiLPM*((0.5*gLPM+phiLPM)*Phi1(screenvar) +
407  0.5*gLPM*Phi2(screenvar) - 0.5*(gLPM+phiLPM)*FZ)/F20;
408  }
409  else {
410  greject = (ScreenFunction2(screenvar) - FZ)/F20;
411  }
412  }
413 
414  } while( greject < G4UniformRand() );
415 
416  } // end of epsil sampling
417 
418  //
419  // fixe charges randomly
420  //
421 
422  G4double ElectTotEnergy, PositTotEnergy;
423  if (G4UniformRand() > 0.5) {
424 
425  ElectTotEnergy = (1.-epsil)*GammaEnergy;
426  PositTotEnergy = epsil*GammaEnergy;
427 
428  } else {
429 
430  PositTotEnergy = (1.-epsil)*GammaEnergy;
431  ElectTotEnergy = epsil*GammaEnergy;
432  }
433 
434  //
435  // scattered electron (positron) angles. ( Z - axis along the parent photon)
436  //
437  // universal distribution suggested by L. Urban
438  // (Geant3 manual (1993) Phys211),
439  // derived from Tsai distribution (Rev Mod Phys 49,421(1977))
440 
441  G4double u;
442  static const G4double a1 = 0.625;
443  static const G4double a2 = 3.*a1;
444  static const G4double d = 27. ;
445 
446  if (9./(9.+d) >G4UniformRand()) u= - G4Log(G4UniformRand()*G4UniformRand())/a1;
447  else u= - G4Log(G4UniformRand()*G4UniformRand())/a2;
448 
449  G4double TetEl = u*electron_mass_c2/ElectTotEnergy;
450  G4double TetPo = u*electron_mass_c2/PositTotEnergy;
451  G4double Phi = twopi * G4UniformRand();
452  G4double dxEl= sin(TetEl)*cos(Phi),dyEl= sin(TetEl)*sin(Phi),dzEl=cos(TetEl);
453  G4double dxPo=-sin(TetPo)*cos(Phi),dyPo=-sin(TetPo)*sin(Phi),dzPo=cos(TetPo);
454 
455  //
456  // kinematic of the created pair
457  //
458  // the electron and positron are assumed to have a symetric
459  // angular distribution with respect to the Z axis along the parent photon.
460 
461  G4double ElectKineEnergy = max(0.,ElectTotEnergy - electron_mass_c2);
462 
463  G4ThreeVector ElectDirection (dxEl, dyEl, dzEl);
464  ElectDirection.rotateUz(GammaDirection);
465 
466  // create G4DynamicParticle object for the particle1
467  G4DynamicParticle* aParticle1= new G4DynamicParticle(
468  theElectron,ElectDirection,ElectKineEnergy);
469 
470  // the e+ is always created (even with Ekine=0) for further annihilation.
471 
472  G4double PositKineEnergy = max(0.,PositTotEnergy - electron_mass_c2);
473 
474  G4ThreeVector PositDirection (dxPo, dyPo, dzPo);
475  PositDirection.rotateUz(GammaDirection);
476 
477  // create G4DynamicParticle object for the particle2
478  G4DynamicParticle* aParticle2= new G4DynamicParticle(
479  thePositron,PositDirection,PositKineEnergy);
480 
481  // Fill output vector
482  fvect->push_back(aParticle1);
483  fvect->push_back(aParticle2);
484 
485  // kill incident photon
488 }
489 
490 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
491 
492 
494  const G4Material* mat, G4double)
495 {
497  // G4cout<<" lpmEnergy="<<lpmEnergy<<G4endl;
498 }
499 
500 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
G4double Phi1(G4double delta) const
G4double ComputeXSectionPerAtom(G4double totalEnergy, G4double Z)
G4double GetKineticEnergy() const
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:135
const char * p
Definition: xmltok.h:285
G4double Phi2(G4double delta) const
G4PairProductionRelModel(const G4ParticleDefinition *p=0, const G4String &nam="BetheHeitlerLPM")
G4double GetZ13(G4double Z)
G4bool IsMaster() const
Definition: G4VEmModel.hh:676
G4double GetfCoulomb() const
Definition: G4Element.hh:190
static const G4double Finel_light[5]
#define F20
G4ParticleDefinition * thePositron
int G4int
Definition: G4Types.hh:78
static G4NistManager * Instance()
static const G4double xgi[8]
G4double ScreenFunction1(G4double ScreenVariable)
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double)
#define G4UniformRand()
Definition: Randomize.hh:87
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
static const G4double wgi[8]
const G4ThreeVector & GetMomentumDirection() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
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
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:760
G4double GetRadlen() const
Definition: G4Material.hh:218
G4double GetlogZ3() const
G4double G4Log(G4double x)
Definition: G4Log.hh:227
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
static G4Positron * Positron()
Definition: G4Positron.cc:94
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:768
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
G4ParticleChangeForGamma * fParticleChange
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4IonisParamElm * GetIonisation() const
Definition: G4Element.hh:198
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
static G4Electron * Electron()
Definition: G4Electron.cc:94
void SetProposedKineticEnergy(G4double proposedKinEnergy)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4double ComputeRelDXSectionPerAtom(G4double eplusEnergy, G4double totalEnergy, G4double Z)
T sqr(const T &x)
Definition: templates.hh:145
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cut=0., G4double emax=DBL_MAX)
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:510
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:121
const G4Material * GetMaterial() const
G4double ComputeDXSectionPerAtom(G4double eplusEnergy, G4double totalEnergy, G4double Z)