Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4LivermoreComptonModel.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: G4LivermoreComptonModel.cc 84216 2014-10-10 14:51:51Z gcosmo $
27 // GEANT4 tag $Name: not supported by cvs2svn $
28 //
29 //
30 // Author: Alexander Bagulya
31 // 11 March 2012
32 // on the base of G4LivermoreComptonModel
33 //
34 // History:
35 // --------
36 // 18 Apr 2009 V Ivanchenko Cleanup initialisation and generation of secondaries:
37 // - apply internal high-energy limit only in constructor
38 // - do not apply low-energy limit (default is 0)
39 // - remove GetMeanFreePath method and table
40 // - added protection against numerical problem in energy sampling
41 // - use G4ElementSelector
42 // 26 Dec 2010 V Ivanchenko Load data tables only once to avoid memory leak
43 
45 #include "G4PhysicalConstants.hh"
46 #include "G4SystemOfUnits.hh"
47 #include "G4Electron.hh"
49 #include "G4LossTableManager.hh"
50 #include "G4VAtomDeexcitation.hh"
51 #include "G4AtomicShell.hh"
52 #include "G4Gamma.hh"
53 #include "G4ShellData.hh"
54 #include "G4DopplerProfile.hh"
55 #include "G4Log.hh"
56 #include "G4Exp.hh"
57 
58 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
59 
60 using namespace std;
61 
62 G4int G4LivermoreComptonModel::maxZ = 99;
63 G4LPhysicsFreeVector* G4LivermoreComptonModel::data[] = {0};
64 G4ShellData* G4LivermoreComptonModel::shellData = 0;
65 G4DopplerProfile* G4LivermoreComptonModel::profileData = 0;
66 
67 static const G4double ln10 = G4Log(10.);
68 
70  const G4String& nam)
71  : G4VEmModel(nam),isInitialised(false)
72 {
73  verboseLevel=1 ;
74  // Verbosity scale:
75  // 0 = nothing
76  // 1 = warning for energy non-conservation
77  // 2 = details of energy budget
78  // 3 = calculation of cross sections, file openings, sampling of atoms
79  // 4 = entering in methods
80 
81  if( verboseLevel>1 ) {
82  G4cout << "Livermore Compton model is constructed " << G4endl;
83  }
84 
85  //Mark this model as "applicable" for atomic deexcitation
86  SetDeexcitationFlag(true);
87 
88  fParticleChange = 0;
89  fAtomDeexcitation = 0;
90 }
91 
92 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
93 
95 {
96  if(IsMaster()) {
97  delete shellData;
98  shellData = 0;
99  delete profileData;
100  profileData = 0;
101  for(G4int i=0; i<maxZ; ++i) {
102  if(data[i]) {
103  delete data[i];
104  data[i] = 0;
105  }
106  }
107  }
108 }
109 
110 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
111 
113  const G4DataVector& cuts)
114 {
115  if (verboseLevel > 1) {
116  G4cout << "Calling G4LivermoreComptonModel::Initialise()" << G4endl;
117  }
118 
119  // Initialise element selector
120 
121  if(IsMaster()) {
122 
123  // Access to elements
124 
125  char* path = getenv("G4LEDATA");
126 
127  G4ProductionCutsTable* theCoupleTable =
129  G4int numOfCouples = theCoupleTable->GetTableSize();
130 
131  for(G4int i=0; i<numOfCouples; ++i) {
132  const G4Material* material =
133  theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
134  const G4ElementVector* theElementVector = material->GetElementVector();
135  G4int nelm = material->GetNumberOfElements();
136 
137  for (G4int j=0; j<nelm; ++j) {
138  G4int Z = G4lrint((*theElementVector)[j]->GetZ());
139  if(Z < 1) { Z = 1; }
140  else if(Z > maxZ){ Z = maxZ; }
141 
142  if( (!data[Z]) ) { ReadData(Z, path); }
143  }
144  }
145 
146  // For Doppler broadening
147  if(!shellData) {
148  shellData = new G4ShellData();
149  shellData->SetOccupancyData();
150  G4String file = "/doppler/shell-doppler";
151  shellData->LoadData(file);
152  }
153  if(!profileData) { profileData = new G4DopplerProfile(); }
154 
155  InitialiseElementSelectors(particle, cuts);
156  }
157 
158  if (verboseLevel > 2) {
159  G4cout << "Loaded cross section files" << G4endl;
160  }
161 
162  if( verboseLevel>1 ) {
163  G4cout << "G4LivermoreComptonModel is initialized " << G4endl
164  << "Energy range: "
165  << LowEnergyLimit() / eV << " eV - "
166  << HighEnergyLimit() / GeV << " GeV"
167  << G4endl;
168  }
169  //
170  if(isInitialised) { return; }
171 
172  fParticleChange = GetParticleChangeForGamma();
173  fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
174  isInitialised = true;
175 }
176 
177 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
178 
180  G4VEmModel* masterModel)
181 {
183 }
184 
185 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
186 
187 void G4LivermoreComptonModel::ReadData(size_t Z, const char* path)
188 {
189  if (verboseLevel > 1)
190  {
191  G4cout << "G4LivermoreComptonModel::ReadData()"
192  << G4endl;
193  }
194  if(data[Z]) { return; }
195  const char* datadir = path;
196  if(!datadir)
197  {
198  datadir = getenv("G4LEDATA");
199  if(!datadir)
200  {
201  G4Exception("G4LivermoreComptonModel::ReadData()",
202  "em0006",FatalException,
203  "Environment variable G4LEDATA not defined");
204  return;
205  }
206  }
207 
208  data[Z] = new G4LPhysicsFreeVector();
209 
210  // Activation of spline interpolation
211  data[Z]->SetSpline(false);
212 
213  std::ostringstream ost;
214  ost << datadir << "/livermore/comp/ce-cs-" << Z <<".dat";
215  std::ifstream fin(ost.str().c_str());
216 
217  if( !fin.is_open())
218  {
220  ed << "G4LivermoreComptonModel data file <" << ost.str().c_str()
221  << "> is not opened!" << G4endl;
222  G4Exception("G4LivermoreComptonModel::ReadData()",
223  "em0003",FatalException,
224  ed,"G4LEDATA version should be G4EMLOW6.34 or later");
225  return;
226  } else {
227  if(verboseLevel > 3) {
228  G4cout << "File " << ost.str()
229  << " is opened by G4LivermoreComptonModel" << G4endl;
230  }
231  data[Z]->Retrieve(fin, true);
232  data[Z]->ScaleVector(MeV, MeV*barn);
233  }
234  fin.close();
235 }
236 
237 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
238 
239 G4double
241  G4double GammaEnergy,
242  G4double Z, G4double,
244 {
245  if (verboseLevel > 3) {
246  G4cout << "G4LivermoreComptonModel::ComputeCrossSectionPerAtom()"
247  << G4endl;
248  }
249  G4double cs = 0.0;
250 
251  if (GammaEnergy < LowEnergyLimit()) { return 0.0; }
252 
253  G4int intZ = G4lrint(Z);
254  if(intZ < 1 || intZ > maxZ) { return cs; }
255 
256  G4LPhysicsFreeVector* pv = data[intZ];
257 
258  // if element was not initialised
259  // do initialisation safely for MT mode
260  if(!pv)
261  {
262  InitialiseForElement(0, intZ);
263  pv = data[intZ];
264  if(!pv) { return cs; }
265  }
266 
267  G4int n = pv->GetVectorLength() - 1;
268  G4double e1 = pv->Energy(0);
269  G4double e2 = pv->Energy(n);
270 
271  if(GammaEnergy <= e1) { cs = GammaEnergy/(e1*e1)*pv->Value(e1); }
272  else if(GammaEnergy <= e2) { cs = pv->Value(GammaEnergy)/GammaEnergy; }
273  else if(GammaEnergy > e2) { cs = pv->Value(e2)/GammaEnergy; }
274 
275  return cs;
276 }
277 
278 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
279 
280 
282  std::vector<G4DynamicParticle*>* fvect,
283  const G4MaterialCutsCouple* couple,
284  const G4DynamicParticle* aDynamicGamma,
286 {
287 
288  // The scattered gamma energy is sampled according to Klein - Nishina
289  // formula then accepted or rejected depending on the Scattering Function
290  // multiplied by factor from Klein - Nishina formula.
291  // Expression of the angular distribution as Klein Nishina
292  // angular and energy distribution and Scattering fuctions is taken from
293  // D. E. Cullen "A simple model of photon transport" Nucl. Instr. Meth.
294  // Phys. Res. B 101 (1995). Method of sampling with form factors is different
295  // data are interpolated while in the article they are fitted.
296  // Reference to the article is from J. Stepanek New Photon, Positron
297  // and Electron Interaction Data for GEANT in Energy Range from 1 eV to 10
298  // TeV (draft).
299  // The random number techniques of Butcher & Messel are used
300  // (Nucl Phys 20(1960),15).
301 
302  G4double photonEnergy0 = aDynamicGamma->GetKineticEnergy();
303 
304  if (verboseLevel > 3) {
305  G4cout << "G4LivermoreComptonModel::SampleSecondaries() E(MeV)= "
306  << photonEnergy0/MeV << " in " << couple->GetMaterial()->GetName()
307  << G4endl;
308  }
309 
310  // do nothing below the threshold
311  // should never get here because the XS is zero below the limit
312  if (photonEnergy0 < LowEnergyLimit())
313  return ;
314 
315  G4double e0m = photonEnergy0 / electron_mass_c2 ;
316  G4ParticleMomentum photonDirection0 = aDynamicGamma->GetMomentumDirection();
317 
318  // Select randomly one element in the current material
319  const G4ParticleDefinition* particle = aDynamicGamma->GetDefinition();
320  const G4Element* elm = SelectRandomAtom(couple,particle,photonEnergy0);
321 
322  G4int Z = G4lrint(elm->GetZ());
323 
324  G4double epsilon0Local = 1. / (1. + 2. * e0m);
325  G4double epsilon0Sq = epsilon0Local * epsilon0Local;
326  G4double alpha1 = -G4Log(epsilon0Local);
327  G4double alpha2 = 0.5 * (1. - epsilon0Sq);
328 
329  G4double wlPhoton = h_Planck*c_light/photonEnergy0;
330 
331  // Sample the energy of the scattered photon
333  G4double epsilonSq;
334  G4double oneCosT;
335  G4double sinT2;
336  G4double gReject;
337 
338  if (verboseLevel > 3) {
339  G4cout << "Started loop to sample gamma energy" << G4endl;
340  }
341 
342  do {
343  if ( alpha1/(alpha1+alpha2) > G4UniformRand())
344  {
345  epsilon = G4Exp(-alpha1 * G4UniformRand());
346  epsilonSq = epsilon * epsilon;
347  }
348  else
349  {
350  epsilonSq = epsilon0Sq + (1. - epsilon0Sq) * G4UniformRand();
351  epsilon = std::sqrt(epsilonSq);
352  }
353 
354  oneCosT = (1. - epsilon) / ( epsilon * e0m);
355  sinT2 = oneCosT * (2. - oneCosT);
356  G4double x = std::sqrt(oneCosT/2.) * cm / wlPhoton;
357  G4double scatteringFunction = ComputeScatteringFunction(x, Z);
358  gReject = (1. - epsilon * sinT2 / (1. + epsilonSq)) * scatteringFunction;
359 
360  } while(gReject < G4UniformRand()*Z);
361 
362  G4double cosTheta = 1. - oneCosT;
363  G4double sinTheta = std::sqrt (sinT2);
364  G4double phi = twopi * G4UniformRand() ;
365  G4double dirx = sinTheta * std::cos(phi);
366  G4double diry = sinTheta * std::sin(phi);
367  G4double dirz = cosTheta ;
368 
369  // Doppler broadening - Method based on:
370  // Y. Namito, S. Ban and H. Hirayama,
371  // "Implementation of the Doppler Broadening of a Compton-Scattered Photon
372  // into the EGS4 Code", NIM A 349, pp. 489-494, 1994
373 
374  // Maximum number of sampling iterations
375  static G4int maxDopplerIterations = 1000;
376  G4double bindingE = 0.;
377  G4double photonEoriginal = epsilon * photonEnergy0;
378  G4double photonE = -1.;
379  G4int iteration = 0;
380  G4double eMax = photonEnergy0;
381 
382  G4int shellIdx = 0;
383 
384  if (verboseLevel > 3) {
385  G4cout << "Started loop to sample broading" << G4endl;
386  }
387 
388  do
389  {
390  ++iteration;
391  // Select shell based on shell occupancy
392  shellIdx = shellData->SelectRandomShell(Z);
393  bindingE = shellData->BindingEnergy(Z,shellIdx);
394 
395  if (verboseLevel > 3) {
396  G4cout << "Shell ID= " << shellIdx
397  << " Ebind(keV)= " << bindingE/keV << G4endl;
398  }
399 
400  eMax = photonEnergy0 - bindingE;
401 
402  // Randomly sample bound electron momentum
403  // (memento: the data set is in Atomic Units)
404  G4double pSample = profileData->RandomSelectMomentum(Z,shellIdx);
405  if (verboseLevel > 3) {
406  G4cout << "pSample= " << pSample << G4endl;
407  }
408  // Rescale from atomic units
409  G4double pDoppler = pSample * fine_structure_const;
410  G4double pDoppler2 = pDoppler * pDoppler;
411  G4double var2 = 1. + oneCosT * e0m;
412  G4double var3 = var2*var2 - pDoppler2;
413  G4double var4 = var2 - pDoppler2 * cosTheta;
414  G4double var = var4*var4 - var3 + pDoppler2 * var3;
415  if (var > 0.)
416  {
417  G4double varSqrt = std::sqrt(var);
418  G4double scale = photonEnergy0 / var3;
419  // Random select either root
420  if (G4UniformRand() < 0.5) { photonE = (var4 - varSqrt) * scale; }
421  else { photonE = (var4 + varSqrt) * scale; }
422  }
423  else
424  {
425  photonE = -1.;
426  }
427  } while (iteration <= maxDopplerIterations && photonE > eMax);
428  // (photonE < 0. || photonE > eMax || photonE < eMax*G4UniformRand()) );
429 
430 
431  // End of recalculation of photon energy with Doppler broadening
432  // Revert to original if maximum number of iterations threshold
433  // has been reached
434 
435  if (iteration >= maxDopplerIterations)
436  {
437  photonE = photonEoriginal;
438  bindingE = 0.;
439  }
440 
441  // Update G4VParticleChange for the scattered photon
442 
443  G4ThreeVector photonDirection1(dirx,diry,dirz);
444  photonDirection1.rotateUz(photonDirection0);
445  fParticleChange->ProposeMomentumDirection(photonDirection1) ;
446 
447  G4double photonEnergy1 = photonE;
448 
449  if (photonEnergy1 > 0.) {
450  fParticleChange->SetProposedKineticEnergy(photonEnergy1) ;
451 
452  } else {
453  // photon absorbed
454  photonEnergy1 = 0.;
455  fParticleChange->SetProposedKineticEnergy(0.) ;
456  fParticleChange->ProposeTrackStatus(fStopAndKill);
457  fParticleChange->ProposeLocalEnergyDeposit(photonEnergy0);
458  return;
459  }
460 
461  // Kinematics of the scattered electron
462  G4double eKineticEnergy = photonEnergy0 - photonEnergy1 - bindingE;
463 
464  // protection against negative final energy: no e- is created
465  if(eKineticEnergy < 0.0) {
466  fParticleChange->ProposeLocalEnergyDeposit(photonEnergy0 - photonEnergy1);
467  return;
468  }
469 
470  G4double eTotalEnergy = eKineticEnergy + electron_mass_c2;
471 
472  G4double electronE = photonEnergy0 * (1. - epsilon) + electron_mass_c2;
473  G4double electronP2 =
474  electronE*electronE - electron_mass_c2*electron_mass_c2;
475  G4double sinThetaE = -1.;
476  G4double cosThetaE = 0.;
477  if (electronP2 > 0.)
478  {
479  cosThetaE = (eTotalEnergy + photonEnergy1 )*
480  (1. - epsilon) / std::sqrt(electronP2);
481  sinThetaE = -1. * sqrt(1. - cosThetaE * cosThetaE);
482  }
483 
484  G4double eDirX = sinThetaE * std::cos(phi);
485  G4double eDirY = sinThetaE * std::sin(phi);
486  G4double eDirZ = cosThetaE;
487 
488  G4ThreeVector eDirection(eDirX,eDirY,eDirZ);
489  eDirection.rotateUz(photonDirection0);
491  eDirection,eKineticEnergy) ;
492  fvect->push_back(dp);
493 
494  // sample deexcitation
495  //
496 
497  if (verboseLevel > 3) {
498  G4cout << "Started atomic de-excitation " << fAtomDeexcitation << G4endl;
499  }
500 
501  if(fAtomDeexcitation && iteration < maxDopplerIterations) {
502  G4int index = couple->GetIndex();
503  if(fAtomDeexcitation->CheckDeexcitationActiveRegion(index)) {
504  size_t nbefore = fvect->size();
506  const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as);
507  fAtomDeexcitation->GenerateParticles(fvect, shell, Z, index);
508  size_t nafter = fvect->size();
509  if(nafter > nbefore) {
510  for (size_t i=nbefore; i<nafter; ++i) {
511  //Check if there is enough residual energy
512  if (bindingE >= ((*fvect)[i])->GetKineticEnergy())
513  {
514  //Ok, this is a valid secondary: keep it
515  bindingE -= ((*fvect)[i])->GetKineticEnergy();
516  }
517  else
518  {
519  //Invalid secondary: not enough energy to create it!
520  //Keep its energy in the local deposit
521  delete (*fvect)[i];
522  (*fvect)[i]=0;
523  }
524  }
525  }
526  }
527  }
528  //This should never happen
529  if(bindingE < 0.0)
530  G4Exception("G4LivermoreComptonModel::SampleSecondaries()",
531  "em2050",FatalException,"Negative local energy deposit");
532 
533  fParticleChange->ProposeLocalEnergyDeposit(bindingE);
534 }
535 
536 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
537 
538 G4double
539 G4LivermoreComptonModel::ComputeScatteringFunction(G4double x, G4int Z)
540 {
541  G4double value = Z;
542  if (x <= ScatFuncFitParam[Z][2]) {
543 
544  G4double lgq = G4Log(x)/ln10;
545 
546  if (lgq < ScatFuncFitParam[Z][1]) {
547  value = ScatFuncFitParam[Z][3] + lgq*ScatFuncFitParam[Z][4];
548  } else {
549  value = ScatFuncFitParam[Z][5] + lgq*ScatFuncFitParam[Z][6] +
550  lgq*lgq*ScatFuncFitParam[Z][7] + lgq*lgq*lgq*ScatFuncFitParam[Z][8];
551  }
552  value = G4Exp(value*ln10);
553  }
554  //G4cout << " value= " << value << G4endl;
555  return value;
556 }
557 
558 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
559 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
560 
561 #include "G4AutoLock.hh"
562 namespace { G4Mutex LivermoreComptonModelMutex = G4MUTEX_INITIALIZER; }
563 
564 void
566  G4int Z)
567 {
568  G4AutoLock l(&LivermoreComptonModelMutex);
569  // G4cout << "G4LivermoreComptonModel::InitialiseForElement Z= "
570  // << Z << G4endl;
571  if(!data[Z]) { ReadData(Z); }
572  l.unlock();
573 }
574 
575 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
576 
577 //Fitting data to compute scattering function
578 const G4double G4LivermoreComptonModel::ScatFuncFitParam[101][9] = {
579 { 0, 0., 0., 0., 0., 0., 0., 0., 0.},
580 { 1, 6.673, 1.49968E+08, -14.352, 1.999, -143.374, 50.787, -5.951, 0.2304418},
581 { 2, 6.500, 2.50035E+08, -14.215, 1.970, -53.649, 13.892, -0.948, 0.006996759},
582 { 3, 6.551, 3.99945E+08, -13.555, 1.993, -62.090, 21.462, -2.453, 0.093416},
583 { 4, 6.500, 5.00035E+08, -13.746, 1.998, -127.906, 46.491, -5.614, 0.2262103},
584 { 5, 6.500, 5.99791E+08, -13.800, 1.998, -131.153, 47.132, -5.619, 0.2233819},
585 { 6, 6.708, 6.99842E+08, -13.885, 1.999, -128.143, 45.379, -5.325, 0.2083009},
586 { 7, 6.685, 7.99834E+08, -13.885, 2.000, -131.048, 46.314, -5.421, 0.2114925},
587 { 8, 6.669, 7.99834E+08, -13.962, 2.001, -128.225, 44.818, -5.183, 0.1997155},
588 { 9, 6.711, 7.99834E+08, -13.999, 2.000, -122.112, 42.103, -4.796, 0.1819099},
589 { 10, 6.702, 7.99834E+08, -14.044, 1.999, -110.143, 37.225, -4.143, 0.1532094},
590 { 11, 6.425, 1.00000E+09, -13.423, 1.993, -41.137, 12.313, -1.152, 0.03384553},
591 { 12, 6.542, 1.00000E+09, -13.389, 1.997, -53.549, 17.420, -1.840, 0.06431849},
592 { 13, 6.570, 1.49968E+09, -13.401, 1.997, -66.243, 22.297, -2.460, 0.09045854},
593 { 14, 6.364, 1.49968E+09, -13.452, 1.999, -78.271, 26.757, -3.008, 0.1128195},
594 { 15, 6.500, 1.49968E+09, -13.488, 1.998, -85.069, 29.164, -3.291, 0.1239113},
595 { 16, 6.500, 1.49968E+09, -13.532, 1.998, -93.640, 32.274, -3.665, 0.1388633},
596 { 17, 6.500, 1.49968E+09, -13.584, 2.000, -98.534, 33.958, -3.857, 0.1461557},
597 { 18, 6.500, 1.49968E+09, -13.618, 1.999, -100.077, 34.379, -3.891, 0.1468902},
598 { 19, 6.500, 1.99986E+09, -13.185, 1.992, -53.819, 17.528, -1.851, 0.0648722},
599 { 20, 6.490, 1.99986E+09, -13.123, 1.993, -52.221, 17.169, -1.832, 0.06502094},
600 { 21, 6.498, 1.99986E+09, -13.157, 1.994, -55.365, 18.276, -1.961, 0.07002778},
601 { 22, 6.495, 1.99986E+09, -13.183, 1.994, -57.412, 18.957, -2.036, 0.07278856},
602 { 23, 6.487, 1.99986E+09, -13.216, 1.995, -58.478, 19.270, -2.065, 0.07362722},
603 { 24, 6.500, 1.99986E+09, -13.330, 1.997, -62.192, 20.358, -2.167, 0.07666583},
604 { 25, 6.488, 1.99986E+09, -13.277, 1.997, -58.007, 18.924, -2.003, 0.0704305},
605 { 26, 6.500, 5.00035E+09, -13.292, 1.997, -61.176, 20.067, -2.141, 0.0760269},
606 { 27, 6.500, 5.00035E+09, -13.321, 1.998, -61.909, 20.271, -2.159, 0.07653559},
607 { 28, 6.500, 5.00035E+09, -13.340, 1.998, -62.402, 20.391, -2.167, 0.07664243},
608 { 29, 6.500, 5.00035E+09, -13.439, 1.998, -67.305, 21.954, -2.331, 0.0823267},
609 { 30, 6.500, 5.00035E+09, -13.383, 1.999, -62.064, 20.136, -2.122, 0.07437589},
610 { 31, 6.500, 5.00035E+09, -13.349, 1.997, -61.068, 19.808, -2.086, 0.07307488},
611 { 32, 6.500, 5.00035E+09, -13.373, 1.999, -63.126, 20.553, -2.175, 0.07660222},
612 { 33, 6.500, 5.00035E+09, -13.395, 1.999, -65.674, 21.445, -2.278, 0.08054694},
613 { 34, 6.500, 5.00035E+09, -13.417, 1.999, -69.457, 22.811, -2.442, 0.08709536},
614 { 35, 6.500, 5.00035E+09, -13.442, 2.000, -72.283, 23.808, -2.558, 0.09156808},
615 { 36, 6.500, 5.00035E+09, -13.451, 1.998, -74.696, 24.641, -2.653, 0.09516597},
616 { 37, 6.500, 5.00035E+09, -13.082, 1.991, -46.235, 14.519, -1.458, 0.04837659},
617 { 38, 6.465, 5.00035E+09, -13.022, 1.993, -41.784, 13.065, -1.300, 0.04267703},
618 { 39, 6.492, 5.00035E+09, -13.043, 1.994, -44.609, 14.114, -1.429, 0.0479348},
619 { 40, 6.499, 5.00035E+09, -13.064, 1.994, -47.142, 15.019, -1.536, 0.0521347},
620 { 41, 6.384, 5.00035E+09, -13.156, 1.996, -53.114, 17.052, -1.766, 0.06079426},
621 { 42, 6.500, 5.00035E+09, -13.176, 1.996, -54.590, 17.550, -1.822, 0.06290335},
622 { 43, 6.500, 5.00035E+09, -13.133, 1.997, -51.272, 16.423, -1.694, 0.05806108},
623 { 44, 6.500, 5.00035E+09, -13.220, 1.996, -58.314, 18.839, -1.969, 0.0684608},
624 { 45, 6.500, 5.00035E+09, -13.246, 1.998, -59.674, 19.295, -2.020, 0.07037294},
625 { 46, 6.500, 5.00035E+09, -13.407, 1.999, -72.228, 23.693, -2.532, 0.09017969},
626 { 47, 6.500, 5.00035E+09, -13.277, 1.998, -60.890, 19.647, -2.053, 0.07138694},
627 { 48, 6.500, 5.00035E+09, -13.222, 1.998, -56.152, 18.002, -1.863, 0.06410123},
628 { 49, 6.500, 5.00035E+09, -13.199, 1.997, -56.208, 18.052, -1.872, 0.06456884},
629 { 50, 6.500, 5.00035E+09, -13.215, 1.998, -58.478, 18.887, -1.973, 0.06860356},
630 { 51, 6.500, 5.00035E+09, -13.230, 1.998, -60.708, 19.676, -2.066, 0.07225841},
631 { 52, 6.500, 7.99834E+09, -13.246, 1.998, -63.341, 20.632, -2.180, 0.0767412},
632 { 53, 6.500, 5.00035E+09, -13.262, 1.998, -66.339, 21.716, -2.310, 0.08191981},
633 { 54, 6.500, 7.99834E+09, -13.279, 1.998, -67.649, 22.151, -2.357, 0.08357825},
634 { 55, 6.500, 5.00035E+09, -12.951, 1.990, -45.302, 14.219, -1.423, 0.04712317},
635 { 56, 6.425, 5.00035E+09, -12.882, 1.992, -39.825, 12.363, -1.214, 0.03931009},
636 { 57, 6.466, 2.82488E+09, -12.903, 1.992, -38.952, 11.982, -1.160, 0.03681554},
637 { 58, 6.451, 5.00035E+09, -12.915, 1.993, -41.959, 13.118, -1.302, 0.04271291},
638 { 59, 6.434, 5.00035E+09, -12.914, 1.993, -40.528, 12.555, -1.230, 0.03971407},
639 { 60, 6.444, 5.00035E+09, -12.922, 1.992, -39.986, 12.329, -1.200, 0.03843737},
640 { 61, 6.414, 7.99834E+09, -12.930, 1.993, -42.756, 13.362, -1.327, 0.0436124},
641 { 62, 6.420, 7.99834E+09, -12.938, 1.992, -42.682, 13.314, -1.319, 0.04322509},
642 { 63, 6.416, 7.99834E+09, -12.946, 1.993, -42.399, 13.185, -1.301, 0.04243861},
643 { 64, 6.443, 7.99834E+09, -12.963, 1.993, -43.226, 13.475, -1.335, 0.04377341},
644 { 65, 6.449, 7.99834E+09, -12.973, 1.993, -43.232, 13.456, -1.330, 0.04347536},
645 { 66, 6.419, 7.99834E+09, -12.966, 1.993, -42.047, 12.990, -1.270, 0.04095499},
646 { 67, 6.406, 7.99834E+09, -12.976, 1.993, -42.405, 13.106, -1.283, 0.04146024},
647 { 68, 6.424, 7.99834E+09, -12.986, 1.993, -41.974, 12.926, -1.259, 0.040435},
648 { 69, 6.417, 7.99834E+09, -12.989, 1.993, -42.132, 12.967, -1.262, 0.04048908},
649 { 70, 6.405, 7.99834E+09, -13.000, 1.994, -42.582, 13.122, -1.280, 0.04119599},
650 { 71, 6.449, 7.99834E+09, -13.015, 1.994, -42.586, 13.115, -1.278, 0.04107587},
651 { 72, 6.465, 7.99834E+09, -13.030, 1.994, -43.708, 13.509, -1.324, 0.04286491},
652 { 73, 6.447, 7.99834E+09, -13.048, 1.996, -44.838, 13.902, -1.369, 0.04457132},
653 { 74, 6.452, 7.99834E+09, -13.073, 1.997, -45.545, 14.137, -1.395, 0.04553459},
654 { 75, 6.432, 7.99834E+09, -13.082, 1.997, -46.426, 14.431, -1.428, 0.04678218},
655 { 76, 6.439, 7.99834E+09, -13.100, 1.997, -47.513, 14.806, -1.471, 0.04842566},
656 { 77, 6.432, 7.99834E+09, -13.110, 1.997, -48.225, 15.042, -1.497, 0.04938364},
657 { 78, 6.500, 7.99834E+09, -13.185, 1.997, -53.256, 16.739, -1.687, 0.05645173},
658 { 79, 6.500, 7.99834E+09, -13.200, 1.997, -53.900, 16.946, -1.709, 0.05723134},
659 { 80, 6.500, 7.99834E+09, -13.156, 1.998, -49.801, 15.536, -1.547, 0.05103522},
660 { 81, 6.500, 7.99834E+09, -13.128, 1.997, -49.651, 15.512, -1.548, 0.05123203},
661 { 82, 6.500, 7.99834E+09, -13.134, 1.997, -51.021, 16.018, -1.609, 0.05364831},
662 { 83, 6.500, 7.99834E+09, -13.148, 1.998, -52.693, 16.612, -1.679, 0.05638698},
663 { 84, 6.500, 7.99834E+09, -13.161, 1.998, -54.415, 17.238, -1.754, 0.05935566},
664 { 85, 6.500, 7.99834E+09, -13.175, 1.998, -56.083, 17.834, -1.824, 0.06206068},
665 { 86, 6.500, 7.99834E+09, -13.189, 1.998, -57.860, 18.463, -1.898, 0.0649633},
666 { 87, 6.500, 7.99834E+09, -12.885, 1.990, -39.973, 12.164, -1.162, 0.0364598},
667 { 88, 6.417, 7.99834E+09, -12.816, 1.991, -34.591, 10.338, -0.956, 0.0287409},
668 { 89, 6.442, 7.99834E+09, -12.831, 1.992, -36.002, 10.867, -1.021, 0.03136835},
669 { 90, 6.463, 7.99834E+09, -12.850, 1.993, -37.660, 11.475, -1.095, 0.03435334},
670 { 91, 6.447, 7.99834E+09, -12.852, 1.993, -37.268, 11.301, -1.071, 0.0330539},
671 { 92, 6.439, 7.99834E+09, -12.858, 1.993, -37.695, 11.438, -1.085, 0.03376669},
672 { 93, 6.437, 1.00000E+10, -12.866, 1.993, -39.010, 11.927, -1.146, 0.03630848},
673 { 94, 6.432, 7.99834E+09, -12.862, 1.993, -37.192, 11.229, -1.057, 0.0325621},
674 { 95, 6.435, 7.99834E+09, -12.869, 1.993, -37.589, 11.363, -1.072, 0.03312393},
675 { 96, 6.449, 1.00000E+10, -12.886, 1.993, -39.573, 12.095, -1.162, 0.03680527},
676 { 97, 6.446, 1.00000E+10, -12.892, 1.993, -40.007, 12.242, -1.178, 0.03737377},
677 { 98, 6.421, 1.00000E+10, -12.887, 1.993, -39.509, 12.041, -1.152, 0.03629023},
678 { 99, 6.414, 1.00000E+10, -12.894, 1.993, -39.939, 12.183, -1.168, 0.03690464},
679 {100, 6.412, 1.00000E+10, -12.900, 1.993, -39.973, 12.180, -1.166, 0.036773}
680  };
void SetOccupancyData()
Definition: G4ShellData.hh:70
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:640
static constexpr double h_Planck
G4bool CheckDeexcitationActiveRegion(G4int coupleIndex)
static G4LossTableManager * Instance()
std::vector< G4Element * > G4ElementVector
G4LivermoreComptonModel(const G4ParticleDefinition *p=0, const G4String &nam="LivermoreCompton")
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4double GetKineticEnergy() const
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:146
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:633
const G4String & GetName() const
Definition: G4Material.hh:178
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
G4double GetZ() const
Definition: G4Element.hh:131
G4bool IsMaster() const
Definition: G4VEmModel.hh:717
G4ParticleDefinition * GetDefinition() const
size_t GetVectorLength() const
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
int G4int
Definition: G4Types.hh:78
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:175
G4int SelectRandomShell(G4int Z) const
Definition: G4ShellData.cc:363
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
static constexpr double twopi
Definition: G4SIunits.hh:76
static constexpr double electron_mass_c2
void LoadData(const G4String &fileName)
Definition: G4ShellData.cc:234
const XML_Char const XML_Char * data
Definition: expat.h:268
static const G4double ln10
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
const XML_Char int const XML_Char * value
Definition: expat.h:331
const G4ThreeVector & GetMomentumDirection() const
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX)
G4double BindingEnergy(G4int Z, G4int shellIndex) const
Definition: G4ShellData.cc:166
static constexpr double cm
Definition: G4SIunits.hh:119
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
static constexpr double eV
Definition: G4SIunits.hh:215
G4double Energy(size_t index) const
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:801
G4double Value(G4double theEnergy, size_t &lastidx) 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:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
static G4ProductionCutsTable * GetProductionCutsTable()
G4int G4Mutex
Definition: G4Threading.hh:173
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:809
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
int G4lrint(double ad)
Definition: templates.hh:163
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
static constexpr double c_light
static constexpr double GeV
Definition: G4SIunits.hh:217
static G4Electron * Electron()
Definition: G4Electron.cc:94
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4VAtomDeexcitation * AtomDeexcitation()
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:780
static constexpr double fine_structure_const
static constexpr double barn
Definition: G4SIunits.hh:105
static constexpr double keV
Definition: G4SIunits.hh:216
double epsilon(double density, double temperature)
G4double RandomSelectMomentum(G4int Z, G4int shellIndex) const
G4AtomicShellEnumerator
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:541
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:132
const G4Material * GetMaterial() const