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

#include <G4PEEffectFluoModel.hh>

+ Inheritance diagram for G4PEEffectFluoModel:

Public Member Functions

 G4PEEffectFluoModel (const G4String &nam="PhotoElectric")
 
 ~G4PEEffectFluoModel () override
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double, G4double) override
 
G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
G4PEEffectFluoModeloperator= (const G4PEEffectFluoModel &right)=delete
 
 G4PEEffectFluoModel (const G4PEEffectFluoModel &)=delete
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
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 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 *, 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 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 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)
 
void SetLPMFlag (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

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 *)
 
- 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
 
std::size_t currentCoupleIndex = 0
 
std::size_t basedCoupleIndex = 0
 
G4bool lossFlucFlag = true
 

Detailed Description

Definition at line 57 of file G4PEEffectFluoModel.hh.

Constructor & Destructor Documentation

◆ G4PEEffectFluoModel() [1/2]

G4PEEffectFluoModel::G4PEEffectFluoModel ( const G4String & nam = "PhotoElectric")
explicit

Definition at line 65 of file G4PEEffectFluoModel.cc.

66 : G4VEmModel(nam)
67{
68 theGamma = G4Gamma::Gamma();
69 theElectron = G4Electron::Electron();
70 fminimalEnergy = 1.0*CLHEP::eV;
72
73 fSandiaCof.resize(4,0.0);
74
75 // default generator
77}
static G4Electron * Electron()
Definition G4Electron.cc:91
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
void SetDeexcitationFlag(G4bool val)
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67
void SetAngularDistribution(G4VEmAngularDistribution *)

◆ ~G4PEEffectFluoModel()

G4PEEffectFluoModel::~G4PEEffectFluoModel ( )
overridedefault

◆ G4PEEffectFluoModel() [2/2]

G4PEEffectFluoModel::G4PEEffectFluoModel ( const G4PEEffectFluoModel & )
delete

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

G4double G4PEEffectFluoModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition * ,
G4double kinEnergy,
G4double Z,
G4double A,
G4double ,
G4double  )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 106 of file G4PEEffectFluoModel.cc.

110{
111 // This method may be used only if G4MaterialCutsCouple pointer
112 // has been set properly
114 ->GetSandiaTable()->GetSandiaCofPerAtom((G4int)Z, energy, fSandiaCof);
115
116 G4double x1 = 1 / energy;
117
118 return x1 * (fSandiaCof[0] + x1 * (fSandiaCof[1] +
119 x1 * (fSandiaCof[2] + x1 * fSandiaCof[3])));
120}
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
const G4Material * GetMaterial() const
G4SandiaTable * GetSandiaTable() const
void GetSandiaCofPerAtom(G4int Z, G4double energy, std::vector< G4double > &coeff) const
const G4MaterialCutsCouple * CurrentCouple() const
G4double energy(const ThreeVector &p, const G4double m)

◆ CrossSectionPerVolume()

G4double G4PEEffectFluoModel::CrossSectionPerVolume ( const G4Material * material,
const G4ParticleDefinition * ,
G4double kineticEnergy,
G4double cutEnergy,
G4double maxEnergy )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 125 of file G4PEEffectFluoModel.cc.

129{
130 // This method may be used only if G4MaterialCutsCouple pointer
131 // has been set properly
132 energy = std::max(energy, fMatEnergyTh[material->GetIndex()]);
133 const G4double* SandiaCof =
134 material->GetSandiaTable()->GetSandiaCofForMaterial(energy);
135
136 G4double x1 = 1 / energy;
137
138 return x1 * (SandiaCof[0] + x1 * (SandiaCof[1] +
139 x1 * (SandiaCof[2] + x1 * SandiaCof[3])));
140}
std::size_t GetIndex() const
G4double GetSandiaCofForMaterial(G4int, G4int) const

◆ Initialise()

void G4PEEffectFluoModel::Initialise ( const G4ParticleDefinition * ,
const G4DataVector &  )
overridevirtual

Implements G4VEmModel.

Reimplemented in G4PolarizedPhotoElectricModel.

Definition at line 85 of file G4PEEffectFluoModel.cc.

87{
88 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
90 if(nullptr == fParticleChange) {
91 fParticleChange = GetParticleChangeForGamma();
92 }
93 std::size_t nmat = G4Material::GetNumberOfMaterials();
94 fMatEnergyTh.resize(nmat, 0.0);
95 for(std::size_t i=0; i<nmat; ++i) {
96 fMatEnergyTh[i] = (*(G4Material::GetMaterialTable()))[i]
97 ->GetSandiaTable()->GetSandiaCofForMaterial(0, 0);
98 //G4cout << "G4PEEffectFluoModel::Initialise Eth(eV)= "
99 // << fMatEnergyTh[i]/eV << G4endl;
100 }
101}
static G4EmParameters * Instance()
G4bool PhotoeffectBelowKShell() const
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
static std::size_t GetNumberOfMaterials()
static G4MaterialTable * GetMaterialTable()
G4ParticleChangeForGamma * GetParticleChangeForGamma()

Referenced by G4PolarizedPhotoElectricModel::Initialise().

◆ operator=()

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Reimplemented in G4PolarizedPhotoElectricModel.

Definition at line 145 of file G4PEEffectFluoModel.cc.

150{
151 SetCurrentCouple(couple);
152 const G4Material* aMaterial = couple->GetMaterial();
153
154 G4double energy = aDynamicPhoton->GetKineticEnergy();
155
156 // select randomly one element constituing the material.
157 const G4Element* anElement = SelectRandomAtom(aMaterial,theGamma,energy);
158
159 //
160 // Photo electron
161 //
162
163 // Select atomic shell
164 G4int nShells = anElement->GetNbOfAtomicShells();
165 G4int i = 0;
166 for(; i<nShells; ++i) {
167 /*
168 G4cout << "i= " << i << " E(eV)= " << energy/eV
169 << " Eb(eV)= " << anElement->GetAtomicShell(i)/eV
170 << " " << anElement->GetName()
171 << G4endl;
172 */
173 if(energy >= anElement->GetAtomicShell(i)) { break; }
174 }
175
176 G4double edep = energy;
177
178 // photo-electron is not sampled if shell is not found or
179 // the flag of photoeffect is "false" and shell is no K
180 if ( (fPEBelowKShell || 0 == i) && i < nShells ) {
181
182 G4double bindingEnergy = anElement->GetAtomicShell(i);
183 edep = bindingEnergy;
184 G4double esec = 0.0;
185
186 // sample deexcitation cascade
187 //
188 if(nullptr != fAtomDeexcitation) {
189 G4int index = couple->GetIndex();
190 if(fAtomDeexcitation->CheckDeexcitationActiveRegion(index)) {
191 G4int Z = G4lrint(anElement->GetZ());
192 auto as = (G4AtomicShellEnumerator)(i);
193 const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as);
194 G4double eshell = shell->BindingEnergy();
195 if(eshell > bindingEnergy && eshell <= energy) {
196 bindingEnergy = eshell;
197 edep = eshell;
198 }
199 std::size_t nbefore = fvect->size();
200 fAtomDeexcitation->GenerateParticles(fvect, shell, Z, index);
201 std::size_t nafter = fvect->size();
202 for (std::size_t j=nbefore; j<nafter; ++j) {
203 G4double e = ((*fvect)[j])->GetKineticEnergy();
204 if(esec + e > edep) {
205 // correct energy in order to have energy balance
206 e = edep - esec;
207 ((*fvect)[j])->SetKineticEnergy(e);
208 esec += e;
209 /*
210 G4cout << "### G4PEffectFluoModel Edep(eV)= " << edep/eV
211 << " Esec(eV)= " << esec/eV
212 << " E["<< j << "](eV)= " << e/eV
213 << " N= " << nafter
214 << " Z= " << Z << " shell= " << i
215 << " Ebind(keV)= " << bindingEnergy/keV
216 << " Eshell(keV)= " << shell->BindingEnergy()/keV
217 << G4endl;
218 */
219 // delete the rest of secondaries (should not happens)
220 for (std::size_t jj=nafter-1; jj>j; --jj) {
221 delete (*fvect)[jj];
222 fvect->pop_back();
223 }
224 break;
225 }
226 esec += e;
227 }
228 edep -= esec;
229 }
230 }
231 // create photo electron
232 //
233 G4double elecKineEnergy = energy - bindingEnergy;
234 if (elecKineEnergy > fminimalEnergy) {
235 auto aParticle = new G4DynamicParticle(theElectron,
236 GetAngularDistribution()->SampleDirection(aDynamicPhoton,
237 elecKineEnergy,
238 i, couple->GetMaterial()),
239 elecKineEnergy);
240 fvect->push_back(aParticle);
241 } else {
242 edep += elecKineEnergy;
243 elecKineEnergy = 0.0;
244 }
245 if(std::abs(energy - elecKineEnergy - esec - edep) > CLHEP::eV) {
246 G4cout << "### G4PEffectFluoModel dE(eV)= "
247 << (energy - elecKineEnergy - esec - edep)/eV
248 << " shell= " << i
249 << " E(keV)= " << energy/keV
250 << " Ebind(keV)= " << bindingEnergy/keV
251 << " Ee(keV)= " << elecKineEnergy/keV
252 << " Esec(keV)= " << esec/keV
253 << " Edep(keV)= " << edep/keV
254 << G4endl;
255 }
256 }
257
258 // kill primary photon
259 fParticleChange->SetProposedKineticEnergy(0.);
260 fParticleChange->ProposeTrackStatus(fStopAndKill);
261 if(edep > 0.0) {
262 fParticleChange->ProposeLocalEnergyDeposit(edep);
263 }
264}
@ fStopAndKill
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4double BindingEnergy() const
G4double GetKineticEnergy() const
G4double GetZ() const
Definition G4Element.hh:119
G4int GetNbOfAtomicShells() const
Definition G4Element.hh:134
G4double GetAtomicShell(G4int index) const
Definition G4Element.cc:361
void SetProposedKineticEnergy(G4double proposedKinEnergy)
G4bool CheckDeexcitationActiveRegion(G4int coupleIndex)
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
G4VEmAngularDistribution * GetAngularDistribution()
void SetCurrentCouple(const G4MaterialCutsCouple *)
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4double bindingEnergy(G4int A, G4int Z)
int G4lrint(double ad)
Definition templates.hh:134

Referenced by G4PolarizedPhotoElectricModel::SampleSecondaries().


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