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