Geant4  9.6.p02
 All Classes 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$
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 = 2.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,
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 
172  // Do it fast if photon energy < 2. MeV
173  if (photonEnergy < smallEnergy )
174  {
175  epsilon = epsilon0Local + (0.5 - epsilon0Local) * G4UniformRand();
176 
177  if (G4int(2*G4UniformRand()))
178  {
179  electronTotEnergy = (1. - epsilon) * photonEnergy;
180  positronTotEnergy = epsilon * photonEnergy;
181  }
182  else
183  {
184  positronTotEnergy = (1. - epsilon) * photonEnergy;
185  electronTotEnergy = epsilon * photonEnergy;
186  }
187  }
188  else
189  {
190  // Select randomly one element in the current material
191  //const G4Element* element = crossSectionHandler->SelectRandomElement(couple,photonEnergy);
192  const G4ParticleDefinition* particle = aDynamicGamma->GetDefinition();
193  const G4Element* element = SelectRandomAtom(couple,particle,photonEnergy);
194  G4cout << "G4LivermoreGammaConversionModelRC::SampleSecondaries" << G4endl;
195 
196  if (element == 0)
197  {
198  G4cout << "G4LivermoreGammaConversionModelRC::SampleSecondaries - element = 0"
199  << G4endl;
200  return;
201  }
202  G4IonisParamElm* ionisation = element->GetIonisation();
203  if (ionisation == 0)
204  {
205  G4cout << "G4LivermoreGammaConversionModelRC::SampleSecondaries - ionisation = 0"
206  << G4endl;
207  return;
208  }
209 
210  // Extract Coulomb factor for this Element
211  G4double fZ = 8. * (ionisation->GetlogZ3());
212  if (photonEnergy > 50. * MeV) fZ += 8. * (element->GetfCoulomb());
213 
214  // Limits of the screening variable
215  G4double screenFactor = 136. * epsilon0Local / (element->GetIonisation()->GetZ3()) ;
216  G4double screenMax = std::exp ((42.24 - fZ)/8.368) - 0.952 ;
217  G4double screenMin = std::min(4.*screenFactor,screenMax) ;
218 
219  // Limits of the energy sampling
220  G4double epsilon1 = 0.5 - 0.5 * std::sqrt(1. - screenMin / screenMax) ;
221  G4double epsilonMin = std::max(epsilon0Local,epsilon1);
222  G4double epsilonRange = 0.5 - epsilonMin ;
223 
224  // Sample the energy rate of the created electron (or positron)
225  G4double screen;
226  G4double gReject ;
227 
228  G4double f10 = ScreenFunction1(screenMin) - fZ;
229  G4double f20 = ScreenFunction2(screenMin) - fZ;
230  G4double normF1 = std::max(f10 * epsilonRange * epsilonRange,0.);
231  G4double normF2 = std::max(1.5 * f20,0.);
232  G4double a=393.3750918, b=115.3070201, c=810.6428451, d=19.96497475, e=1016.874592, f=1.936685510,
233  gLocal=751.2140962, h=0.099751048, i=299.9466339, j=0.002057250, k=49.81034926;
234  G4double aa=-18.6371131, bb=-1729.95248, cc=9450.971186, dd=106336.0145, ee=55143.09287, ff=-117602.840,
235  gg=-721455.467, hh=693957.8635, ii=156266.1085, jj=533209.9347;
236  G4double Rechazo = 0.;
237  G4double logepsMin = log(epsilonMin);
238  G4double NormaRC = a + b*logepsMin + c/logepsMin + d*pow(logepsMin,2.) + e/pow(logepsMin,2.) + f*pow(logepsMin,3.) +
239  gLocal/pow(logepsMin,3.) + h*pow(logepsMin,4.) + i/pow(logepsMin,4.) + j*pow(logepsMin,5.) +
240  k/pow(logepsMin,5.);
241 
242  do {
243  do {
244  if (normF1 / (normF1 + normF2) > G4UniformRand() )
245  {
246  epsilon = 0.5 - epsilonRange * std::pow(G4UniformRand(), 0.3333) ;
247  screen = screenFactor / (epsilon * (1. - epsilon));
248  gReject = (ScreenFunction1(screen) - fZ) / f10 ;
249  }
250  else
251  {
252  epsilon = epsilonMin + epsilonRange * G4UniformRand();
253  screen = screenFactor / (epsilon * (1 - epsilon));
254  gReject = (ScreenFunction2(screen) - fZ) / f20 ;
255  }
256  } while ( gReject < G4UniformRand() );
257 
258  if (G4int(2*G4UniformRand())) epsilon = (1. - epsilon); // Extención de Epsilon hasta 1.
259 
260  G4double logepsilon = log(epsilon);
261  G4double deltaP_R1 = 1. + (a + b*logepsilon + c/logepsilon + d*pow(logepsilon,2.) + e/pow(logepsilon,2.) +
262  f*pow(logepsilon,3.) + gLocal/pow(logepsilon,3.) + h*pow(logepsilon,4.) + i/pow(logepsilon,4.) +
263  j*pow(logepsilon,5.) + k/pow(logepsilon,5.))/100.;
264  G4double deltaP_R2 = 1.+((aa + cc*logepsilon + ee*pow(logepsilon,2.) + gg*pow(logepsilon,3.) + ii*pow(logepsilon,4.))
265  / (1. + bb*logepsilon + dd*pow(logepsilon,2.) + ff*pow(logepsilon,3.) + hh*pow(logepsilon,4.)
266  + jj*pow(logepsilon,5.) ))/100.;
267 
268  if (epsilon <= 0.5)
269  {
270  Rechazo = deltaP_R1/NormaRC;
271  }
272  else
273  {
274  Rechazo = deltaP_R2/NormaRC;
275  }
276  G4cout << Rechazo << " " << NormaRC << " " << epsilon << G4endl;
277  } while (Rechazo < G4UniformRand() );
278 
279  electronTotEnergy = (1. - epsilon) * photonEnergy;
280  positronTotEnergy = epsilon * photonEnergy;
281 
282  } // End of epsilon sampling
283 
284  // Fix charges randomly
285 
286  // Scattered electron (positron) angles. ( Z - axis along the parent photon)
287  // Universal distribution suggested by L. Urban (Geant3 manual (1993) Phys211),
288  // derived from Tsai distribution (Rev. Mod. Phys. 49, 421 (1977)
289 
290  G4double u;
291  const G4double a1 = 0.625;
292  G4double a2 = 3. * a1;
293  // G4double d = 27. ;
294 
295  // if (9. / (9. + d) > G4UniformRand())
296  if (0.25 > G4UniformRand())
297  {
298  u = - std::log(G4UniformRand() * G4UniformRand()) / a1 ;
299  }
300  else
301  {
302  u = - std::log(G4UniformRand() * G4UniformRand()) / a2 ;
303  }
304 
305  G4double thetaEle = u*electron_mass_c2/electronTotEnergy;
306  G4double thetaPos = u*electron_mass_c2/positronTotEnergy;
307  G4double phi = twopi * G4UniformRand();
308 
309  G4double dxEle= std::sin(thetaEle)*std::cos(phi),dyEle= std::sin(thetaEle)*std::sin(phi),dzEle=std::cos(thetaEle);
310  G4double dxPos=-std::sin(thetaPos)*std::cos(phi),dyPos=-std::sin(thetaPos)*std::sin(phi),dzPos=std::cos(thetaPos);
311 
312 
313  // Kinematics of the created pair:
314  // the electron and positron are assumed to have a symetric angular
315  // distribution with respect to the Z axis along the parent photon
316 
317  G4double electronKineEnergy = std::max(0.,electronTotEnergy - electron_mass_c2) ;
318 
319  // SI - The range test has been removed wrt original G4LowEnergyGammaconversion class
320 
321  G4ThreeVector electronDirection (dxEle, dyEle, dzEle);
322  electronDirection.rotateUz(photonDirection);
323 
325  electronDirection,
326  electronKineEnergy);
327 
328  // The e+ is always created (even with kinetic energy = 0) for further annihilation
329  G4double positronKineEnergy = std::max(0.,positronTotEnergy - electron_mass_c2) ;
330 
331  // SI - The range test has been removed wrt original G4LowEnergyGammaconversion class
332 
333  G4ThreeVector positronDirection (dxPos, dyPos, dzPos);
334  positronDirection.rotateUz(photonDirection);
335 
336  // Create G4DynamicParticle object for the particle2
338  positronDirection, positronKineEnergy);
339  // Fill output vector
340 // G4cout << "Cree el e+ " << epsilon << G4endl;
341  fvect->push_back(particle1);
342  fvect->push_back(particle2);
343 
344  // kill incident photon
347 
348 }
349 
350 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
351 
352 G4double G4LivermoreGammaConversionModelRC::ScreenFunction1(G4double screenVariable)
353 {
354  // Compute the value of the screening function 3*phi1 - phi2
355 
356  G4double value;
357 
358  if (screenVariable > 1.)
359  value = 42.24 - 8.368 * std::log(screenVariable + 0.952);
360  else
361  value = 42.392 - screenVariable * (7.796 - 1.961 * screenVariable);
362 
363  return value;
364 }
365 
366 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
367 
368 G4double G4LivermoreGammaConversionModelRC::ScreenFunction2(G4double screenVariable)
369 {
370  // Compute the value of the screening function 1.5*phi1 - 0.5*phi2
371 
372  G4double value;
373 
374  if (screenVariable > 1.)
375  value = 42.24 - 8.368 * std::log(screenVariable + 0.952);
376  else
377  value = 41.405 - screenVariable * (5.828 - 0.8945 * screenVariable);
378 
379  return value;
380 }
381