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

#include <G4PEEffectModel.hh>

+ Inheritance diagram for G4PEEffectModel:

Public Member Functions

 G4PEEffectModel (const G4ParticleDefinition *p=0, const G4String &nam="PhotoElectric")
 
virtual ~G4PEEffectModel ()
 
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 61 of file G4PEEffectModel.hh.

Constructor & Destructor Documentation

◆ G4PEEffectModel()

G4PEEffectModel::G4PEEffectModel ( const G4ParticleDefinition p = 0,
const G4String nam = "PhotoElectric" 
)

Definition at line 66 of file G4PEEffectModel.cc.

68 : G4VEmModel(nam)
69{
70 G4cout << "### G4PEEffectModel is obsolete "
71 << "and will be removed for the next release." << G4endl;
72
73 theGamma = G4Gamma::Gamma();
74 theElectron = G4Electron::Electron();
75 fminimalEnergy = 1.0*eV;
76 fParticleChange = 0;
77
78 // default generator
80}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
static G4Electron * Electron()
Definition: G4Electron.cc:94
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:515

◆ ~G4PEEffectModel()

G4PEEffectModel::~G4PEEffectModel ( )
virtual

Definition at line 84 of file G4PEEffectModel.cc.

85{}

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 99 of file G4PEEffectModel.cc.

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

Referenced by G4AdjointPhotoElectricModel::AdjointCrossSectionPerAtom().

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 116 of file G4PEEffectModel.cc.

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

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 89 of file G4PEEffectModel.cc.

91{
92 // always false before the run
94 if(!fParticleChange) { fParticleChange = GetParticleChangeForGamma(); }
95}
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:109
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:641

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 134 of file G4PEEffectModel.cc.

139{
140 const G4Material* aMaterial = couple->GetMaterial();
141
142 G4double energy = aDynamicPhoton->GetKineticEnergy();
143 G4ParticleMomentum PhotonDirection = aDynamicPhoton->GetMomentumDirection();
144
145 // select randomly one element constituing the material.
146 const G4Element* anElement = SelectRandomAtom(aMaterial,theGamma,energy);
147
148 //
149 // Photo electron
150 //
151
152 // Select atomic shell
153 G4int nShells = anElement->GetNbOfAtomicShells();
154 G4int i = 0;
155 while ((i<nShells) && (energy<anElement->GetAtomicShell(i))) { ++i; }
156
157 G4double edep = energy;
158
159 // no shell available
160 if (i < nShells) {
161
162 G4double bindingEnergy = anElement->GetAtomicShell(i);
163 G4double elecKineEnergy = energy - bindingEnergy;
164
165 if (elecKineEnergy > fminimalEnergy) {
166
167 edep = bindingEnergy;
168 G4ThreeVector elecDirection =
169 GetAngularDistribution()->SampleDirection(aDynamicPhoton,
170 elecKineEnergy,
171 i,
172 couple->GetMaterial());
173
174 G4DynamicParticle* aParticle =
175 new G4DynamicParticle(theElectron, elecDirection, elecKineEnergy);
176 fvect->push_back(aParticle);
177
178 }
179 }
180
181 fParticleChange->SetProposedKineticEnergy(0.);
182 fParticleChange->ProposeTrackStatus(fStopAndKill);
183 if(edep > 0.0) {
184 fParticleChange->ProposeLocalEnergyDeposit(edep);
185 }
186}
@ fStopAndKill
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4int GetNbOfAtomicShells() const
Definition: G4Element.hh:146
G4double GetAtomicShell(G4int index) const
Definition: G4Element.cc:367
const G4Material * GetMaterial() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
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)

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