Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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// $Id$
27//
28// Author: Luciano Pandola
29//
30// History:
31// --------
32// 08 Jan 2010 L Pandola First implementation
33// 01 Feb 2011 L Pandola Suppress fake energy-violation warning when Auger is active.
34// Make sure that fluorescence/Auger is generated only if
35// 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//
40
43#include "G4SystemOfUnits.hh"
46#include "G4DynamicParticle.hh"
47#include "G4PhysicsTable.hh"
49#include "G4ElementTable.hh"
50#include "G4Element.hh"
52#include "G4AtomicShell.hh"
53#include "G4Gamma.hh"
54#include "G4Electron.hh"
55#include "G4LossTableManager.hh"
56
57//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
58
59
61 const G4String& nam)
62 :G4VEmModel(nam),fParticleChange(0),isInitialised(false),fAtomDeexcitation(0),
63 logAtomicShellXS(0)
64{
65 fIntrinsicLowEnergyLimit = 100.0*eV;
66 fIntrinsicHighEnergyLimit = 100.0*GeV;
67 // SetLowEnergyLimit(fIntrinsicLowEnergyLimit);
68 SetHighEnergyLimit(fIntrinsicHighEnergyLimit);
69 //
70 verboseLevel= 0;
71 // Verbosity scale:
72 // 0 = nothing
73 // 1 = warning for energy non-conservation
74 // 2 = details of energy budget
75 // 3 = calculation of cross sections, file openings, sampling of atoms
76 // 4 = entering in methods
77
78 //Mark this model as "applicable" for atomic deexcitation
80
81 fTransitionManager = G4AtomicTransitionManager::Instance();
82}
83
84//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
85
87{
88 std::map <G4int,G4PhysicsTable*>::iterator i;
89 if (logAtomicShellXS)
90 {
91 for (i=logAtomicShellXS->begin();i != logAtomicShellXS->end();i++)
92 {
93 G4PhysicsTable* tab = i->second;
94 tab->clearAndDestroy();
95 delete tab;
96 }
97 }
98 delete logAtomicShellXS;
99}
100
101//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
102
104 const G4DataVector& cuts)
105{
106 if (verboseLevel > 3)
107 G4cout << "Calling G4PenelopePhotoElectricModel::Initialise()" << G4endl;
108
109 // logAtomicShellXS is created only once, since it is never cleared
110 if (!logAtomicShellXS)
111 logAtomicShellXS = new std::map<G4int,G4PhysicsTable*>;
112
113 InitialiseElementSelectors(particle,cuts);
114
115 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
116 //Issue warning if the AtomicDeexcitation has not been declared
117 if (!fAtomDeexcitation)
118 {
119 G4cout << G4endl;
120 G4cout << "WARNING from G4PenelopePhotoElectricModel " << G4endl;
121 G4cout << "Atomic de-excitation module is not instantiated, so there will not be ";
122 G4cout << "any fluorescence/Auger emission." << G4endl;
123 G4cout << "Please make sure this is intended" << G4endl;
124 }
125
126 if (verboseLevel > 0) {
127 G4cout << "Penelope Photo-Electric model v2008 is initialized " << G4endl
128 << "Energy range: "
129 << LowEnergyLimit() / MeV << " MeV - "
130 << HighEnergyLimit() / GeV << " GeV";
131 }
132
133 if(isInitialised) return;
135 isInitialised = true;
136
137}
138
139//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
140
143 G4double energy,
146{
147 //
148 // Penelope model v2008
149 //
150
151 if (verboseLevel > 3)
152 G4cout << "Calling ComputeCrossSectionPerAtom() of G4PenelopePhotoElectricModel" << G4endl;
153
154 G4int iZ = (G4int) Z;
155
156 //read data files
157 if (!logAtomicShellXS->count(iZ))
158 ReadDataFile(iZ);
159 //now it should be ok
160 if (!logAtomicShellXS->count(iZ))
161 G4Exception("G4PenelopePhotoElectricModel::ComputeCrossSectionPerAtom()",
162 "em2038",FatalException,
163 "Unable to retrieve the shell cross section table");
164
165 G4double cross = 0;
166
167 G4PhysicsTable* theTable = logAtomicShellXS->find(iZ)->second;
168 G4PhysicsFreeVector* totalXSLog = (G4PhysicsFreeVector*) (*theTable)[0];
169
170 if (!totalXSLog)
171 {
172 G4Exception("G4PenelopePhotoElectricModel::ComputeCrossSectionPerAtom()",
173 "em2039",FatalException,
174 "Unable to retrieve the total cross section table");
175 return 0;
176 }
177 G4double logene = std::log(energy);
178 G4double logXS = totalXSLog->Value(logene);
179 cross = std::exp(logXS);
180
181 if (verboseLevel > 2)
182 G4cout << "Photoelectric cross section at " << energy/MeV << " MeV for Z=" << Z <<
183 " = " << cross/barn << " barn" << G4endl;
184 return cross;
185}
186
187//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
188
189void G4PenelopePhotoElectricModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
190 const G4MaterialCutsCouple* couple,
191 const G4DynamicParticle* aDynamicGamma,
192 G4double,
193 G4double)
194{
195 //
196 // Photoelectric effect, Penelope model v2008
197 //
198 // The target atom and the target shell are sampled according to the Livermore
199 // database
200 // D.E. Cullen et al., Report UCRL-50400 (1989)
201 // The angular distribution of the electron in the final state is sampled
202 // according to the Sauter distribution from
203 // F. Sauter, Ann. Phys. 11 (1931) 454
204 // The energy of the final electron is given by the initial photon energy minus
205 // the binding energy. Fluorescence de-excitation is subsequently produced
206 // (to fill the vacancy) according to the general Geant4 G4DeexcitationManager:
207 // J. Stepanek, Comp. Phys. Comm. 1206 pp 1-1-9 (1997)
208
209 if (verboseLevel > 3)
210 G4cout << "Calling SamplingSecondaries() of G4PenelopePhotoElectricModel" << G4endl;
211
212 G4double photonEnergy = aDynamicGamma->GetKineticEnergy();
213
214 // always kill primary
217
218 if (photonEnergy <= fIntrinsicLowEnergyLimit)
219 {
221 return ;
222 }
223
224 G4ParticleMomentum photonDirection = aDynamicGamma->GetMomentumDirection();
225
226 // Select randomly one element in the current material
227 if (verboseLevel > 2)
228 G4cout << "Going to select element in " << couple->GetMaterial()->GetName() << G4endl;
229
230 // atom can be selected efficiently if element selectors are initialised
231 const G4Element* anElement =
232 SelectRandomAtom(couple,G4Gamma::GammaDefinition(),photonEnergy);
233 G4int Z = (G4int) anElement->GetZ();
234 if (verboseLevel > 2)
235 G4cout << "Selected " << anElement->GetName() << G4endl;
236
237 // Select the ionised shell in the current atom according to shell cross sections
238 //shellIndex = 0 --> K shell
239 // 1-3 --> L shells
240 // 4-8 --> M shells
241 // 9 --> outer shells cumulatively
242 //
243 size_t shellIndex = SelectRandomShell(Z,photonEnergy);
244
245 if (verboseLevel > 2)
246 G4cout << "Selected shell " << shellIndex << " of element " << anElement->GetName() << G4endl;
247
248 // Retrieve the corresponding identifier and binding energy of the selected shell
250
251 //The number of shell cross section possibly reported in the Penelope database
252 //might be different from the number of shells in the G4AtomicTransitionManager
253 //(namely, Penelope may contain more shell, especially for very light elements).
254 //In order to avoid a warning message from the G4AtomicTransitionManager, I
255 //add this protection. Results are anyway changed, because when G4AtomicTransitionManager
256 //has a shellID>maxID, it sets the shellID to the last valid shell.
257 size_t numberOfShells = (size_t) transitionManager->NumberOfShells(Z);
258 if (shellIndex >= numberOfShells)
259 shellIndex = numberOfShells-1;
260
261 const G4AtomicShell* shell = fTransitionManager->Shell(Z,shellIndex);
262 G4double bindingEnergy = shell->BindingEnergy();
263 //G4int shellId = shell->ShellId();
264
265 //Penelope considers only K, L and M shells. Cross sections of outer shells are
266 //not included in the Penelope database. If SelectRandomShell() returns
267 //shellIndex = 9, it means that an outer shell was ionized. In this case the
268 //Penelope recipe is to set bindingEnergy = 0 (the energy is entirely assigned
269 //to the electron) and to disregard fluorescence.
270 if (shellIndex == 9)
271 bindingEnergy = 0.*eV;
272
273
274 G4double localEnergyDeposit = 0.0;
275 G4double cosTheta = 1.0;
276
277 // Primary outcoming electron
278 G4double eKineticEnergy = photonEnergy - bindingEnergy;
279
280 // There may be cases where the binding energy of the selected shell is > photon energy
281 // In such cases do not generate secondaries
282 if (eKineticEnergy > 0.)
283 {
284 // The electron is created
285 // Direction sampled from the Sauter distribution
286 cosTheta = SampleElectronDirection(eKineticEnergy);
287 G4double sinTheta = std::sqrt(1-cosTheta*cosTheta);
288 G4double phi = twopi * G4UniformRand() ;
289 G4double dirx = sinTheta * std::cos(phi);
290 G4double diry = sinTheta * std::sin(phi);
291 G4double dirz = cosTheta ;
292 G4ThreeVector electronDirection(dirx,diry,dirz); //electron direction
293 electronDirection.rotateUz(photonDirection);
295 electronDirection,
296 eKineticEnergy);
297 fvect->push_back(electron);
298 }
299 else
300 bindingEnergy = photonEnergy;
301
302
303 G4double energyInFluorescence = 0; //testing purposes
304 G4double energyInAuger = 0; //testing purposes
305
306 //Now, take care of fluorescence, if required. According to the Penelope
307 //recipe, I have to skip fluoresence completely if shellIndex == 9
308 //(= sampling of a shell outer than K,L,M)
309 if (fAtomDeexcitation && shellIndex<9)
310 {
311 G4int index = couple->GetIndex();
312 if (fAtomDeexcitation->CheckDeexcitationActiveRegion(index))
313 {
314 size_t nBefore = fvect->size();
315 fAtomDeexcitation->GenerateParticles(fvect,shell,Z,index);
316 size_t nAfter = fvect->size();
317
318 if (nAfter > nBefore) //actual production of fluorescence
319 {
320 for (size_t j=nBefore;j<nAfter;j++) //loop on products
321 {
322 G4double itsEnergy = ((*fvect)[j])->GetKineticEnergy();
323 bindingEnergy -= itsEnergy;
324 if (((*fvect)[j])->GetParticleDefinition() == G4Gamma::Definition())
325 energyInFluorescence += itsEnergy;
326 else if (((*fvect)[j])->GetParticleDefinition() == G4Electron::Definition())
327 energyInAuger += itsEnergy;
328
329 }
330 }
331
332 }
333 }
334
335 //Residual energy is deposited locally
336 localEnergyDeposit += bindingEnergy;
337
338 if (localEnergyDeposit < 0)
339 {
340 G4cout << "WARNING - "
341 << "G4PenelopePhotoElectricModel::SampleSecondaries() - Negative energy deposit"
342 << G4endl;
343 localEnergyDeposit = 0;
344 }
345
346 fParticleChange->ProposeLocalEnergyDeposit(localEnergyDeposit);
347
348 if (verboseLevel > 1)
349 {
350 G4cout << "-----------------------------------------------------------" << G4endl;
351 G4cout << "Energy balance from G4PenelopePhotoElectric" << G4endl;
352 G4cout << "Selected shell: " << WriteTargetShell(shellIndex) << " of element " <<
353 anElement->GetName() << G4endl;
354 G4cout << "Incoming photon energy: " << photonEnergy/keV << " keV" << G4endl;
355 G4cout << "-----------------------------------------------------------" << G4endl;
356 if (eKineticEnergy)
357 G4cout << "Outgoing electron " << eKineticEnergy/keV << " keV" << G4endl;
358 if (energyInFluorescence)
359 G4cout << "Fluorescence x-rays: " << energyInFluorescence/keV << " keV" << G4endl;
360 if (energyInAuger)
361 G4cout << "Auger electrons: " << energyInAuger/keV << " keV" << G4endl;
362 G4cout << "Local energy deposit " << localEnergyDeposit/keV << " keV" << G4endl;
363 G4cout << "Total final state: " <<
364 (eKineticEnergy+energyInFluorescence+localEnergyDeposit+energyInAuger)/keV <<
365 " keV" << G4endl;
366 G4cout << "-----------------------------------------------------------" << G4endl;
367 }
368 if (verboseLevel > 0)
369 {
370 G4double energyDiff =
371 std::fabs(eKineticEnergy+energyInFluorescence+localEnergyDeposit+energyInAuger-photonEnergy);
372 if (energyDiff > 0.05*keV)
373 {
374 G4cout << "Warning from G4PenelopePhotoElectric: problem with energy conservation: " <<
375 (eKineticEnergy+energyInFluorescence+localEnergyDeposit+energyInAuger)/keV
376 << " keV (final) vs. " <<
377 photonEnergy/keV << " keV (initial)" << G4endl;
378 G4cout << "-----------------------------------------------------------" << G4endl;
379 G4cout << "Energy balance from G4PenelopePhotoElectric" << G4endl;
380 G4cout << "Selected shell: " << WriteTargetShell(shellIndex) << " of element " <<
381 anElement->GetName() << G4endl;
382 G4cout << "Incoming photon energy: " << photonEnergy/keV << " keV" << G4endl;
383 G4cout << "-----------------------------------------------------------" << G4endl;
384 if (eKineticEnergy)
385 G4cout << "Outgoing electron " << eKineticEnergy/keV << " keV" << G4endl;
386 if (energyInFluorescence)
387 G4cout << "Fluorescence x-rays: " << energyInFluorescence/keV << " keV" << G4endl;
388 if (energyInAuger)
389 G4cout << "Auger electrons: " << energyInAuger/keV << " keV" << G4endl;
390 G4cout << "Local energy deposit " << localEnergyDeposit/keV << " keV" << G4endl;
391 G4cout << "Total final state: " <<
392 (eKineticEnergy+energyInFluorescence+localEnergyDeposit+energyInAuger)/keV <<
393 " keV" << G4endl;
394 G4cout << "-----------------------------------------------------------" << G4endl;
395 }
396 }
397}
398
399//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
400
401G4double G4PenelopePhotoElectricModel::SampleElectronDirection(G4double energy)
402{
403 G4double costheta = 1.0;
404 if (energy>1*GeV) return costheta;
405
406 //1) initialize energy-dependent variables
407 // Variable naming according to Eq. (2.24) of Penelope Manual
408 // (pag. 44)
409 G4double gamma = 1.0 + energy/electron_mass_c2;
410 G4double gamma2 = gamma*gamma;
411 G4double beta = std::sqrt((gamma2-1.0)/gamma2);
412
413 // ac corresponds to "A" of Eq. (2.31)
414 //
415 G4double ac = (1.0/beta) - 1.0;
416 G4double a1 = 0.5*beta*gamma*(gamma-1.0)*(gamma-2.0);
417 G4double a2 = ac + 2.0;
418 G4double gtmax = 2.0*(a1 + 1.0/ac);
419
420 G4double tsam = 0;
421 G4double gtr = 0;
422
423 //2) sampling. Eq. (2.31) of Penelope Manual
424 // tsam = 1-std::cos(theta)
425 // gtr = rejection function according to Eq. (2.28)
426 do{
427 G4double rand = G4UniformRand();
428 tsam = 2.0*ac * (2.0*rand + a2*std::sqrt(rand)) / (a2*a2 - 4.0*rand);
429 gtr = (2.0 - tsam) * (a1 + 1.0/(ac+tsam));
430 }while(G4UniformRand()*gtmax > gtr);
431 costheta = 1.0-tsam;
432
433
434 return costheta;
435}
436
437//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
438
439void G4PenelopePhotoElectricModel::ReadDataFile(G4int Z)
440{
441 if (verboseLevel > 2)
442 {
443 G4cout << "G4PenelopePhotoElectricModel::ReadDataFile()" << G4endl;
444 G4cout << "Going to read PhotoElectric data files for Z=" << Z << G4endl;
445 }
446
447 char* path = getenv("G4LEDATA");
448 if (!path)
449 {
450 G4String excep = "G4PenelopePhotoElectricModel - G4LEDATA environment variable not set!";
451 G4Exception("G4PenelopePhotoElectricModel::ReadDataFile()",
452 "em0006",FatalException,excep);
453 return;
454 }
455
456 /*
457 Read the cross section file
458 */
459 std::ostringstream ost;
460 if (Z>9)
461 ost << path << "/penelope/photoelectric/pdgph" << Z << ".p08";
462 else
463 ost << path << "/penelope/photoelectric/pdgph0" << Z << ".p08";
464 std::ifstream file(ost.str().c_str());
465 if (!file.is_open())
466 {
467 G4String excep = "G4PenelopePhotoElectricModel - data file " + G4String(ost.str()) + " not found!";
468 G4Exception("G4PenelopePhotoElectricModel::ReadDataFile()",
469 "em0003",FatalException,excep);
470 }
471 //I have to know in advance how many points are in the data list
472 //to initialize the G4PhysicsFreeVector()
473 size_t ndata=0;
474 G4String line;
475 while( getline(file, line) )
476 ndata++;
477 ndata -= 1;
478 //G4cout << "Found: " << ndata << " lines" << G4endl;
479
480 file.clear();
481 file.close();
482 file.open(ost.str().c_str());
483
484 G4int readZ =0;
485 size_t nShells= 0;
486 file >> readZ >> nShells;
487
488 if (verboseLevel > 3)
489 G4cout << "Element Z=" << Z << " , nShells = " << nShells << G4endl;
490
491 //check the right file is opened.
492 if (readZ != Z || nShells <= 0 || nShells > 50) //protect nShell against large values
493 {
495 ed << "Corrupted data file for Z=" << Z << G4endl;
496 G4Exception("G4PenelopePhotoElectricModel::ReadDataFile()",
497 "em0005",FatalException,ed);
498 return;
499 }
500 G4PhysicsTable* thePhysicsTable = new G4PhysicsTable();
501
502 //the table has to contain nShell+1 G4PhysicsFreeVectors,
503 //(theTable)[0] --> total cross section
504 //(theTable)[ishell] --> cross section for shell (ishell-1)
505
506 //reserve space for the vectors
507 //everything is log-log
508 for (size_t i=0;i<nShells+1;i++)
509 thePhysicsTable->push_back(new G4PhysicsFreeVector(ndata));
510
511 size_t k =0;
512 for (k=0;k<ndata && !file.eof();k++)
513 {
514 G4double energy = 0;
515 G4double aValue = 0;
516 file >> energy ;
517 energy *= eV;
518 G4double logene = std::log(energy);
519 //loop on the columns
520 for (size_t i=0;i<nShells+1;i++)
521 {
522 file >> aValue;
523 aValue *= barn;
524 G4PhysicsFreeVector* theVec = (G4PhysicsFreeVector*) ((*thePhysicsTable)[i]);
525 if (aValue < 1e-40*cm2) //protection against log(0)
526 aValue = 1e-40*cm2;
527 theVec->PutValue(k,logene,std::log(aValue));
528 }
529 }
530
531 if (verboseLevel > 2)
532 {
533 G4cout << "G4PenelopePhotoElectricModel: read " << k << " points for element Z = "
534 << Z << G4endl;
535 }
536
537 logAtomicShellXS->insert(std::make_pair(Z,thePhysicsTable));
538
539 file.close();
540 return;
541}
542
543//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
544
545size_t G4PenelopePhotoElectricModel::SelectRandomShell(G4int Z,G4double energy)
546{
547 G4double logEnergy = std::log(energy);
548
549 //Check if data have been read (it should be!)
550 if (!logAtomicShellXS->count(Z))
551 {
553 ed << "Cannot find shell cross section data for Z=" << Z << G4endl;
554 G4Exception("G4PenelopePhotoElectricModel::SelectRandomShell()",
555 "em2038",FatalException,ed);
556 }
557
558 size_t shellIndex = 0;
559
560 G4PhysicsTable* theTable = logAtomicShellXS->find(Z)->second;
561
562 G4DataVector* tempVector = new G4DataVector();
563
564 G4double sum = 0;
565 //loop on shell partial XS, retrieve the value for the given energy and store on
566 //a temporary vector
567 tempVector->push_back(sum); //first element is zero
568
569 G4PhysicsFreeVector* totalXSLog = (G4PhysicsFreeVector*) (*theTable)[0];
570 G4double logXS = totalXSLog->Value(logEnergy);
571 G4double totalXS = std::exp(logXS);
572
573 //Notice: totalXS is the total cross section and it does *not* correspond to
574 //the sum of partialXS's, since these include only K, L and M shells.
575 //
576 // Therefore, here one have to consider the possibility of ionisation of
577 // an outer shell. Conventionally, it is indicated with id=10 in Penelope
578 //
579
580 for (size_t k=1;k<theTable->entries();k++)
581 {
582 G4PhysicsFreeVector* partialXSLog = (G4PhysicsFreeVector*) (*theTable)[k];
583 G4double logXSLocal = partialXSLog->Value(logEnergy);
584 G4double partialXS = std::exp(logXSLocal);
585 sum += partialXS;
586 tempVector->push_back(sum);
587 }
588
589 tempVector->push_back(totalXS); //last element
590
591 G4double random = G4UniformRand()*totalXS;
592
593 /*
594 for (size_t i=0;i<tempVector->size(); i++)
595 G4cout << i << " " << (*tempVector)[i]/totalXS << G4endl;
596 */
597
598 //locate bin of tempVector
599 //Now one has to sample according to the elements in tempVector
600 //This gives the left edge of the interval...
601 size_t lowerBound = 0;
602 size_t upperBound = tempVector->size()-1;
603 while (lowerBound <= upperBound)
604 {
605 size_t midBin = (lowerBound + upperBound)/2;
606 if( random < (*tempVector)[midBin])
607 upperBound = midBin-1;
608 else
609 lowerBound = midBin+1;
610 }
611
612 shellIndex = upperBound;
613
614 delete tempVector;
615 return shellIndex;
616}
617
618//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
619
621{
622 //read data files
623 if (!logAtomicShellXS->count(Z))
624 ReadDataFile(Z);
625 //now it should be ok
626 if (!logAtomicShellXS->count(Z))
627 {
629 ed << "Cannot find shell cross section data for Z=" << Z << G4endl;
630 G4Exception("G4PenelopePhotoElectricModel::GetNumberOfShellXS()",
631 "em2038",FatalException,ed);
632 }
633 //one vector is allocated for the _total_ cross section
634 size_t nEntries = logAtomicShellXS->find(Z)->second->entries();
635 return (nEntries-1);
636}
637
638//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
639
641{
642 //this forces also the loading of the data
643 size_t entries = GetNumberOfShellXS(Z);
644
645 if (shellID >= entries)
646 {
647 G4cout << "Element Z=" << Z << " has data for " << entries << " shells only" << G4endl;
648 G4cout << "so shellID should be from 0 to " << entries-1 << G4endl;
649 return 0;
650 }
651
652 G4PhysicsTable* theTable = logAtomicShellXS->find(Z)->second;
653 //[0] is the total XS, shellID is in the element [shellID+1]
654 G4PhysicsFreeVector* totalXSLog = (G4PhysicsFreeVector*) (*theTable)[shellID+1];
655
656 if (!totalXSLog)
657 {
658 G4Exception("G4PenelopePhotoElectricModel::GetShellCrossSection()",
659 "em2039",FatalException,
660 "Unable to retrieve the total cross section table");
661 return 0;
662 }
663 G4double logene = std::log(energy);
664 G4double logXS = totalXSLog->Value(logene);
665 G4double cross = std::exp(logXS);
666 if (cross < 2e-40*cm2) cross = 0;
667 return cross;
668}
669
670//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
671
672G4String G4PenelopePhotoElectricModel::WriteTargetShell(size_t shellID)
673{
674 G4String theShell = "outer shell";
675 if (shellID == 0)
676 theShell = "K";
677 else if (shellID == 1)
678 theShell = "L1";
679 else if (shellID == 2)
680 theShell = "L2";
681 else if (shellID == 3)
682 theShell = "L3";
683 else if (shellID == 4)
684 theShell = "M1";
685 else if (shellID == 5)
686 theShell = "M2";
687 else if (shellID == 6)
688 theShell = "M3";
689 else if (shellID == 7)
690 theShell = "M4";
691 else if (shellID == 8)
692 theShell = "M5";
693
694 return theShell;
695}
@ FatalException
@ fStopAndKill
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
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:49
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4double GetZ() const
Definition: G4Element.hh:131
const G4String & GetName() const
Definition: G4Element.hh:127
static G4Gamma * GammaDefinition()
Definition: G4Gamma.cc:81
static G4Gamma * Definition()
Definition: G4Gamma.cc:49
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
const G4Material * GetMaterial() const
const G4String & GetName() const
Definition: G4Material.hh:177
void SetProposedKineticEnergy(G4double proposedKinEnergy)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4PenelopePhotoElectricModel(const G4ParticleDefinition *p=0, const G4String &processName="PenPhotoElec")
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX)
G4double GetShellCrossSection(G4int Z, size_t shellID, G4double energy)
G4ParticleChangeForGamma * fParticleChange
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
void PutValue(size_t binNumber, G4double binValue, G4double dataValue)
void push_back(G4PhysicsVector *)
void clearAndDestroy()
size_t entries() const
G4double Value(G4double theEnergy)
G4bool CheckDeexcitationActiveRegion(G4int coupleIndex)
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:109
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:529
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:522
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:459
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:641
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:123
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76