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
G4PEEffectModel.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: G4PEEffectModel
34//
35// Author: Vladimir Ivanchenko on base of Michel Maire code
36//
37// Creation date: 21.03.2005
38//
39// Modifications:
40//
41// 04.12.05 : SetProposedKineticEnergy(0.) for the killed photon (mma)
42// 20.02.09 : Added initialisation of deexcitation flag and method
43// CrossSectionPerVolume instead of mfp (V.Ivanchenko)
44//
45// Class Description:
46//
47// -------------------------------------------------------------------
48//
49//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
50//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
51
52#include "G4PEEffectModel.hh"
54#include "G4SystemOfUnits.hh"
55#include "G4Electron.hh"
56#include "G4Gamma.hh"
57#include "Randomize.hh"
58#include "G4DataVector.hh"
61
62//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
63
64using namespace std;
65
67 const G4String& nam)
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}
81
82//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
83
85{}
86
87//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
88
90 const G4DataVector&)
91{
92 // always false before the run
94 if(!fParticleChange) { fParticleChange = GetParticleChangeForGamma(); }
95}
96
97//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
98
100 G4double energy,
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}
113
114//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
115
118 G4double energy,
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}
131
132//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
133
134void G4PEEffectModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
135 const G4MaterialCutsCouple* couple,
136 const G4DynamicParticle* aDynamicPhoton,
137 G4double,
138 G4double)
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}
187
188//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
@ fStopAndKill
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4int GetNbOfAtomicShells() const
Definition: G4Element.hh:146
G4double GetAtomicShell(G4int index) const
Definition: G4Element.cc:367
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
const G4Material * GetMaterial() const
G4SandiaTable * GetSandiaTable() const
Definition: G4Material.hh:228
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
G4PEEffectModel(const G4ParticleDefinition *p=0, const G4String &nam="PhotoElectric")
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double, G4double)
virtual ~G4PEEffectModel()
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)
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)