Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PenelopePhotoElectricModel Class Reference

#include <G4PenelopePhotoElectricModel.hh>

+ Inheritance 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 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 Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
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 ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., 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 *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
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)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *)
 
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
 
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 ForceBuildTable (G4bool val)
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 

Protected Attributes

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

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 *)
 

Detailed Description

Definition at line 58 of file G4PenelopePhotoElectricModel.hh.

Constructor & Destructor Documentation

◆ G4PenelopePhotoElectricModel()

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

Definition at line 60 of file G4PenelopePhotoElectricModel.cc.

62 :G4VEmModel(nam),fParticleChange(0),isInitialised(false),fAtomDeexcitation(0),
63 logAtomicShellXS(0)
64{
65 fIntrinsicLowEnergyLimit = 100.0*eV;
66 fIntrinsicHighEnergyLimit = 100.0*GeV;
67 // SetLowEnergyLimit(fIntrinsicLowEnergyLimit);
68 SetHighEnergyLimit(fIntrinsicHighEnergyLimit);
69 //
70 verboseLevel= 0;
71 // Verbosity scale:
72 // 0 = nothing
73 // 1 = warning for energy non-conservation
74 // 2 = details of energy budget
75 // 3 = calculation of cross sections, file openings, sampling of atoms
76 // 4 = entering in methods
77
78 //Mark this model as "applicable" for atomic deexcitation
80
81 fTransitionManager = G4AtomicTransitionManager::Instance();
82}
static G4AtomicTransitionManager * Instance()
G4ParticleChangeForGamma * fParticleChange
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:641

◆ ~G4PenelopePhotoElectricModel()

G4PenelopePhotoElectricModel::~G4PenelopePhotoElectricModel ( )
virtual

Definition at line 86 of file G4PenelopePhotoElectricModel.cc.

87{
88 std::map <G4int,G4PhysicsTable*>::iterator i;
89 if (logAtomicShellXS)
90 {
91 for (i=logAtomicShellXS->begin();i != logAtomicShellXS->end();i++)
92 {
93 G4PhysicsTable* tab = i->second;
94 tab->clearAndDestroy();
95 delete tab;
96 }
97 }
98 delete logAtomicShellXS;
99}
void clearAndDestroy()

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

146{
147 //
148 // Penelope model v2008
149 //
150
151 if (verboseLevel > 3)
152 G4cout << "Calling ComputeCrossSectionPerAtom() of G4PenelopePhotoElectricModel" << G4endl;
153
154 G4int iZ = (G4int) Z;
155
156 //read data files
157 if (!logAtomicShellXS->count(iZ))
158 ReadDataFile(iZ);
159 //now it should be ok
160 if (!logAtomicShellXS->count(iZ))
161 G4Exception("G4PenelopePhotoElectricModel::ComputeCrossSectionPerAtom()",
162 "em2038",FatalException,
163 "Unable to retrieve the shell cross section table");
164
165 G4double cross = 0;
166
167 G4PhysicsTable* theTable = logAtomicShellXS->find(iZ)->second;
168 G4PhysicsFreeVector* totalXSLog = (G4PhysicsFreeVector*) (*theTable)[0];
169
170 if (!totalXSLog)
171 {
172 G4Exception("G4PenelopePhotoElectricModel::ComputeCrossSectionPerAtom()",
173 "em2039",FatalException,
174 "Unable to retrieve the total cross section table");
175 return 0;
176 }
177 G4double logene = std::log(energy);
178 G4double logXS = totalXSLog->Value(logene);
179 cross = std::exp(logXS);
180
181 if (verboseLevel > 2)
182 G4cout << "Photoelectric cross section at " << energy/MeV << " MeV for Z=" << Z <<
183 " = " << cross/barn << " barn" << G4endl;
184 return cross;
185}
@ FatalException
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
G4double Value(G4double theEnergy)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

◆ GetNumberOfShellXS()

size_t G4PenelopePhotoElectricModel::GetNumberOfShellXS ( G4int  Z)

Definition at line 620 of file G4PenelopePhotoElectricModel.cc.

621{
622 //read data files
623 if (!logAtomicShellXS->count(Z))
624 ReadDataFile(Z);
625 //now it should be ok
626 if (!logAtomicShellXS->count(Z))
627 {
629 ed << "Cannot find shell cross section data for Z=" << Z << G4endl;
630 G4Exception("G4PenelopePhotoElectricModel::GetNumberOfShellXS()",
631 "em2038",FatalException,ed);
632 }
633 //one vector is allocated for the _total_ cross section
634 size_t nEntries = logAtomicShellXS->find(Z)->second->entries();
635 return (nEntries-1);
636}
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76

Referenced by GetShellCrossSection().

◆ GetShellCrossSection()

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

Definition at line 640 of file G4PenelopePhotoElectricModel.cc.

641{
642 //this forces also the loading of the data
643 size_t entries = GetNumberOfShellXS(Z);
644
645 if (shellID >= entries)
646 {
647 G4cout << "Element Z=" << Z << " has data for " << entries << " shells only" << G4endl;
648 G4cout << "so shellID should be from 0 to " << entries-1 << G4endl;
649 return 0;
650 }
651
652 G4PhysicsTable* theTable = logAtomicShellXS->find(Z)->second;
653 //[0] is the total XS, shellID is in the element [shellID+1]
654 G4PhysicsFreeVector* totalXSLog = (G4PhysicsFreeVector*) (*theTable)[shellID+1];
655
656 if (!totalXSLog)
657 {
658 G4Exception("G4PenelopePhotoElectricModel::GetShellCrossSection()",
659 "em2039",FatalException,
660 "Unable to retrieve the total cross section table");
661 return 0;
662 }
663 G4double logene = std::log(energy);
664 G4double logXS = totalXSLog->Value(logene);
665 G4double cross = std::exp(logXS);
666 if (cross < 2e-40*cm2) cross = 0;
667 return cross;
668}

◆ GetVerbosityLevel()

G4int G4PenelopePhotoElectricModel::GetVerbosityLevel ( )
inline

Definition at line 85 of file G4PenelopePhotoElectricModel.hh.

85{return verboseLevel;};

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 103 of file G4PenelopePhotoElectricModel.cc.

105{
106 if (verboseLevel > 3)
107 G4cout << "Calling G4PenelopePhotoElectricModel::Initialise()" << G4endl;
108
109 // logAtomicShellXS is created only once, since it is never cleared
110 if (!logAtomicShellXS)
111 logAtomicShellXS = new std::map<G4int,G4PhysicsTable*>;
112
113 InitialiseElementSelectors(particle,cuts);
114
115 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
116 //Issue warning if the AtomicDeexcitation has not been declared
117 if (!fAtomDeexcitation)
118 {
119 G4cout << G4endl;
120 G4cout << "WARNING from G4PenelopePhotoElectricModel " << G4endl;
121 G4cout << "Atomic de-excitation module is not instantiated, so there will not be ";
122 G4cout << "any fluorescence/Auger emission." << G4endl;
123 G4cout << "Please make sure this is intended" << G4endl;
124 }
125
126 if (verboseLevel > 0) {
127 G4cout << "Penelope Photo-Electric model v2008 is initialized " << G4endl
128 << "Energy range: "
129 << LowEnergyLimit() / MeV << " MeV - "
130 << HighEnergyLimit() / GeV << " GeV";
131 }
132
133 if(isInitialised) return;
135 isInitialised = true;
136
137}
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:109
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:529
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:522
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:123

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 189 of file G4PenelopePhotoElectricModel.cc.

194{
195 //
196 // Photoelectric effect, Penelope model v2008
197 //
198 // The target atom and the target shell are sampled according to the Livermore
199 // database
200 // D.E. Cullen et al., Report UCRL-50400 (1989)
201 // The angular distribution of the electron in the final state is sampled
202 // according to the Sauter distribution from
203 // F. Sauter, Ann. Phys. 11 (1931) 454
204 // The energy of the final electron is given by the initial photon energy minus
205 // the binding energy. Fluorescence de-excitation is subsequently produced
206 // (to fill the vacancy) according to the general Geant4 G4DeexcitationManager:
207 // J. Stepanek, Comp. Phys. Comm. 1206 pp 1-1-9 (1997)
208
209 if (verboseLevel > 3)
210 G4cout << "Calling SamplingSecondaries() of G4PenelopePhotoElectricModel" << G4endl;
211
212 G4double photonEnergy = aDynamicGamma->GetKineticEnergy();
213
214 // always kill primary
217
218 if (photonEnergy <= fIntrinsicLowEnergyLimit)
219 {
221 return ;
222 }
223
224 G4ParticleMomentum photonDirection = aDynamicGamma->GetMomentumDirection();
225
226 // Select randomly one element in the current material
227 if (verboseLevel > 2)
228 G4cout << "Going to select element in " << couple->GetMaterial()->GetName() << G4endl;
229
230 // atom can be selected efficiently if element selectors are initialised
231 const G4Element* anElement =
232 SelectRandomAtom(couple,G4Gamma::GammaDefinition(),photonEnergy);
233 G4int Z = (G4int) anElement->GetZ();
234 if (verboseLevel > 2)
235 G4cout << "Selected " << anElement->GetName() << G4endl;
236
237 // Select the ionised shell in the current atom according to shell cross sections
238 //shellIndex = 0 --> K shell
239 // 1-3 --> L shells
240 // 4-8 --> M shells
241 // 9 --> outer shells cumulatively
242 //
243 size_t shellIndex = SelectRandomShell(Z,photonEnergy);
244
245 if (verboseLevel > 2)
246 G4cout << "Selected shell " << shellIndex << " of element " << anElement->GetName() << G4endl;
247
248 // Retrieve the corresponding identifier and binding energy of the selected shell
250
251 //The number of shell cross section possibly reported in the Penelope database
252 //might be different from the number of shells in the G4AtomicTransitionManager
253 //(namely, Penelope may contain more shell, especially for very light elements).
254 //In order to avoid a warning message from the G4AtomicTransitionManager, I
255 //add this protection. Results are anyway changed, because when G4AtomicTransitionManager
256 //has a shellID>maxID, it sets the shellID to the last valid shell.
257 size_t numberOfShells = (size_t) transitionManager->NumberOfShells(Z);
258 if (shellIndex >= numberOfShells)
259 shellIndex = numberOfShells-1;
260
261 const G4AtomicShell* shell = fTransitionManager->Shell(Z,shellIndex);
263 //G4int shellId = shell->ShellId();
264
265 //Penelope considers only K, L and M shells. Cross sections of outer shells are
266 //not included in the Penelope database. If SelectRandomShell() returns
267 //shellIndex = 9, it means that an outer shell was ionized. In this case the
268 //Penelope recipe is to set bindingEnergy = 0 (the energy is entirely assigned
269 //to the electron) and to disregard fluorescence.
270 if (shellIndex == 9)
271 bindingEnergy = 0.*eV;
272
273
274 G4double localEnergyDeposit = 0.0;
275 G4double cosTheta = 1.0;
276
277 // Primary outcoming electron
278 G4double eKineticEnergy = photonEnergy - bindingEnergy;
279
280 // There may be cases where the binding energy of the selected shell is > photon energy
281 // In such cases do not generate secondaries
282 if (eKineticEnergy > 0.)
283 {
284 // The electron is created
285 // Direction sampled from the Sauter distribution
286 cosTheta = SampleElectronDirection(eKineticEnergy);
287 G4double sinTheta = std::sqrt(1-cosTheta*cosTheta);
288 G4double phi = twopi * G4UniformRand() ;
289 G4double dirx = sinTheta * std::cos(phi);
290 G4double diry = sinTheta * std::sin(phi);
291 G4double dirz = cosTheta ;
292 G4ThreeVector electronDirection(dirx,diry,dirz); //electron direction
293 electronDirection.rotateUz(photonDirection);
295 electronDirection,
296 eKineticEnergy);
297 fvect->push_back(electron);
298 }
299 else
300 bindingEnergy = photonEnergy;
301
302
303 G4double energyInFluorescence = 0; //testing purposes
304 G4double energyInAuger = 0; //testing purposes
305
306 //Now, take care of fluorescence, if required. According to the Penelope
307 //recipe, I have to skip fluoresence completely if shellIndex == 9
308 //(= sampling of a shell outer than K,L,M)
309 if (fAtomDeexcitation && shellIndex<9)
310 {
311 G4int index = couple->GetIndex();
312 if (fAtomDeexcitation->CheckDeexcitationActiveRegion(index))
313 {
314 size_t nBefore = fvect->size();
315 fAtomDeexcitation->GenerateParticles(fvect,shell,Z,index);
316 size_t nAfter = fvect->size();
317
318 if (nAfter > nBefore) //actual production of fluorescence
319 {
320 for (size_t j=nBefore;j<nAfter;j++) //loop on products
321 {
322 G4double itsEnergy = ((*fvect)[j])->GetKineticEnergy();
323 bindingEnergy -= itsEnergy;
324 if (((*fvect)[j])->GetParticleDefinition() == G4Gamma::Definition())
325 energyInFluorescence += itsEnergy;
326 else if (((*fvect)[j])->GetParticleDefinition() == G4Electron::Definition())
327 energyInAuger += itsEnergy;
328
329 }
330 }
331
332 }
333 }
334
335 //Residual energy is deposited locally
336 localEnergyDeposit += bindingEnergy;
337
338 if (localEnergyDeposit < 0)
339 {
340 G4cout << "WARNING - "
341 << "G4PenelopePhotoElectricModel::SampleSecondaries() - Negative energy deposit"
342 << G4endl;
343 localEnergyDeposit = 0;
344 }
345
346 fParticleChange->ProposeLocalEnergyDeposit(localEnergyDeposit);
347
348 if (verboseLevel > 1)
349 {
350 G4cout << "-----------------------------------------------------------" << G4endl;
351 G4cout << "Energy balance from G4PenelopePhotoElectric" << G4endl;
352 G4cout << "Selected shell: " << WriteTargetShell(shellIndex) << " of element " <<
353 anElement->GetName() << G4endl;
354 G4cout << "Incoming photon energy: " << photonEnergy/keV << " keV" << G4endl;
355 G4cout << "-----------------------------------------------------------" << G4endl;
356 if (eKineticEnergy)
357 G4cout << "Outgoing electron " << eKineticEnergy/keV << " keV" << G4endl;
358 if (energyInFluorescence)
359 G4cout << "Fluorescence x-rays: " << energyInFluorescence/keV << " keV" << G4endl;
360 if (energyInAuger)
361 G4cout << "Auger electrons: " << energyInAuger/keV << " keV" << G4endl;
362 G4cout << "Local energy deposit " << localEnergyDeposit/keV << " keV" << G4endl;
363 G4cout << "Total final state: " <<
364 (eKineticEnergy+energyInFluorescence+localEnergyDeposit+energyInAuger)/keV <<
365 " keV" << G4endl;
366 G4cout << "-----------------------------------------------------------" << G4endl;
367 }
368 if (verboseLevel > 0)
369 {
370 G4double energyDiff =
371 std::fabs(eKineticEnergy+energyInFluorescence+localEnergyDeposit+energyInAuger-photonEnergy);
372 if (energyDiff > 0.05*keV)
373 {
374 G4cout << "Warning from G4PenelopePhotoElectric: problem with energy conservation: " <<
375 (eKineticEnergy+energyInFluorescence+localEnergyDeposit+energyInAuger)/keV
376 << " keV (final) vs. " <<
377 photonEnergy/keV << " keV (initial)" << G4endl;
378 G4cout << "-----------------------------------------------------------" << G4endl;
379 G4cout << "Energy balance from G4PenelopePhotoElectric" << G4endl;
380 G4cout << "Selected shell: " << WriteTargetShell(shellIndex) << " of element " <<
381 anElement->GetName() << G4endl;
382 G4cout << "Incoming photon energy: " << photonEnergy/keV << " keV" << G4endl;
383 G4cout << "-----------------------------------------------------------" << G4endl;
384 if (eKineticEnergy)
385 G4cout << "Outgoing electron " << eKineticEnergy/keV << " keV" << G4endl;
386 if (energyInFluorescence)
387 G4cout << "Fluorescence x-rays: " << energyInFluorescence/keV << " keV" << G4endl;
388 if (energyInAuger)
389 G4cout << "Auger electrons: " << energyInAuger/keV << " keV" << G4endl;
390 G4cout << "Local energy deposit " << localEnergyDeposit/keV << " keV" << G4endl;
391 G4cout << "Total final state: " <<
392 (eKineticEnergy+energyInFluorescence+localEnergyDeposit+energyInAuger)/keV <<
393 " keV" << G4endl;
394 G4cout << "-----------------------------------------------------------" << G4endl;
395 }
396 }
397}
@ fStopAndKill
#define G4UniformRand()
Definition: Randomize.hh:53
G4double BindingEnergy() const
G4AtomicShell * Shell(G4int Z, size_t shellIndex) const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
static G4Electron * Definition()
Definition: G4Electron.cc:49
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4double GetZ() const
Definition: G4Element.hh:131
const G4String & GetName() const
Definition: G4Element.hh:127
static G4Gamma * GammaDefinition()
Definition: G4Gamma.cc:81
static G4Gamma * Definition()
Definition: G4Gamma.cc:49
const G4Material * GetMaterial() const
const G4String & GetName() const
Definition: G4Material.hh:177
void SetProposedKineticEnergy(G4double proposedKinEnergy)
G4bool CheckDeexcitationActiveRegion(G4int coupleIndex)
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:459
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4double bindingEnergy(G4int A, G4int Z)

◆ SetVerbosityLevel()

void G4PenelopePhotoElectricModel::SetVerbosityLevel ( G4int  lev)
inline

Definition at line 84 of file G4PenelopePhotoElectricModel.hh.

84{verboseLevel = lev;};

Member Data Documentation

◆ fParticleChange

G4ParticleChangeForGamma* G4PenelopePhotoElectricModel::fParticleChange
protected

Definition at line 93 of file G4PenelopePhotoElectricModel.hh.

Referenced by Initialise(), and SampleSecondaries().


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