Geant4  10.02.p03
G4PenelopePhotoElectricModel Class Reference

#include <G4PenelopePhotoElectricModel.hh>

Inheritance diagram for G4PenelopePhotoElectricModel:
Collaboration diagram for G4PenelopePhotoElectricModel:

Public Member Functions

 G4PenelopePhotoElectricModel (const G4ParticleDefinition *p=0, const G4String &processName="PenPhotoElec")
 
virtual ~G4PenelopePhotoElectricModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle *> *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
void SetVerbosityLevel (G4int lev)
 
G4int GetVerbosityLevel ()
 
size_t GetNumberOfShellXS (G4int)
 
G4double GetShellCrossSection (G4int Z, size_t shellID, G4double energy)
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int, const G4ParticleDefinition *, G4double)
 
virtual G4double ComputeCrossSectionPerShell (const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ChargeSquareRatio (const G4Track &)
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void StartTracking (G4Track *)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector *> *)
 
G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectIsotopeNumber (const G4Element *)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool LPMFlag () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy)
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 

Protected Attributes

G4ParticleChangeForGamma * fParticleChange
 
const G4ParticleDefinitionfParticle
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChange * pParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 

Private Member Functions

G4PenelopePhotoElectricModeloperator= (const G4PenelopePhotoElectricModel &right)
 
 G4PenelopePhotoElectricModel (const G4PenelopePhotoElectricModel &)
 
void SetParticle (const G4ParticleDefinition *)
 
G4double SampleElectronDirection (G4double energy)
 
void ReadDataFile (G4int Z)
 
size_t SelectRandomShell (G4int Z, G4double energy)
 
G4String WriteTargetShell (size_t shellID)
 

Private Attributes

G4double fIntrinsicLowEnergyLimit
 
G4double fIntrinsicHighEnergyLimit
 
G4int verboseLevel
 
G4bool isInitialised
 
G4VAtomDeexcitationfAtomDeexcitation
 
const G4AtomicTransitionManagerfTransitionManager
 
std::map< G4int, G4PhysicsTable * > * logAtomicShellXS
 
G4bool fLocalTable
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLoss * GetParticleChangeForLoss ()
 
G4ParticleChangeForGamma * GetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 
- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

Definition at line 59 of file G4PenelopePhotoElectricModel.hh.

Constructor & Destructor Documentation

◆ G4PenelopePhotoElectricModel() [1/2]

G4PenelopePhotoElectricModel::G4PenelopePhotoElectricModel ( const G4ParticleDefinition p = 0,
const G4String processName = "PenPhotoElec" 
)

Definition at line 66 of file G4PenelopePhotoElectricModel.cc.

70 {
73  // SetLowEnergyLimit(fIntrinsicLowEnergyLimit);
75  //
76 
77  if (part)
79 
80  verboseLevel= 0;
81  // Verbosity scale:
82  // 0 = nothing
83  // 1 = warning for energy non-conservation
84  // 2 = details of energy budget
85  // 3 = calculation of cross sections, file openings, sampling of atoms
86  // 4 = entering in methods
87 
88  //Mark this model as "applicable" for atomic deexcitation
89  SetDeexcitationFlag(true);
90 
92 }
G4ParticleChangeForGamma * fParticleChange
std::map< G4int, G4PhysicsTable * > * logAtomicShellXS
const G4AtomicTransitionManager * fTransitionManager
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:69
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:725
TString part[npart]
static const double GeV
Definition: G4SIunits.hh:214
static const double eV
Definition: G4SIunits.hh:212
void SetParticle(const G4ParticleDefinition *)
const G4ParticleDefinition * fParticle
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:781
static G4AtomicTransitionManager * Instance()
Here is the call graph for this function:

◆ ~G4PenelopePhotoElectricModel()

G4PenelopePhotoElectricModel::~G4PenelopePhotoElectricModel ( )
virtual

Definition at line 96 of file G4PenelopePhotoElectricModel.cc.

97 {
98  if (IsMaster() || fLocalTable)
99  {
100  std::map <G4int,G4PhysicsTable*>::iterator i;
101  if (logAtomicShellXS)
102  {
103  for (i=logAtomicShellXS->begin();i != logAtomicShellXS->end();i++)
104  {
105  G4PhysicsTable* tab = i->second;
106  //tab->clearAndDestroy();
107  delete tab;
108  }
109  }
110  delete logAtomicShellXS;
111  }
112 }
G4bool IsMaster() const
Definition: G4VEmModel.hh:718
std::map< G4int, G4PhysicsTable * > * logAtomicShellXS
std::ostream & tab(std::ostream &)
Definition: CCalutils.cc:89
Here is the call graph for this function:

◆ G4PenelopePhotoElectricModel() [2/2]

G4PenelopePhotoElectricModel::G4PenelopePhotoElectricModel ( const G4PenelopePhotoElectricModel )
private

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

G4double G4PenelopePhotoElectricModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition ,
G4double  kinEnergy,
G4double  Z,
G4double  A = 0,
G4double  cut = 0,
G4double  emax = DBL_MAX 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 207 of file G4PenelopePhotoElectricModel.cc.

212 {
213  //
214  // Penelope model v2008
215  //
216 
217  if (verboseLevel > 3)
218  G4cout << "Calling ComputeCrossSectionPerAtom() of G4PenelopePhotoElectricModel" << G4endl;
219 
220  G4int iZ = (G4int) Z;
221 
222  //Either Initialize() was not called, or we are in a slave and InitializeLocal() was
223  //not invoked
224  if (!logAtomicShellXS)
225  {
226  //create a **thread-local** version of the table. Used only for G4EmCalculator and
227  //Unit Tests
228  fLocalTable = true;
229  logAtomicShellXS = new std::map<G4int,G4PhysicsTable*>;
230  }
231 
232  //now it should be ok
233  if (!logAtomicShellXS->count(iZ))
234  {
235  //If we are here, it means that Initialize() was inkoved, but the MaterialTable was
236  //not filled up. This can happen in a UnitTest or via G4EmCalculator
237  if (verboseLevel > 0)
238  {
239  //Issue a G4Exception (warning) only in verbose mode
241  ed << "Unable to retrieve the shell cross section table for Z=" << iZ << G4endl;
242  ed << "This can happen only in Unit Tests or via G4EmCalculator" << G4endl;
243  G4Exception("G4PenelopePhotoElectricModel::ComputeCrossSectionPerAtom()",
244  "em2038",JustWarning,ed);
245  }
246  //protect file reading via autolock
247  G4AutoLock lock(&PenelopePhotoElectricModelMutex);
248  ReadDataFile(iZ);
249  lock.unlock();
250 
251  }
252 
253  G4double cross = 0;
254 
255  G4PhysicsTable* theTable = logAtomicShellXS->find(iZ)->second;
256  G4PhysicsFreeVector* totalXSLog = (G4PhysicsFreeVector*) (*theTable)[0];
257 
258  if (!totalXSLog)
259  {
260  G4Exception("G4PenelopePhotoElectricModel::ComputeCrossSectionPerAtom()",
261  "em2039",FatalException,
262  "Unable to retrieve the total cross section table");
263  return 0;
264  }
265  G4double logene = std::log(energy);
266  G4double logXS = totalXSLog->Value(logene);
267  cross = G4Exp(logXS);
268 
269  if (verboseLevel > 2)
270  G4cout << "Photoelectric cross section at " << energy/MeV << " MeV for Z=" << Z <<
271  " = " << cross/barn << " barn" << G4endl;
272  return cross;
273 }
static const double MeV
Definition: G4SIunits.hh:211
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
int G4int
Definition: G4Types.hh:78
std::map< G4int, G4PhysicsTable * > * logAtomicShellXS
G4GLOB_DLL std::ostream G4cout
double energy
Definition: plottest35.C:25
Float_t Z
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
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
#define G4endl
Definition: G4ios.hh:61
static const double barn
Definition: G4SIunits.hh:104
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ GetNumberOfShellXS()

size_t G4PenelopePhotoElectricModel::GetNumberOfShellXS ( G4int  Z)

Definition at line 636 of file G4PenelopePhotoElectricModel.cc.

637 {
638  if (!IsMaster())
639  //Should not be here!
640  G4Exception("G4PenelopePhotoElectricModel::GetNumberOfShellXS()",
641  "em0100",FatalException,"Worker thread in this method");
642 
643  //read data files
644  if (!logAtomicShellXS->count(Z))
645  ReadDataFile(Z);
646  //now it should be ok
647  if (!logAtomicShellXS->count(Z))
648  {
650  ed << "Cannot find shell cross section data for Z=" << Z << G4endl;
651  G4Exception("G4PenelopePhotoElectricModel::GetNumberOfShellXS()",
652  "em2038",FatalException,ed);
653  }
654  //one vector is allocated for the _total_ cross section
655  size_t nEntries = logAtomicShellXS->find(Z)->second->entries();
656  return (nEntries-1);
657 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4bool IsMaster() const
Definition: G4VEmModel.hh:718
std::map< G4int, G4PhysicsTable * > * logAtomicShellXS
Float_t Z
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetShellCrossSection()

G4double G4PenelopePhotoElectricModel::GetShellCrossSection ( G4int  Z,
size_t  shellID,
G4double  energy 
)

Definition at line 661 of file G4PenelopePhotoElectricModel.cc.

662 {
663  //this forces also the loading of the data
664  size_t entries = GetNumberOfShellXS(Z);
665 
666  if (shellID >= entries)
667  {
668  G4cout << "Element Z=" << Z << " has data for " << entries << " shells only" << G4endl;
669  G4cout << "so shellID should be from 0 to " << entries-1 << G4endl;
670  return 0;
671  }
672 
673  G4PhysicsTable* theTable = logAtomicShellXS->find(Z)->second;
674  //[0] is the total XS, shellID is in the element [shellID+1]
675  G4PhysicsFreeVector* totalXSLog = (G4PhysicsFreeVector*) (*theTable)[shellID+1];
676 
677  if (!totalXSLog)
678  {
679  G4Exception("G4PenelopePhotoElectricModel::GetShellCrossSection()",
680  "em2039",FatalException,
681  "Unable to retrieve the total cross section table");
682  return 0;
683  }
684  G4double logene = std::log(energy);
685  G4double logXS = totalXSLog->Value(logene);
686  G4double cross = G4Exp(logXS);
687  if (cross < 2e-40*cm2) cross = 0;
688  return cross;
689 }
static const double cm2
Definition: G4SIunits.hh:119
std::map< G4int, G4PhysicsTable * > * logAtomicShellXS
G4GLOB_DLL std::ostream G4cout
double energy
Definition: plottest35.C:25
Float_t Z
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
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetVerbosityLevel()

G4int G4PenelopePhotoElectricModel::GetVerbosityLevel ( )
inline

Definition at line 88 of file G4PenelopePhotoElectricModel.hh.

Here is the call graph for this function:

◆ Initialise()

void G4PenelopePhotoElectricModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector cuts 
)
virtual

Implements G4VEmModel.

Definition at line 116 of file G4PenelopePhotoElectricModel.cc.

118 {
119  if (verboseLevel > 3)
120  G4cout << "Calling G4PenelopePhotoElectricModel::Initialise()" << G4endl;
121 
123  //Issue warning if the AtomicDeexcitation has not been declared
124  if (!fAtomDeexcitation)
125  {
126  G4cout << G4endl;
127  G4cout << "WARNING from G4PenelopePhotoElectricModel " << G4endl;
128  G4cout << "Atomic de-excitation module is not instantiated, so there will not be ";
129  G4cout << "any fluorescence/Auger emission." << G4endl;
130  G4cout << "Please make sure this is intended" << G4endl;
131  }
132 
133  SetParticle(particle);
134 
135  //Only the master model creates/fills/destroys the tables
136  if (IsMaster() && particle == fParticle)
137  {
138 
139  // logAtomicShellXS is created only once, since it is never cleared
140  if (!logAtomicShellXS)
141  logAtomicShellXS = new std::map<G4int,G4PhysicsTable*>;
142 
143  G4ProductionCutsTable* theCoupleTable =
145 
146  for (size_t i=0;i<theCoupleTable->GetTableSize();i++)
147  {
148  const G4Material* material =
149  theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
150  const G4ElementVector* theElementVector = material->GetElementVector();
151 
152  for (size_t j=0;j<material->GetNumberOfElements();j++)
153  {
154  G4int iZ = (G4int) theElementVector->at(j)->GetZ();
155  //read data files only in the master
156  if (!logAtomicShellXS->count(iZ))
157  ReadDataFile(iZ);
158  }
159  }
160 
161 
162  InitialiseElementSelectors(particle,cuts);
163 
164  if (verboseLevel > 0) {
165  G4cout << "Penelope Photo-Electric model v2008 is initialized " << G4endl
166  << "Energy range: "
167  << LowEnergyLimit() / MeV << " MeV - "
168  << HighEnergyLimit() / GeV << " GeV";
169  }
170  }
171 
172  if(isInitialised) return;
174  isInitialised = true;
175 
176 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
static const double MeV
Definition: G4SIunits.hh:211
static G4LossTableManager * Instance()
std::vector< G4Element * > G4ElementVector
const G4Material * GetMaterial() const
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:148
G4ParticleChangeForGamma * fParticleChange
G4bool IsMaster() const
Definition: G4VEmModel.hh:718
int G4int
Definition: G4Types.hh:78
std::map< G4int, G4PhysicsTable * > * logAtomicShellXS
string material
Definition: eplot.py:19
G4GLOB_DLL std::ostream G4cout
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
static const double GeV
Definition: G4SIunits.hh:214
static G4ProductionCutsTable * GetProductionCutsTable()
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
void SetParticle(const G4ParticleDefinition *)
const G4ParticleDefinition * fParticle
#define G4endl
Definition: G4ios.hh:61
G4VAtomDeexcitation * AtomDeexcitation()
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:134
Here is the call graph for this function:

◆ InitialiseLocal()

void G4PenelopePhotoElectricModel::InitialiseLocal ( const G4ParticleDefinition part,
G4VEmModel masterModel 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 178 of file G4PenelopePhotoElectricModel.cc.

180 {
181  if (verboseLevel > 3)
182  G4cout << "Calling G4PenelopePhotoElectricModel::InitialiseLocal()" << G4endl;
183 
184  //
185  //Check that particle matches: one might have multiple master models (e.g.
186  //for e+ and e-).
187  //
188  if (part == fParticle)
189  {
191 
192  //Get the const table pointers from the master to the workers
193  const G4PenelopePhotoElectricModel* theModel =
194  static_cast<G4PenelopePhotoElectricModel*> (masterModel);
195 
196  logAtomicShellXS = theModel->logAtomicShellXS;
197 
198  //Same verbosity for all workers, as the master
199  verboseLevel = theModel->verboseLevel;
200  }
201 
202  return;
203 }
void SetElementSelectors(std::vector< G4EmElementSelector *> *)
Definition: G4VEmModel.hh:810
std::map< G4int, G4PhysicsTable * > * logAtomicShellXS
G4GLOB_DLL std::ostream G4cout
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:802
const G4ParticleDefinition * fParticle
#define G4endl
Definition: G4ios.hh:61
Here is the call graph for this function:

◆ operator=()

G4PenelopePhotoElectricModel& G4PenelopePhotoElectricModel::operator= ( const G4PenelopePhotoElectricModel right)
private

◆ ReadDataFile()

void G4PenelopePhotoElectricModel::ReadDataFile ( G4int  Z)
private

Definition at line 525 of file G4PenelopePhotoElectricModel.cc.

526 {
527  if (!IsMaster())
528  //Should not be here!
529  G4Exception("G4PenelopePhotoElectricModel::ReadDataFile()",
530  "em0100",FatalException,"Worker thread in this method");
531 
532  if (verboseLevel > 2)
533  {
534  G4cout << "G4PenelopePhotoElectricModel::ReadDataFile()" << G4endl;
535  G4cout << "Going to read PhotoElectric data files for Z=" << Z << G4endl;
536  }
537 
538  char* path = getenv("G4LEDATA");
539  if (!path)
540  {
541  G4String excep = "G4PenelopePhotoElectricModel - G4LEDATA environment variable not set!";
542  G4Exception("G4PenelopePhotoElectricModel::ReadDataFile()",
543  "em0006",FatalException,excep);
544  return;
545  }
546 
547  /*
548  Read the cross section file
549  */
550  std::ostringstream ost;
551  if (Z>9)
552  ost << path << "/penelope/photoelectric/pdgph" << Z << ".p08";
553  else
554  ost << path << "/penelope/photoelectric/pdgph0" << Z << ".p08";
555  std::ifstream file(ost.str().c_str());
556  if (!file.is_open())
557  {
558  G4String excep = "G4PenelopePhotoElectricModel - data file " + G4String(ost.str()) + " not found!";
559  G4Exception("G4PenelopePhotoElectricModel::ReadDataFile()",
560  "em0003",FatalException,excep);
561  }
562  //I have to know in advance how many points are in the data list
563  //to initialize the G4PhysicsFreeVector()
564  size_t ndata=0;
565  G4String line;
566  while( getline(file, line) )
567  ndata++;
568  ndata -= 1;
569  //G4cout << "Found: " << ndata << " lines" << G4endl;
570 
571  file.clear();
572  file.close();
573  file.open(ost.str().c_str());
574 
575  G4int readZ =0;
576  size_t nShells= 0;
577  file >> readZ >> nShells;
578 
579  if (verboseLevel > 3)
580  G4cout << "Element Z=" << Z << " , nShells = " << nShells << G4endl;
581 
582  //check the right file is opened.
583  if (readZ != Z || nShells <= 0 || nShells > 50) //protect nShell against large values
584  {
586  ed << "Corrupted data file for Z=" << Z << G4endl;
587  G4Exception("G4PenelopePhotoElectricModel::ReadDataFile()",
588  "em0005",FatalException,ed);
589  return;
590  }
591  G4PhysicsTable* thePhysicsTable = new G4PhysicsTable();
592 
593  //the table has to contain nShell+1 G4PhysicsFreeVectors,
594  //(theTable)[0] --> total cross section
595  //(theTable)[ishell] --> cross section for shell (ishell-1)
596 
597  //reserve space for the vectors
598  //everything is log-log
599  for (size_t i=0;i<nShells+1;i++)
600  thePhysicsTable->push_back(new G4PhysicsFreeVector(ndata));
601 
602  size_t k =0;
603  for (k=0;k<ndata && !file.eof();k++)
604  {
605  G4double energy = 0;
606  G4double aValue = 0;
607  file >> energy ;
608  energy *= eV;
609  G4double logene = std::log(energy);
610  //loop on the columns
611  for (size_t i=0;i<nShells+1;i++)
612  {
613  file >> aValue;
614  aValue *= barn;
615  G4PhysicsFreeVector* theVec = (G4PhysicsFreeVector*) ((*thePhysicsTable)[i]);
616  if (aValue < 1e-40*cm2) //protection against log(0)
617  aValue = 1e-40*cm2;
618  theVec->PutValue(k,logene,std::log(aValue));
619  }
620  }
621 
622  if (verboseLevel > 2)
623  {
624  G4cout << "G4PenelopePhotoElectricModel: read " << k << " points for element Z = "
625  << Z << G4endl;
626  }
627 
628  logAtomicShellXS->insert(std::make_pair(Z,thePhysicsTable));
629 
630  file.close();
631  return;
632 }
static const double cm2
Definition: G4SIunits.hh:119
void PutValue(size_t binNumber, G4double binValue, G4double dataValue)
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
void push_back(G4PhysicsVector *)
G4bool IsMaster() const
Definition: G4VEmModel.hh:718
TFile * file
int G4int
Definition: G4Types.hh:78
std::map< G4int, G4PhysicsTable * > * logAtomicShellXS
G4GLOB_DLL std::ostream G4cout
double energy
Definition: plottest35.C:25
Float_t Z
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static const double eV
Definition: G4SIunits.hh:212
#define G4endl
Definition: G4ios.hh:61
static const double barn
Definition: G4SIunits.hh:104
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SampleElectronDirection()

G4double G4PenelopePhotoElectricModel::SampleElectronDirection ( G4double  energy)
private

Definition at line 487 of file G4PenelopePhotoElectricModel.cc.

488 {
489  G4double costheta = 1.0;
490  if (energy>1*GeV) return costheta;
491 
492  //1) initialize energy-dependent variables
493  // Variable naming according to Eq. (2.24) of Penelope Manual
494  // (pag. 44)
495  G4double gamma = 1.0 + energy/electron_mass_c2;
496  G4double gamma2 = gamma*gamma;
497  G4double beta = std::sqrt((gamma2-1.0)/gamma2);
498 
499  // ac corresponds to "A" of Eq. (2.31)
500  //
501  G4double ac = (1.0/beta) - 1.0;
502  G4double a1 = 0.5*beta*gamma*(gamma-1.0)*(gamma-2.0);
503  G4double a2 = ac + 2.0;
504  G4double gtmax = 2.0*(a1 + 1.0/ac);
505 
506  G4double tsam = 0;
507  G4double gtr = 0;
508 
509  //2) sampling. Eq. (2.31) of Penelope Manual
510  // tsam = 1-std::cos(theta)
511  // gtr = rejection function according to Eq. (2.28)
512  do{
513  G4double rand = G4UniformRand();
514  tsam = 2.0*ac * (2.0*rand + a2*std::sqrt(rand)) / (a2*a2 - 4.0*rand);
515  gtr = (2.0 - tsam) * (a1 + 1.0/(ac+tsam));
516  }while(G4UniformRand()*gtmax > gtr);
517  costheta = 1.0-tsam;
518 
519 
520  return costheta;
521 }
static const G4double a1
#define G4UniformRand()
Definition: Randomize.hh:97
double energy
Definition: plottest35.C:25
static const double GeV
Definition: G4SIunits.hh:214
float electron_mass_c2
Definition: hepunit.py:274
double G4double
Definition: G4Types.hh:76
static const G4double a2
Here is the caller graph for this function:

◆ SampleSecondaries()

void G4PenelopePhotoElectricModel::SampleSecondaries ( std::vector< G4DynamicParticle *> *  fvect,
const G4MaterialCutsCouple couple,
const G4DynamicParticle aDynamicGamma,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 277 of file G4PenelopePhotoElectricModel.cc.

282 {
283  //
284  // Photoelectric effect, Penelope model v2008
285  //
286  // The target atom and the target shell are sampled according to the Livermore
287  // database
288  // D.E. Cullen et al., Report UCRL-50400 (1989)
289  // The angular distribution of the electron in the final state is sampled
290  // according to the Sauter distribution from
291  // F. Sauter, Ann. Phys. 11 (1931) 454
292  // The energy of the final electron is given by the initial photon energy minus
293  // the binding energy. Fluorescence de-excitation is subsequently produced
294  // (to fill the vacancy) according to the general Geant4 G4DeexcitationManager:
295  // J. Stepanek, Comp. Phys. Comm. 1206 pp 1-1-9 (1997)
296 
297  if (verboseLevel > 3)
298  G4cout << "Calling SamplingSecondaries() of G4PenelopePhotoElectricModel" << G4endl;
299 
300  G4double photonEnergy = aDynamicGamma->GetKineticEnergy();
301 
302  // always kill primary
303  fParticleChange->ProposeTrackStatus(fStopAndKill);
304  fParticleChange->SetProposedKineticEnergy(0.);
305 
306  if (photonEnergy <= fIntrinsicLowEnergyLimit)
307  {
308  fParticleChange->ProposeLocalEnergyDeposit(photonEnergy);
309  return ;
310  }
311 
312  G4ParticleMomentum photonDirection = aDynamicGamma->GetMomentumDirection();
313 
314  // Select randomly one element in the current material
315  if (verboseLevel > 2)
316  G4cout << "Going to select element in " << couple->GetMaterial()->GetName() << G4endl;
317 
318  // atom can be selected efficiently if element selectors are initialised
319  const G4Element* anElement =
320  SelectRandomAtom(couple,G4Gamma::GammaDefinition(),photonEnergy);
321  G4int Z = (G4int) anElement->GetZ();
322  if (verboseLevel > 2)
323  G4cout << "Selected " << anElement->GetName() << G4endl;
324 
325  // Select the ionised shell in the current atom according to shell cross sections
326  //shellIndex = 0 --> K shell
327  // 1-3 --> L shells
328  // 4-8 --> M shells
329  // 9 --> outer shells cumulatively
330  //
331  size_t shellIndex = SelectRandomShell(Z,photonEnergy);
332 
333  if (verboseLevel > 2)
334  G4cout << "Selected shell " << shellIndex << " of element " << anElement->GetName() << G4endl;
335 
336  // Retrieve the corresponding identifier and binding energy of the selected shell
338 
339  //The number of shell cross section possibly reported in the Penelope database
340  //might be different from the number of shells in the G4AtomicTransitionManager
341  //(namely, Penelope may contain more shell, especially for very light elements).
342  //In order to avoid a warning message from the G4AtomicTransitionManager, I
343  //add this protection. Results are anyway changed, because when G4AtomicTransitionManager
344  //has a shellID>maxID, it sets the shellID to the last valid shell.
345  size_t numberOfShells = (size_t) transitionManager->NumberOfShells(Z);
346  if (shellIndex >= numberOfShells)
347  shellIndex = numberOfShells-1;
348 
349  const G4AtomicShell* shell = fTransitionManager->Shell(Z,shellIndex);
350  G4double bindingEnergy = shell->BindingEnergy();
351  //G4int shellId = shell->ShellId();
352 
353  //Penelope considers only K, L and M shells. Cross sections of outer shells are
354  //not included in the Penelope database. If SelectRandomShell() returns
355  //shellIndex = 9, it means that an outer shell was ionized. In this case the
356  //Penelope recipe is to set bindingEnergy = 0 (the energy is entirely assigned
357  //to the electron) and to disregard fluorescence.
358  if (shellIndex == 9)
359  bindingEnergy = 0.*eV;
360 
361 
362  G4double localEnergyDeposit = 0.0;
363  G4double cosTheta = 1.0;
364 
365  // Primary outcoming electron
366  G4double eKineticEnergy = photonEnergy - bindingEnergy;
367 
368  // There may be cases where the binding energy of the selected shell is > photon energy
369  // In such cases do not generate secondaries
370  if (eKineticEnergy > 0.)
371  {
372  // The electron is created
373  // Direction sampled from the Sauter distribution
374  cosTheta = SampleElectronDirection(eKineticEnergy);
375  G4double sinTheta = std::sqrt(1-cosTheta*cosTheta);
376  G4double phi = twopi * G4UniformRand() ;
377  G4double dirx = sinTheta * std::cos(phi);
378  G4double diry = sinTheta * std::sin(phi);
379  G4double dirz = cosTheta ;
380  G4ThreeVector electronDirection(dirx,diry,dirz); //electron direction
381  electronDirection.rotateUz(photonDirection);
383  electronDirection,
384  eKineticEnergy);
385  fvect->push_back(electron);
386  }
387  else
388  bindingEnergy = photonEnergy;
389 
390 
391  G4double energyInFluorescence = 0; //testing purposes
392  G4double energyInAuger = 0; //testing purposes
393 
394  //Now, take care of fluorescence, if required. According to the Penelope
395  //recipe, I have to skip fluoresence completely if shellIndex == 9
396  //(= sampling of a shell outer than K,L,M)
397  if (fAtomDeexcitation && shellIndex<9)
398  {
399  G4int index = couple->GetIndex();
401  {
402  size_t nBefore = fvect->size();
403  fAtomDeexcitation->GenerateParticles(fvect,shell,Z,index);
404  size_t nAfter = fvect->size();
405 
406  if (nAfter > nBefore) //actual production of fluorescence
407  {
408  for (size_t j=nBefore;j<nAfter;j++) //loop on products
409  {
410  G4double itsEnergy = ((*fvect)[j])->GetKineticEnergy();
411  bindingEnergy -= itsEnergy;
412  if (((*fvect)[j])->GetParticleDefinition() == G4Gamma::Definition())
413  energyInFluorescence += itsEnergy;
414  else if (((*fvect)[j])->GetParticleDefinition() == G4Electron::Definition())
415  energyInAuger += itsEnergy;
416  }
417  }
418  }
419  }
420 
421  //Residual energy is deposited locally
422  localEnergyDeposit += bindingEnergy;
423 
424  if (localEnergyDeposit < 0)
425  {
426  G4cout << "WARNING - "
427  << "G4PenelopePhotoElectricModel::SampleSecondaries() - Negative energy deposit"
428  << G4endl;
429  localEnergyDeposit = 0;
430  }
431 
432  fParticleChange->ProposeLocalEnergyDeposit(localEnergyDeposit);
433 
434  if (verboseLevel > 1)
435  {
436  G4cout << "-----------------------------------------------------------" << G4endl;
437  G4cout << "Energy balance from G4PenelopePhotoElectric" << G4endl;
438  G4cout << "Selected shell: " << WriteTargetShell(shellIndex) << " of element " <<
439  anElement->GetName() << G4endl;
440  G4cout << "Incoming photon energy: " << photonEnergy/keV << " keV" << G4endl;
441  G4cout << "-----------------------------------------------------------" << G4endl;
442  if (eKineticEnergy)
443  G4cout << "Outgoing electron " << eKineticEnergy/keV << " keV" << G4endl;
444  if (energyInFluorescence)
445  G4cout << "Fluorescence x-rays: " << energyInFluorescence/keV << " keV" << G4endl;
446  if (energyInAuger)
447  G4cout << "Auger electrons: " << energyInAuger/keV << " keV" << G4endl;
448  G4cout << "Local energy deposit " << localEnergyDeposit/keV << " keV" << G4endl;
449  G4cout << "Total final state: " <<
450  (eKineticEnergy+energyInFluorescence+localEnergyDeposit+energyInAuger)/keV <<
451  " keV" << G4endl;
452  G4cout << "-----------------------------------------------------------" << G4endl;
453  }
454  if (verboseLevel > 0)
455  {
456  G4double energyDiff =
457  std::fabs(eKineticEnergy+energyInFluorescence+localEnergyDeposit+energyInAuger-photonEnergy);
458  if (energyDiff > 0.05*keV)
459  {
460  G4cout << "Warning from G4PenelopePhotoElectric: problem with energy conservation: " <<
461  (eKineticEnergy+energyInFluorescence+localEnergyDeposit+energyInAuger)/keV
462  << " keV (final) vs. " <<
463  photonEnergy/keV << " keV (initial)" << G4endl;
464  G4cout << "-----------------------------------------------------------" << G4endl;
465  G4cout << "Energy balance from G4PenelopePhotoElectric" << G4endl;
466  G4cout << "Selected shell: " << WriteTargetShell(shellIndex) << " of element " <<
467  anElement->GetName() << G4endl;
468  G4cout << "Incoming photon energy: " << photonEnergy/keV << " keV" << G4endl;
469  G4cout << "-----------------------------------------------------------" << G4endl;
470  if (eKineticEnergy)
471  G4cout << "Outgoing electron " << eKineticEnergy/keV << " keV" << G4endl;
472  if (energyInFluorescence)
473  G4cout << "Fluorescence x-rays: " << energyInFluorescence/keV << " keV" << G4endl;
474  if (energyInAuger)
475  G4cout << "Auger electrons: " << energyInAuger/keV << " keV" << G4endl;
476  G4cout << "Local energy deposit " << localEnergyDeposit/keV << " keV" << G4endl;
477  G4cout << "Total final state: " <<
478  (eKineticEnergy+energyInFluorescence+localEnergyDeposit+energyInAuger)/keV <<
479  " keV" << G4endl;
480  G4cout << "-----------------------------------------------------------" << G4endl;
481  }
482  }
483 }
G4bool CheckDeexcitationActiveRegion(G4int coupleIndex)
G4AtomicShell * Shell(G4int Z, size_t shellIndex) const
const G4Material * GetMaterial() const
Int_t index
G4ParticleChangeForGamma * fParticleChange
static G4Electron * Definition()
Definition: G4Electron.cc:49
int G4int
Definition: G4Types.hh:78
const G4AtomicTransitionManager * fTransitionManager
G4double GetKineticEnergy() const
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
Float_t Z
static const double twopi
Definition: G4SIunits.hh:75
const G4String & GetName() const
Definition: G4Element.hh:127
static const double eV
Definition: G4SIunits.hh:212
const G4ThreeVector & GetMomentumDirection() const
void GenerateParticles(std::vector< G4DynamicParticle *> *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
size_t SelectRandomShell(G4int Z, G4double energy)
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
static const double keV
Definition: G4SIunits.hh:213
double G4double
Definition: G4Types.hh:76
G4double BindingEnergy() const
const G4String & GetName() const
Definition: G4Material.hh:178
static G4AtomicTransitionManager * Instance()
G4double SampleElectronDirection(G4double energy)
G4double GetZ() const
Definition: G4Element.hh:131
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:544
static G4Gamma * Definition()
Definition: G4Gamma.cc:49
static G4Gamma * GammaDefinition()
Definition: G4Gamma.cc:81
Here is the call graph for this function:

◆ SelectRandomShell()

size_t G4PenelopePhotoElectricModel::SelectRandomShell ( G4int  Z,
G4double  energy 
)
private

Definition at line 729 of file G4PenelopePhotoElectricModel.cc.

730 {
731  G4double logEnergy = std::log(energy);
732 
733  //Check if data have been read (it should be!)
734  if (!logAtomicShellXS->count(Z))
735  {
737  ed << "Cannot find shell cross section data for Z=" << Z << G4endl;
738  G4Exception("G4PenelopePhotoElectricModel::SelectRandomShell()",
739  "em2038",FatalException,ed);
740  }
741 
742  // size_t shellIndex = 0;
743 
744  G4PhysicsTable* theTable = logAtomicShellXS->find(Z)->second;
745 
746  G4double sum = 0;
747 
748  G4PhysicsFreeVector* totalXSLog = (G4PhysicsFreeVector*) (*theTable)[0];
749  G4double logXS = totalXSLog->Value(logEnergy);
750  G4double totalXS = G4Exp(logXS);
751 
752  //Notice: totalXS is the total cross section and it does *not* correspond to
753  //the sum of partialXS's, since these include only K, L and M shells.
754  //
755  // Therefore, here one have to consider the possibility of ionisation of
756  // an outer shell. Conventionally, it is indicated with id=10 in Penelope
757  //
758  G4double random = G4UniformRand()*totalXS;
759 
760  for (size_t k=1;k<theTable->entries();k++)
761  {
762  //Add one shell
763  G4PhysicsFreeVector* partialXSLog = (G4PhysicsFreeVector*) (*theTable)[k];
764  G4double logXSLocal = partialXSLog->Value(logEnergy);
765  G4double partialXS = G4Exp(logXSLocal);
766  sum += partialXS;
767  if (random <= sum)
768  return k-1;
769  }
770  //none of the shells K, L, M: return outer shell
771  return 9;
772 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
std::map< G4int, G4PhysicsTable * > * logAtomicShellXS
#define G4UniformRand()
Definition: Randomize.hh:97
double energy
Definition: plottest35.C:25
Float_t Z
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
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
size_t entries() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetParticle()

void G4PenelopePhotoElectricModel::SetParticle ( const G4ParticleDefinition p)
private

Definition at line 720 of file G4PenelopePhotoElectricModel.cc.

721 {
722  if(!fParticle) {
723  fParticle = p;
724  }
725 }
const G4ParticleDefinition * fParticle
Here is the caller graph for this function:

◆ SetVerbosityLevel()

void G4PenelopePhotoElectricModel::SetVerbosityLevel ( G4int  lev)
inline

◆ WriteTargetShell()

G4String G4PenelopePhotoElectricModel::WriteTargetShell ( size_t  shellID)
private

Definition at line 693 of file G4PenelopePhotoElectricModel.cc.

694 {
695  G4String theShell = "outer shell";
696  if (shellID == 0)
697  theShell = "K";
698  else if (shellID == 1)
699  theShell = "L1";
700  else if (shellID == 2)
701  theShell = "L2";
702  else if (shellID == 3)
703  theShell = "L3";
704  else if (shellID == 4)
705  theShell = "M1";
706  else if (shellID == 5)
707  theShell = "M2";
708  else if (shellID == 6)
709  theShell = "M3";
710  else if (shellID == 7)
711  theShell = "M4";
712  else if (shellID == 8)
713  theShell = "M5";
714 
715  return theShell;
716 }
Here is the caller graph for this function:

Member Data Documentation

◆ fAtomDeexcitation

G4VAtomDeexcitation* G4PenelopePhotoElectricModel::fAtomDeexcitation
private

Definition at line 114 of file G4PenelopePhotoElectricModel.hh.

◆ fIntrinsicHighEnergyLimit

G4double G4PenelopePhotoElectricModel::fIntrinsicHighEnergyLimit
private

Definition at line 109 of file G4PenelopePhotoElectricModel.hh.

◆ fIntrinsicLowEnergyLimit

G4double G4PenelopePhotoElectricModel::fIntrinsicLowEnergyLimit
private

Definition at line 108 of file G4PenelopePhotoElectricModel.hh.

◆ fLocalTable

G4bool G4PenelopePhotoElectricModel::fLocalTable
private

Definition at line 129 of file G4PenelopePhotoElectricModel.hh.

◆ fParticle

const G4ParticleDefinition* G4PenelopePhotoElectricModel::fParticle
protected

Definition at line 97 of file G4PenelopePhotoElectricModel.hh.

◆ fParticleChange

G4ParticleChangeForGamma* G4PenelopePhotoElectricModel::fParticleChange
protected

Definition at line 96 of file G4PenelopePhotoElectricModel.hh.

◆ fTransitionManager

const G4AtomicTransitionManager* G4PenelopePhotoElectricModel::fTransitionManager
private

Definition at line 115 of file G4PenelopePhotoElectricModel.hh.

◆ isInitialised

G4bool G4PenelopePhotoElectricModel::isInitialised
private

Definition at line 112 of file G4PenelopePhotoElectricModel.hh.

◆ logAtomicShellXS

std::map<G4int,G4PhysicsTable*>* G4PenelopePhotoElectricModel::logAtomicShellXS
private

Definition at line 123 of file G4PenelopePhotoElectricModel.hh.

◆ verboseLevel

G4int G4PenelopePhotoElectricModel::verboseLevel
private

Definition at line 111 of file G4PenelopePhotoElectricModel.hh.


The documentation for this class was generated from the following files: