Geant4 10.7.0
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")
 
virtual ~G4PEEffectFluoModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double, G4double) override
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
- 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 *, 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)
 
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)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
G4int SelectIsotopeNumber (const G4Element *)
 
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 *)
 
const G4ElementGetCurrentElement () const
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (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
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const G4MaterialpBaseMaterial
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
G4bool lossFlucFlag
 
G4double inveplus
 
G4double pFactor
 

Detailed Description

Definition at line 57 of file G4PEEffectFluoModel.hh.

Constructor & Destructor Documentation

◆ G4PEEffectFluoModel()

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*eV;
72 fParticleChange = nullptr;
73 fAtomDeexcitation = nullptr;
74
75 fSandiaCof.resize(4,0.0);
76
77 // default generator
79}
static G4Electron * Electron()
Definition: G4Electron.cc:93
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:813
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:618

◆ ~G4PEEffectFluoModel()

G4PEEffectFluoModel::~G4PEEffectFluoModel ( )
virtual

Definition at line 83 of file G4PEEffectFluoModel.cc.

84{}

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 energy2 = energy*energy;
117 G4double energy3 = energy*energy2;
118 G4double energy4 = energy2*energy2;
119
120 return fSandiaCof[0]/energy + fSandiaCof[1]/energy2 +
121 fSandiaCof[2]/energy3 + fSandiaCof[3]/energy4;
122}
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4Material * GetMaterial() const
G4SandiaTable * GetSandiaTable() const
Definition: G4Material.hh:227
void GetSandiaCofPerAtom(G4int Z, G4double energy, std::vector< G4double > &coeff) const
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:480
G4double energy(const ThreeVector &p, const G4double m)

Referenced by G4AdjointPhotoElectricModel::AdjointCrossSectionPerAtom().

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 127 of file G4PEEffectFluoModel.cc.

131{
132 // This method may be used only if G4MaterialCutsCouple pointer
133 // has been set properly
134 energy = std::max(energy, fMatEnergyTh[material->GetIndex()]);
135 const G4double* SandiaCof =
136 material->GetSandiaTable()->GetSandiaCofForMaterial(energy);
137
138 G4double energy2 = energy*energy;
139 G4double energy3 = energy*energy2;
140 G4double energy4 = energy2*energy2;
141
142 return SandiaCof[0]/energy + SandiaCof[1]/energy2 +
143 SandiaCof[2]/energy3 + SandiaCof[3]/energy4;
144}
size_t GetIndex() const
Definition: G4Material.hh:258
G4double GetSandiaCofForMaterial(G4int, G4int) const

◆ Initialise()

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

Implements G4VEmModel.

Reimplemented in G4PolarizedPEEffectModel.

Definition at line 88 of file G4PEEffectFluoModel.cc.

90{
91 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
92 if(nullptr == fParticleChange) { fParticleChange = GetParticleChangeForGamma(); }
94 fMatEnergyTh.resize(nmat, 0.0);
95 for(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 G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:644
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:637
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:133

Referenced by G4PolarizedPEEffectModel::Initialise().

◆ SampleSecondaries()

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

Implements G4VEmModel.

Reimplemented in G4PolarizedPEEffectModel.

Definition at line 149 of file G4PEEffectFluoModel.cc.

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

Referenced by G4PolarizedPEEffectModel::SampleSecondaries().


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