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
G4KleinNishinaModel.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: G4KleinNishinaModel
34//
35// Author: Vladimir Ivanchenko on base of G4KleinNishinaCompton
36//
37// Creation date: 13.06.2010
38//
39// Modifications:
40//
41// Class Description:
42//
43// -------------------------------------------------------------------
44//
45//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
46//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
47
50#include "G4SystemOfUnits.hh"
51#include "G4Electron.hh"
52#include "G4Gamma.hh"
53#include "Randomize.hh"
54#include "G4RandomDirection.hh"
55#include "G4DataVector.hh"
58#include "G4AtomicShells.hh"
59#include "G4LossTableManager.hh"
60
61//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
62
63using namespace std;
64
66 : G4VEmModel(nam)
67{
70 lowestGammaEnergy = 1.0*eV;
71 limitFactor = 4;
72 fProbabilities.resize(9,0.0);
75 fAtomDeexcitation = 0;
76}
77
78//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
79
81{}
82
83//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
84
86 const G4DataVector& cuts)
87{
88 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
91}
92
93//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
94
97 G4double GammaEnergy,
100{
101 G4double xSection = 0.0 ;
102 if ( Z < 0.9999 || GammaEnergy < 0.1*keV) { return xSection; }
103
104 static const G4double a = 20.0 , b = 230.0 , c = 440.0;
105
106 static const G4double
107 d1= 2.7965e-1*barn, d2=-1.8300e-1*barn, d3= 6.7527 *barn, d4=-1.9798e+1*barn,
108 e1= 1.9756e-5*barn, e2=-1.0205e-2*barn, e3=-7.3913e-2*barn, e4= 2.7079e-2*barn,
109 f1=-3.9178e-7*barn, f2= 6.8241e-5*barn, f3= 6.0480e-5*barn, f4= 3.0274e-4*barn;
110
111 G4double p1Z = Z*(d1 + e1*Z + f1*Z*Z), p2Z = Z*(d2 + e2*Z + f2*Z*Z),
112 p3Z = Z*(d3 + e3*Z + f3*Z*Z), p4Z = Z*(d4 + e4*Z + f4*Z*Z);
113
114 G4double T0 = 15.0*keV;
115 if (Z < 1.5) { T0 = 40.0*keV; }
116
117 G4double X = max(GammaEnergy, T0) / electron_mass_c2;
118 xSection = p1Z*std::log(1.+2.*X)/X
119 + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
120
121 // modification for low energy. (special case for Hydrogen)
122 if (GammaEnergy < T0) {
123 G4double dT0 = keV;
124 X = (T0+dT0) / electron_mass_c2 ;
125 G4double sigma = p1Z*log(1.+2*X)/X
126 + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
127 G4double c1 = -T0*(sigma-xSection)/(xSection*dT0);
128 G4double c2 = 0.150;
129 if (Z > 1.5) { c2 = 0.375-0.0556*log(Z); }
130 G4double y = log(GammaEnergy/T0);
131 xSection *= exp(-y*(c1+c2*y));
132 }
133
134 if(xSection < 0.0) { xSection = 0.0; }
135 // G4cout << "e= " << GammaEnergy << " Z= " << Z
136 // << " cross= " << xSection << G4endl;
137 return xSection;
138}
139
140//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
141
143 std::vector<G4DynamicParticle*>* fvect,
144 const G4MaterialCutsCouple* couple,
145 const G4DynamicParticle* aDynamicGamma,
146 G4double,
147 G4double)
148{
149 // primary gamma
150 G4double energy = aDynamicGamma->GetKineticEnergy();
151 G4ThreeVector direction = aDynamicGamma->GetMomentumDirection();
152
153 // select atom
154 const G4Element* elm = SelectRandomAtom(couple, theGamma, energy);
155
156 // select shell first
157 G4int nShells = elm->GetNbOfAtomicShells();
158 if(nShells > (G4int)fProbabilities.size()) { fProbabilities.resize(nShells); }
159 G4double totprob = 0.0;
160 G4int i;
161 for(i=0; i<nShells; ++i) {
162 //G4double bindingEnergy = elm->GetAtomicShell(i);
163 totprob += elm->GetNbOfShellElectrons(i);
164 //totprob += elm->GetNbOfShellElectrons(i)/(bindingEnergy*bindingEnergy);
165 fProbabilities[i] = totprob;
166 }
167 //if(totprob == 0.0) { return; }
168
169 // Loop on sampling
170 // const G4int nlooplim = 100;
171 //G4int nloop = 0;
172
173 G4double bindingEnergy, ePotEnergy, eKinEnergy;
174 G4double gamEnergy0, gamEnergy1;
175
176 //static const G4double eminus2 = 1.0 - exp(-2.0);
177
178 do {
179 //++nloop;
180 G4double xprob = totprob*G4UniformRand();
181
182 // select shell
183 for(i=0; i<nShells; ++i) { if(xprob <= fProbabilities[i]) { break; } }
184
185 bindingEnergy = elm->GetAtomicShell(i);
186 // ePotEnergy = bindingEnergy;
187 // gamEnergy0 = energy;
188 lv1.set(0.0,0.0,energy,energy);
189
190 //G4cout << "nShells= " << nShells << " i= " << i
191 // << " Egamma= " << energy << " Ebind= " << bindingEnergy
192 // << " Elim= " << limitEnergy
193 // << G4endl;
194
195 // for rest frame of the electron
196 G4double x = -log(G4UniformRand());
197 eKinEnergy = bindingEnergy*x;
198 ePotEnergy = bindingEnergy*(1.0 + x);
199
200 // for rest frame of the electron
201 G4double eTotMomentum = sqrt(eKinEnergy*(eKinEnergy + 2*electron_mass_c2));
202 G4double phi = G4UniformRand()*twopi;
203 G4double costet = 2*G4UniformRand() - 1;
204 G4double sintet = sqrt((1 - costet)*(1 + costet));
205 lv2.set(eTotMomentum*sintet*cos(phi),eTotMomentum*sintet*sin(phi),
206 eTotMomentum*costet,eKinEnergy + electron_mass_c2);
207 bst = lv2.boostVector();
208 lv1.boost(-bst);
209
210 gamEnergy0 = lv1.e();
211
212 // In the rest frame of the electron
213 // The scattered gamma energy is sampled according to Klein-Nishina formula
214 // The random number techniques of Butcher & Messel are used
215 // (Nuc Phys 20(1960),15).
216 G4double E0_m = gamEnergy0/electron_mass_c2;
217
218 //
219 // sample the energy rate of the scattered gamma
220 //
221
222 G4double epsilon, epsilonsq, onecost, sint2, greject ;
223
224 G4double eps0 = 1./(1 + 2*E0_m);
225 G4double epsilon0sq = eps0*eps0;
226 G4double alpha1 = - log(eps0);
227 G4double alpha2 = 0.5*(1 - epsilon0sq);
228
229 do {
230 if ( alpha1/(alpha1+alpha2) > G4UniformRand() ) {
231 epsilon = exp(-alpha1*G4UniformRand()); // epsilon0**r
232 epsilonsq = epsilon*epsilon;
233
234 } else {
235 epsilonsq = epsilon0sq + (1.- epsilon0sq)*G4UniformRand();
236 epsilon = sqrt(epsilonsq);
237 }
238
239 onecost = (1.- epsilon)/(epsilon*E0_m);
240 sint2 = onecost*(2.-onecost);
241 greject = 1. - epsilon*sint2/(1.+ epsilonsq);
242
243 } while (greject < G4UniformRand());
244 gamEnergy1 = epsilon*gamEnergy0;
245
246 // before scattering total 4-momentum in e- system
247 lv2.set(0.0,0.0,0.0,electron_mass_c2);
248 lv2 += lv1;
249
250 //
251 // scattered gamma angles. ( Z - axis along the parent gamma)
252 //
253 if(sint2 < 0.0) { sint2 = 0.0; }
254 costet = 1. - onecost;
255 sintet = sqrt(sint2);
256 phi = twopi * G4UniformRand();
257
258 // e- recoil
259 //
260 // in rest frame of the electron
261 G4ThreeVector gamDir = lv1.vect().unit();
262 G4ThreeVector v = G4ThreeVector(sintet*cos(phi),sintet*sin(phi),costet);
263 v.rotateUz(gamDir);
264 lv1.set(gamEnergy1*v.x(),gamEnergy1*v.y(),gamEnergy1*v.z(),gamEnergy1);
265 lv2 -= lv1;
266 //G4cout<<"Egam= "<<lv1.e()<<" Ee= "<< lv2.e()-electron_mass_c2 << G4endl;
267 lv2.boost(bst);
268 eKinEnergy = lv2.e() - electron_mass_c2 - ePotEnergy;
269 //G4cout << "eKinEnergy= " << eKinEnergy << G4endl;
270
271 } while ( eKinEnergy < 0.0 );
272
273 //
274 // update G4VParticleChange for the scattered gamma
275 //
276
277 lv1.boost(bst);
278 gamEnergy1 = lv1.e();
279 if(gamEnergy1 > lowestGammaEnergy) {
280 G4ThreeVector gamDirection1 = lv1.vect().unit();
281 gamDirection1.rotateUz(direction);
283 } else {
285 gamEnergy1 = 0.0;
286 }
288
289 //
290 // kinematic of the scattered electron
291 //
292
293 if(eKinEnergy > lowestGammaEnergy) {
294 G4ThreeVector eDirection = lv2.vect().unit();
295 eDirection.rotateUz(direction);
296 G4DynamicParticle* dp =
297 new G4DynamicParticle(theElectron,eDirection,eKinEnergy);
298 fvect->push_back(dp);
299 } else { eKinEnergy = 0.0; }
300
301 G4double edep = energy - gamEnergy1 - eKinEnergy;
302
303 // sample deexcitation
304 //
305 if(fAtomDeexcitation) {
306 G4int index = couple->GetIndex();
307 if(fAtomDeexcitation->CheckDeexcitationActiveRegion(index)) {
308 G4int Z = G4lrint(elm->GetZ());
310 const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as);
311 size_t nbefore = fvect->size();
312 fAtomDeexcitation->GenerateParticles(fvect, shell, Z, index);
313 size_t nafter = fvect->size();
314 if(nafter > nbefore) {
315 for (size_t j=nbefore; j<nafter; ++j) {
316 edep -= ((*fvect)[j])->GetKineticEnergy();
317 }
318 }
319 }
320 }
321 // energy balance
322 if(edep < 0.0) { edep = 0.0; }
324}
325
326//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
327
G4AtomicShellEnumerator
CLHEP::Hep3Vector G4ThreeVector
@ fStopAndKill
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4UniformRand()
Definition: Randomize.hh:53
double z() const
Hep3Vector unit() const
double x() const
double y() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
Hep3Vector boostVector() const
HepLorentzVector & boost(double, double, double)
Hep3Vector vect() const
void set(double x, double y, double z, double t)
const G4ThreeVector & GetMomentumDirection() 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
G4int GetNbOfShellElectrons(G4int index) const
Definition: G4Element.cc:383
G4double GetAtomicShell(G4int index) const
Definition: G4Element.cc:367
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
G4ParticleDefinition * theElectron
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, G4double cut, G4double emax)
G4KleinNishinaModel(const G4String &nam="KleinNishina")
G4ParticleChangeForGamma * fParticleChange
G4ParticleDefinition * theGamma
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
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)
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 InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:123
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
int G4lrint(double ad)
Definition: templates.hh:163