Geant4  10.01.p02
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 90583 2015-06-04 11:16:41Z 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*CLHEP::barn, d2=-1.8300e-1*CLHEP::barn,
67  d3= 6.7527 *CLHEP::barn, d4=-1.9798e+1*CLHEP::barn,
68  e1= 1.9756e-5*CLHEP::barn, e2=-1.0205e-2*CLHEP::barn,
69  e3=-7.3913e-2*CLHEP::barn, e4= 2.7079e-2*CLHEP::barn,
70  f1=-3.9178e-7*CLHEP::barn, f2= 6.8241e-5*CLHEP::barn,
71  f3= 6.0480e-5*CLHEP::barn, f4= 3.0274e-4*CLHEP::barn;
72 static const G4int nlooplim = 1000;
73 
75  const G4String& nam)
76  : G4VEmModel(nam)
77 {
80  lowestSecondaryEnergy = 100.0*eV;
81  fParticleChange = 0;
82 }
83 
84 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
85 
87 {}
88 
89 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
90 
92  const G4DataVector& cuts)
93 {
94  if(IsMaster()) { InitialiseElementSelectors(p, cuts); }
96 }
97 
98 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
99 
101  G4VEmModel* masterModel)
102 {
104 }
105 
106 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
107 
109  const G4ParticleDefinition*,
110  G4double GammaEnergy,
111  G4double Z, G4double,
113 {
114  G4double xSection = 0.0 ;
115  if (GammaEnergy <= LowEnergyLimit()) { return xSection; }
116 
117  static const G4double a = 20.0 , b = 230.0 , c = 440.0;
118 
119  G4double p1Z = Z*(d1 + e1*Z + f1*Z*Z), p2Z = Z*(d2 + e2*Z + f2*Z*Z),
120  p3Z = Z*(d3 + e3*Z + f3*Z*Z), p4Z = Z*(d4 + e4*Z + f4*Z*Z);
121 
122  G4double T0 = 15.0*keV;
123  if (Z < 1.5) { T0 = 40.0*keV; }
124 
125  G4double X = max(GammaEnergy, T0) / electron_mass_c2;
126  xSection = p1Z*G4Log(1.+2.*X)/X
127  + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
128 
129  // modification for low energy. (special case for Hydrogen)
130  if (GammaEnergy < T0) {
131  static const G4double dT0 = keV;
132  X = (T0+dT0) / electron_mass_c2 ;
133  G4double sigma = p1Z*G4Log(1.+2*X)/X
134  + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
135  G4double c1 = -T0*(sigma-xSection)/(xSection*dT0);
136  G4double c2 = 0.150;
137  if (Z > 1.5) { c2 = 0.375-0.0556*G4Log(Z); }
138  G4double y = G4Log(GammaEnergy/T0);
139  xSection *= G4Exp(-y*(c1+c2*y));
140  }
141  // G4cout<<"e= "<< GammaEnergy<<" Z= "<<Z<<" cross= " << xSection << G4endl;
142  return xSection;
143 }
144 
145 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
146 
148  std::vector<G4DynamicParticle*>* fvect,
149  const G4MaterialCutsCouple*,
150  const G4DynamicParticle* aDynamicGamma,
151  G4double,
152  G4double)
153 {
154  // The scattered gamma energy is sampled according to Klein - Nishina formula.
155  // The random number techniques of Butcher & Messel are used
156  // (Nuc Phys 20(1960),15).
157  // Note : Effects due to binding of atomic electrons are negliged.
158 
159  G4double gamEnergy0 = aDynamicGamma->GetKineticEnergy();
160 
161  // do nothing below the threshold
162  if(gamEnergy0 <= LowEnergyLimit()) { return; }
163 
164  G4double E0_m = gamEnergy0 / electron_mass_c2 ;
165 
166  G4ThreeVector gamDirection0 = aDynamicGamma->GetMomentumDirection();
167 
168  //
169  // sample the energy rate of the scattered gamma
170  //
171 
172  G4double epsilon, epsilonsq, onecost, sint2, greject ;
173 
174  G4double eps0 = 1./(1. + 2.*E0_m);
175  G4double epsilon0sq = eps0*eps0;
176  G4double alpha1 = - G4Log(eps0);
177  G4double alpha2 = 0.5*(1.- epsilon0sq);
178 
179  rndmEngineMod = G4Random::getTheEngine();
180 
181  G4int nloop = 0;
182  do {
183  ++nloop;
184  // false interaction if too many iterations
185  if(nloop > nlooplim) { return; }
186 
187  if ( alpha1/(alpha1+alpha2) > rndmEngineMod->flat() ) {
188  epsilon = G4Exp(-alpha1*rndmEngineMod->flat()); // eps0**r
189  epsilonsq = epsilon*epsilon;
190 
191  } else {
192  epsilonsq = epsilon0sq + (1.- epsilon0sq)*rndmEngineMod->flat();
193  epsilon = sqrt(epsilonsq);
194  };
195 
196  onecost = (1.- epsilon)/(epsilon*E0_m);
197  sint2 = onecost*(2.-onecost);
198  greject = 1. - epsilon*sint2/(1.+ epsilonsq);
199 
200  } while (greject < rndmEngineMod->flat());
201 
202  //
203  // scattered gamma angles. ( Z - axis along the parent gamma)
204  //
205 
206  if(sint2 < 0.0) { sint2 = 0.0; }
207  G4double cosTeta = 1. - onecost;
208  G4double sinTeta = sqrt (sint2);
209  G4double Phi = twopi * rndmEngineMod->flat();
210 
211  //
212  // update G4VParticleChange for the scattered gamma
213  //
214 
215  G4ThreeVector gamDirection1(sinTeta*cos(Phi), sinTeta*sin(Phi), cosTeta);
216  gamDirection1.rotateUz(gamDirection0);
217  G4double gamEnergy1 = epsilon*gamEnergy0;
218  G4double edep = 0.0;
219  if(gamEnergy1 > lowestSecondaryEnergy) {
222  } else {
225  edep = gamEnergy1;
226  }
227 
228  //
229  // kinematic of the scattered electron
230  //
231 
232  G4double eKinEnergy = gamEnergy0 - gamEnergy1;
233 
234  if(eKinEnergy > lowestSecondaryEnergy) {
235  G4ThreeVector eDirection = gamEnergy0*gamDirection0 - gamEnergy1*gamDirection1;
236  eDirection = eDirection.unit();
237 
238  // create G4DynamicParticle object for the electron.
239  G4DynamicParticle* dp = new G4DynamicParticle(theElectron,eDirection,eKinEnergy);
240  fvect->push_back(dp);
241  } else {
242  edep += eKinEnergy;
243  }
244  // energy balance
245  if(edep > 0.0) {
247  }
248 }
249 
250 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
251 
252 
static const G4double d3
static c2_factory< G4double > c2
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:623
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:148
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
static const G4double e2
G4bool IsMaster() const
Definition: G4VEmModel.hh:700
G4double a
Definition: TRTMaterials.hh:39
static const G4int nlooplim
static const G4double e4
int G4int
Definition: G4Types.hh:78
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)
CLHEP::HepRandomEngine * rndmEngineMod
Definition: G4VEmModel.hh:400
const G4ThreeVector & GetMomentumDirection() const
static const G4double f4
double flat()
Definition: G4AblaRandom.cc:47
static const G4double c1
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:784
static const G4double f1
static const G4double e1
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
static const double eV
Definition: G4SIunits.hh:194
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:792
static const G4double f3
T max(const T t1, const T t2)
brief Return the largest of the two arguments
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 dT0
static const G4double d4
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:134
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)