source: trunk/source/processes/electromagnetic/standard/src/G4PairProductionRelModel.cc

Last change on this file was 1340, checked in by garnier, 14 years ago

update ti head

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