Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4LivermoreGammaConversionModel.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 // Author: Sebastien Incerti
27 // 22 January 2012
28 // on base of G4LivermoreGammaConversionModel (original version)
29 // and G4LivermoreRayleighModel (MT version)
30 
32 #include "G4PhysicalConstants.hh"
33 #include "G4SystemOfUnits.hh"
34 #include "G4Log.hh"
35 #include "G4Exp.hh"
36 
37 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
38 
39 using namespace std;
40 
41 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
42 
43 G4int G4LivermoreGammaConversionModel::maxZ = 99;
44 G4LPhysicsFreeVector* G4LivermoreGammaConversionModel::data[] = {0};
45 
47 (const G4ParticleDefinition*, const G4String& nam)
48 :G4VEmModel(nam),isInitialised(false),smallEnergy(2.*MeV)
49 {
50  fParticleChange = 0;
51 
52  lowEnergyLimit = 2.0*electron_mass_c2;
53 
54  verboseLevel= 0;
55  // Verbosity scale for debugging purposes:
56  // 0 = nothing
57  // 1 = calculation of cross sections, file openings...
58  // 2 = entering in methods
59 
60  if(verboseLevel > 0)
61  {
62  G4cout << "G4LivermoreGammaConversionModel is constructed " << G4endl;
63  }
64 }
65 
66 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
67 
69 {}
70 
71 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
72 
74  const G4ParticleDefinition* particle,
75  const G4DataVector& cuts)
76 {
77  if (verboseLevel > 1)
78  {
79  G4cout << "Calling Initialise() of G4LivermoreGammaConversionModel."
80  << G4endl
81  << "Energy range: "
82  << LowEnergyLimit() / MeV << " MeV - "
83  << HighEnergyLimit() / GeV << " GeV"
84  << G4endl;
85  }
86 
87  if(IsMaster())
88  {
89 
90  // Initialise element selector
91 
92  InitialiseElementSelectors(particle, cuts);
93 
94  // Access to elements
95 
96  char* path = getenv("G4LEDATA");
97 
98  G4ProductionCutsTable* theCoupleTable =
100 
101  G4int numOfCouples = theCoupleTable->GetTableSize();
102 
103  for(G4int i=0; i<numOfCouples; ++i)
104  {
105  const G4Material* material =
106  theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
107  const G4ElementVector* theElementVector = material->GetElementVector();
108  G4int nelm = material->GetNumberOfElements();
109 
110  for (G4int j=0; j<nelm; ++j)
111  {
112  G4int Z = (G4int)(*theElementVector)[j]->GetZ();
113  if(Z < 1) { Z = 1; }
114  else if(Z > maxZ) { Z = maxZ; }
115  if(!data[Z]) { ReadData(Z, path); }
116  }
117  }
118  }
119 
120  //
121 
122  if(isInitialised) { return; }
123  fParticleChange = GetParticleChangeForGamma();
124  isInitialised = true;
125 }
126 
127 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
128 
130  G4VEmModel* masterModel)
131 {
132  SetElementSelectors(masterModel->GetElementSelectors());
133 }
134 
135 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
136 
137 G4double
139  const G4ParticleDefinition*,
140  G4double)
141 {
142  return lowEnergyLimit;
143 }
144 
145 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
146 
147 void G4LivermoreGammaConversionModel::ReadData(size_t Z, const char* path)
148 {
149  if (verboseLevel > 1)
150  {
151  G4cout << "Calling ReadData() of G4LivermoreGammaConversionModel"
152  << G4endl;
153  }
154 
155  if(data[Z]) { return; }
156 
157  const char* datadir = path;
158 
159  if(!datadir)
160  {
161  datadir = getenv("G4LEDATA");
162  if(!datadir)
163  {
164  G4Exception("G4LivermoreGammaConversionModel::ReadData()",
165  "em0006",FatalException,
166  "Environment variable G4LEDATA not defined");
167  return;
168  }
169  }
170 
171  //
172 
173  data[Z] = new G4LPhysicsFreeVector();
174 
175  //
176 
177  std::ostringstream ost;
178  ost << datadir << "/livermore/pair/pp-cs-" << Z <<".dat";
179  std::ifstream fin(ost.str().c_str());
180 
181  if( !fin.is_open())
182  {
184  ed << "G4LivermoreGammaConversionModel data file <" << ost.str().c_str()
185  << "> is not opened!" << G4endl;
186  G4Exception("G4LivermoreGammaConversionModel::ReadData()",
187  "em0003",FatalException,
188  ed,"G4LEDATA version should be G4EMLOW6.27 or later.");
189  return;
190  }
191 
192  else
193  {
194 
195  if(verboseLevel > 3) { G4cout << "File " << ost.str()
196  << " is opened by G4LivermoreGammaConversionModel" << G4endl;}
197 
198  data[Z]->Retrieve(fin, true);
199  }
200 
201  // Activation of spline interpolation
202  data[Z] ->SetSpline(true);
203 
204 }
205 
206 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
207 
208 G4double
210  G4double GammaEnergy,
211  G4double Z, G4double,
213 {
214  if (verboseLevel > 1)
215  {
216  G4cout << "Calling ComputeCrossSectionPerAtom() of G4LivermoreGammaConversionModel"
217  << G4endl;
218  }
219 
220  if (GammaEnergy < lowEnergyLimit) { return 0.0; }
221 
222  G4double xs = 0.0;
223 
224  G4int intZ=G4int(Z);
225 
226  if(intZ < 1 || intZ > maxZ) { return xs; }
227 
228  G4LPhysicsFreeVector* pv = data[intZ];
229 
230  // if element was not initialised
231  // do initialisation safely for MT mode
232  if(!pv)
233  {
234  InitialiseForElement(0, intZ);
235  pv = data[intZ];
236  if(!pv) { return xs; }
237  }
238  // x-section is taken from the table
239  xs = pv->Value(GammaEnergy);
240 
241  if(verboseLevel > 0)
242  {
243  G4int n = pv->GetVectorLength() - 1;
244  G4cout << "****** DEBUG: tcs value for Z=" << Z << " at energy (MeV)="
245  << GammaEnergy/MeV << G4endl;
246  G4cout << " cs (Geant4 internal unit)=" << xs << G4endl;
247  G4cout << " -> first cs value in EADL data file (iu) =" << (*pv)[0] << G4endl;
248  G4cout << " -> last cs value in EADL data file (iu) =" << (*pv)[n] << G4endl;
249  G4cout << "*********************************************************" << G4endl;
250  }
251 
252  return xs;
253 
254 }
255 
256 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
257 
259  std::vector<G4DynamicParticle*>* fvect,
260  const G4MaterialCutsCouple* couple,
261  const G4DynamicParticle* aDynamicGamma,
263 {
264 
265 // The energies of the e+ e- secondaries are sampled using the Bethe - Heitler
266 // cross sections with Coulomb correction. A modified version of the random
267 // number techniques of Butcher & Messel is used (Nuc Phys 20(1960),15).
268 
269 // Note 1 : Effects due to the breakdown of the Born approximation at low
270 // energy are ignored.
271 // Note 2 : The differential cross section implicitly takes account of
272 // pair creation in both nuclear and atomic electron fields. However triplet
273 // prodution is not generated.
274 
275  if (verboseLevel > 1) {
276  G4cout << "Calling SampleSecondaries() of G4LivermoreGammaConversionModel"
277  << G4endl;
278  }
279 
280  G4double photonEnergy = aDynamicGamma->GetKineticEnergy();
281  G4ParticleMomentum photonDirection = aDynamicGamma->GetMomentumDirection();
282 
283  G4double epsilon ;
284  G4double epsilon0Local = electron_mass_c2 / photonEnergy ;
285 
286  // Do it fast if photon energy < 2. MeV
287  if (photonEnergy < smallEnergy )
288  {
289  epsilon = epsilon0Local + (0.5 - epsilon0Local) * G4UniformRand();
290  }
291  else
292  {
293  // Select randomly one element in the current material
294 
295  const G4ParticleDefinition* particle = aDynamicGamma->GetDefinition();
296  const G4Element* element = SelectRandomAtom(couple,particle,photonEnergy);
297 
298  if (element == 0)
299  {
300  G4cout << "G4LivermoreGammaConversionModel::SampleSecondaries - element = 0"
301  << G4endl;
302  return;
303  }
304  G4IonisParamElm* ionisation = element->GetIonisation();
305  if (ionisation == 0)
306  {
307  G4cout << "G4LivermoreGammaConversionModel::SampleSecondaries - ionisation = 0"
308  << G4endl;
309  return;
310  }
311 
312  // Extract Coulomb factor for this Element
313  G4double fZ = 8. * (ionisation->GetlogZ3());
314  if (photonEnergy > 50. * MeV) fZ += 8. * (element->GetfCoulomb());
315 
316  // Limits of the screening variable
317  G4double screenFactor = 136. * epsilon0Local / (element->GetIonisation()->GetZ3()) ;
318  G4double screenMax = G4Exp ((42.24 - fZ)/8.368) - 0.952 ;
319  G4double screenMin = std::min(4.*screenFactor,screenMax) ;
320 
321  // Limits of the energy sampling
322  G4double epsilon1 = 0.5 - 0.5 * std::sqrt(1. - screenMin / screenMax) ;
323  G4double epsilonMin = std::max(epsilon0Local,epsilon1);
324  G4double epsilonRange = 0.5 - epsilonMin ;
325 
326  // Sample the energy rate of the created electron (or positron)
327  G4double screen;
328  G4double gReject ;
329 
330  G4double f10 = ScreenFunction1(screenMin) - fZ;
331  G4double f20 = ScreenFunction2(screenMin) - fZ;
332  G4double normF1 = std::max(f10 * epsilonRange * epsilonRange,0.);
333  G4double normF2 = std::max(1.5 * f20,0.);
334 
335  do
336  {
337  if (normF1 / (normF1 + normF2) > G4UniformRand() )
338  {
339  epsilon = 0.5 - epsilonRange * std::pow(G4UniformRand(), 0.333333) ;
340  screen = screenFactor / (epsilon * (1. - epsilon));
341  gReject = (ScreenFunction1(screen) - fZ) / f10 ;
342  }
343  else
344  {
345  epsilon = epsilonMin + epsilonRange * G4UniformRand();
346  screen = screenFactor / (epsilon * (1 - epsilon));
347  gReject = (ScreenFunction2(screen) - fZ) / f20 ;
348  }
349  } while ( gReject < G4UniformRand() );
350 
351  } // End of epsilon sampling
352 
353  // Fix charges randomly
354 
355  G4double electronTotEnergy;
356  G4double positronTotEnergy;
357 
358  if (G4UniformRand() > 0.5)
359  {
360  electronTotEnergy = (1. - epsilon) * photonEnergy;
361  positronTotEnergy = epsilon * photonEnergy;
362  }
363  else
364  {
365  positronTotEnergy = (1. - epsilon) * photonEnergy;
366  electronTotEnergy = epsilon * photonEnergy;
367  }
368 
369  // Scattered electron (positron) angles. ( Z - axis along the parent photon)
370  // Universal distribution suggested by L. Urban (Geant3 manual (1993) Phys211),
371  // derived from Tsai distribution (Rev. Mod. Phys. 49, 421 (1977)
372 
373  G4double u;
374  const G4double a1 = 0.625;
375  G4double a2 = 3. * a1;
376  // G4double d = 27. ;
377 
378  // if (9. / (9. + d) > G4UniformRand())
379  if (0.25 > G4UniformRand())
380  {
381  u = - G4Log(G4UniformRand() * G4UniformRand()) / a1 ;
382  }
383  else
384  {
385  u = - G4Log(G4UniformRand() * G4UniformRand()) / a2 ;
386  }
387 
388  G4double thetaEle = u*electron_mass_c2/electronTotEnergy;
389  G4double thetaPos = u*electron_mass_c2/positronTotEnergy;
390  G4double phi = twopi * G4UniformRand();
391 
392  G4double dxEle= std::sin(thetaEle)*std::cos(phi),dyEle= std::sin(thetaEle)*std::sin(phi),dzEle=std::cos(thetaEle);
393  G4double dxPos=-std::sin(thetaPos)*std::cos(phi),dyPos=-std::sin(thetaPos)*std::sin(phi),dzPos=std::cos(thetaPos);
394 
395 
396  // Kinematics of the created pair:
397  // the electron and positron are assumed to have a symetric angular
398  // distribution with respect to the Z axis along the parent photon
399 
400  G4double electronKineEnergy = std::max(0.,electronTotEnergy - electron_mass_c2) ;
401 
402  G4ThreeVector electronDirection (dxEle, dyEle, dzEle);
403  electronDirection.rotateUz(photonDirection);
404 
406  electronDirection,
407  electronKineEnergy);
408 
409  // The e+ is always created
410  G4double positronKineEnergy = std::max(0.,positronTotEnergy - electron_mass_c2) ;
411 
412  G4ThreeVector positronDirection (dxPos, dyPos, dzPos);
413  positronDirection.rotateUz(photonDirection);
414 
415  // Create G4DynamicParticle object for the particle2
417  positronDirection,
418  positronKineEnergy);
419  // Fill output vector
420  fvect->push_back(particle1);
421  fvect->push_back(particle2);
422 
423  // kill incident photon
424  fParticleChange->SetProposedKineticEnergy(0.);
425  fParticleChange->ProposeTrackStatus(fStopAndKill);
426 
427 }
428 
429 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
430 
431 G4double
432 G4LivermoreGammaConversionModel::ScreenFunction1(G4double screenVariable)
433 {
434  // Compute the value of the screening function 3*phi1 - phi2
435 
436  G4double value;
437 
438  if (screenVariable > 1.)
439  value = 42.24 - 8.368 * G4Log(screenVariable + 0.952);
440  else
441  value = 42.392 - screenVariable * (7.796 - 1.961 * screenVariable);
442 
443  return value;
444 }
445 
446 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
447 
448 G4double
449 G4LivermoreGammaConversionModel::ScreenFunction2(G4double screenVariable)
450 {
451  // Compute the value of the screening function 1.5*phi1 - 0.5*phi2
452 
453  G4double value;
454 
455  if (screenVariable > 1.)
456  value = 42.24 - 8.368 * G4Log(screenVariable + 0.952);
457  else
458  value = 41.405 - screenVariable * (5.828 - 0.8945 * screenVariable);
459 
460  return value;
461 }
462 
463 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
464 
465 #include "G4AutoLock.hh"
466 namespace { G4Mutex LivermoreGammaConversionModelMutex = G4MUTEX_INITIALIZER; }
467 
469  const G4ParticleDefinition*,
470  G4int Z)
471 {
472  G4AutoLock l(&LivermoreGammaConversionModelMutex);
473  // G4cout << "G4LivermoreGammaConversionModel::InitialiseForElement Z= "
474  // << Z << G4endl;
475  if(!data[Z]) { ReadData(Z); }
476  l.unlock();
477 }
478 
479 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double)
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
std::vector< G4Element * > G4ElementVector
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4double GetKineticEnergy() const
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX)
G4double GetfCoulomb() const
Definition: G4Element.hh:190
G4ParticleDefinition * GetDefinition() const
size_t GetVectorLength() const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
int G4int
Definition: G4Types.hh:78
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:158
string material
Definition: eplot.py:19
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
G4LivermoreGammaConversionModel(const G4ParticleDefinition *p=0, const G4String &nam="LivermoreConversion")
const G4ThreeVector & GetMomentumDirection() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
float electron_mass_c2
Definition: hepunit.py:274
const G4int n
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:760
G4double Value(G4double theEnergy, size_t &lastidx) const
G4double GetlogZ3() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4double G4Log(G4double x)
Definition: G4Log.hh:227
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
static G4ProductionCutsTable * GetProductionCutsTable()
G4int G4Mutex
Definition: G4Threading.hh:156
static G4Positron * Positron()
Definition: G4Positron.cc:94
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4IonisParamElm * GetIonisation() const
Definition: G4Element.hh:198
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
#define G4endl
Definition: G4ios.hh:61
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
double G4double
Definition: G4Types.hh:76
G4double GetZ3() const
const XML_Char const XML_Char * data
const G4Material * GetMaterial() const