Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4LivermorePolarizedPhotoElectricModel.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: G4LivermorePolarizedPhotoElectricModel.hh,v 1.2 2010-11-23 16:42:15 flongo Exp $
27// GEANT4 tag $Name: not supported by cvs2svn $
28//
29// Authors: G.Depaola & F.Longo
30//
31
34#include "G4SystemOfUnits.hh"
35#include "G4LossTableManager.hh"
39#include "G4AtomicShell.hh"
41#include "G4Gamma.hh"
42
43//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
44
45using namespace std;
46
47//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
48
50 const G4String& nam)
51 :G4VEmModel(nam),fParticleChange(0),
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}
85
86//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
87
89{
90 delete crossSectionHandler;
91 delete shellCrossSectionHandler;
92}
93
94//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
95
96void
99 const G4DataVector&)
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}
143
144//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
145
148 G4double GammaEnergy,
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}
163
164//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
165
167 std::vector<G4DynamicParticle*>* fvect,
168 const G4MaterialCutsCouple* couple,
169 const G4DynamicParticle* aDynamicGamma,
170 G4double,
171 G4double)
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;
203 G4double bindingEnergy;
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}
264
265//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
G4AtomicShellEnumerator
@ fStopAndKill
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
G4double BindingEnergy() const
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
G4LivermorePolarizedPhotoElectricModel(const G4String &nam="LivermorePolarizedPhotoElectric")
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX)
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
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)
void LoadShellData(const G4String &dataFile)
G4double FindValue(G4int Z, G4double e) const
void LoadData(const G4String &dataFile)
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
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:109
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:529
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:522
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