Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4KleinNishinaCompton.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//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class file
30//
31//
32// File name: G4KleinNishinaCompton
33//
34// Author: Vladimir Ivanchenko on base of Michel Maire code
35//
36// Creation date: 15.03.2005
37//
38// Modifications:
39// 18-04-05 Use G4ParticleChangeForGamma (V.Ivantchenko)
40// 27-03-06 Remove upper limit of cross section (V.Ivantchenko)
41//
42// Class Description:
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"
57#include "G4Log.hh"
58#include "G4Exp.hh"
59
60//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
61
62using namespace std;
63
65 const G4String& nam)
66 : G4VEmModel(nam)
67{
70 lowestSecondaryEnergy = 100.0*eV;
71 fParticleChange = nullptr;
72}
73
74//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
75
77
78//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
79
81 const G4DataVector& cuts)
82{
83 if(IsMaster()) { InitialiseElementSelectors(p, cuts); }
84 if(nullptr == fParticleChange) {
86 }
87}
88
89//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
90
92 G4VEmModel* masterModel)
93{
95}
96
97//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
98
101 G4double GammaEnergy,
104{
105 G4double xSection = 0.0 ;
106 if (GammaEnergy <= LowEnergyLimit()) { return xSection; }
107
108 static const G4double a = 20.0 , b = 230.0 , c = 440.0;
109
110 static const G4double
111 d1= 2.7965e-1*CLHEP::barn, d2=-1.8300e-1*CLHEP::barn,
112 d3= 6.7527 *CLHEP::barn, d4=-1.9798e+1*CLHEP::barn,
113 e1= 1.9756e-5*CLHEP::barn, e2=-1.0205e-2*CLHEP::barn,
114 e3=-7.3913e-2*CLHEP::barn, e4= 2.7079e-2*CLHEP::barn,
115 f1=-3.9178e-7*CLHEP::barn, f2= 6.8241e-5*CLHEP::barn,
116 f3= 6.0480e-5*CLHEP::barn, f4= 3.0274e-4*CLHEP::barn;
117
118 G4double p1Z = Z*(d1 + e1*Z + f1*Z*Z), p2Z = Z*(d2 + e2*Z + f2*Z*Z),
119 p3Z = Z*(d3 + e3*Z + f3*Z*Z), p4Z = Z*(d4 + e4*Z + f4*Z*Z);
120
121 G4double T0 = 15.0*keV;
122 if (Z < 1.5) { T0 = 40.0*keV; }
123
124 G4double X = max(GammaEnergy, T0) / electron_mass_c2;
125 xSection = p1Z*G4Log(1.+2.*X)/X
126 + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
127
128 // modification for low energy. (special case for Hydrogen)
129 if (GammaEnergy < T0) {
130 static const G4double dT0 = keV;
131 X = (T0+dT0) / electron_mass_c2 ;
132 G4double sigma = p1Z*G4Log(1.+2*X)/X
133 + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
134 G4double c1 = -T0*(sigma-xSection)/(xSection*dT0);
135 G4double c2 = 0.150;
136 if (Z > 1.5) { c2 = 0.375-0.0556*G4Log(Z); }
137 G4double y = G4Log(GammaEnergy/T0);
138 xSection *= G4Exp(-y*(c1+c2*y));
139 }
140 // G4cout<<"e= "<< GammaEnergy<<" Z= "<<Z<<" cross= " << xSection << G4endl;
141 return xSection;
142}
143
144//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
145
147 std::vector<G4DynamicParticle*>* fvect,
149 const G4DynamicParticle* aDynamicGamma,
150 G4double,
151 G4double)
152{
153 // The scattered gamma energy is sampled according to Klein - Nishina formula.
154 // The random number techniques of Butcher & Messel are used
155 // (Nuc Phys 20(1960),15).
156 // Note : Effects due to binding of atomic electrons are negliged.
157
158 G4double gamEnergy0 = aDynamicGamma->GetKineticEnergy();
159
160 // do nothing below the threshold
161 if(gamEnergy0 <= LowEnergyLimit()) { return; }
162
163 G4double E0_m = gamEnergy0 / electron_mass_c2 ;
164
165 G4ThreeVector gamDirection0 = aDynamicGamma->GetMomentumDirection();
166
167 //
168 // sample the energy rate of the scattered gamma
169 //
170
171 G4double epsilon, epsilonsq, onecost, sint2, greject ;
172
173 G4double eps0 = 1./(1. + 2.*E0_m);
174 G4double epsilon0sq = eps0*eps0;
175 G4double alpha1 = - G4Log(eps0);
176 G4double alpha2 = alpha1 + 0.5*(1.- epsilon0sq);
177
178 CLHEP::HepRandomEngine* rndmEngineMod = G4Random::getTheEngine();
179 G4double rndm[3];
180
181 static const G4int nlooplim = 1000;
182 G4int nloop = 0;
183 do {
184 ++nloop;
185 // false interaction if too many iterations
186 if(nloop > nlooplim) { return; }
187
188 // 3 random numbers to sample scattering
189 rndmEngineMod->flatArray(3, rndm);
190
191 if ( alpha1 > alpha2*rndm[0] ) {
192 epsilon = G4Exp(-alpha1*rndm[1]); // eps0**r
193 epsilonsq = epsilon*epsilon;
194
195 } else {
196 epsilonsq = epsilon0sq + (1.- epsilon0sq)*rndm[1];
197 epsilon = sqrt(epsilonsq);
198 };
199
200 onecost = (1.- epsilon)/(epsilon*E0_m);
201 sint2 = onecost*(2.-onecost);
202 greject = 1. - epsilon*sint2/(1.+ epsilonsq);
203
204 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
205 } while (greject < rndm[2]);
206
207 //
208 // scattered gamma angles. ( Z - axis along the parent gamma)
209 //
210
211 if(sint2 < 0.0) { sint2 = 0.0; }
212 G4double cosTeta = 1. - onecost;
213 G4double sinTeta = sqrt (sint2);
214 G4double Phi = twopi * rndmEngineMod->flat();
215
216 //
217 // update G4VParticleChange for the scattered gamma
218 //
219
220 G4ThreeVector gamDirection1(sinTeta*cos(Phi), sinTeta*sin(Phi), cosTeta);
221 gamDirection1.rotateUz(gamDirection0);
222 G4double gamEnergy1 = epsilon*gamEnergy0;
223 G4double edep = 0.0;
224 if(gamEnergy1 > lowestSecondaryEnergy) {
227 } else {
230 edep = gamEnergy1;
231 }
232
233 //
234 // kinematic of the scattered electron
235 //
236
237 G4double eKinEnergy = gamEnergy0 - gamEnergy1;
238
239 if(eKinEnergy > lowestSecondaryEnergy) {
240 G4ThreeVector eDirection = gamEnergy0*gamDirection0 - gamEnergy1*gamDirection1;
241 eDirection = eDirection.unit();
242
243 // create G4DynamicParticle object for the electron.
244 auto dp = new G4DynamicParticle(theElectron,eDirection,eKinEnergy);
245 fvect->push_back(dp);
246 } else {
247 edep += eKinEnergy;
248 }
249 // energy balance
250 if(edep > 0.0) {
252 }
253}
254
255//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
256
257
G4double epsilon(G4double density, G4double temperature)
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
G4double G4Log(G4double x)
Definition: G4Log.hh:227
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double alpha2
Hep3Vector unit() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
Definition: G4Electron.cc:93
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax) override
void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
G4KleinNishinaCompton(const G4ParticleDefinition *p=nullptr, const G4String &nam="Klein-Nishina")
G4ParticleDefinition * theGamma
G4ParticleChangeForGamma * fParticleChange
G4ParticleDefinition * theElectron
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
~G4KleinNishinaCompton() override
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:831
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:124
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:823
G4bool IsMaster() const
Definition: G4VEmModel.hh:725
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:139
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)