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