Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4LivermorePhotoElectricModel.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// Author: Sebastien Incerti
28// 30 October 2008
29// on base of G4LowEnergyPhotoElectric developed by A.Forti and M.G.Pia
30//
31// 22 Oct 2012 A & V Ivanchenko Migration data structure to G4PhysicsVector
32// 1 June 2017 M Bandieramonte: New model based on livermore/epics2014
33// evaluated data - parameterization fits in two ranges
34
36#include "G4SystemOfUnits.hh"
38#include "G4LossTableManager.hh"
39#include "G4Electron.hh"
40#include "G4Gamma.hh"
46#include "G4AtomicShell.hh"
47
48//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
49
50G4PhysicsFreeVector* G4LivermorePhotoElectricModel::fCrossSection[] = {nullptr};
51G4PhysicsFreeVector* G4LivermorePhotoElectricModel::fCrossSectionLE[] = {nullptr};
52std::vector<G4double>* G4LivermorePhotoElectricModel::fParamHigh[] = {nullptr};
53std::vector<G4double>* G4LivermorePhotoElectricModel::fParamLow[] = {nullptr};
54G4int G4LivermorePhotoElectricModel::fNShells[] = {0};
55G4int G4LivermorePhotoElectricModel::fNShellsUsed[] = {0};
56G4ElementData* G4LivermorePhotoElectricModel::fShellCrossSection = nullptr;
57G4Material* G4LivermorePhotoElectricModel::fWater = nullptr;
58G4double G4LivermorePhotoElectricModel::fWaterEnergyLimit = 0.0;
59G4String G4LivermorePhotoElectricModel::fDataDirectory = "";
60
61#ifdef G4MULTITHREADED
62 G4Mutex G4LivermorePhotoElectricModel::livPhotoeffMutex = G4MUTEX_INITIALIZER;
63#endif
64
65using namespace std;
66
67//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
68
70 : G4VEmModel(nam),fParticleChange(nullptr), fAtomDeexcitation(nullptr),
71 maxZ(100),nShellLimit(100),fDeexcitationActive(false),isInitialised(false)
72{
73 verboseLevel= 0;
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 theGamma = G4Gamma::Gamma();
82 theElectron = G4Electron::Electron();
83
84 // default generator
86
87 if(verboseLevel>0) {
88 G4cout << "Livermore PhotoElectric is constructed "
89 << " nShellLimit= " << nShellLimit << G4endl;
90 }
91
92 //Mark this model as "applicable" for atomic deexcitation
94 fSandiaCof.resize(4,0.0);
95 fCurrSection = 0.0;
96}
97
98//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
99
101{
102 if(IsMaster())
103 {
104 delete fShellCrossSection;
105 fShellCrossSection = nullptr;
106 for(G4int i = 0; i <= maxZ; ++i)
107 {
108 if(fParamHigh[i]){
109 delete fParamHigh[i];
110 fParamHigh[i] = nullptr;
111 }
112 if(fParamLow[i]){
113 delete fParamLow[i];
114 fParamLow[i] = nullptr;
115 }
116 if(fCrossSection[i]){
117 delete fCrossSection[i];
118 fCrossSection[i] = nullptr;
119 }
120 if(fCrossSectionLE[i]){
121 delete fCrossSectionLE[i];
122 fCrossSectionLE[i] = nullptr;
123 }
124
125 }
126 }
127}
128
129//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
130
131void
133 const G4DataVector&)
134{
135 if (verboseLevel > 2) {
136 G4cout << "Calling G4LivermorePhotoElectricModel::Initialise() " << G4endl;
137 }
138
139 if(IsMaster()) {
140
141 if(fWater == nullptr) {
142 fWater = G4Material::GetMaterial("G4_WATER", false);
143 if(fWater == nullptr) { fWater = G4Material::GetMaterial("Water", false); }
144 if(fWater) { fWaterEnergyLimit = 13.6*eV; }
145 }
146
147 if(fShellCrossSection == nullptr) { fShellCrossSection = new G4ElementData(); }
148
149 const G4ElementTable* elemTable = G4Element::GetElementTable();
150 std::size_t numElems = (*elemTable).size();
151 for(std::size_t ie = 0; ie < numElems; ++ie)
152 {
153 const G4Element* elem = (*elemTable)[ie];
154 const G4int Z = std::min(maxZ, elem->GetZasInt());
155 if(fCrossSection[Z] == nullptr)
156 {
157 ReadData(Z);
158 }
159 }
160 }
161
162 if (verboseLevel > 2) {
163 G4cout << "Loaded cross section files for new LivermorePhotoElectric model"
164 << G4endl;
165 }
166 if(!isInitialised) {
167 isInitialised = true;
169 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
170 }
171
172 fDeexcitationActive = false;
173 if(fAtomDeexcitation) {
174 fDeexcitationActive = fAtomDeexcitation->IsFluoActive();
175 }
176
177 if (verboseLevel > 0) {
178 G4cout << "LivermorePhotoElectric model is initialized " << G4endl
179 << G4endl;
180 }
181}
182
183//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
184
186 const G4Material* material,
187 const G4ParticleDefinition* p,
188 G4double energy,
190{
191 fCurrSection = 0.0;
192 if(fWater && (material == fWater ||
193 material->GetBaseMaterial() == fWater)) {
194 if(energy <= fWaterEnergyLimit) {
195 fWater->GetSandiaTable()->GetSandiaCofWater(energy, fSandiaCof);
196
197 G4double energy2 = energy*energy;
198 G4double energy3 = energy*energy2;
199 G4double energy4 = energy2*energy2;
200
201 fCurrSection = material->GetDensity()*
202 (fSandiaCof[0]/energy + fSandiaCof[1]/energy2 +
203 fSandiaCof[2]/energy3 + fSandiaCof[3]/energy4);
204 }
205 }
206 if(0.0 == fCurrSection) {
207 fCurrSection = G4VEmModel::CrossSectionPerVolume(material, p, energy);
208 }
209 return fCurrSection;
210}
211
212//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
213
216 G4double energy,
217 G4double ZZ, G4double,
219{
220 if (verboseLevel > 3) {
221 G4cout << "\n G4LivermorePhotoElectricModel::ComputeCrossSectionPerAtom():"
222 << " Z= " << ZZ << " R(keV)= " << energy/keV << G4endl;
223 }
224 G4double cs = 0.0;
225 G4int Z = G4lrint(ZZ);
226 if(Z > maxZ) { return cs; }
227 // if element was not initialised
228
229 // do initialisation safely for MT mode
230 if(fCrossSection[Z] == nullptr) { InitialiseForElement(theGamma, Z); }
231
232 //7: rows in the parameterization file; 5: number of parameters
233 G4int idx = fNShells[Z]*7 - 5;
234
235 energy = std::max(energy, (*(fParamHigh[Z]))[idx-1]);
236
237 G4double x1 = 1.0/energy;
238 G4double x2 = x1*x1;
239 G4double x3 = x2*x1;
240
241 // high energy parameterisation
242 if(energy >= (*(fParamHigh[Z]))[0]) {
243
244 G4double x4 = x2*x2;
245 G4double x5 = x4*x1;
246
247 cs = x1*((*(fParamHigh[Z]))[idx] + x1*(*(fParamHigh[Z]))[idx+1]
248 + x2*(*(fParamHigh[Z]))[idx+2] + x3*(*(fParamHigh[Z]))[idx+3]
249 + x4*(*(fParamHigh[Z]))[idx+4]+ x5*(*(fParamHigh[Z]))[idx+5]);
250
251 }
252 // low energy parameterisation
253 else if(energy >= (*(fParamLow[Z]))[0]) {
254
255 G4double x4 = x2*x2;
256 G4double x5 = x4*x1; //this variable usage can probably be optimized
257 cs = x1*((*(fParamLow[Z]))[idx] + x1*(*(fParamLow[Z]))[idx+1]
258 + x2*(*(fParamLow[Z]))[idx+2] + x3*(*(fParamLow[Z]))[idx+3]
259 + x4*(*(fParamLow[Z]))[idx+4]+ x5*(*(fParamLow[Z]))[idx+5]);
260
261 }
262 // Tabulated values above k-shell ionization energy
263 else if(energy >= (*(fParamHigh[Z]))[1]) {
264 cs = x3*(fCrossSection[Z])->Value(energy);
265 }
266 // Tabulated values below k-shell ionization energy
267 else
268 {
269 cs = x3*(fCrossSectionLE[Z])->Value(energy);
270 }
271 if (verboseLevel > 1) {
272 G4cout << "G4LivermorePhotoElectricModel: E(keV)= " << energy/keV
273 << " Z= " << Z << " cross(barn)= " << cs/barn << G4endl;
274 }
275 return cs;
276}
277
278//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
279
280void
282 std::vector<G4DynamicParticle*>* fvect,
283 const G4MaterialCutsCouple* couple,
284 const G4DynamicParticle* aDynamicGamma,
286{
287 G4double gammaEnergy = aDynamicGamma->GetKineticEnergy();
288 if (verboseLevel > 3) {
289 G4cout << "G4LivermorePhotoElectricModel::SampleSecondaries() Egamma(keV)= "
290 << gammaEnergy/keV << G4endl;
291 }
292
293 // kill incident photon
296
297 // low-energy photo-effect in water - full absorption
298 const G4Material* material = couple->GetMaterial();
299 if(fWater && (material == fWater ||
300 material->GetBaseMaterial() == fWater)) {
301 if(gammaEnergy <= fWaterEnergyLimit) {
303 return;
304 }
305 }
306
307 // Returns the normalized direction of the momentum
308 G4ThreeVector photonDirection = aDynamicGamma->GetMomentumDirection();
309
310 // Select randomly one element in the current material
311 const G4Element* elm = SelectRandomAtom(material, theGamma, gammaEnergy);
312 G4int Z = elm->GetZasInt();
313
314 // Select the ionised shell in the current atom according to shell
315 // cross sections
316 // G4cout << "Select random shell Z= " << Z << G4endl;
317 if(Z > maxZ) { Z = maxZ; }
318
319 // element was not initialised gamma should be absorbed
320 if(fCrossSection[Z] == nullptr) {
322 return;
323 }
324
325 // SAMPLING OF THE SHELL INDEX
326 std::size_t shellIdx = 0;
327 std::size_t nn = fNShellsUsed[Z];
328 if(nn > 1)
329 {
330 if(gammaEnergy >= (*(fParamHigh[Z]))[0])
331 {
332 G4double x1 = 1.0/gammaEnergy;
333 G4double x2 = x1*x1;
334 G4double x3 = x2*x1;
335 G4double x4 = x3*x1;
336 G4double x5 = x4*x1;
337 std::size_t idx = nn*7 - 5;
338 // when do sampling common factors are not taken into account
339 // so cross section is not real
340
341 G4double rand=G4UniformRand();
342 G4double cs0 = rand*( (*(fParamHigh[Z]))[idx]
343 + x1*(*(fParamHigh[Z]))[idx+1]
344 + x2*(*(fParamHigh[Z]))[idx+2]
345 + x3*(*(fParamHigh[Z]))[idx+3]
346 + x4*(*(fParamHigh[Z]))[idx+4]
347 + x5*(*(fParamHigh[Z]))[idx+5]);
348
349 for(shellIdx=0; shellIdx<nn; ++shellIdx)
350 {
351 idx = shellIdx*7 + 2;
352 if(gammaEnergy > (*(fParamHigh[Z]))[idx-1])
353 {
354 G4double cs =
355 (*(fParamHigh[Z]))[idx]
356 + x1*(*(fParamHigh[Z]))[idx+1]
357 + x2*(*(fParamHigh[Z]))[idx+2]
358 + x3*(*(fParamHigh[Z]))[idx+3]
359 + x4*(*(fParamHigh[Z]))[idx+4]
360 + x5*(*(fParamHigh[Z]))[idx+5];
361
362 if(cs >= cs0) { break; }
363 }
364 }
365 if(shellIdx >= nn) { shellIdx = nn-1; }
366 }
367 else if(gammaEnergy >= (*(fParamLow[Z]))[0])
368 {
369 G4double x1 = 1.0/gammaEnergy;
370 G4double x2 = x1*x1;
371 G4double x3 = x2*x1;
372 G4double x4 = x3*x1;
373 G4double x5 = x4*x1;
374 std::size_t idx = nn*7 - 5;
375 // when do sampling common factors are not taken into account
376 // so cross section is not real
377 G4double cs0 = G4UniformRand()*((*(fParamLow[Z]))[idx]
378 + x1*(*(fParamLow[Z]))[idx+1]
379 + x2*(*(fParamLow[Z]))[idx+2]
380 + x3*(*(fParamLow[Z]))[idx+3]
381 + x4*(*(fParamLow[Z]))[idx+4]
382 + x5*(*(fParamLow[Z]))[idx+5]);
383 for(shellIdx=0; shellIdx<nn; ++shellIdx)
384 {
385 idx = shellIdx*7 + 2;
386 if(gammaEnergy > (*(fParamLow[Z]))[idx-1])
387 {
388 G4double cs = (*(fParamLow[Z]))[idx] + x1*(*(fParamLow[Z]))[idx+1]
389 + x2*(*(fParamLow[Z]))[idx+2] + x3*(*(fParamLow[Z]))[idx+3]
390 + x4*(*(fParamLow[Z]))[idx+4]+ x5*(*(fParamLow[Z]))[idx+5];
391 if(cs >= cs0) { break; }
392 }
393 }
394 if(shellIdx >= nn) {shellIdx = nn-1;}
395 }
396 else
397 {
398 // when do sampling common factors are not taken into account
399 // so cross section is not real
400 G4double cs = G4UniformRand();
401
402 if(gammaEnergy >= (*(fParamHigh[Z]))[1]) {
403 //above K-shell binding energy
404 cs*= (fCrossSection[Z])->Value(gammaEnergy);
405 }
406 else
407 {
408 //below K-shell binding energy
409 cs *= (fCrossSectionLE[Z])->Value(gammaEnergy);
410 }
411
412 for(G4int j=0; j<(G4int)nn; ++j) {
413
414 shellIdx = (std::size_t)fShellCrossSection->GetComponentID(Z, j);
415 if(gammaEnergy > (*(fParamLow[Z]))[7*shellIdx+1]) {
416 cs -= fShellCrossSection->GetValueForComponent(Z, j, gammaEnergy);
417 }
418 if(cs <= 0.0 || j+1 == (G4int)nn) {break;}
419 }
420 }
421 }
422 // END: SAMPLING OF THE SHELL
423
424 G4double bindingEnergy = (*(fParamHigh[Z]))[shellIdx*7 + 1];
425 const G4AtomicShell* shell = nullptr;
426
427 // no de-excitation from the last shell
428 if(fDeexcitationActive && shellIdx + 1 < nn) {
430 shell = fAtomDeexcitation->GetAtomicShell(Z, as);
431 }
432
433 // If binding energy of the selected shell is larger than photon energy
434 // do not generate secondaries
435 if(gammaEnergy < bindingEnergy) {
437 return;
438 }
439
440 // Primary outcoming electron
441 G4double eKineticEnergy = gammaEnergy - bindingEnergy;
442 G4double edep = bindingEnergy;
443
444 // Calculate direction of the photoelectron
445 G4ThreeVector electronDirection =
447 eKineticEnergy,
448 (G4int)shellIdx,
449 couple->GetMaterial());
450
451 // The electron is created
452 G4DynamicParticle* electron = new G4DynamicParticle (theElectron,
453 electronDirection,
454 eKineticEnergy);
455 fvect->push_back(electron);
456
457 // Sample deexcitation
458 if(shell) {
459 G4int index = couple->GetIndex();
460 if(fAtomDeexcitation->CheckDeexcitationActiveRegion(index)) {
461 std::size_t nbefore = fvect->size();
462
463 fAtomDeexcitation->GenerateParticles(fvect, shell, Z, index);
464 std::size_t nafter = fvect->size();
465 if(nafter > nbefore) {
466 G4double esec = 0.0;
467 for (std::size_t j=nbefore; j<nafter; ++j) {
468
469 G4double e = ((*fvect)[j])->GetKineticEnergy();
470 if(esec + e > edep) {
471 // correct energy in order to have energy balance
472 e = edep - esec;
473 ((*fvect)[j])->SetKineticEnergy(e);
474 esec += e;
475 // delete the rest of secondaries (should not happens)
476 for (std::size_t jj=nafter-1; jj>j; --jj) {
477 delete (*fvect)[jj];
478 fvect->pop_back();
479 }
480 break;
481 }
482 esec += e;
483 }
484 edep -= esec;
485 }
486 }
487 }
488 // energy balance - excitation energy left
489 if(edep > 0.0) {
491 }
492}
493
494//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
495
496const G4String& G4LivermorePhotoElectricModel::FindDirectoryPath()
497{
498 // check environment variable
499 // build the complete string identifying the file with the data set
500 if(fDataDirectory.empty()) {
501 const char* path = G4FindDataDir("G4LEDATA");
502 if (path) {
503 std::ostringstream ost;
504 if(G4EmParameters::Instance()->LivermoreDataDir() =="livermore"){
505 ost << path << "/livermore/phot_epics2014/";
506 }else{
507 ost << path << "/epics2017/phot/";
508 }
509
510 fDataDirectory = ost.str();
511 } else {
512 G4Exception("G4SeltzerBergerModel::FindDirectoryPath()","em0006",
514 "Environment variable G4LEDATA not defined");
515 }
516 }
517 return fDataDirectory;
518}
519
520//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
521
522void G4LivermorePhotoElectricModel::ReadData(G4int Z)
523{
524 if (verboseLevel > 1)
525 {
526 G4cout << "Calling ReadData() of G4LivermorePhotoElectricModel"
527 << G4endl;
528 }
529
530 if(fCrossSection[Z]!= nullptr) { return; }
531
532 // spline for photoeffect total x-section above K-shell when using EPDL97
533 // but below the parameterized ones
534
535 if(G4EmParameters::Instance()->LivermoreDataDir() =="livermore"){
536 fCrossSection[Z] = new G4PhysicsFreeVector(true);
537 }else{
538 fCrossSection[Z] = new G4PhysicsFreeVector();
539 }
540 std::ostringstream ost;
541 ost << FindDirectoryPath() << "pe-cs-" << Z <<".dat";
542 std::ifstream fin(ost.str().c_str());
543 if( !fin.is_open()) {
545 ed << "G4LivermorePhotoElectricModel data file <" << ost.str().c_str()
546 << "> is not opened!" << G4endl;
547 G4Exception("G4LivermorePhotoElectricModel::ReadData()",
548 "em0003",FatalException,
549 ed,"G4LEDATA version should be G4EMLOW8.0 or later.");
550 return;
551 }
552 if(verboseLevel > 3) {
553 G4cout << "File " << ost.str().c_str()
554 << " is opened by G4LivermorePhotoElectricModel" << G4endl;
555 }
556 fCrossSection[Z]->Retrieve(fin, true);
557 fCrossSection[Z]->ScaleVector(MeV, barn);
558 fCrossSection[Z]->FillSecondDerivatives();
559 fin.close();
560
561 // read high-energy fit parameters
562 fParamHigh[Z] = new std::vector<G4double>;
563 G4int n1 = 0;
564 G4int n2 = 0;
565 G4double x;
566 std::ostringstream ost1;
567 ost1 << fDataDirectory << "pe-high-" << Z <<".dat";
568 std::ifstream fin1(ost1.str().c_str());
569 if( !fin1.is_open()) {
571 ed << "G4LivermorePhotoElectricModel data file <" << ost1.str().c_str()
572 << "> is not opened!" << G4endl;
573 G4Exception("G4LivermorePhotoElectricModel::ReadData()",
574 "em0003",FatalException,
575 ed,"G4LEDATA version should be G4EMLOW7.2 or later.");
576 return;
577 }
578 if(verboseLevel > 3) {
579 G4cout << "File " << ost1.str().c_str()
580 << " is opened by G4LivermorePhotoElectricModel" << G4endl;
581 }
582 fin1 >> n1;
583 if(fin1.fail()) { return; }
584 if(0 > n1 || n1 >= INT_MAX) { n1 = 0; }
585
586 fin1 >> n2;
587 if(fin1.fail()) { return; }
588 if(0 > n2 || n2 >= INT_MAX) { n2 = 0; }
589
590 fin1 >> x;
591 if(fin1.fail()) { return; }
592
593 fNShells[Z] = n1;
594 fParamHigh[Z]->reserve(7*n1+1);
595 fParamHigh[Z]->push_back(x*MeV);
596 for(G4int i=0; i<n1; ++i) {
597 for(G4int j=0; j<7; ++j) {
598 fin1 >> x;
599 if(0 == j) { x *= MeV; }
600 else { x *= barn; }
601 fParamHigh[Z]->push_back(x);
602 }
603 }
604 fin1.close();
605
606 // read low-energy fit parameters
607 fParamLow[Z] = new std::vector<G4double>;
608 G4int n1_low = 0;
609 G4int n2_low = 0;
610 G4double x_low;
611 std::ostringstream ost1_low;
612 ost1_low << fDataDirectory << "pe-low-" << Z <<".dat";
613 std::ifstream fin1_low(ost1_low.str().c_str());
614 if( !fin1_low.is_open()) {
616 ed << "G4LivermorePhotoElectricModel data file <" << ost1_low.str().c_str()
617 << "> is not opened!" << G4endl;
618 G4Exception("G4LivermorePhotoElectricModel::ReadData()",
619 "em0003",FatalException,
620 ed,"G4LEDATA version should be G4EMLOW8.0 or later.");
621 return;
622 }
623 if(verboseLevel > 3) {
624 G4cout << "File " << ost1_low.str().c_str()
625 << " is opened by G4LivermorePhotoElectricModel" << G4endl;
626 }
627 fin1_low >> n1_low;
628 if(fin1_low.fail()) { return; }
629 if(0 > n1_low || n1_low >= INT_MAX) { n1_low = 0; }
630
631 fin1_low >> n2_low;
632 if(fin1_low.fail()) { return; }
633 if(0 > n2_low || n2_low >= INT_MAX) { n2_low = 0; }
634
635 fin1_low >> x_low;
636 if(fin1_low.fail()) { return; }
637
638 fNShells[Z] = n1_low;
639 fParamLow[Z]->reserve(7*n1_low+1);
640 fParamLow[Z]->push_back(x_low*MeV);
641 for(G4int i=0; i<n1_low; ++i) {
642 for(G4int j=0; j<7; ++j) {
643 fin1_low >> x_low;
644 if(0 == j) { x_low *= MeV; }
645 else { x_low *= barn; }
646 fParamLow[Z]->push_back(x_low);
647 }
648 }
649 fin1_low.close();
650
651 // there is a possibility to use only main shells
652 if(nShellLimit < n2) { n2 = nShellLimit; }
653 fShellCrossSection->InitialiseForComponent(Z, n2);//number of shells
654 fNShellsUsed[Z] = n2;
655
656 if(1 < n2) {
657 std::ostringstream ost2;
658 ost2 << fDataDirectory << "pe-ss-cs-" << Z <<".dat";
659 std::ifstream fin2(ost2.str().c_str());
660 if( !fin2.is_open()) {
662 ed << "G4LivermorePhotoElectricModel data file <" << ost2.str().c_str()
663 << "> is not opened!" << G4endl;
664 G4Exception("G4LivermorePhotoElectricModel::ReadData()",
665 "em0003",FatalException,
666 ed,"G4LEDATA version should be G4EMLOW7.2 or later.");
667 return;
668 }
669 if(verboseLevel > 3) {
670 G4cout << "File " << ost2.str().c_str()
671 << " is opened by G4LivermorePhotoElectricModel" << G4endl;
672 }
673
674 G4int n3, n4;
675 G4double y;
676 for(G4int i=0; i<n2; ++i) {
677 fin2 >> x >> y >> n3 >> n4;
678 G4PhysicsFreeVector* v = new G4PhysicsFreeVector(n3, x, y);
679 for(G4int j=0; j<n3; ++j) {
680 fin2 >> x >> y;
681 v->PutValues(j, x*MeV, y*barn);
682 }
683 fShellCrossSection->AddComponent(Z, n4, v);
684 }
685 fin2.close();
686 }
687
688 // no spline for photoeffect total x-section below K-shell
689 if(1 < fNShells[Z]) {
690 fCrossSectionLE[Z] = new G4PhysicsFreeVector();
691 std::ostringstream ost3;
692 ost3 << fDataDirectory << "pe-le-cs-" << Z <<".dat";
693 std::ifstream fin3(ost3.str().c_str());
694 if( !fin3.is_open()) {
696 ed << "G4LivermorePhotoElectricModel data file <" << ost3.str().c_str()
697 << "> is not opened!" << G4endl;
698 G4Exception("G4LivermorePhotoElectricModel::ReadData()",
699 "em0003",FatalException,
700 ed,"G4LEDATA version should be G4EMLOW8.0 or later.");
701 return;
702 }
703 if(verboseLevel > 3) {
704 G4cout << "File " << ost3.str().c_str()
705 << " is opened by G4LivermorePhotoElectricModel" << G4endl;
706 }
707 fCrossSectionLE[Z]->Retrieve(fin3, true);
708 fCrossSectionLE[Z]->ScaleVector(MeV, barn);
709 fin3.close();
710 }
711}
712
713//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
714
716{
717 if(Z < 1 || Z > maxZ) { return -1;} //If Z is out of the supported return 0
718 //If necessary load data for Z
719 InitialiseForElement(theGamma, Z);
720 if(fCrossSection[Z] == nullptr || shell < 0 || shell >= fNShellsUsed[Z]) { return -1; }
721
722 if(Z>2)
723 return fShellCrossSection->GetComponentDataByIndex(Z, shell)->Energy(0);
724 else
725 return fCrossSection[Z]->Energy(0);
726}
727
728//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
729
732{
733 if (fCrossSection[Z] == nullptr) {
734#ifdef G4MULTITHREADED
735 G4MUTEXLOCK(&livPhotoeffMutex);
736 if (fCrossSection[Z] == nullptr) {
737#endif
738 ReadData(Z);
739#ifdef G4MULTITHREADED
740 }
741 G4MUTEXUNLOCK(&livPhotoeffMutex);
742#endif
743 }
744}
745
746//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
G4AtomicShellEnumerator
std::vector< G4Element * > G4ElementTable
const char * G4FindDataDir(const char *)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
#define elem(i, j)
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:85
#define G4MUTEXLOCK(mutex)
Definition: G4Threading.hh:251
#define G4MUTEXUNLOCK(mutex)
Definition: G4Threading.hh:254
std::mutex G4Mutex
Definition: G4Threading.hh:81
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
Definition: G4Electron.cc:93
G4PhysicsVector * GetComponentDataByIndex(G4int Z, G4int idx)
void InitialiseForComponent(G4int Z, G4int nComponents=0)
G4double GetValueForComponent(G4int Z, G4int idx, G4double kinEnergy)
void AddComponent(G4int Z, G4int id, G4PhysicsVector *v)
G4int GetComponentID(G4int Z, G4int idx)
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:403
G4int GetZasInt() const
Definition: G4Element.hh:132
static G4EmParameters * Instance()
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double energy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX) override
G4ParticleChangeForGamma * fParticleChange
G4double GetBindingEnergy(G4int Z, G4int shell)
void InitialiseForElement(const G4ParticleDefinition *, G4int Z) override
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4LivermorePhotoElectricModel(const G4String &nam="LivermorePhElectric")
G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double energy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX) override
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
const G4Material * GetMaterial() const
G4double GetDensity() const
Definition: G4Material.hh:175
const G4Material * GetBaseMaterial() const
Definition: G4Material.hh:228
G4SandiaTable * GetSandiaTable() const
Definition: G4Material.hh:224
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:691
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void PutValues(const std::size_t index, const G4double energy, const G4double value)
void ScaleVector(const G4double factorE, const G4double factorV)
G4double Energy(const std::size_t index) const
G4bool Retrieve(std::ifstream &fIn, G4bool ascii=false)
void FillSecondDerivatives(const G4SplineType=G4SplineType::Base, const G4double dir1=0.0, const G4double dir2=0.0)
void GetSandiaCofWater(G4double energy, std::vector< G4double > &coeff) const
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)
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:600
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:124
G4bool IsMaster() const
Definition: G4VEmModel.hh:725
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:181
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:349
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:561
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:802
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:607
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
int G4lrint(double ad)
Definition: templates.hh:134
#define INT_MAX
Definition: templates.hh:90