Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PenelopeAnnihilationModel.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// Author: Luciano Pandola
29//
30// History:
31// --------
32// 29 Oct 2008 L Pandola Migration from process to model
33// 15 Apr 2009 V Ivanchenko Cleanup initialisation and generation of secondaries:
34// - apply internal high-energy limit only in constructor
35// - do not apply low-energy limit (default is 0)
36// - do not use G4ElementSelector
37
40#include "G4SystemOfUnits.hh"
44#include "G4DynamicParticle.hh"
45#include "G4Gamma.hh"
46
47//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
48
49
51 const G4String& nam)
52 :G4VEmModel(nam),fParticleChange(0),isInitialised(false)
53{
54 fIntrinsicLowEnergyLimit = 0.0;
55 fIntrinsicHighEnergyLimit = 100.0*GeV;
56 // SetLowEnergyLimit(fIntrinsicLowEnergyLimit);
57 SetHighEnergyLimit(fIntrinsicHighEnergyLimit);
58
59 //Calculate variable that will be used later on
60 fPielr2 = pi*classic_electr_radius*classic_electr_radius;
61
62 verboseLevel= 0;
63 // Verbosity scale:
64 // 0 = nothing
65 // 1 = warning for energy non-conservation
66 // 2 = details of energy budget
67 // 3 = calculation of cross sections, file openings, sampling of atoms
68 // 4 = entering in methods
69
70}
71
72//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
73
75{;}
76
77//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
78
80 const G4DataVector&)
81{
82 if (verboseLevel > 3)
83 G4cout << "Calling G4PenelopeAnnihilationModel::Initialise()" << G4endl;
84
85 if(verboseLevel > 0) {
86 G4cout << "Penelope Annihilation model is initialized " << G4endl
87 << "Energy range: "
88 << LowEnergyLimit() / keV << " keV - "
89 << HighEnergyLimit() / GeV << " GeV"
90 << G4endl;
91 }
92
93 if(isInitialised) return;
95 isInitialised = true;
96}
97
98//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
99
102 G4double energy,
105{
106 if (verboseLevel > 3)
107 G4cout << "Calling ComputeCrossSectionPerAtom() of G4PenelopeAnnihilationModel" <<
108 G4endl;
109
110 G4double cs = Z*ComputeCrossSectionPerElectron(energy);
111
112 if (verboseLevel > 2)
113 G4cout << "Annihilation cross Section at " << energy/keV << " keV for Z=" << Z <<
114 " = " << cs/barn << " barn" << G4endl;
115 return cs;
116}
117
118//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
119
120void G4PenelopeAnnihilationModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
122 const G4DynamicParticle* aDynamicPositron,
123 G4double,
124 G4double)
125{
126 //
127 // Penelope model to sample final state for positron annihilation.
128 // Target eletrons are assumed to be free and at rest. Binding effects enabling
129 // one-photon annihilation are neglected.
130 // For annihilation at rest, two back-to-back photons are emitted, having energy of 511 keV
131 // and isotropic angular distribution.
132 // For annihilation in flight, it is used the theory from
133 // W. Heitler, The quantum theory of radiation, Oxford University Press (1954)
134 // The two photons can have different energy. The efficiency of the sampling algorithm
135 // of the photon energy from the dSigma/dE distribution is practically 100% for
136 // positrons of kinetic energy < 10 keV. It reaches a minimum (about 80%) at energy
137 // of about 10 MeV.
138 // The angle theta is kinematically linked to the photon energy, to ensure momentum
139 // conservation. The angle phi is sampled isotropically for the first gamma.
140 //
141 if (verboseLevel > 3)
142 G4cout << "Calling SamplingSecondaries() of G4PenelopeAnnihilationModel" << G4endl;
143
144 G4double kineticEnergy = aDynamicPositron->GetKineticEnergy();
145
146 // kill primary
149
150 if (kineticEnergy == 0.0)
151 {
152 //Old AtRestDoIt
153 G4double cosTheta = -1.0+2.0*G4UniformRand();
154 G4double sinTheta = std::sqrt(1.0-cosTheta*cosTheta);
155 G4double phi = twopi*G4UniformRand();
156 G4ThreeVector direction (sinTheta*std::cos(phi),sinTheta*std::sin(phi),cosTheta);
158 direction, electron_mass_c2);
160 -direction, electron_mass_c2);
161
162 fvect->push_back(firstGamma);
163 fvect->push_back(secondGamma);
164 return;
165 }
166
167 //This is the "PostStep" case (annihilation in flight)
168 G4ParticleMomentum positronDirection =
169 aDynamicPositron->GetMomentumDirection();
170 G4double gamma = 1.0 + std::max(kineticEnergy,1.0*eV)/electron_mass_c2;
171 G4double gamma21 = std::sqrt(gamma*gamma-1);
172 G4double ani = 1.0+gamma;
173 G4double chimin = 1.0/(ani+gamma21);
174 G4double rchi = (1.0-chimin)/chimin;
175 G4double gt0 = ani*ani-2.0;
176 G4double test=0.0;
177 G4double epsilon = 0;
178 do{
179 epsilon = chimin*std::pow(rchi,G4UniformRand());
180 G4double reject = ani*ani*(1.0-epsilon)+2.0*gamma-(1.0/epsilon);
181 test = G4UniformRand()*gt0-reject;
182 }while(test>0);
183
184 G4double totalAvailableEnergy = kineticEnergy + 2.0*electron_mass_c2;
185 G4double photon1Energy = epsilon*totalAvailableEnergy;
186 G4double photon2Energy = (1.0-epsilon)*totalAvailableEnergy;
187 G4double cosTheta1 = (ani-1.0/epsilon)/gamma21;
188 G4double cosTheta2 = (ani-1.0/(1.0-epsilon))/gamma21;
189
190 //G4double localEnergyDeposit = 0.;
191
192 G4double sinTheta1 = std::sqrt(1.-cosTheta1*cosTheta1);
193 G4double phi1 = twopi * G4UniformRand();
194 G4double dirx1 = sinTheta1 * std::cos(phi1);
195 G4double diry1 = sinTheta1 * std::sin(phi1);
196 G4double dirz1 = cosTheta1;
197
198 G4double sinTheta2 = std::sqrt(1.-cosTheta2*cosTheta2);
199 G4double phi2 = phi1+pi;
200 G4double dirx2 = sinTheta2 * std::cos(phi2);
201 G4double diry2 = sinTheta2 * std::sin(phi2);
202 G4double dirz2 = cosTheta2;
203
204 G4ThreeVector photon1Direction (dirx1,diry1,dirz1);
205 photon1Direction.rotateUz(positronDirection);
206 // create G4DynamicParticle object for the particle1
208 photon1Direction,
209 photon1Energy);
210 fvect->push_back(aParticle1);
211
212 G4ThreeVector photon2Direction(dirx2,diry2,dirz2);
213 photon2Direction.rotateUz(positronDirection);
214 // create G4DynamicParticle object for the particle2
216 photon2Direction,
217 photon2Energy);
218 fvect->push_back(aParticle2);
219
220 if (verboseLevel > 1)
221 {
222 G4cout << "-----------------------------------------------------------" << G4endl;
223 G4cout << "Energy balance from G4PenelopeAnnihilation" << G4endl;
224 G4cout << "Kinetic positron energy: " << kineticEnergy/keV << " keV" << G4endl;
225 G4cout << "Total available energy: " << totalAvailableEnergy/keV << " keV " << G4endl;
226 G4cout << "-----------------------------------------------------------" << G4endl;
227 G4cout << "Photon energy 1: " << photon1Energy/keV << " keV" << G4endl;
228 G4cout << "Photon energy 2: " << photon2Energy/keV << " keV" << G4endl;
229 G4cout << "Total final state: " << (photon1Energy+photon2Energy)/keV <<
230 " keV" << G4endl;
231 G4cout << "-----------------------------------------------------------" << G4endl;
232 }
233 if (verboseLevel > 0)
234 {
235 G4double energyDiff = std::fabs(totalAvailableEnergy-photon1Energy-photon2Energy);
236 if (energyDiff > 0.05*keV)
237 G4cout << "Warning from G4PenelopeAnnihilation: problem with energy conservation: " <<
238 (photon1Energy+photon2Energy)/keV <<
239 " keV (final) vs. " <<
240 totalAvailableEnergy/keV << " keV (initial)" << G4endl;
241 }
242 return;
243}
244
245//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
246
247G4double G4PenelopeAnnihilationModel:: ComputeCrossSectionPerElectron(G4double energy)
248{
249 //
250 // Penelope model to calculate cross section for positron annihilation.
251 // The annihilation cross section per electron is calculated according
252 // to the Heitler formula
253 // W. Heitler, The quantum theory of radiation, Oxford University Press (1954)
254 // in the assumptions of electrons free and at rest.
255 //
256 G4double gamma = 1.0+std::max(energy,1.0*eV)/electron_mass_c2;
257 G4double gamma2 = gamma*gamma;
258 G4double f2 = gamma2-1.0;
259 G4double f1 = std::sqrt(f2);
260 G4double crossSection = fPielr2*((gamma2+4.0*gamma+1.0)*std::log(gamma+f1)/f2
261 - (gamma+3.0)/f1)/(gamma+1.0);
262 return crossSection;
263}
@ fStopAndKill
double G4double
Definition: G4Types.hh:64
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
void SetProposedKineticEnergy(G4double proposedKinEnergy)
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)
G4PenelopeAnnihilationModel(const G4ParticleDefinition *p=0, const G4String &processName="PenAnnih")
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4ParticleChangeForGamma * fParticleChange
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:109
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:529
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:522
void ProposeTrackStatus(G4TrackStatus status)