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

Last change on this file since 1103 was 1058, checked in by garnier, 17 years ago

file release beta

File size: 17.0 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.3 2009/05/15 17:12:33 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-03-beta-cand-01 $
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 theCrossSectionTable(0),
88 nbins(10),
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{
94 fParticleChange = 0;
95 theGamma = G4Gamma::Gamma();
96 thePositron = G4Positron::Positron();
97 theElectron = G4Electron::Electron();
98
99 nist = G4NistManager::Instance();
100
101}
102
103//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
104
105G4PairProductionRelModel::~G4PairProductionRelModel()
106{
107 if(theCrossSectionTable) {
108 theCrossSectionTable->clearAndDestroy();
109 delete theCrossSectionTable;
110 }
111}
112
113//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
114
115void G4PairProductionRelModel::Initialise(const G4ParticleDefinition*,
116 const G4DataVector&)
117{
118 fParticleChange = GetParticleChangeForGamma();
119
120 if(theCrossSectionTable) {
121 theCrossSectionTable->clearAndDestroy();
122 delete theCrossSectionTable;
123 }
124
125 const G4ElementTable* theElementTable = G4Element::GetElementTable();
126 size_t nvect = G4Element::GetNumberOfElements();
127 theCrossSectionTable = new G4PhysicsTable(nvect);
128 G4PhysicsLogVector* ptrVector;
129 G4double emin = LowEnergyLimit();
130 G4double emax = HighEnergyLimit();
131 G4int n = nbins*G4int(log10(emax/emin));
132 G4bool spline = G4LossTableManager::Instance()->SplineFlag();
133 G4double e, value;
134
135 for(size_t j=0; j<nvect ; j++) {
136
137 ptrVector = new G4PhysicsLogVector(emin, emax, n);
138 ptrVector->SetSpline(spline);
139 G4double Z = (*theElementTable)[j]->GetZ();
140 G4VEmModel::SetCurrentElement((*theElementTable)[j]);
141 G4int iz = G4int(Z);
142 indexZ[iz] = j;
143
144 for(G4int i=0; i<nbins; i++) {
145 e = ptrVector->GetLowEdgeEnergy( i ) ;
146 value = ComputeCrossSectionPerAtom(theGamma, e, Z);
147 ptrVector->PutValue( i, value );
148 }
149
150 theCrossSectionTable->insert(ptrVector);
151 }
152}
153
154//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
155/*
156G4double G4PairProductionRelModel::ComputeRelXSectionPerAtom(G4double k, G4double Z)
157{
158
159 G4double cross = 0.0;
160
161
162
163
164}
165*/
166
167//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
168
169G4double G4PairProductionRelModel::ComputeXSectionPerAtom(G4double totalEnergy, G4double Z)
170{
171 G4double cross = 0.0;
172
173 // number of intervals and integration step
174 G4double vcut = electron_mass_c2/totalEnergy ;
175
176 // limits by the screening variable
177 G4double dmax = DeltaMax();
178 G4double dmin = min(DeltaMin(totalEnergy),dmax);
179 G4double vcut1 = 0.5 - 0.5*sqrt(1. - dmin/dmax) ;
180 vcut = max(vcut, vcut1);
181
182
183 G4double vmax = 0.5;
184 G4int n = 1; // needs optimisation
185
186 G4double delta = (vmax - vcut)*totalEnergy/G4double(n);
187
188 G4double e0 = vcut*totalEnergy;
189 G4double xs;
190
191 // simple integration
192 for(G4int l=0; l<n; l++,e0 += delta) {
193 for(G4int i=0; i<8; i++) {
194
195 G4double eg = (e0 + xgi[i]*delta);
196 if (fLPMflag && totalEnergy>100.*GeV)
197 xs = ComputeRelDXSectionPerAtom(eg,totalEnergy,Z);
198 else
199 xs = ComputeDXSectionPerAtom(eg,totalEnergy,Z);
200 cross += wgi[i]*xs;
201
202 }
203 }
204
205 cross *= delta*2.;
206
207 return cross;
208}
209
210//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
211
212G4double
213G4PairProductionRelModel::ComputeDXSectionPerAtom(G4double eplusEnergy,
214 G4double totalEnergy,
215 G4double /*Z*/)
216{
217 // most simple case - complete screening:
218
219 // dsig/dE+ = 4 * alpha * Z**2 * r0**2 / k
220 // * [ (y**2 + (1-y**2) + 2/3*y*(1-y) ) * ( log (183 * Z**-1/3) + 1/9 * y*(1-y) ]
221 // y = E+/k
222 G4double yp=eplusEnergy/totalEnergy;
223 G4double ym=1.-yp;
224
225 G4double cross = 0.;
226 if (use_completescreening)
227 cross = (yp*yp + ym*ym + 2./3.*ym*yp)*(Fel - fCoulomb) + 1./9.*yp*ym;
228 else {
229 G4double delta = 0.25*DeltaMin(totalEnergy)/(yp*ym);
230 cross = (yp*yp + ym*ym)*(0.25*Phi1(delta) - lnZ/3. - fCoulomb)
231 + 2./3.*ym*yp*(0.25*Phi2(delta) - lnZ/3. - fCoulomb);
232 }
233 return cross/totalEnergy;
234
235}
236//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
237
238G4double
239G4PairProductionRelModel::ComputeRelDXSectionPerAtom(G4double eplusEnergy,
240 G4double totalEnergy,
241 G4double /*Z*/)
242{
243 // most simple case - complete screening:
244
245 // dsig/dE+ = 4 * alpha * Z**2 * r0**2 / k
246 // * [ (y**2 + (1-y**2) + 2/3*y*(1-y) ) * ( log (183 * Z**-1/3) + 1/9 * y*(1-y) ]
247 // y = E+/k
248 G4double yp=eplusEnergy/totalEnergy;
249 G4double ym=1.-yp;
250
251 CalcLPMFunctions(totalEnergy,eplusEnergy); // gamma
252
253 G4double cross = 0.;
254 if (use_completescreening)
255 cross = xiLPM*(2./3.*phiLPM*(yp*yp + ym*ym) + gLPM)*(Fel - fCoulomb);
256 else {
257 G4double delta = 0.25*DeltaMin(totalEnergy)/(yp*ym);
258 cross = (1./3.*gLPM + 2./3.*phiLPM)*(yp*yp + ym*ym)
259 *(0.25*Phi1(delta) - lnZ/3. - fCoulomb)
260 + 2./3.*gLPM*ym*yp*(0.25*Phi2(delta) - lnZ/3. - fCoulomb);
261 cross *= xiLPM;
262 }
263 return cross/totalEnergy;
264
265}
266
267//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
268
269void G4PairProductionRelModel::CalcLPMFunctions(G4double k, G4double eplus)
270{
271 // *** calculate lpm variable s & sprime ***
272 // Klein eqs. (78) & (79)
273 G4double sprime = sqrt(0.125*k*lpmEnergy/(eplus*(k-eplus)));
274
275 G4double s1 = preS1*z23;
276 G4double logS1 = 2./3.*lnZ-2.*facFel;
277 G4double logTS1 = logTwo+logS1;
278
279 xiLPM = 2.;
280
281 if (sprime>1)
282 xiLPM = 1.;
283 else if (sprime>sqrt(2.)*s1) {
284 G4double h = log(sprime)/logTS1;
285 xiLPM = 1+h-0.08*(1-h)*(1-sqr(1-h))/logTS1;
286 }
287
288 G4double s = sprime/sqrt(xiLPM);
289// G4cout<<"k="<<k<<" y="<<eplus/k<<G4endl;
290// G4cout<<"s="<<s<<G4endl;
291
292 // *** calculate supression functions phi and G ***
293 // Klein eqs. (77)
294 G4double s2=s*s;
295 G4double s3=s*s2;
296 G4double s4=s2*s2;
297
298 if (s<0.1) {
299 // high suppression limit
300 phiLPM = 6.*s - 18.84955592153876*s2 + 39.47841760435743*s3
301 - 57.69873135166053*s4;
302 gLPM = 37.69911184307752*s2 - 236.8705056261446*s3 + 807.7822389*s4;
303 }
304 else if (s<1.9516) {
305 // intermediate suppression
306 // using eq.77 approxim. valid s<2.
307 phiLPM = 1.-exp(-6.*s*(1.+(3.-pi)*s)
308 +s3/(0.623+0.795*s+0.658*s2));
309 if (s<0.415827397755) {
310 // using eq.77 approxim. valid 0.07<s<2
311 G4double psiLPM = 1-exp(-4*s-8*s2/(1+3.936*s+4.97*s2-0.05*s3+7.50*s4));
312 gLPM = 3*psiLPM-2*phiLPM;
313 }
314 else {
315 // using alternative parametrisiation
316 G4double pre = -0.16072300849123999 + s*3.7550300067531581 + s2*-1.7981383069010097
317 + s3*0.67282686077812381 + s4*-0.1207722909879257;
318 gLPM = tanh(pre);
319 }
320 }
321 else {
322 // low suppression limit valid s>2.
323 phiLPM = 1. - 0.0119048/s4;
324 gLPM = 1. - 0.0230655/s4;
325 }
326
327 // *** make sure suppression is smaller than 1 ***
328 // *** caused by Migdal approximation in xi ***
329 if (xiLPM*phiLPM>1. || s>0.57) xiLPM=1./phiLPM;
330}
331
332
333//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
334
335G4double
336G4PairProductionRelModel::ComputeCrossSectionPerAtom(const G4ParticleDefinition*,
337 G4double gammaEnergy, G4double Z,
338 G4double, G4double, G4double)
339{
340 // static const G4double gammaEnergyLimit = 1.5*MeV;
341 G4double crossSection = 0.0 ;
342 if ( Z < 1. ) return crossSection;
343 if ( gammaEnergy <= 2.0*electron_mass_c2 ) return crossSection;
344
345 SetCurrentElement(Z);
346 // choose calculator according to parameters and switches
347 // in the moment only one calculator:
348 crossSection=ComputeXSectionPerAtom(gammaEnergy,Z);
349
350 G4double xi = Finel/(Fel - fCoulomb); // inelastic contribution
351 crossSection*=4.*Z*(Z+xi)*fine_structure_const*classic_electr_radius*classic_electr_radius;
352
353 return crossSection;
354}
355
356//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
357
358void
359G4PairProductionRelModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
360 const G4MaterialCutsCouple* couple,
361 const G4DynamicParticle* aDynamicGamma,
362 G4double,
363 G4double)
364// The secondaries e+e- energies are sampled using the Bethe - Heitler
365// cross sections with Coulomb correction.
366// A modified version of the random number techniques of Butcher & Messel
367// is used (Nuc Phys 20(1960),15).
368//
369// GEANT4 internal units.
370//
371// Note 1 : Effects due to the breakdown of the Born approximation at
372// low energy are ignored.
373// Note 2 : The differential cross section implicitly takes account of
374// pair creation in both nuclear and atomic electron fields.
375// However triplet prodution is not generated.
376{
377 const G4Material* aMaterial = couple->GetMaterial();
378
379 G4double GammaEnergy = aDynamicGamma->GetKineticEnergy();
380 G4ParticleMomentum GammaDirection = aDynamicGamma->GetMomentumDirection();
381
382 G4double epsil ;
383 G4double epsil0 = electron_mass_c2/GammaEnergy ;
384 if(epsil0 > 1.0) return;
385
386 // do it fast if GammaEnergy < 2. MeV
387 static const G4double Egsmall=2.*MeV;
388
389 // select randomly one element constituing the material
390 const G4Element* anElement = SelectRandomAtom(aMaterial, theGamma, GammaEnergy);
391
392 if (GammaEnergy < Egsmall) {
393
394 epsil = epsil0 + (0.5-epsil0)*G4UniformRand();
395
396 } else {
397 // now comes the case with GammaEnergy >= 2. MeV
398
399 // Extract Coulomb factor for this Element
400 G4double FZ = 8.*(anElement->GetIonisation()->GetlogZ3());
401 if (GammaEnergy > 50.*MeV) FZ += 8.*(anElement->GetfCoulomb());
402
403 // limits of the screening variable
404 G4double screenfac = 136.*epsil0/(anElement->GetIonisation()->GetZ3());
405 G4double screenmax = exp ((42.24 - FZ)/8.368) - 0.952 ;
406 G4double screenmin = min(4.*screenfac,screenmax);
407
408 // limits of the energy sampling
409 G4double epsil1 = 0.5 - 0.5*sqrt(1. - screenmin/screenmax) ;
410 G4double epsilmin = max(epsil0,epsil1) , epsilrange = 0.5 - epsilmin;
411
412 //
413 // sample the energy rate of the created electron (or positron)
414 //
415 //G4double epsil, screenvar, greject ;
416 G4double screenvar, greject ;
417
418 G4double F10 = ScreenFunction1(screenmin) - FZ;
419 G4double F20 = ScreenFunction2(screenmin) - FZ;
420 G4double NormF1 = max(F10*epsilrange*epsilrange,0.);
421 G4double NormF2 = max(1.5*F20,0.);
422
423 do {
424 if ( NormF1/(NormF1+NormF2) > G4UniformRand() ) {
425 epsil = 0.5 - epsilrange*pow(G4UniformRand(), 0.333333);
426 screenvar = screenfac/(epsil*(1-epsil));
427 if (fLPMflag && GammaEnergy>100.*GeV) {
428 CalcLPMFunctions(GammaEnergy,GammaEnergy*epsil);
429 greject = xiLPM*((gLPM+2.*phiLPM)*Phi1(screenvar) - gLPM*Phi2(screenvar) - phiLPM*FZ)/F10;
430 }
431 else {
432 greject = (ScreenFunction1(screenvar) - FZ)/F10;
433 }
434
435 } else {
436 epsil = epsilmin + epsilrange*G4UniformRand();
437 screenvar = screenfac/(epsil*(1-epsil));
438 if (fLPMflag && GammaEnergy>100.*GeV) {
439 CalcLPMFunctions(GammaEnergy,GammaEnergy*epsil);
440 greject = xiLPM*((0.5*gLPM+phiLPM)*Phi1(screenvar) + 0.5*gLPM*Phi2(screenvar) - 0.5*(gLPM+phiLPM)*FZ)/F20;
441 }
442 else {
443 greject = (ScreenFunction2(screenvar) - FZ)/F20;
444 }
445 }
446
447 } while( greject < G4UniformRand() );
448
449 } // end of epsil sampling
450
451 //
452 // fixe charges randomly
453 //
454
455 G4double ElectTotEnergy, PositTotEnergy;
456 if (G4UniformRand() > 0.5) {
457
458 ElectTotEnergy = (1.-epsil)*GammaEnergy;
459 PositTotEnergy = epsil*GammaEnergy;
460
461 } else {
462
463 PositTotEnergy = (1.-epsil)*GammaEnergy;
464 ElectTotEnergy = epsil*GammaEnergy;
465 }
466
467 //
468 // scattered electron (positron) angles. ( Z - axis along the parent photon)
469 //
470 // universal distribution suggested by L. Urban
471 // (Geant3 manual (1993) Phys211),
472 // derived from Tsai distribution (Rev Mod Phys 49,421(1977))
473
474 G4double u;
475 const G4double a1 = 0.625 , a2 = 3.*a1 , d = 27. ;
476
477 if (9./(9.+d) >G4UniformRand()) u= - log(G4UniformRand()*G4UniformRand())/a1;
478 else u= - log(G4UniformRand()*G4UniformRand())/a2;
479
480 G4double TetEl = u*electron_mass_c2/ElectTotEnergy;
481 G4double TetPo = u*electron_mass_c2/PositTotEnergy;
482 G4double Phi = twopi * G4UniformRand();
483 G4double dxEl= sin(TetEl)*cos(Phi),dyEl= sin(TetEl)*sin(Phi),dzEl=cos(TetEl);
484 G4double dxPo=-sin(TetPo)*cos(Phi),dyPo=-sin(TetPo)*sin(Phi),dzPo=cos(TetPo);
485
486 //
487 // kinematic of the created pair
488 //
489 // the electron and positron are assumed to have a symetric
490 // angular distribution with respect to the Z axis along the parent photon.
491
492 G4double ElectKineEnergy = max(0.,ElectTotEnergy - electron_mass_c2);
493
494 G4ThreeVector ElectDirection (dxEl, dyEl, dzEl);
495 ElectDirection.rotateUz(GammaDirection);
496
497 // create G4DynamicParticle object for the particle1
498 G4DynamicParticle* aParticle1= new G4DynamicParticle(
499 theElectron,ElectDirection,ElectKineEnergy);
500
501 // the e+ is always created (even with Ekine=0) for further annihilation.
502
503 G4double PositKineEnergy = max(0.,PositTotEnergy - electron_mass_c2);
504
505 G4ThreeVector PositDirection (dxPo, dyPo, dzPo);
506 PositDirection.rotateUz(GammaDirection);
507
508 // create G4DynamicParticle object for the particle2
509 G4DynamicParticle* aParticle2= new G4DynamicParticle(
510 thePositron,PositDirection,PositKineEnergy);
511
512 // Fill output vector
513 fvect->push_back(aParticle1);
514 fvect->push_back(aParticle2);
515
516 // kill incident photon
517 fParticleChange->SetProposedKineticEnergy(0.);
518 fParticleChange->ProposeTrackStatus(fStopAndKill);
519}
520
521//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
522
523
524void G4PairProductionRelModel::SetupForMaterial(const G4ParticleDefinition*,
525 const G4Material* mat, G4double)
526{
527 lpmEnergy = mat->GetRadlen()*fLPMconstant;
528 // G4cout<<" lpmEnergy="<<lpmEnergy<<G4endl;
529}
530
531//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
Note: See TracBrowser for help on using the repository browser.