Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PEEffectFluoModel.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// $Id$
27//
28// -------------------------------------------------------------------
29//
30// GEANT4 Class file
31//
32//
33// File name: G4PEEffectFluoModel
34//
35// Author: Vladimir Ivanchenko on base of G4PEEffectModel
36//
37// Creation date: 13.06.2010
38//
39// Modifications:
40//
41// Class Description:
42// Implementation of the photo-electric effect with deexcitation
43//
44// -------------------------------------------------------------------
45//
46//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
47//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
48
51#include "G4SystemOfUnits.hh"
52#include "G4Electron.hh"
53#include "G4Gamma.hh"
54#include "Randomize.hh"
55#include "G4DataVector.hh"
58#include "G4LossTableManager.hh"
60
61//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
62
63using namespace std;
64
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}
78
79//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
80
82{}
83
84//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
85
87 const G4DataVector&)
88{
89 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
90 if(!fParticleChange) { fParticleChange = GetParticleChangeForGamma(); }
91}
92
93//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
94
97 G4double energy,
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}
110
111//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
112
116 G4double energy,
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}
129
130//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
131
132void
133G4PEEffectFluoModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
134 const G4MaterialCutsCouple* couple,
135 const G4DynamicParticle* aDynamicPhoton,
136 G4double,
137 G4double)
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}
213
214//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4AtomicShellEnumerator
@ fStopAndKill
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
G4double GetKineticEnergy() const
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4double GetZ() const
Definition: G4Element.hh:131
G4int GetNbOfAtomicShells() const
Definition: G4Element.hh:146
G4double GetAtomicShell(G4int index) const
Definition: G4Element.cc:367
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
const G4Material * GetMaterial() const
G4SandiaTable * GetSandiaTable() const
Definition: G4Material.hh:228
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double, G4double)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
G4PEEffectFluoModel(const G4String &nam="PhotoElectric")
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
void SetProposedKineticEnergy(G4double proposedKinEnergy)
static G4double * GetSandiaCofPerAtom(G4int Z, G4double energy)
G4double GetSandiaCofForMaterial(G4int, G4int)
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
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:109
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:459
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:641
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:515
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
int G4lrint(double ad)
Definition: templates.hh:163