Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 level, const G4ParticleDefinition *, G4double kineticEnergy)
 
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 * > *)
 
virtual 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=nullptr)
 
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

G4ParticleChangeForGammafParticleChange
 
const G4ParticleDefinitionfParticle
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
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::G4PenelopePhotoElectricModel ( const G4ParticleDefinition p = 0,
const G4String processName = "PenPhotoElec" 
)

Definition at line 66 of file G4PenelopePhotoElectricModel.cc.

69  isInitialised(false),fAtomDeexcitation(0),logAtomicShellXS(0),fLocalTable(false)
70 {
71  fIntrinsicLowEnergyLimit = 100.0*eV;
72  fIntrinsicHighEnergyLimit = 100.0*GeV;
73  // SetLowEnergyLimit(fIntrinsicLowEnergyLimit);
74  SetHighEnergyLimit(fIntrinsicHighEnergyLimit);
75  //
76 
77  if (part)
78  SetParticle(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 
91  fTransitionManager = G4AtomicTransitionManager::Instance();
92 }
G4ParticleChangeForGamma * fParticleChange
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:68
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:724
static constexpr double eV
Definition: G4SIunits.hh:215
static constexpr double GeV
Definition: G4SIunits.hh:217
const G4ParticleDefinition * fParticle
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:780
static G4AtomicTransitionManager * Instance()

Here is the call graph for this function:

G4PenelopePhotoElectricModel::~G4PenelopePhotoElectricModel ( )
virtual

Definition at line 96 of file G4PenelopePhotoElectricModel.cc.

97 {
98  if (IsMaster() || fLocalTable)
99  {
100  if (logAtomicShellXS)
101  {
102  for (auto& item : (*logAtomicShellXS))
103  {
104  //G4PhysicsTable* tab = item.second;
105  //tab->clearAndDestroy();
106  delete item.second;
107  }
108  }
109  delete logAtomicShellXS;
110  }
111 }
G4bool IsMaster() const
Definition: G4VEmModel.hh:717

Here is the call graph for this function:

Member Function Documentation

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 206 of file G4PenelopePhotoElectricModel.cc.

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
236  if (verboseLevel > 0)
237  {
238  //Issue a G4Exception (warning) only in verbose mode
240  ed << "Unable to retrieve the shell cross section table for Z=" << iZ << G4endl;
241  ed << "This can happen only in Unit Tests or via G4EmCalculator" << G4endl;
242  G4Exception("G4PenelopePhotoElectricModel::ComputeCrossSectionPerAtom()",
243  "em2038",JustWarning,ed);
244  }
245  //protect file reading via autolock
246  G4AutoLock lock(&PenelopePhotoElectricModelMutex);
247  ReadDataFile(iZ);
248  lock.unlock();
249 
250  }
251 
252  G4double cross = 0;
253 
254  G4PhysicsTable* theTable = logAtomicShellXS->find(iZ)->second;
255  G4PhysicsFreeVector* totalXSLog = (G4PhysicsFreeVector*) (*theTable)[0];
256 
257  if (!totalXSLog)
258  {
259  G4Exception("G4PenelopePhotoElectricModel::ComputeCrossSectionPerAtom()",
260  "em2039",FatalException,
261  "Unable to retrieve the total cross section table");
262  return 0;
263  }
264  G4double logene = std::log(energy);
265  G4double logXS = totalXSLog->Value(logene);
266  cross = G4Exp(logXS);
267 
268  if (verboseLevel > 2)
269  G4cout << "Photoelectric cross section at " << energy/MeV << " MeV for Z=" << Z <<
270  " = " << cross/barn << " barn" << G4endl;
271  return cross;
272 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
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
G4double energy(const ThreeVector &p, const G4double m)
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
double G4double
Definition: G4Types.hh:76
static constexpr double barn
Definition: G4SIunits.hh:105

Here is the call graph for this function:

size_t G4PenelopePhotoElectricModel::GetNumberOfShellXS ( G4int  Z)

Definition at line 635 of file G4PenelopePhotoElectricModel.cc.

636 {
637  if (!IsMaster())
638  //Should not be here!
639  G4Exception("G4PenelopePhotoElectricModel::GetNumberOfShellXS()",
640  "em0100",FatalException,"Worker thread in this method");
641 
642  //read data files
643  if (!logAtomicShellXS->count(Z))
644  ReadDataFile(Z);
645  //now it should be ok
646  if (!logAtomicShellXS->count(Z))
647  {
649  ed << "Cannot find shell cross section data for Z=" << Z << G4endl;
650  G4Exception("G4PenelopePhotoElectricModel::GetNumberOfShellXS()",
651  "em2038",FatalException,ed);
652  }
653  //one vector is allocated for the _total_ cross section
654  size_t nEntries = logAtomicShellXS->find(Z)->second->entries();
655  return (nEntries-1);
656 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4bool IsMaster() const
Definition: G4VEmModel.hh:717
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:

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

Definition at line 660 of file G4PenelopePhotoElectricModel.cc.

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

Here is the call graph for this function:

G4int G4PenelopePhotoElectricModel::GetVerbosityLevel ( )
inline

Definition at line 88 of file G4PenelopePhotoElectricModel.hh.

88 {return verboseLevel;};
void G4PenelopePhotoElectricModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector cuts 
)
virtual

Implements G4VEmModel.

Definition at line 115 of file G4PenelopePhotoElectricModel.cc.

117 {
118  if (verboseLevel > 3)
119  G4cout << "Calling G4PenelopePhotoElectricModel::Initialise()" << G4endl;
120 
121  fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
122  //Issue warning if the AtomicDeexcitation has not been declared
123  if (!fAtomDeexcitation)
124  {
125  G4cout << G4endl;
126  G4cout << "WARNING from G4PenelopePhotoElectricModel " << G4endl;
127  G4cout << "Atomic de-excitation module is not instantiated, so there will not be ";
128  G4cout << "any fluorescence/Auger emission." << G4endl;
129  G4cout << "Please make sure this is intended" << G4endl;
130  }
131 
132  SetParticle(particle);
133 
134  //Only the master model creates/fills/destroys the tables
135  if (IsMaster() && particle == fParticle)
136  {
137 
138  // logAtomicShellXS is created only once, since it is never cleared
139  if (!logAtomicShellXS)
140  logAtomicShellXS = new std::map<G4int,G4PhysicsTable*>;
141 
142  G4ProductionCutsTable* theCoupleTable =
144 
145  for (size_t i=0;i<theCoupleTable->GetTableSize();i++)
146  {
147  const G4Material* material =
148  theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
149  const G4ElementVector* theElementVector = material->GetElementVector();
150 
151  for (size_t j=0;j<material->GetNumberOfElements();j++)
152  {
153  G4int iZ = (G4int) theElementVector->at(j)->GetZ();
154  //read data files only in the master
155  if (!logAtomicShellXS->count(iZ))
156  ReadDataFile(iZ);
157  }
158  }
159 
160 
161  InitialiseElementSelectors(particle,cuts);
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 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:640
static G4LossTableManager * Instance()
std::vector< G4Element * > G4ElementVector
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:146
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:633
G4ParticleChangeForGamma * fParticleChange
G4bool IsMaster() const
Definition: G4VEmModel.hh:717
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
static G4ProductionCutsTable * GetProductionCutsTable()
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
static constexpr double GeV
Definition: G4SIunits.hh:217
const G4ParticleDefinition * fParticle
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
G4VAtomDeexcitation * AtomDeexcitation()
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:132
const G4Material * GetMaterial() const

Here is the call graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 177 of file G4PenelopePhotoElectricModel.cc.

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 }
G4GLOB_DLL std::ostream G4cout
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:801
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:809
const G4ParticleDefinition * fParticle
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

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

Implements G4VEmModel.

Definition at line 276 of file G4PenelopePhotoElectricModel.cc.

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

Here is the call graph for this function:

void G4PenelopePhotoElectricModel::SetVerbosityLevel ( G4int  lev)
inline

Definition at line 87 of file G4PenelopePhotoElectricModel.hh.

87 {verboseLevel = lev;};

Member Data Documentation

const G4ParticleDefinition* G4PenelopePhotoElectricModel::fParticle
protected

Definition at line 97 of file G4PenelopePhotoElectricModel.hh.

G4ParticleChangeForGamma* G4PenelopePhotoElectricModel::fParticleChange
protected

Definition at line 96 of file G4PenelopePhotoElectricModel.hh.


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