Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4LivermoreGammaConversionModelRC.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: G4LivermoreGammaConversionModelRC.cc 74822 2013-10-22 14:42:13Z gcosmo $
27 //
28 // Author: Francesco Longo & Gerardo Depaola
29 // on base of G4LivermoreGammaConversionModel
30 //
31 // History:
32 // --------
33 // 12 Apr 2009 V Ivanchenko Cleanup initialisation and generation of secondaries:
34 // - apply internal high-energy limit only in constructor
35 // - do not apply low-energy limit (default is 0)
36 // - use CLHEP electron mass for low-enegry limit
37 // - remove MeanFreePath method and table
38 
39 
41 #include "G4PhysicalConstants.hh"
42 #include "G4SystemOfUnits.hh"
43 
44 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
45 
46 using namespace std;
47 
48 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
49 
51  const G4String& nam)
52  :G4VEmModel(nam),fParticleChange(0),smallEnergy(2.*MeV),isInitialised(false),
53  crossSectionHandler(0),meanFreePathTable(0)
54 {
55  lowEnergyLimit = 4.0*electron_mass_c2;
56  highEnergyLimit = 100 * GeV;
57  SetHighEnergyLimit(highEnergyLimit);
58 
59  verboseLevel= 0;
60  // Verbosity scale:
61  // 0 = nothing
62  // 1 = warning for energy non-conservation
63  // 2 = details of energy budget
64  // 3 = calculation of cross sections, file openings, sampling of atoms
65  // 4 = entering in methods
66 
67  if(verboseLevel > 0) {
68  G4cout << "Livermore Gamma conversion is constructed " << G4endl
69  << "Energy range: "
70  << lowEnergyLimit / MeV << " MeV - "
71  << highEnergyLimit / GeV << " GeV"
72  << G4endl;
73  }
74 }
75 
76 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
77 
79 {
80  if (crossSectionHandler) delete crossSectionHandler;
81 }
82 
83 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
84 
85 void
87  const G4DataVector&)
88 {
89  if (verboseLevel > 3)
90  G4cout << "Calling G4LivermoreGammaConversionModelRC::Initialise()" << G4endl;
91 
92  if (crossSectionHandler)
93  {
94  crossSectionHandler->Clear();
95  delete crossSectionHandler;
96  }
97 
98  // Read data tables for all materials
99 
100  crossSectionHandler = new G4CrossSectionHandler();
101  crossSectionHandler->Initialise(0,lowEnergyLimit,100.*GeV,400);
102  G4String crossSectionFile = "pair/pp-cs-";
103  crossSectionHandler->LoadData(crossSectionFile);
104 
105  //
106 
107  if (verboseLevel > 2)
108  G4cout << "Loaded cross section files for Livermore Gamma Conversion model RC" << G4endl;
109 
110  if (verboseLevel > 0) {
111  G4cout << "Livermore Gamma Conversion model is initialized " << G4endl
112  << "Energy range: "
113  << LowEnergyLimit() / MeV << " MeV - "
114  << HighEnergyLimit() / GeV << " GeV"
115  << G4endl;
116  }
117 
118  if(isInitialised) return;
120  isInitialised = true;
121 }
122 
123 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
124 
125 G4double
127  G4double GammaEnergy,
128  G4double Z, G4double,
130 {
131  if (verboseLevel > 3) {
132  G4cout << "Calling ComputeCrossSectionPerAtom() of G4LivermoreGammaConversionModelRC"
133  << G4endl;
134  }
135  if (GammaEnergy < lowEnergyLimit || GammaEnergy > highEnergyLimit) return 0;
136 
137  G4double cs = crossSectionHandler->FindValue(G4int(Z), GammaEnergy);
138  return cs;
139 }
140 
141 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
142 
143 void G4LivermoreGammaConversionModelRC::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
144  const G4MaterialCutsCouple* couple,
145  const G4DynamicParticle* aDynamicGamma,
146  G4double,
147  G4double)
148 {
149 
150 // The energies of the e+ e- secondaries are sampled using the Bethe - Heitler
151 // cross sections with Coulomb correction. A modified version of the random
152 // number techniques of Butcher & Messel is used (Nuc Phys 20(1960),15).
153 
154 // Note 1 : Effects due to the breakdown of the Born approximation at low
155 // energy are ignored.
156 // Note 2 : The differential cross section implicitly takes account of
157 // pair creation in both nuclear and atomic electron fields. However triplet
158 // prodution is not generated.
159 
160  if (verboseLevel > 3)
161  G4cout << "Calling SampleSecondaries() of G4LivermoreGammaConversionModelRC" << G4endl;
162 
163  G4double photonEnergy = aDynamicGamma->GetKineticEnergy();
164  G4ParticleMomentum photonDirection = aDynamicGamma->GetMomentumDirection();
165 
166  G4double epsilon ;
167  G4double epsilon0Local = electron_mass_c2 / photonEnergy ;
168  G4double electronTotEnergy;
169  G4double positronTotEnergy;
170 
171  G4double HardPhotonEnergy = 0.0;
172 
173 
174  // Do it fast if photon energy < 2. MeV
175  if (photonEnergy < smallEnergy )
176  {
177  epsilon = epsilon0Local + (0.5 - epsilon0Local) * G4UniformRand();
178 
179  if (G4int(2*G4UniformRand()))
180  {
181  electronTotEnergy = (1. - epsilon) * photonEnergy;
182  positronTotEnergy = epsilon * photonEnergy;
183  }
184  else
185  {
186  positronTotEnergy = (1. - epsilon) * photonEnergy;
187  electronTotEnergy = epsilon * photonEnergy;
188  }
189  }
190  else
191  {
192  // Select randomly one element in the current material
193  //const G4Element* element = crossSectionHandler->SelectRandomElement(couple,photonEnergy);
194  const G4ParticleDefinition* particle = aDynamicGamma->GetDefinition();
195  const G4Element* element = SelectRandomAtom(couple,particle,photonEnergy);
196  G4cout << "G4LivermoreGammaConversionModelRC::SampleSecondaries" << G4endl;
197 
198  if (element == 0)
199  {
200  G4cout << "G4LivermoreGammaConversionModelRC::SampleSecondaries - element = 0"
201  << G4endl;
202  return;
203  }
204  G4IonisParamElm* ionisation = element->GetIonisation();
205  if (ionisation == 0)
206  {
207  G4cout << "G4LivermoreGammaConversionModelRC::SampleSecondaries - ionisation = 0"
208  << G4endl;
209  return;
210  }
211 
212  // Extract Coulomb factor for this Element
213  G4double fZ = 8. * (ionisation->GetlogZ3());
214  if (photonEnergy > 50. * MeV) fZ += 8. * (element->GetfCoulomb());
215 
216  // Limits of the screening variable
217  G4double screenFactor = 136. * epsilon0Local / (element->GetIonisation()->GetZ3()) ;
218  G4double screenMax = std::exp ((42.24 - fZ)/8.368) - 0.952 ;
219  G4double screenMin = std::min(4.*screenFactor,screenMax) ;
220 
221  // Limits of the energy sampling
222  G4double epsilon1 = 0.5 - 0.5 * std::sqrt(1. - screenMin / screenMax) ;
223  G4double epsilonMin = std::max(epsilon0Local,epsilon1);
224  G4double epsilonRange = 0.5 - epsilonMin ;
225 
226  // Sample the energy rate of the created electron (or positron)
227  G4double screen;
228  G4double gReject ;
229 
230  G4double f10 = ScreenFunction1(screenMin) - fZ;
231  G4double f20 = ScreenFunction2(screenMin) - fZ;
232  G4double normF1 = std::max(f10 * epsilonRange * epsilonRange,0.);
233  G4double normF2 = std::max(1.5 * f20,0.);
234  G4double a=393.3750918, b=115.3070201, c=810.6428451, d=19.96497475, e=1016.874592, f=1.936685510,
235  gLocal=751.2140962, h=0.099751048, i=299.9466339, j=0.002057250, k=49.81034926;
236  G4double aa=-18.6371131, bb=-1729.95248, cc=9450.971186, dd=106336.0145, ee=55143.09287, ff=-117602.840,
237  gg=-721455.467, hh=693957.8635, ii=156266.1085, jj=533209.9347;
238  G4double Rechazo = 0.;
239  G4double logepsMin = log(epsilonMin);
240  G4double NormaRC = a + b*logepsMin + c/logepsMin + d*pow(logepsMin,2.) + e/pow(logepsMin,2.) + f*pow(logepsMin,3.) +
241  gLocal/pow(logepsMin,3.) + h*pow(logepsMin,4.) + i/pow(logepsMin,4.) + j*pow(logepsMin,5.) +
242  k/pow(logepsMin,5.);
243 
244 
245  G4double HardPhotonThreshold = 0.08;
246  G4double r1, r2, r3, beta=0, gbeta, sigt = 582.068, sigh, rejet;
247  // , Pi = 2.*acos(0.);
248  G4double cg = (11./2.)/(exp(-11.*HardPhotonThreshold/2.)-exp(-11./2.));
249 
250 
251  r1 = G4UniformRand();
252  sigh = 1028.58*exp(-HardPhotonThreshold/0.09033) + 136.63; // sigma hard
253  if (r1 > 1.- sigh/sigt) {
254  r2 = G4UniformRand();
255  rejet = 0.;
256  while (r2 > rejet) {
257  r3 = G4UniformRand();
258  beta = (-2./11.)*log(exp(-0.08*11./2.)-r3*11./(2.*cg));
259  gbeta = exp(-11.*beta/2.);
260  rejet = fbeta(beta)/(8000.*gbeta);
261  }
262  HardPhotonEnergy = beta * photonEnergy;
263  }
264  else{
265  HardPhotonEnergy = 0.;
266  }
267 
268  photonEnergy -= HardPhotonEnergy;
269 
270  do {
271  do {
272  if (normF1 / (normF1 + normF2) > G4UniformRand() )
273  {
274  epsilon = 0.5 - epsilonRange * std::pow(G4UniformRand(), 0.3333) ;
275  screen = screenFactor / (epsilon * (1. - epsilon));
276  gReject = (ScreenFunction1(screen) - fZ) / f10 ;
277  }
278  else
279  {
280  epsilon = epsilonMin + epsilonRange * G4UniformRand();
281  screen = screenFactor / (epsilon * (1 - epsilon));
282  gReject = (ScreenFunction2(screen) - fZ) / f20 ;
283  }
284  } while ( gReject < G4UniformRand() );
285 
286  if (G4int(2*G4UniformRand())) epsilon = (1. - epsilon); // Extención de Epsilon hasta 1.
287 
288  G4double logepsilon = log(epsilon);
289  G4double deltaP_R1 = 1. + (a + b*logepsilon + c/logepsilon + d*pow(logepsilon,2.) + e/pow(logepsilon,2.) +
290  f*pow(logepsilon,3.) + gLocal/pow(logepsilon,3.) + h*pow(logepsilon,4.) + i/pow(logepsilon,4.) +
291  j*pow(logepsilon,5.) + k/pow(logepsilon,5.))/100.;
292  G4double deltaP_R2 = 1.+((aa + cc*logepsilon + ee*pow(logepsilon,2.) + gg*pow(logepsilon,3.) + ii*pow(logepsilon,4.))
293  / (1. + bb*logepsilon + dd*pow(logepsilon,2.) + ff*pow(logepsilon,3.) + hh*pow(logepsilon,4.)
294  + jj*pow(logepsilon,5.) ))/100.;
295 
296  if (epsilon <= 0.5)
297  {
298  Rechazo = deltaP_R1/NormaRC;
299  }
300  else
301  {
302  Rechazo = deltaP_R2/NormaRC;
303  }
304  G4cout << Rechazo << " " << NormaRC << " " << epsilon << G4endl;
305  } while (Rechazo < G4UniformRand() );
306 
307  electronTotEnergy = (1. - epsilon) * photonEnergy;
308  positronTotEnergy = epsilon * photonEnergy;
309 
310  } // End of epsilon sampling
311 
312  // Fix charges randomly
313 
314  // Scattered electron (positron) angles. ( Z - axis along the parent photon)
315  // Universal distribution suggested by L. Urban (Geant3 manual (1993) Phys211),
316  // derived from Tsai distribution (Rev. Mod. Phys. 49, 421 (1977)
317 
318  G4double u;
319  const G4double a1 = 0.625;
320  G4double a2 = 3. * a1;
321  // G4double d = 27. ;
322 
323  // if (9. / (9. + d) > G4UniformRand())
324  if (0.25 > G4UniformRand())
325  {
326  u = - std::log(G4UniformRand() * G4UniformRand()) / a1 ;
327  }
328  else
329  {
330  u = - std::log(G4UniformRand() * G4UniformRand()) / a2 ;
331  }
332 
333  G4double thetaEle = u*electron_mass_c2/electronTotEnergy;
334  G4double thetaPos = u*electron_mass_c2/positronTotEnergy;
335  G4double phi = twopi * G4UniformRand();
336 
337  G4double dxEle= std::sin(thetaEle)*std::cos(phi),dyEle= std::sin(thetaEle)*std::sin(phi),dzEle=std::cos(thetaEle);
338  G4double dxPos=-std::sin(thetaPos)*std::cos(phi),dyPos=-std::sin(thetaPos)*std::sin(phi),dzPos=std::cos(thetaPos);
339 
340 
341  // Kinematics of the created pair:
342  // the electron and positron are assumed to have a symetric angular
343  // distribution with respect to the Z axis along the parent photon
344 
345  G4double electronKineEnergy = std::max(0.,electronTotEnergy - electron_mass_c2) ;
346 
347  // SI - The range test has been removed wrt original G4LowEnergyGammaconversion class
348 
349  G4ThreeVector electronDirection (dxEle, dyEle, dzEle);
350  electronDirection.rotateUz(photonDirection);
351 
353  electronDirection,
354  electronKineEnergy);
355 
356  // The e+ is always created (even with kinetic energy = 0) for further annihilation
357  G4double positronKineEnergy = std::max(0.,positronTotEnergy - electron_mass_c2) ;
358 
359  // SI - The range test has been removed wrt original G4LowEnergyGammaconversion class
360 
361  G4ThreeVector positronDirection (dxPos, dyPos, dzPos);
362  positronDirection.rotateUz(photonDirection);
363 
364  // Create G4DynamicParticle object for the particle2
366  positronDirection,
367  positronKineEnergy);
368  // Fill output vector
369 // G4cout << "Cree el e+ " << epsilon << G4endl;
370  fvect->push_back(particle1);
371  fvect->push_back(particle2);
372 
373  if (HardPhotonEnergy > 0.)
374  {
375  G4double thetaHardPhoton = u*electron_mass_c2/HardPhotonEnergy;
376  phi = twopi * G4UniformRand();
377  G4double dxHardP= std::sin(thetaHardPhoton)*std::cos(phi);
378  G4double dyHardP= std::sin(thetaHardPhoton)*std::sin(phi);
379  G4double dzHardP =std::cos(thetaHardPhoton);
380 
381  G4ThreeVector hardPhotonDirection (dxHardP, dyHardP, dzHardP);
382  hardPhotonDirection.rotateUz(photonDirection);
384  hardPhotonDirection,
385  HardPhotonEnergy);
386  fvect->push_back(particle3);
387  }
388 
389 
390  // kill incident photon
393 
394 }
395 
396 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
397 
398 G4double G4LivermoreGammaConversionModelRC::ScreenFunction1(G4double screenVariable)
399 {
400  // Compute the value of the screening function 3*phi1 - phi2
401 
402  G4double value;
403 
404  if (screenVariable > 1.)
405  value = 42.24 - 8.368 * std::log(screenVariable + 0.952);
406  else
407  value = 42.392 - screenVariable * (7.796 - 1.961 * screenVariable);
408 
409  return value;
410 }
411 
412 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
413 
414 G4double G4LivermoreGammaConversionModelRC::ScreenFunction2(G4double screenVariable)
415 {
416  // Compute the value of the screening function 1.5*phi1 - 0.5*phi2
417 
418  G4double value;
419 
420  if (screenVariable > 1.)
421  value = 42.24 - 8.368 * std::log(screenVariable + 0.952);
422  else
423  value = 41.405 - screenVariable * (5.828 - 0.8945 * screenVariable);
424 
425  return value;
426 }
427 
428 
429 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
430 
431 G4double G4LivermoreGammaConversionModelRC::fbeta(G4double x)
432 {
433  // compute the probabililty distribution for hard photon
434  G4double Pi, gamma, eta, d, p1, p2, p3, p4, p5, p6, p7, ffbeta;
435  gamma = (1.-x)*(1.-x)/x;
436  eta = (1.-x)/(1.+x);
437  d = Dilog(1./x)-Dilog(x);
438  Pi = 2.*acos(0.);
439  p1 = -1.*(25528.*pow(gamma,2) + 116044.* gamma +151556.)/105.;
440  p2 = 256.* pow(gamma,3) + 1092.* pow(gamma,2) +1260.*gamma + 420.;
441  p3 = (676.*pow(gamma,3) + 9877.*pow(gamma,2) + 58415.*gamma + 62160.)/105.;
442  p4 = 64.*pow(gamma,3) + 305.*pow(gamma,2) + 475.*gamma + 269. - 276./gamma;
443  p5 = (676.*pow(gamma,3) + 38109.*pow(gamma,2) + 211637.*gamma + 266660. - 53632./gamma)/105.;
444  p6 = 32.*pow(gamma,2) + 416.*gamma + 1310. +1184./gamma;
445  p7 = 128.*pow(gamma,3) + 802.*pow(gamma,2) + 1028.*gamma - 470. - 1184./gamma;
446  ffbeta = (1.-x) * (p1 + p2*Pi*Pi/6. + p3*log(gamma) +
447  p4*pow(log(x),2) + (p5 + p6*log(gamma))*eta*log(x) + p7*d*eta);
448 
449  return ffbeta;
450 }
451 
452 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
453 
454 G4double G4LivermoreGammaConversionModelRC::Dilog(G4double y)
455 {
456  G4double fdilog = 0.0;
457  G4double Pi = 2.*acos(0.); // serve?
458  if (y <= 0.5) {
459  fdilog = pow(Pi,2)/6. + (1.-y)*(log(1-y)-1.)+pow((1.-y),2)*((1./2.)*log(1.-y)-1./4.)
460  +pow((1.-y),3)*((1./3.)*log(1.-y)-1./9.)+pow((1.-y),4)*((1./4.)*log(1.-y)-1./16.);
461  }
462  if (0.5 < y && y < 2.) {
463  fdilog = 1.-y+pow((1.-y),2)/4.+pow((1.-y),3)/9.+pow((1.-y),4)/16.+
464  pow((1.-y),5)/25.+pow((1.-y),6)/36.+pow((1.-y),7)/49.;
465  }
466  if (y >= 2.) {
467  fdilog = -pow(log(y),2)/2. - pow(Pi,2)/6. + (log(y)+1.)/y +
468  (log(y)/2.+1./4.)/pow(y,2) + (log(y)/3.+1./9.)/pow(y,3);
469  }
470  return fdilog;
471 
472 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:599
G4double GetKineticEnergy() const
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:592
G4double GetfCoulomb() const
Definition: G4Element.hh:190
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:683
G4double FindValue(G4int Z, G4double e) const
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
float electron_mass_c2
Definition: hepunit.py:274
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX)
G4LivermoreGammaConversionModelRC(const G4ParticleDefinition *p=0, const G4String &nam="LivermoreConversion")
G4double GetlogZ3() const
void Initialise(G4VDataSetAlgorithm *interpolation=0, G4double minE=250 *CLHEP::eV, G4double maxE=100 *CLHEP::GeV, G4int numberOfBins=200, G4double unitE=CLHEP::MeV, G4double unitData=CLHEP::barn, G4int minZ=1, G4int maxZ=99)
static G4Positron * Positron()
Definition: G4Positron.cc:94
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4IonisParamElm * GetIonisation() const
Definition: G4Element.hh:198
void LoadData(const G4String &dataFile)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
const XML_Char int const XML_Char * value
static G4Electron * Electron()
Definition: G4Electron.cc:94
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
G4double GetZ3() const
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