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