Geant4  10.00.p01
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 // $Id: G4KleinNishinaCompton.cc 74309 2013-10-03 06:42:30Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4KleinNishinaCompton
34 //
35 // Author: Vladimir Ivanchenko on base of Michel Maire code
36 //
37 // Creation date: 15.03.2005
38 //
39 // Modifications:
40 // 18-04-05 Use G4ParticleChangeForGamma (V.Ivantchenko)
41 // 27-03-06 Remove upper limit of cross section (V.Ivantchenko)
42 //
43 // Class Description:
44 //
45 // -------------------------------------------------------------------
46 //
47 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
48 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
49 
50 #include "G4KleinNishinaCompton.hh"
51 #include "G4PhysicalConstants.hh"
52 #include "G4SystemOfUnits.hh"
53 #include "G4Electron.hh"
54 #include "G4Gamma.hh"
55 #include "Randomize.hh"
56 #include "G4DataVector.hh"
58 #include "G4Log.hh"
59 #include "G4Exp.hh"
60 
61 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
62 
63 using namespace std;
64 
65 static const G4double
66  d1= 2.7965e-1*barn, d2=-1.8300e-1*barn, d3= 6.7527 *barn, d4=-1.9798e+1*barn,
67  e1= 1.9756e-5*barn, e2=-1.0205e-2*barn, e3=-7.3913e-2*barn, e4= 2.7079e-2*barn,
68  f1=-3.9178e-7*barn, f2= 6.8241e-5*barn, f3= 6.0480e-5*barn, f4= 3.0274e-4*barn;
69 
71  const G4String& nam)
72  : G4VEmModel(nam)
73 {
76  lowestGammaEnergy = 1.0*eV;
77  fParticleChange = 0;
78 }
79 
80 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
81 
83 {}
84 
85 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
86 
88  const G4DataVector& cuts)
89 {
90  if(IsMaster()) { InitialiseElementSelectors(p, cuts); }
92 }
93 
94 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
95 
97  G4VEmModel* masterModel)
98 {
100 }
101 
102 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
103 
105  const G4ParticleDefinition*,
106  G4double GammaEnergy,
107  G4double Z, G4double,
109 {
110  G4double xSection = 0.0 ;
111  if ( Z < 0.9999 ) return xSection;
112  if ( GammaEnergy < 0.1*keV ) return xSection;
113  // if ( GammaEnergy > (100.*GeV/Z) ) return xSection;
114 
115  static const G4double a = 20.0 , b = 230.0 , c = 440.0;
116 
117  G4double p1Z = Z*(d1 + e1*Z + f1*Z*Z), p2Z = Z*(d2 + e2*Z + f2*Z*Z),
118  p3Z = Z*(d3 + e3*Z + f3*Z*Z), p4Z = Z*(d4 + e4*Z + f4*Z*Z);
119 
120  G4double T0 = 15.0*keV;
121  if (Z < 1.5) T0 = 40.0*keV;
122 
123  G4double X = max(GammaEnergy, T0) / electron_mass_c2;
124  xSection = p1Z*G4Log(1.+2.*X)/X
125  + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
126 
127  // modification for low energy. (special case for Hydrogen)
128  if (GammaEnergy < T0) {
129  G4double dT0 = 1.*keV;
130  X = (T0+dT0) / electron_mass_c2 ;
131  G4double sigma = p1Z*G4Log(1.+2*X)/X
132  + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
133  G4double c1 = -T0*(sigma-xSection)/(xSection*dT0);
134  G4double c2 = 0.150;
135  if (Z > 1.5) c2 = 0.375-0.0556*G4Log(Z);
136  G4double y = G4Log(GammaEnergy/T0);
137  xSection *= G4Exp(-y*(c1+c2*y));
138  }
139  // G4cout << "e= " << GammaEnergy << " Z= " << Z << " cross= " << xSection << G4endl;
140  return xSection;
141 }
142 
143 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
144 
145 void G4KleinNishinaCompton::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
146  const G4MaterialCutsCouple*,
147  const G4DynamicParticle* aDynamicGamma,
148  G4double,
149  G4double)
150 {
151  // The scattered gamma energy is sampled according to Klein - Nishina formula.
152  // The random number techniques of Butcher & Messel are used
153  // (Nuc Phys 20(1960),15).
154  // Note : Effects due to binding of atomic electrons are negliged.
155 
156  G4double gamEnergy0 = aDynamicGamma->GetKineticEnergy();
157 
158  // extra protection
159  if(gamEnergy0 < lowestGammaEnergy) {
163  return;
164  }
165 
166  G4double E0_m = gamEnergy0 / electron_mass_c2 ;
167 
168  G4ThreeVector gamDirection0 = aDynamicGamma->GetMomentumDirection();
169 
170  //
171  // sample the energy rate of the scattered gamma
172  //
173 
174  G4double epsilon, epsilonsq, onecost, sint2, greject ;
175 
176  G4double eps0 = 1./(1. + 2.*E0_m);
177  G4double epsilon0sq = eps0*eps0;
178  G4double alpha1 = - G4Log(eps0);
179  G4double alpha2 = 0.5*(1.- epsilon0sq);
180 
181  do {
182  if ( alpha1/(alpha1+alpha2) > G4UniformRand() ) {
183  epsilon = G4Exp(-alpha1*G4UniformRand()); // eps0**r
184  epsilonsq = epsilon*epsilon;
185 
186  } else {
187  epsilonsq = epsilon0sq + (1.- epsilon0sq)*G4UniformRand();
188  epsilon = sqrt(epsilonsq);
189  };
190 
191  onecost = (1.- epsilon)/(epsilon*E0_m);
192  sint2 = onecost*(2.-onecost);
193  greject = 1. - epsilon*sint2/(1.+ epsilonsq);
194 
195  } while (greject < G4UniformRand());
196 
197  //
198  // scattered gamma angles. ( Z - axis along the parent gamma)
199  //
200 
201  if(sint2 < 0.0) { sint2 = 0.0; }
202  G4double cosTeta = 1. - onecost;
203  G4double sinTeta = sqrt (sint2);
204  G4double Phi = twopi * G4UniformRand();
205 
206  //
207  // update G4VParticleChange for the scattered gamma
208  //
209 
210  G4ThreeVector gamDirection1(sinTeta*cos(Phi), sinTeta*sin(Phi), cosTeta);
211  gamDirection1.rotateUz(gamDirection0);
212  G4double gamEnergy1 = epsilon*gamEnergy0;
213  if(gamEnergy1 > lowestGammaEnergy) {
216  } else {
220  }
221 
222  //
223  // kinematic of the scattered electron
224  //
225 
226  G4double eKinEnergy = gamEnergy0 - gamEnergy1;
227 
228  if(eKinEnergy > DBL_MIN) {
229  G4ThreeVector eDirection = gamEnergy0*gamDirection0 - gamEnergy1*gamDirection1;
230  eDirection = eDirection.unit();
231 
232  // create G4DynamicParticle object for the electron.
233  G4DynamicParticle* dp = new G4DynamicParticle(theElectron,eDirection,eKinEnergy);
234  fvect->push_back(dp);
235  }
236 }
237 
238 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
239 
240 
static const G4double d3
static c2_factory< G4double > c2
G4ParticleChangeForGamma * fParticleChange
static const G4double d1
G4KleinNishinaCompton(const G4ParticleDefinition *p=0, const G4String &nam="Klein-Nishina")
static const G4double f2
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:135
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
static const G4double e2
G4bool IsMaster() const
Definition: G4VEmModel.hh:676
G4double a
Definition: TRTMaterials.hh:39
static const G4double e4
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax)
#define G4UniformRand()
Definition: Randomize.hh:87
const G4ThreeVector & GetMomentumDirection() const
static const G4double f4
static const G4double c1
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:760
static const G4double f1
static const G4double e1
G4double G4Log(G4double x)
Definition: G4Log.hh:227
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
static const double eV
Definition: G4SIunits.hh:194
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:768
static const G4double f3
T max(const T t1, const T t2)
brief Return the largest of the two arguments
#define DBL_MIN
Definition: templates.hh:75
static G4Electron * Electron()
Definition: G4Electron.cc:94
void SetProposedKineticEnergy(G4double proposedKinEnergy)
G4ParticleDefinition * theElectron
static const double keV
Definition: G4SIunits.hh:195
static const double barn
Definition: G4SIunits.hh:95
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
static const G4double e3
static const G4double d2
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
G4ParticleDefinition * theGamma
static const G4double d4
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:121
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)