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

#include <G4LivermorePolarizedPhotoElectricModel.hh>

+ Inheritance diagram for G4LivermorePolarizedPhotoElectricModel:

Public Member Functions

 G4LivermorePolarizedPhotoElectricModel (const G4String &nam="LivermorePolarizedPhotoElectric")
 
virtual ~G4LivermorePolarizedPhotoElectricModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX)
 
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
 

Protected Attributes

G4ParticleChangeForGammafParticleChange
 
- Protected Attributes inherited from G4VEmModel
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 

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 41 of file G4LivermorePolarizedPhotoElectricModel.hh.

Constructor & Destructor Documentation

◆ G4LivermorePolarizedPhotoElectricModel()

G4LivermorePolarizedPhotoElectricModel::G4LivermorePolarizedPhotoElectricModel ( const G4String nam = "LivermorePolarizedPhotoElectric")

Definition at line 49 of file G4LivermorePolarizedPhotoElectricModel.cc.

52 crossSectionHandler(0), shellCrossSectionHandler(0)
53{
54 lowEnergyLimit = 10 * eV; // SI - Could be 10 eV ?
55 highEnergyLimit = 100 * GeV;
56
57 verboseLevel= 0;
58 // Verbosity scale:
59 // 0 = nothing
60 // 1 = warning for energy non-conservation
61 // 2 = details of energy budget
62 // 3 = calculation of cross sections, file openings, sampling of atoms
63 // 4 = entering in methods
64
65 theGamma = G4Gamma::Gamma();
66 theElectron = G4Electron::Electron();
67
68 // use default generator
70
71 // polarized generator needs fix
72 //SetAngularDistribution(new G4PhotoElectricAngularGeneratorPolarized());
74 fAtomDeexcitation = 0;
75 fDeexcitationActive = false;
76
77 if (verboseLevel > 1) {
78 G4cout << "Livermore Polarized PhotoElectric is constructed " << G4endl
79 << "Energy range: "
80 << lowEnergyLimit / keV << " keV - "
81 << highEnergyLimit / GeV << " GeV"
82 << G4endl;
83 }
84}
#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 SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:641
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:515

◆ ~G4LivermorePolarizedPhotoElectricModel()

G4LivermorePolarizedPhotoElectricModel::~G4LivermorePolarizedPhotoElectricModel ( )
virtual

Definition at line 88 of file G4LivermorePolarizedPhotoElectricModel.cc.

89{
90 delete crossSectionHandler;
91 delete shellCrossSectionHandler;
92}

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 146 of file G4LivermorePolarizedPhotoElectricModel.cc.

151{
152 G4double cs = crossSectionHandler->FindValue(G4lrint(Z), GammaEnergy);
153
154 if (verboseLevel > 3) {
155 G4cout << "G4LivermorePolarizedPhotoElectricModel::ComputeCrossSectionPerAtom()"
156 << G4endl;
157 G4cout << " E(keV)= " << GammaEnergy/keV << " Z= " << Z
158 << " CrossSection(barn)= "
159 << cs/barn << G4endl;
160 }
161 return cs;
162}
double G4double
Definition: G4Types.hh:64
G4double FindValue(G4int Z, G4double e) const
int G4lrint(double ad)
Definition: templates.hh:163

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 97 of file G4LivermorePolarizedPhotoElectricModel.cc.

100{
101 if (verboseLevel > 3) {
102 G4cout << "Calling G4LivermorePolarizedPhotoElectricModel::Initialise()" << G4endl;
103 }
104 if (crossSectionHandler)
105 {
106 crossSectionHandler->Clear();
107 delete crossSectionHandler;
108 }
109
110 if (shellCrossSectionHandler)
111 {
112 shellCrossSectionHandler->Clear();
113 delete shellCrossSectionHandler;
114 }
115
116
117 // Reading of data files - all materials are read
118 crossSectionHandler = new G4CrossSectionHandler;
119 crossSectionHandler->Clear();
120 G4String crossSectionFile = "phot/pe-cs-";
121 crossSectionHandler->LoadData(crossSectionFile);
122
123 shellCrossSectionHandler = new G4CrossSectionHandler();
124 shellCrossSectionHandler->Clear();
125 G4String shellCrossSectionFile = "phot/pe-ss-cs-";
126 shellCrossSectionHandler->LoadShellData(shellCrossSectionFile);
127
129
130 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
131 if(fAtomDeexcitation) {
132 fDeexcitationActive = fAtomDeexcitation->IsFluoActive();
133 }
134 if (verboseLevel > 1) {
135 G4cout << "Livermore Polarized PhotoElectric model is initialized "
136 << G4endl
137 << "Energy range: "
138 << LowEnergyLimit() / keV << " keV - "
139 << HighEnergyLimit() / GeV << " GeV"
140 << G4endl;
141 }
142}
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
void LoadShellData(const G4String &dataFile)
void LoadData(const G4String &dataFile)
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:109
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:529
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:522

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 166 of file G4LivermorePolarizedPhotoElectricModel.cc.

172{
173 if (verboseLevel > 3) {
174 G4cout << "Calling SampleSecondaries() of G4LivermorePolarizedPhotoElectricModel"
175 << G4endl;
176 }
177
178 G4double photonEnergy = aDynamicGamma->GetKineticEnergy();
179
180 // kill incident photon
183
184 // low-energy gamma is absorpted by this process
185
186 if (photonEnergy <= lowEnergyLimit)
187 {
189 return;
190 }
191
192 // Select randomly one element in the current material
193 //G4cout << "Select random atom Egamma(keV)= " << photonEnergy/keV << G4endl;
194 const G4Element* elm =
195 SelectRandomAtom(couple->GetMaterial(),theGamma,photonEnergy);
196 G4int Z = G4lrint(elm->GetZ());
197
198 // Select the ionised shell in the current atom according to shell cross sections
199 //G4cout << "Select random shell Z= " << Z << G4endl;
200 size_t shellIndex = shellCrossSectionHandler->SelectRandomShell(Z,photonEnergy);
201 //G4cout << "Shell index= " << shellIndex
202 // << " nShells= " << elm->GetNbOfAtomicShells() << G4endl;
204 const G4AtomicShell* shell = 0;
205 if(fDeexcitationActive) {
207 shell = fAtomDeexcitation->GetAtomicShell(Z, as);
208 bindingEnergy = shell->BindingEnergy();
209 } else {
210 G4int nshells = elm->GetNbOfAtomicShells() - 1;
211 if(G4int(shellIndex) > nshells) { shellIndex = std::max(0, nshells); }
212 bindingEnergy = elm->GetAtomicShell(shellIndex);
213 }
214
215 // There may be cases where the binding energy of the selected
216 // shell is > photon energy
217 // In such cases do not generate secondaries
218 if(photonEnergy < bindingEnergy) {
220 return;
221 }
222 //G4cout << "Z= " << Z << " shellIndex= " << shellIndex
223 // << " Ebind(keV)= " << bindingEnergy/keV << G4endl;
224
225
226 // Primary outcoming electron
227 G4double eKineticEnergy = photonEnergy - bindingEnergy;
228 G4double edep = bindingEnergy;
229
230 // Calculate direction of the photoelectron
231 G4ThreeVector electronDirection =
232 GetAngularDistribution()->SampleDirection(aDynamicGamma,
233 eKineticEnergy,
234 shellIndex,
235 couple->GetMaterial());
236
237 // The electron is created ...
238 G4DynamicParticle* electron = new G4DynamicParticle (theElectron,
239 electronDirection,
240 eKineticEnergy);
241 fvect->push_back(electron);
242
243 // Sample deexcitation
244 if(fDeexcitationActive) {
245 G4int index = couple->GetIndex();
246 if(fAtomDeexcitation->CheckDeexcitationActiveRegion(index)) {
247 size_t nbefore = fvect->size();
248
249 fAtomDeexcitation->GenerateParticles(fvect, shell, Z, index);
250 size_t nafter = fvect->size();
251 if(nafter > nbefore) {
252 for (size_t j=nbefore; j<nafter; ++j) {
253 edep -= ((*fvect)[j])->GetKineticEnergy();
254 }
255 }
256 }
257 }
258
259 // energy balance - excitation energy left
260 if(edep > 0.0) {
262 }
263}
G4AtomicShellEnumerator
@ fStopAndKill
int G4int
Definition: G4Types.hh:66
G4double BindingEnergy() const
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)
G4int SelectRandomShell(G4int Z, G4double e) const
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)

Member Data Documentation

◆ fParticleChange

G4ParticleChangeForGamma* G4LivermorePolarizedPhotoElectricModel::fParticleChange
protected

Definition at line 69 of file G4LivermorePolarizedPhotoElectricModel.hh.

Referenced by Initialise(), and SampleSecondaries().


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