Geant4 9.6.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 &)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double, G4double)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
- 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
 

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
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 

Detailed Description

Definition at line 57 of file G4PEEffectFluoModel.hh.

Constructor & Destructor Documentation

◆ G4PEEffectFluoModel()

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

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 = 0;
73 fAtomDeexcitation = 0;
74
75 // default generator
77}
static G4Electron * Electron()
Definition: G4Electron.cc:94
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:641
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:515

◆ ~G4PEEffectFluoModel()

G4PEEffectFluoModel::~G4PEEffectFluoModel ( )
virtual

Definition at line 81 of file G4PEEffectFluoModel.cc.

82{}

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 96 of file G4PEEffectFluoModel.cc.

100{
101 G4double* SandiaCof = G4SandiaTable::GetSandiaCofPerAtom((G4int)Z, energy);
102
103 G4double energy2 = energy*energy;
104 G4double energy3 = energy*energy2;
105 G4double energy4 = energy2*energy2;
106
107 return SandiaCof[0]/energy + SandiaCof[1]/energy2 +
108 SandiaCof[2]/energy3 + SandiaCof[3]/energy4;
109}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
static G4double * GetSandiaCofPerAtom(G4int Z, G4double energy)

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 114 of file G4PEEffectFluoModel.cc.

118{
119 G4double* SandiaCof =
120 material->GetSandiaTable()->GetSandiaCofForMaterial(energy);
121
122 G4double energy2 = energy*energy;
123 G4double energy3 = energy*energy2;
124 G4double energy4 = energy2*energy2;
125
126 return SandiaCof[0]/energy + SandiaCof[1]/energy2 +
127 SandiaCof[2]/energy3 + SandiaCof[3]/energy4;
128}
G4SandiaTable * GetSandiaTable() const
Definition: G4Material.hh:228
G4double GetSandiaCofForMaterial(G4int, G4int)

◆ Initialise()

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

Implements G4VEmModel.

Reimplemented in G4PolarizedPEEffectModel.

Definition at line 86 of file G4PEEffectFluoModel.cc.

88{
89 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
90 if(!fParticleChange) { fParticleChange = GetParticleChangeForGamma(); }
91}
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:109

Referenced by G4PolarizedPEEffectModel::Initialise().

◆ SampleSecondaries()

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

Implements G4VEmModel.

Reimplemented in G4PolarizedPEEffectModel.

Definition at line 133 of file G4PEEffectFluoModel.cc.

138{
139 const G4Material* aMaterial = couple->GetMaterial();
140
141 G4double energy = aDynamicPhoton->GetKineticEnergy();
142
143 // select randomly one element constituing the material.
144 const G4Element* anElement = SelectRandomAtom(aMaterial,theGamma,energy);
145
146 //
147 // Photo electron
148 //
149
150 // Select atomic shell
151 G4int nShells = anElement->GetNbOfAtomicShells();
152 G4int i = 0;
153 for(; i<nShells; ++i) {
154 /*
155 G4cout << "i= " << i << " E(eV)= " << energy/eV
156 << " Eb(eV)= " << anElement->GetAtomicShell(i)/eV
157 << " " << anElement->GetName()
158 << G4endl;
159 */
160 if(energy >= anElement->GetAtomicShell(i)) { break; }
161 }
162
163 G4double edep = energy;
164
165 // Normally one shell is available
166 if (i < nShells) {
167
168 G4double bindingEnergy = anElement->GetAtomicShell(i);
169 G4double elecKineEnergy = energy - bindingEnergy;
170
171 // create photo electron
172 //
173 if (elecKineEnergy > fminimalEnergy) {
174 edep = bindingEnergy;
175 G4ThreeVector elecDirection =
176 GetAngularDistribution()->SampleDirection(aDynamicPhoton,
177 elecKineEnergy,
178 i,
179 couple->GetMaterial());
180
181 G4DynamicParticle* aParticle =
182 new G4DynamicParticle(theElectron, elecDirection, elecKineEnergy);
183 fvect->push_back(aParticle);
184 }
185
186 // sample deexcitation
187 //
188 if(fAtomDeexcitation) {
189 G4int index = couple->GetIndex();
190 if(fAtomDeexcitation->CheckDeexcitationActiveRegion(index)) {
191 G4int Z = G4lrint(anElement->GetZ());
193 const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as);
194 size_t nbefore = fvect->size();
195 fAtomDeexcitation->GenerateParticles(fvect, shell, Z, index);
196 size_t nafter = fvect->size();
197 if(nafter > nbefore) {
198 for (size_t j=nbefore; j<nafter; ++j) {
199 edep -= ((*fvect)[j])->GetKineticEnergy();
200 }
201 }
202 }
203 }
204 }
205
206 // kill primary photon
207 fParticleChange->SetProposedKineticEnergy(0.);
208 fParticleChange->ProposeTrackStatus(fStopAndKill);
209 if(edep > 0.0) {
210 fParticleChange->ProposeLocalEnergyDeposit(edep);
211 }
212}
G4AtomicShellEnumerator
@ fStopAndKill
G4double GetKineticEnergy() const
G4double GetZ() const
Definition: G4Element.hh:131
G4int GetNbOfAtomicShells() const
Definition: G4Element.hh:146
G4double GetAtomicShell(G4int index) const
Definition: G4Element.cc:367
const G4Material * GetMaterial() const
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)
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:508
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)
int G4lrint(double ad)
Definition: templates.hh:163

Referenced by G4PolarizedPEEffectModel::SampleSecondaries().


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