Geant4-11
G4PenelopePhotoElectricModel.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: Luciano Pandola
28//
29// History:
30// --------
31// 08 Jan 2010 L Pandola First implementation
32// 01 Feb 2011 L Pandola Suppress fake energy-violation warning when Auger
33// is active.
34// Make sure that fluorescence/Auger is generated
35// only if above threshold
36// 25 May 2011 L Pandola Renamed (make v2008 as default Penelope)
37// 10 Jun 2011 L Pandola Migrate atomic deexcitation interface
38// 07 Oct 2011 L Pandola Bug fix (potential violation of energy conservation)
39// 27 Sep 2013 L Pandola Migrate to MT paradigm, only master model manages
40// tables.
41// 02 Oct 2013 L Pandola Rewrite sampling algorithm of SelectRandomShell()
42// to improve CPU performances
43//
44
47#include "G4SystemOfUnits.hh"
50#include "G4DynamicParticle.hh"
51#include "G4PhysicsTable.hh"
53#include "G4ElementTable.hh"
54#include "G4Element.hh"
56#include "G4AtomicShell.hh"
57#include "G4Gamma.hh"
58#include "G4Electron.hh"
59#include "G4AutoLock.hh"
60#include "G4LossTableManager.hh"
61#include "G4Exp.hh"
62
63//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
64
67
68//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
69
71 const G4String& nam)
72 :G4VEmModel(nam),fParticleChange(nullptr),fParticle(nullptr),
73 fAtomDeexcitation(nullptr),fIsInitialised(false),fLocalTable(false)
74{
77 // SetLowEnergyLimit(fIntrinsicLowEnergyLimit);
79 //
80
81 if (part)
82 SetParticle(part);
83
85 // Verbosity scale:
86 // 0 = nothing
87 // 1 = warning for energy non-conservation
88 // 2 = details of energy budget
89 // 3 = calculation of cross sections, file openings, sampling of atoms
90 // 4 = entering in methods
91
92 //Mark this model as "applicable" for atomic deexcitation
94
96}
97
98//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
99
101{
102 if (IsMaster() || fLocalTable)
103 {
104 for(G4int i=0; i<=fMaxZ; ++i)
105 {
106 if(fLogAtomicShellXS[i]) {
108 delete fLogAtomicShellXS[i];
109 fLogAtomicShellXS[i] = nullptr;
110 }
111 }
112 }
113}
114
115//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
116
118 const G4DataVector& cuts)
119{
120 if (fVerboseLevel > 3)
121 G4cout << "Calling G4PenelopePhotoElectricModel::Initialise()" << G4endl;
122
124 //Issue warning if the AtomicDeexcitation has not been declared
126 {
127 G4cout << G4endl;
128 G4cout << "WARNING from G4PenelopePhotoElectricModel " << G4endl;
129 G4cout << "Atomic de-excitation module is not instantiated, so there will not be ";
130 G4cout << "any fluorescence/Auger emission." << G4endl;
131 G4cout << "Please make sure this is intended" << G4endl;
132 }
133
134 SetParticle(particle);
135
136 //Only the master model creates/fills/destroys the tables
137 if (IsMaster() && particle == fParticle)
138 {
139 G4ProductionCutsTable* theCoupleTable =
141
142 for (size_t i=0;i<theCoupleTable->GetTableSize();i++)
143 {
144 const G4Material* material =
145 theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
146 const G4ElementVector* theElementVector = material->GetElementVector();
147
148 for (size_t j=0;j<material->GetNumberOfElements();j++)
149 {
150 G4int iZ = theElementVector->at(j)->GetZasInt();
151 //read data files only in the master
152 if (!fLogAtomicShellXS[iZ])
153 ReadDataFile(iZ);
154 }
155 }
156
157 InitialiseElementSelectors(particle,cuts);
158
159 if (fVerboseLevel > 0) {
160 G4cout << "Penelope Photo-Electric model v2008 is initialized " << G4endl
161 << "Energy range: "
162 << LowEnergyLimit() / MeV << " MeV - "
163 << HighEnergyLimit() / GeV << " GeV";
164 }
165 }
166
167 if(fIsInitialised) return;
169 fIsInitialised = true;
170
171}
172
174 G4VEmModel *masterModel)
175{
176 if (fVerboseLevel > 3)
177 G4cout << "Calling G4PenelopePhotoElectricModel::InitialiseLocal()" << G4endl;
178 //
179 //Check that particle matches: one might have multiple master models (e.g.
180 //for e+ and e-).
181 //
182 if (part == fParticle)
183 {
185
186 //Get the const table pointers from the master to the workers
187 const G4PenelopePhotoElectricModel* theModel =
188 static_cast<G4PenelopePhotoElectricModel*> (masterModel);
189 for(G4int i=0; i<=fMaxZ; ++i)
190 fLogAtomicShellXS[i] = theModel->fLogAtomicShellXS[i];
191 //Same verbosity for all workers, as the master
192 fVerboseLevel = theModel->fVerboseLevel;
193 }
194
195 return;
196}
197
198//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
205{
206 //
207 // Penelope model v2008
208 //
209 if (fVerboseLevel > 3)
210 G4cout << "Calling ComputeCrossSectionPerAtom() of G4PenelopePhotoElectricModel" << G4endl;
211
212 G4int iZ = G4int(Z);
213
214 if (!fLogAtomicShellXS[iZ])
215 {
216 //If we are here, it means that Initialize() was inkoved, but the MaterialTable was
217 //not filled up. This can happen in a UnitTest or via G4EmCalculator
218 if (fVerboseLevel > 0)
219 {
220 //Issue a G4Exception (warning) only in verbose mode
222 ed << "Unable to retrieve the shell cross section table for Z=" << iZ << G4endl;
223 ed << "This can happen only in Unit Tests or via G4EmCalculator" << G4endl;
224 G4Exception("G4PenelopePhotoElectricModel::ComputeCrossSectionPerAtom()",
225 "em2038",JustWarning,ed);
226 }
227 //protect file reading via autolock
229 ReadDataFile(iZ);
230 lock.unlock();
231 }
232
233 G4double cross = 0;
234 G4PhysicsTable* theTable = fLogAtomicShellXS[iZ];
235 G4PhysicsFreeVector* totalXSLog = (G4PhysicsFreeVector*) (*theTable)[0];
236
237 if (!totalXSLog)
238 {
239 G4Exception("G4PenelopePhotoElectricModel::ComputeCrossSectionPerAtom()",
240 "em2039",FatalException,
241 "Unable to retrieve the total cross section table");
242 return 0;
243 }
244 G4double logene = G4Log(energy);
245 G4double logXS = totalXSLog->Value(logene);
246 cross = G4Exp(logXS);
247
248 if (fVerboseLevel > 2)
249 G4cout << "Photoelectric cross section at " << energy/MeV << " MeV for Z=" << Z <<
250 " = " << cross/barn << " barn" << G4endl;
251 return cross;
252}
253
254//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
255
256void G4PenelopePhotoElectricModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
257 const G4MaterialCutsCouple* couple,
258 const G4DynamicParticle* aDynamicGamma,
259 G4double,
260 G4double)
261{
262 //
263 // Photoelectric effect, Penelope model v2008
264 //
265 // The target atom and the target shell are sampled according to the Livermore
266 // database
267 // D.E. Cullen et al., Report UCRL-50400 (1989)
268 // The angular distribution of the electron in the final state is sampled
269 // according to the Sauter distribution from
270 // F. Sauter, Ann. Phys. 11 (1931) 454
271 // The energy of the final electron is given by the initial photon energy minus
272 // the binding energy. Fluorescence de-excitation is subsequently produced
273 // (to fill the vacancy) according to the general Geant4 G4DeexcitationManager:
274 // J. Stepanek, Comp. Phys. Comm. 1206 pp 1-1-9 (1997)
275
276 if (fVerboseLevel > 3)
277 G4cout << "Calling SamplingSecondaries() of G4PenelopePhotoElectricModel" << G4endl;
278
279 G4double photonEnergy = aDynamicGamma->GetKineticEnergy();
280
281 // always kill primary
284
285 if (photonEnergy <= fIntrinsicLowEnergyLimit)
286 {
288 return ;
289 }
290
291 G4ParticleMomentum photonDirection = aDynamicGamma->GetMomentumDirection();
292
293 // Select randomly one element in the current material
294 if (fVerboseLevel > 2)
295 G4cout << "Going to select element in " << couple->GetMaterial()->GetName() << G4endl;
296
297 // atom can be selected efficiently if element selectors are initialised
298 const G4Element* anElement =
299 SelectRandomAtom(couple,G4Gamma::GammaDefinition(),photonEnergy);
300 G4int Z = anElement->GetZasInt();
301 if (fVerboseLevel > 2)
302 G4cout << "Selected " << anElement->GetName() << G4endl;
303
304 // Select the ionised shell in the current atom according to shell cross sections
305 //shellIndex = 0 --> K shell
306 // 1-3 --> L shells
307 // 4-8 --> M shells
308 // 9 --> outer shells cumulatively
309 //
310 size_t shellIndex = SelectRandomShell(Z,photonEnergy);
311
312 if (fVerboseLevel > 2)
313 G4cout << "Selected shell " << shellIndex << " of element " << anElement->GetName() << G4endl;
314
315 // Retrieve the corresponding identifier and binding energy of the selected shell
317
318 //The number of shell cross section possibly reported in the Penelope database
319 //might be different from the number of shells in the G4AtomicTransitionManager
320 //(namely, Penelope may contain more shell, especially for very light elements).
321 //In order to avoid a warning message from the G4AtomicTransitionManager, I
322 //add this protection. Results are anyway changed, because when G4AtomicTransitionManager
323 //has a shellID>maxID, it sets the shellID to the last valid shell.
324 size_t numberOfShells = (size_t) transitionManager->NumberOfShells(Z);
325 if (shellIndex >= numberOfShells)
326 shellIndex = numberOfShells-1;
327
328 const G4AtomicShell* shell = fTransitionManager->Shell(Z,shellIndex);
330
331 //Penelope considers only K, L and M shells. Cross sections of outer shells are
332 //not included in the Penelope database. If SelectRandomShell() returns
333 //shellIndex = 9, it means that an outer shell was ionized. In this case the
334 //Penelope recipe is to set bindingEnergy = 0 (the energy is entirely assigned
335 //to the electron) and to disregard fluorescence.
336 if (shellIndex == 9)
337 bindingEnergy = 0.*eV;
338
339 G4double localEnergyDeposit = 0.0;
340 G4double cosTheta = 1.0;
341
342 // Primary outcoming electron
343 G4double eKineticEnergy = photonEnergy - bindingEnergy;
344
345 // There may be cases where the binding energy of the selected shell is > photon energy
346 // In such cases do not generate secondaries
347 if (eKineticEnergy > 0.)
348 {
349 // The electron is created
350 // Direction sampled from the Sauter distribution
351 cosTheta = SampleElectronDirection(eKineticEnergy);
352 G4double sinTheta = std::sqrt(1-cosTheta*cosTheta);
353 G4double phi = twopi * G4UniformRand() ;
354 G4double dirx = sinTheta * std::cos(phi);
355 G4double diry = sinTheta * std::sin(phi);
356 G4double dirz = cosTheta ;
357 G4ThreeVector electronDirection(dirx,diry,dirz); //electron direction
358 electronDirection.rotateUz(photonDirection);
360 electronDirection,
361 eKineticEnergy);
362 fvect->push_back(electron);
363 }
364 else
365 bindingEnergy = photonEnergy;
366
367 G4double energyInFluorescence = 0; //testing purposes
368 G4double energyInAuger = 0; //testing purposes
369
370 //Now, take care of fluorescence, if required. According to the Penelope
371 //recipe, I have to skip fluoresence completely if shellIndex == 9
372 //(= sampling of a shell outer than K,L,M)
373 if (fAtomDeexcitation && shellIndex<9)
374 {
375 G4int index = couple->GetIndex();
377 {
378 size_t nBefore = fvect->size();
379 fAtomDeexcitation->GenerateParticles(fvect,shell,Z,index);
380 size_t nAfter = fvect->size();
381
382 if (nAfter > nBefore) //actual production of fluorescence
383 {
384 for (size_t j=nBefore;j<nAfter;j++) //loop on products
385 {
386 G4double itsEnergy = ((*fvect)[j])->GetKineticEnergy();
387 if (itsEnergy < bindingEnergy) // valid secondary, generate it
388 {
389 bindingEnergy -= itsEnergy;
390 if (((*fvect)[j])->GetParticleDefinition() == G4Gamma::Definition())
391 energyInFluorescence += itsEnergy;
392 else if (((*fvect)[j])->GetParticleDefinition() == G4Electron::Definition())
393 energyInAuger += itsEnergy;
394 }
395 else //invalid secondary: takes more than the available energy: delete it
396 {
397 delete (*fvect)[j];
398 (*fvect)[j] = nullptr;
399 }
400 }
401 }
402 }
403 }
404
405 //Residual energy is deposited locally
406 localEnergyDeposit += bindingEnergy;
407
408 if (localEnergyDeposit < 0) //Should not be: issue a G4Exception (warning)
409 {
410 G4Exception("G4PenelopePhotoElectricModel::SampleSecondaries()",
411 "em2099",JustWarning,"WARNING: Negative local energy deposit");
412 localEnergyDeposit = 0;
413 }
414
415 fParticleChange->ProposeLocalEnergyDeposit(localEnergyDeposit);
416
417 if (fVerboseLevel > 1)
418 {
419 G4cout << "-----------------------------------------------------------" << G4endl;
420 G4cout << "Energy balance from G4PenelopePhotoElectric" << G4endl;
421 G4cout << "Selected shell: " << WriteTargetShell(shellIndex) << " of element " <<
422 anElement->GetName() << G4endl;
423 G4cout << "Incoming photon energy: " << photonEnergy/keV << " keV" << G4endl;
424 G4cout << "-----------------------------------------------------------" << G4endl;
425 if (eKineticEnergy)
426 G4cout << "Outgoing electron " << eKineticEnergy/keV << " keV" << G4endl;
427 if (energyInFluorescence)
428 G4cout << "Fluorescence x-rays: " << energyInFluorescence/keV << " keV" << G4endl;
429 if (energyInAuger)
430 G4cout << "Auger electrons: " << energyInAuger/keV << " keV" << G4endl;
431 G4cout << "Local energy deposit " << localEnergyDeposit/keV << " keV" << G4endl;
432 G4cout << "Total final state: " <<
433 (eKineticEnergy+energyInFluorescence+localEnergyDeposit+energyInAuger)/keV <<
434 " keV" << G4endl;
435 G4cout << "-----------------------------------------------------------" << G4endl;
436 }
437 if (fVerboseLevel > 0)
438 {
439 G4double energyDiff =
440 std::fabs(eKineticEnergy+energyInFluorescence+localEnergyDeposit+energyInAuger-photonEnergy);
441 if (energyDiff > 0.05*keV)
442 {
443 G4cout << "Warning from G4PenelopePhotoElectric: problem with energy conservation: " <<
444 (eKineticEnergy+energyInFluorescence+localEnergyDeposit+energyInAuger)/keV
445 << " keV (final) vs. " <<
446 photonEnergy/keV << " keV (initial)" << G4endl;
447 G4cout << "-----------------------------------------------------------" << G4endl;
448 G4cout << "Energy balance from G4PenelopePhotoElectric" << G4endl;
449 G4cout << "Selected shell: " << WriteTargetShell(shellIndex) << " of element " <<
450 anElement->GetName() << G4endl;
451 G4cout << "Incoming photon energy: " << photonEnergy/keV << " keV" << G4endl;
452 G4cout << "-----------------------------------------------------------" << G4endl;
453 if (eKineticEnergy)
454 G4cout << "Outgoing electron " << eKineticEnergy/keV << " keV" << G4endl;
455 if (energyInFluorescence)
456 G4cout << "Fluorescence x-rays: " << energyInFluorescence/keV << " keV" << G4endl;
457 if (energyInAuger)
458 G4cout << "Auger electrons: " << energyInAuger/keV << " keV" << G4endl;
459 G4cout << "Local energy deposit " << localEnergyDeposit/keV << " keV" << G4endl;
460 G4cout << "Total final state: " <<
461 (eKineticEnergy+energyInFluorescence+localEnergyDeposit+energyInAuger)/keV <<
462 " keV" << G4endl;
463 G4cout << "-----------------------------------------------------------" << G4endl;
464 }
465 }
466}
467
468//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
469
471{
472 G4double costheta = 1.0;
473 if (energy>1*GeV) return costheta;
474
475 //1) initialize energy-dependent variables
476 // Variable naming according to Eq. (2.24) of Penelope Manual
477 // (pag. 44)
478 G4double gamma = 1.0 + energy/electron_mass_c2;
479 G4double gamma2 = gamma*gamma;
480 G4double beta = std::sqrt((gamma2-1.0)/gamma2);
481
482 // ac corresponds to "A" of Eq. (2.31)
483 //
484 G4double ac = (1.0/beta) - 1.0;
485 G4double a1 = 0.5*beta*gamma*(gamma-1.0)*(gamma-2.0);
486 G4double a2 = ac + 2.0;
487 G4double gtmax = 2.0*(a1 + 1.0/ac);
488
489 G4double tsam = 0;
490 G4double gtr = 0;
491
492 //2) sampling. Eq. (2.31) of Penelope Manual
493 // tsam = 1-std::cos(theta)
494 // gtr = rejection function according to Eq. (2.28)
495 do{
496 G4double rand = G4UniformRand();
497 tsam = 2.0*ac * (2.0*rand + a2*std::sqrt(rand)) / (a2*a2 - 4.0*rand);
498 gtr = (2.0 - tsam) * (a1 + 1.0/(ac+tsam));
499 }while(G4UniformRand()*gtmax > gtr);
500 costheta = 1.0-tsam;
501
502 return costheta;
503}
504
505//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
506
508{
509 if (!IsMaster())
510 //Should not be here!
511 G4Exception("G4PenelopePhotoElectricModel::ReadDataFile()",
512 "em0100",FatalException,"Worker thread in this method");
513
514 if (fVerboseLevel > 2)
515 {
516 G4cout << "G4PenelopePhotoElectricModel::ReadDataFile()" << G4endl;
517 G4cout << "Going to read PhotoElectric data files for Z=" << Z << G4endl;
518 }
519
520 char* path = std::getenv("G4LEDATA");
521 if (!path)
522 {
523 G4String excep = "G4PenelopePhotoElectricModel - G4LEDATA environment variable not set!";
524 G4Exception("G4PenelopePhotoElectricModel::ReadDataFile()",
525 "em0006",FatalException,excep);
526 return;
527 }
528
529 /*
530 Read the cross section file
531 */
532 std::ostringstream ost;
533 if (Z>9)
534 ost << path << "/penelope/photoelectric/pdgph" << Z << ".p08";
535 else
536 ost << path << "/penelope/photoelectric/pdgph0" << Z << ".p08";
537 std::ifstream file(ost.str().c_str());
538 if (!file.is_open())
539 {
540 G4String excep = "G4PenelopePhotoElectricModel - data file " + G4String(ost.str()) + " not found!";
541 G4Exception("G4PenelopePhotoElectricModel::ReadDataFile()",
542 "em0003",FatalException,excep);
543 }
544 //I have to know in advance how many points are in the data list
545 //to initialize the G4PhysicsFreeVector()
546 size_t ndata=0;
547 G4String line;
548 while( getline(file, line) )
549 ndata++;
550 ndata -= 1;
551 //G4cout << "Found: " << ndata << " lines" << G4endl;
552
553 file.clear();
554 file.close();
555 file.open(ost.str().c_str());
556
557 G4int readZ =0;
558 size_t nShells= 0;
559 file >> readZ >> nShells;
560
561 if (fVerboseLevel > 3)
562 G4cout << "Element Z=" << Z << " , nShells = " << nShells << G4endl;
563
564 //check the right file is opened.
565 if (readZ != Z || nShells <= 0 || nShells > 50) //protect nShell against large values
566 {
568 ed << "Corrupted data file for Z=" << Z << G4endl;
569 G4Exception("G4PenelopePhotoElectricModel::ReadDataFile()",
570 "em0005",FatalException,ed);
571 return;
572 }
573 G4PhysicsTable* thePhysicsTable = new G4PhysicsTable();
574
575 //the table has to contain nShell+1 G4PhysicsFreeVectors,
576 //(theTable)[0] --> total cross section
577 //(theTable)[ishell] --> cross section for shell (ishell-1)
578
579 //reserve space for the vectors
580 //everything is log-log
581 for (size_t i=0;i<nShells+1;i++)
582 thePhysicsTable->push_back(new G4PhysicsFreeVector(ndata));
583
584 size_t k =0;
585 for (k=0;k<ndata && !file.eof();k++)
586 {
587 G4double energy = 0;
588 G4double aValue = 0;
589 file >> energy ;
590 energy *= eV;
591 G4double logene = G4Log(energy);
592 //loop on the columns
593 for (size_t i=0;i<nShells+1;i++)
594 {
595 file >> aValue;
596 aValue *= barn;
597 G4PhysicsFreeVector* theVec = (G4PhysicsFreeVector*) ((*thePhysicsTable)[i]);
598 if (aValue < 1e-40*cm2) //protection against log(0)
599 aValue = 1e-40*cm2;
600 theVec->PutValue(k,logene,G4Log(aValue));
601 }
602 }
603
604 if (fVerboseLevel > 2)
605 {
606 G4cout << "G4PenelopePhotoElectricModel: read " << k << " points for element Z = "
607 << Z << G4endl;
608 }
609
610 fLogAtomicShellXS[Z] = thePhysicsTable;
611
612 file.close();
613 return;
614}
615
616//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
617
619{
620 if (!IsMaster())
621 //Should not be here!
622 G4Exception("G4PenelopePhotoElectricModel::GetNumberOfShellXS()",
623 "em0100",FatalException,"Worker thread in this method");
624
625 //read data files
626 if (!fLogAtomicShellXS[Z])
628 //now it should be ok
629 if (!fLogAtomicShellXS[Z])
630 {
632 ed << "Cannot find shell cross section data for Z=" << Z << G4endl;
633 G4Exception("G4PenelopePhotoElectricModel::GetNumberOfShellXS()",
634 "em2038",FatalException,ed);
635 }
636 //one vector is allocated for the _total_ cross section
637 size_t nEntries = fLogAtomicShellXS[Z]->entries();
638 return (nEntries-1);
639}
640
641//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
642
644{
645 //this forces also the loading of the data
646 size_t entries = GetNumberOfShellXS(Z);
647
648 if (shellID >= entries)
649 {
650 G4cout << "Element Z=" << Z << " has data for " << entries << " shells only" << G4endl;
651 G4cout << "so shellID should be from 0 to " << entries-1 << G4endl;
652 return 0;
653 }
654
656 //[0] is the total XS, shellID is in the element [shellID+1]
657 G4PhysicsFreeVector* totalXSLog = (G4PhysicsFreeVector*) (*theTable)[shellID+1];
658
659 if (!totalXSLog)
660 {
661 G4Exception("G4PenelopePhotoElectricModel::GetShellCrossSection()",
662 "em2039",FatalException,
663 "Unable to retrieve the total cross section table");
664 return 0;
665 }
666 G4double logene = G4Log(energy);
667 G4double logXS = totalXSLog->Value(logene);
668 G4double cross = G4Exp(logXS);
669 if (cross < 2e-40*cm2) cross = 0;
670 return cross;
671}
672
673//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
674
676{
677 G4String theShell = "outer shell";
678 if (shellID == 0)
679 theShell = "K";
680 else if (shellID == 1)
681 theShell = "L1";
682 else if (shellID == 2)
683 theShell = "L2";
684 else if (shellID == 3)
685 theShell = "L3";
686 else if (shellID == 4)
687 theShell = "M1";
688 else if (shellID == 5)
689 theShell = "M2";
690 else if (shellID == 6)
691 theShell = "M3";
692 else if (shellID == 7)
693 theShell = "M4";
694 else if (shellID == 8)
695 theShell = "M5";
696
697 return theShell;
698}
699
700//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...
701
703{
704 if(!fParticle) {
705 fParticle = p;
706 }
707}
708
709//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...
710
712{
713 G4double logEnergy = G4Log(energy);
714
715 //Check if data have been read (it should be!)
716 if (!fLogAtomicShellXS[Z])
717 {
719 ed << "Cannot find shell cross section data for Z=" << Z << G4endl;
720 G4Exception("G4PenelopePhotoElectricModel::SelectRandomShell()",
721 "em2038",FatalException,ed);
722 }
723
725
726 G4double sum = 0;
727 G4PhysicsFreeVector* totalXSLog = (G4PhysicsFreeVector*) (*theTable)[0];
728 G4double logXS = totalXSLog->Value(logEnergy);
729 G4double totalXS = G4Exp(logXS);
730
731 //Notice: totalXS is the total cross section and it does *not* correspond to
732 //the sum of partialXS's, since these include only K, L and M shells.
733 //
734 // Therefore, here one have to consider the possibility of ionisation of
735 // an outer shell. Conventionally, it is indicated with id=10 in Penelope
736 //
737 G4double random = G4UniformRand()*totalXS;
738
739 for (size_t k=1;k<theTable->entries();k++)
740 {
741 //Add one shell
742 G4PhysicsFreeVector* partialXSLog = (G4PhysicsFreeVector*) (*theTable)[k];
743 G4double logXSLocal = partialXSLog->Value(logEnergy);
744 G4double partialXS = G4Exp(logXSLocal);
745 sum += partialXS;
746 if (random <= sum)
747 return k-1;
748 }
749 //none of the shells K, L, M: return outer shell
750 return 9;
751}
std::vector< const G4Element * > G4ElementVector
@ JustWarning
@ 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
static constexpr double twopi
Definition: G4SIunits.hh:56
static constexpr double barn
Definition: G4SIunits.hh:85
static constexpr double cm2
Definition: G4SIunits.hh:100
static constexpr double keV
Definition: G4SIunits.hh:202
static constexpr double eV
Definition: G4SIunits.hh:201
static constexpr double GeV
Definition: G4SIunits.hh:203
static constexpr double MeV
Definition: G4SIunits.hh:200
#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
const G4int Z[17]
#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 BindingEnergy() const
G4AtomicShell * Shell(G4int Z, size_t shellIndex) const
static G4AtomicTransitionManager * Instance()
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
static G4Electron * Definition()
Definition: G4Electron.cc:48
static G4Electron * Electron()
Definition: G4Electron.cc:93
const G4String & GetName() const
Definition: G4Element.hh:127
G4int GetZasInt() const
Definition: G4Element.hh:132
static G4Gamma * GammaDefinition()
Definition: G4Gamma.cc:80
static G4Gamma * Definition()
Definition: G4Gamma.cc:48
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
const G4Material * GetMaterial() const
const G4String & GetName() const
Definition: G4Material.hh:173
void SetProposedKineticEnergy(G4double proposedKinEnergy)
const G4AtomicTransitionManager * fTransitionManager
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4PenelopePhotoElectricModel(const G4ParticleDefinition *p=nullptr, const G4String &processName="PenPhotoElec")
G4double SampleElectronDirection(G4double energy)
G4double GetShellCrossSection(G4int Z, size_t shellID, G4double energy)
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX) override
G4ParticleChangeForGamma * fParticleChange
size_t SelectRandomShell(G4int Z, G4double energy)
static G4PhysicsTable * fLogAtomicShellXS[fMaxZ+1]
void SetParticle(const G4ParticleDefinition *)
const G4ParticleDefinition * fParticle
void PutValue(const std::size_t index, const G4double e, const G4double value)
void push_back(G4PhysicsVector *)
void clearAndDestroy()
std::size_t entries() const
G4double Value(const G4double energy, std::size_t &lastidx) const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
G4bool CheckDeexcitationActiveRegion(G4int coupleIndex)
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:767
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:852
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:123
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:662
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:844
G4bool IsMaster() const
Definition: G4VEmModel.hh:746
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:655
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:580
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:823
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:138
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4double energy(const ThreeVector &p, const G4double m)
G4double bindingEnergy(G4int A, G4int Z)
string material
Definition: eplot.py:19
float electron_mass_c2
Definition: hepunit.py:273