Geant4 11.1.1
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=nullptr, const G4String &processName="PenPhotoElec")
 
virtual ~G4PenelopePhotoElectricModel ()
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel) override
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
void SetVerbosityLevel (G4int lev)
 
G4int GetVerbosityLevel ()
 
std::size_t GetNumberOfShellXS (G4int)
 
G4double GetShellCrossSection (G4int Z, std::size_t shellID, G4double energy)
 
G4PenelopePhotoElectricModeloperator= (const G4PenelopePhotoElectricModel &right)=delete
 
 G4PenelopePhotoElectricModel (const G4PenelopePhotoElectricModel &)=delete
 
- 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 void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
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 ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
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 *, const G4double &length, G4double &eloss)
 
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 FillNumberOfSecondaries (G4int &numberOfTriplets, G4int &numberOfRecoil)
 
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)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectTargetAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementGetCurrentElement (const G4Material *mat=nullptr) const
 
G4int SelectRandomAtomNumber (const G4Material *) const
 
const G4IsotopeGetCurrentIsotope (const G4Element *elm=nullptr) const
 
G4int SelectIsotopeNumber (const G4Element *) const
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
G4VEmModelGetTripletModel ()
 
void SetTripletModel (G4VEmModel *)
 
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) const
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetFluctuationFlag (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
void SetUseBaseMaterials (G4bool val)
 
G4bool UseBaseMaterials () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Attributes

G4ParticleChangeForGammafParticleChange
 
const G4ParticleDefinitionfParticle
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4PhysicsTablexSectionTable = nullptr
 
const G4MaterialpBaseMaterial = nullptr
 
const std::vector< G4double > * theDensityFactor = nullptr
 
const std::vector< G4int > * theDensityIdx = nullptr
 
G4double inveplus
 
G4double pFactor = 1.0
 
size_t currentCoupleIndex = 0
 
size_t basedCoupleIndex = 0
 
G4bool lossFlucFlag = true
 

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() [1/2]

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

Definition at line 70 of file G4PenelopePhotoElectricModel.cc.

72 :G4VEmModel(nam),fParticleChange(nullptr),fParticle(nullptr),
73 fAtomDeexcitation(nullptr),fIsInitialised(false),fLocalTable(false)
74{
75 fIntrinsicLowEnergyLimit = 100.0*eV;
76 fIntrinsicHighEnergyLimit = 100.0*GeV;
77 // SetLowEnergyLimit(fIntrinsicLowEnergyLimit);
78 SetHighEnergyLimit(fIntrinsicHighEnergyLimit);
79 //
80
81 if (part)
82 SetParticle(part);
83
84 fVerboseLevel= 0;
85 // Verbosity scale:
86 // 0 = nothing
87 // 1 = warning for energy non-conservation
88 // 2 = details of energy budget
89 // 3 = calculation of cross sections, file openings, sampling of atoms
90 // 4 = entering in methods
91
92 //Mark this model as "applicable" for atomic deexcitation
94
95 fTransitionManager = G4AtomicTransitionManager::Instance();
96}
static G4AtomicTransitionManager * Instance()
G4ParticleChangeForGamma * fParticleChange
const G4ParticleDefinition * fParticle
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:746
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:802

◆ ~G4PenelopePhotoElectricModel()

G4PenelopePhotoElectricModel::~G4PenelopePhotoElectricModel ( )
virtual

Definition at line 100 of file G4PenelopePhotoElectricModel.cc.

101{
102 if (IsMaster() || fLocalTable)
103 {
104 for(G4int i=0; i<=fMaxZ; ++i)
105 {
106 if(fLogAtomicShellXS[i]) {
107 fLogAtomicShellXS[i]->clearAndDestroy();
108 delete fLogAtomicShellXS[i];
109 fLogAtomicShellXS[i] = nullptr;
110 }
111 }
112 }
113}
int G4int
Definition: G4Types.hh:85
void clearAndDestroy()
G4bool IsMaster() const
Definition: G4VEmModel.hh:725

◆ G4PenelopePhotoElectricModel() [2/2]

G4PenelopePhotoElectricModel::G4PenelopePhotoElectricModel ( const G4PenelopePhotoElectricModel )
delete

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 200 of file G4PenelopePhotoElectricModel.cc.

205{
206 //
207 // Penelope model v2008
208 //
209 if (fVerboseLevel > 3)
210 G4cout << "Calling ComputeCrossSectionPerAtom() of G4PenelopePhotoElectricModel" << G4endl;
211
212 G4int iZ = G4int(Z);
213
214 if (!fLogAtomicShellXS[iZ])
215 {
216 //If we are here, it means that Initialize() was inkoved, but the MaterialTable was
217 //not filled up. This can happen in a UnitTest or via G4EmCalculator
218 if (fVerboseLevel > 0)
219 {
220 //Issue a G4Exception (warning) only in verbose mode
222 ed << "Unable to retrieve the shell cross section table for Z=" << iZ << G4endl;
223 ed << "This can happen only in Unit Tests or via G4EmCalculator" << G4endl;
224 G4Exception("G4PenelopePhotoElectricModel::ComputeCrossSectionPerAtom()",
225 "em2038",JustWarning,ed);
226 }
227 //protect file reading via autolock
228 G4AutoLock lock(&PenelopePhotoElectricModelMutex);
229 ReadDataFile(iZ);
230 lock.unlock();
231 }
232
233 G4double cross = 0;
234 G4PhysicsTable* theTable = fLogAtomicShellXS[iZ];
235 G4PhysicsFreeVector* totalXSLog = (G4PhysicsFreeVector*) (*theTable)[0];
236
237 if (!totalXSLog)
238 {
239 G4Exception("G4PenelopePhotoElectricModel::ComputeCrossSectionPerAtom()",
240 "em2039",FatalException,
241 "Unable to retrieve the total cross section table");
242 return 0;
243 }
244 G4double logene = G4Log(energy);
245 G4double logXS = totalXSLog->Value(logene);
246 cross = G4Exp(logXS);
247
248 if (fVerboseLevel > 2)
249 G4cout << "Photoelectric cross section at " << energy/MeV << " MeV for Z=" << Z <<
250 " = " << cross/barn << " barn" << G4endl;
251 return cross;
252}
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
G4double G4Log(G4double x)
Definition: G4Log.hh:227
double G4double
Definition: G4Types.hh:83
const G4int Z[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4double Value(const G4double energy, std::size_t &lastidx) const
G4double energy(const ThreeVector &p, const G4double m)

◆ GetNumberOfShellXS()

std::size_t G4PenelopePhotoElectricModel::GetNumberOfShellXS ( G4int  Z)

Definition at line 618 of file G4PenelopePhotoElectricModel.cc.

619{
620 if (!IsMaster())
621 //Should not be here!
622 G4Exception("G4PenelopePhotoElectricModel::GetNumberOfShellXS()",
623 "em0100",FatalException,"Worker thread in this method");
624
625 //read data files
626 if (!fLogAtomicShellXS[Z])
627 ReadDataFile(Z);
628 //now it should be ok
629 if (!fLogAtomicShellXS[Z])
630 {
632 ed << "Cannot find shell cross section data for Z=" << Z << G4endl;
633 G4Exception("G4PenelopePhotoElectricModel::GetNumberOfShellXS()",
634 "em2038",FatalException,ed);
635 }
636 //one vector is allocated for the _total_ cross section
637 std::size_t nEntries = fLogAtomicShellXS[Z]->entries();
638 return (nEntries-1);
639}
std::size_t entries() const

Referenced by GetShellCrossSection().

◆ GetShellCrossSection()

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

Definition at line 643 of file G4PenelopePhotoElectricModel.cc.

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

◆ GetVerbosityLevel()

G4int G4PenelopePhotoElectricModel::GetVerbosityLevel ( )
inline

Definition at line 84 of file G4PenelopePhotoElectricModel.hh.

84{return fVerboseLevel;};

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 117 of file G4PenelopePhotoElectricModel.cc.

119{
120 if (fVerboseLevel > 3)
121 G4cout << "Calling G4PenelopePhotoElectricModel::Initialise()" << G4endl;
122
123 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
124 //Issue warning if the AtomicDeexcitation has not been declared
125 if (!fAtomDeexcitation)
126 {
127 G4cout << G4endl;
128 G4cout << "WARNING from G4PenelopePhotoElectricModel " << G4endl;
129 G4cout << "Atomic de-excitation module is not instantiated, so there will not be ";
130 G4cout << "any fluorescence/Auger emission." << G4endl;
131 G4cout << "Please make sure this is intended" << G4endl;
132 }
133
134 SetParticle(particle);
135
136 //Only the master model creates/fills/destroys the tables
137 if (IsMaster() && particle == fParticle)
138 {
139 G4ProductionCutsTable* theCoupleTable =
141
142 for (G4int i=0;i<(G4int)theCoupleTable->GetTableSize();++i)
143 {
144 const G4Material* material =
145 theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
146 const G4ElementVector* theElementVector = material->GetElementVector();
147
148 for (std::size_t j=0;j<material->GetNumberOfElements();++j)
149 {
150 G4int iZ = theElementVector->at(j)->GetZasInt();
151 //read data files only in the master
152 if (!fLogAtomicShellXS[iZ])
153 ReadDataFile(iZ);
154 }
155 }
156
157 InitialiseElementSelectors(particle,cuts);
158
159 if (fVerboseLevel > 0) {
160 G4cout << "Penelope Photo-Electric model v2008 is initialized " << G4endl
161 << "Energy range: "
162 << LowEnergyLimit() / MeV << " MeV - "
163 << HighEnergyLimit() / GeV << " GeV";
164 }
165 }
166
167 if(fIsInitialised) return;
169 fIsInitialised = true;
170
171}
std::vector< const G4Element * > G4ElementVector
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
const G4Material * GetMaterial() const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:185
size_t GetNumberOfElements() const
Definition: G4Material.hh:181
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:124
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:139

◆ InitialiseLocal()

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

Reimplemented from G4VEmModel.

Definition at line 173 of file G4PenelopePhotoElectricModel.cc.

175{
176 if (fVerboseLevel > 3)
177 G4cout << "Calling G4PenelopePhotoElectricModel::InitialiseLocal()" << G4endl;
178 //
179 //Check that particle matches: one might have multiple master models (e.g.
180 //for e+ and e-).
181 //
182 if (part == fParticle)
183 {
185
186 //Get the const table pointers from the master to the workers
187 const G4PenelopePhotoElectricModel* theModel =
188 static_cast<G4PenelopePhotoElectricModel*> (masterModel);
189 for(G4int i=0; i<=fMaxZ; ++i)
190 fLogAtomicShellXS[i] = theModel->fLogAtomicShellXS[i];
191 //Same verbosity for all workers, as the master
192 fVerboseLevel = theModel->fVerboseLevel;
193 }
194
195 return;
196}
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:831
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:823

◆ operator=()

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 256 of file G4PenelopePhotoElectricModel.cc.

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

◆ SetVerbosityLevel()

void G4PenelopePhotoElectricModel::SetVerbosityLevel ( G4int  lev)
inline

Definition at line 83 of file G4PenelopePhotoElectricModel.hh.

83{fVerboseLevel = lev;};

Member Data Documentation

◆ fParticle

const G4ParticleDefinition* G4PenelopePhotoElectricModel::fParticle
protected

Definition at line 95 of file G4PenelopePhotoElectricModel.hh.

Referenced by Initialise(), and InitialiseLocal().

◆ fParticleChange

G4ParticleChangeForGamma* G4PenelopePhotoElectricModel::fParticleChange
protected

Definition at line 94 of file G4PenelopePhotoElectricModel.hh.

Referenced by Initialise(), and SampleSecondaries().


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