Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups 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 // $Id$
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 
59 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
60 
61 using namespace std;
62 
64  const G4String& nam)
65  : G4VEmModel(nam)
66 {
69  lowestGammaEnergy = 1.0*eV;
70  fParticleChange = 0;
71 }
72 
73 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
74 
76 {}
77 
78 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
79 
81  const G4DataVector&)
82 {
84 }
85 
86 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
87 
89  const G4ParticleDefinition*,
90  G4double GammaEnergy,
93 {
94  G4double xSection = 0.0 ;
95  if ( Z < 0.9999 ) return xSection;
96  if ( GammaEnergy < 0.1*keV ) return xSection;
97  // if ( GammaEnergy > (100.*GeV/Z) ) return xSection;
98 
99  static const G4double a = 20.0 , b = 230.0 , c = 440.0;
100 
101  static const G4double
102  d1= 2.7965e-1*barn, d2=-1.8300e-1*barn, d3= 6.7527 *barn, d4=-1.9798e+1*barn,
103  e1= 1.9756e-5*barn, e2=-1.0205e-2*barn, e3=-7.3913e-2*barn, e4= 2.7079e-2*barn,
104  f1=-3.9178e-7*barn, f2= 6.8241e-5*barn, f3= 6.0480e-5*barn, f4= 3.0274e-4*barn;
105 
106  G4double p1Z = Z*(d1 + e1*Z + f1*Z*Z), p2Z = Z*(d2 + e2*Z + f2*Z*Z),
107  p3Z = Z*(d3 + e3*Z + f3*Z*Z), p4Z = Z*(d4 + e4*Z + f4*Z*Z);
108 
109  G4double T0 = 15.0*keV;
110  if (Z < 1.5) T0 = 40.0*keV;
111 
112  G4double X = max(GammaEnergy, T0) / electron_mass_c2;
113  xSection = p1Z*std::log(1.+2.*X)/X
114  + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
115 
116  // modification for low energy. (special case for Hydrogen)
117  if (GammaEnergy < T0) {
118  G4double dT0 = 1.*keV;
119  X = (T0+dT0) / electron_mass_c2 ;
120  G4double sigma = p1Z*log(1.+2*X)/X
121  + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
122  G4double c1 = -T0*(sigma-xSection)/(xSection*dT0);
123  G4double c2 = 0.150;
124  if (Z > 1.5) c2 = 0.375-0.0556*log(Z);
125  G4double y = log(GammaEnergy/T0);
126  xSection *= exp(-y*(c1+c2*y));
127  }
128  // G4cout << "e= " << GammaEnergy << " Z= " << Z << " cross= " << xSection << G4endl;
129  return xSection;
130 }
131 
132 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
133 
134 void G4KleinNishinaCompton::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
135  const G4MaterialCutsCouple*,
136  const G4DynamicParticle* aDynamicGamma,
137  G4double,
138  G4double)
139 {
140  // The scattered gamma energy is sampled according to Klein - Nishina formula.
141  // The random number techniques of Butcher & Messel are used
142  // (Nuc Phys 20(1960),15).
143  // Note : Effects due to binding of atomic electrons are negliged.
144 
145  G4double gamEnergy0 = aDynamicGamma->GetKineticEnergy();
146 
147  // extra protection
148  if(gamEnergy0 < lowestGammaEnergy) {
152  return;
153  }
154 
155  G4double E0_m = gamEnergy0 / electron_mass_c2 ;
156 
157  G4ThreeVector gamDirection0 = aDynamicGamma->GetMomentumDirection();
158 
159  //
160  // sample the energy rate of the scattered gamma
161  //
162 
163  G4double epsilon, epsilonsq, onecost, sint2, greject ;
164 
165  G4double eps0 = 1./(1. + 2.*E0_m);
166  G4double epsilon0sq = eps0*eps0;
167  G4double alpha1 = - log(eps0);
168  G4double alpha2 = 0.5*(1.- epsilon0sq);
169 
170  do {
171  if ( alpha1/(alpha1+alpha2) > G4UniformRand() ) {
172  epsilon = exp(-alpha1*G4UniformRand()); // eps0**r
173  epsilonsq = epsilon*epsilon;
174 
175  } else {
176  epsilonsq = epsilon0sq + (1.- epsilon0sq)*G4UniformRand();
177  epsilon = sqrt(epsilonsq);
178  };
179 
180  onecost = (1.- epsilon)/(epsilon*E0_m);
181  sint2 = onecost*(2.-onecost);
182  greject = 1. - epsilon*sint2/(1.+ epsilonsq);
183 
184  } while (greject < G4UniformRand());
185 
186  //
187  // scattered gamma angles. ( Z - axis along the parent gamma)
188  //
189 
190  if(sint2 < 0.0) { sint2 = 0.0; }
191  G4double cosTeta = 1. - onecost;
192  G4double sinTeta = sqrt (sint2);
193  G4double Phi = twopi * G4UniformRand();
194 
195  //
196  // update G4VParticleChange for the scattered gamma
197  //
198 
199  G4ThreeVector gamDirection1(sinTeta*cos(Phi), sinTeta*sin(Phi), cosTeta);
200  gamDirection1.rotateUz(gamDirection0);
201  G4double gamEnergy1 = epsilon*gamEnergy0;
202  if(gamEnergy1 > lowestGammaEnergy) {
205  } else {
209  }
210 
211  //
212  // kinematic of the scattered electron
213  //
214 
215  G4double eKinEnergy = gamEnergy0 - gamEnergy1;
216 
217  if(eKinEnergy > DBL_MIN) {
218  G4ThreeVector eDirection = gamEnergy0*gamDirection0 - gamEnergy1*gamDirection1;
219  eDirection = eDirection.unit();
220 
221  // create G4DynamicParticle object for the electron.
222  G4DynamicParticle* dp = new G4DynamicParticle(theElectron,eDirection,eKinEnergy);
223  fvect->push_back(dp);
224  }
225 }
226 
227 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
228 
229