source: trunk/source/processes/electromagnetic/standard/src/G4KleinNishinaModel.cc @ 1350

Last change on this file since 1350 was 1350, checked in by garnier, 13 years ago

update to last version 4.9.4

File size: 9.7 KB
Line 
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: G4KleinNishinaModel.cc,v 1.1 2010/09/03 14:11:16 vnivanch Exp $
27// GEANT4 tag $Name: emstand-V09-03-24 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class file
32//
33//
34// File name:     G4KleinNishinaModel
35//
36// Author:        Vladimir Ivanchenko on base of G4KleinNishinaCompton
37//
38// Creation date: 13.06.2010
39//
40// Modifications:
41//
42// Class Description:
43//
44// -------------------------------------------------------------------
45//
46//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
47//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
48
49#include "G4KleinNishinaModel.hh"
50#include "G4Electron.hh"
51#include "G4Gamma.hh"
52#include "Randomize.hh"
53#include "G4RandomDirection.hh"
54#include "G4DataVector.hh"
55#include "G4ParticleChangeForGamma.hh"
56#include "G4VAtomDeexcitation.hh"
57#include "G4LossTableManager.hh"
58
59//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
60
61using namespace std;
62
63G4KleinNishinaModel::G4KleinNishinaModel(const G4String& nam)
64  : G4VEmModel(nam),isInitialized(false)
65{
66  theGamma = G4Gamma::Gamma();
67  theElectron = G4Electron::Electron();
68  lowestGammaEnergy = 1.0*eV;
69  fProbabilities.resize(9,0.0);
70}
71
72//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
73
74G4KleinNishinaModel::~G4KleinNishinaModel()
75{}
76
77//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
78
79void G4KleinNishinaModel::Initialise(const G4ParticleDefinition* p,
80                                     const G4DataVector& cuts)
81{
82  fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
83  InitialiseElementSelectors(p, cuts);
84
85  if (isInitialized) { return; }
86  fParticleChange = GetParticleChangeForGamma();
87  isInitialized = true;
88}
89
90//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
91
92G4double
93G4KleinNishinaModel::ComputeCrossSectionPerAtom(const G4ParticleDefinition*,
94                                                G4double GammaEnergy,
95                                                G4double Z, G4double,
96                                                G4double, G4double)
97{
98  G4double CrossSection = 0.0 ;
99  if ( Z < 0.9999 || GammaEnergy < 0.1*keV) { return CrossSection; }
100
101  static const G4double a = 20.0 , b = 230.0 , c = 440.0;
102 
103  static const G4double
104    d1= 2.7965e-1*barn, d2=-1.8300e-1*barn, d3= 6.7527   *barn, d4=-1.9798e+1*barn,
105    e1= 1.9756e-5*barn, e2=-1.0205e-2*barn, e3=-7.3913e-2*barn, e4= 2.7079e-2*barn,
106    f1=-3.9178e-7*barn, f2= 6.8241e-5*barn, f3= 6.0480e-5*barn, f4= 3.0274e-4*barn;
107     
108  G4double p1Z = Z*(d1 + e1*Z + f1*Z*Z), p2Z = Z*(d2 + e2*Z + f2*Z*Z),
109           p3Z = Z*(d3 + e3*Z + f3*Z*Z), p4Z = Z*(d4 + e4*Z + f4*Z*Z);
110
111  G4double T0  = 15.0*keV; 
112  if (Z < 1.5) { T0 = 40.0*keV; } 
113
114  G4double X   = max(GammaEnergy, T0) / electron_mass_c2;
115  CrossSection = p1Z*std::log(1.+2.*X)/X
116               + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
117               
118  //  modification for low energy. (special case for Hydrogen)
119  if (GammaEnergy < T0) {
120    G4double dT0 = keV;
121    X = (T0+dT0) / electron_mass_c2 ;
122    G4double sigma = p1Z*log(1.+2*X)/X
123                    + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
124    G4double   c1 = -T0*(sigma-CrossSection)/(CrossSection*dT0);             
125    G4double   c2 = 0.150; 
126    if (Z > 1.5) { c2 = 0.375-0.0556*log(Z); }
127    G4double    y = log(GammaEnergy/T0);
128    CrossSection *= exp(-y*(c1+c2*y));         
129  }
130  //  G4cout << "e= " << GammaEnergy << " Z= " << Z
131  //  << " cross= " << CrossSection << G4endl;
132  return CrossSection;
133}
134
135//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
136
137void G4KleinNishinaModel::SampleSecondaries(
138                             std::vector<G4DynamicParticle*>* fvect,
139                             const G4MaterialCutsCouple* couple,
140                             const G4DynamicParticle* aDynamicGamma,
141                             G4double,
142                             G4double)
143{
144  G4double energy = aDynamicGamma->GetKineticEnergy();
145  G4ThreeVector direction = aDynamicGamma->GetMomentumDirection();
146
147  // select atom
148  const G4Element* elm = SelectRandomAtom(couple, theGamma, energy);
149
150  // select shell first
151  G4int Z = (G4int)elm->GetZ();
152  G4int nShells = elm->GetNbOfAtomicShells();
153  if(nShells > (G4int)fProbabilities.size()) { fProbabilities.resize(nShells); }
154  G4double totprob = 0.0;
155  G4int i = 0;
156  for(; i<nShells; ++i) {
157    G4double prob = 0.0;
158    if(energy > elm->GetAtomicShell(i)) { 
159      prob = (G4double)elm->GetNbOfShellElectrons(i);
160    }
161    totprob += prob;
162    fProbabilities[i] = totprob; 
163  }
164  if(totprob == 0.0) { return; }
165
166  G4LorentzVector lv1, lv2, lv3;
167  G4LorentzVector lv0(energy*direction.x(),energy*direction.y(),
168                      energy*direction.z(),energy);
169  G4double eKinEnergy = 0.0;
170  G4double gamEnergy1 = 0.0;
171
172  // Loop on sampling
173  do {
174    G4double xprob = totprob*G4UniformRand();
175
176    for(i=0; i<nShells; ++i) { if(xprob <= fProbabilities[i]) {break;} }
177    if( i == nShells ) { return; }
178   
179    G4double bindingEnergy  = elm->GetAtomicShell(i);
180    G4double tkin =  bindingEnergy*0.5;
181    G4double eEnergy = tkin + electron_mass_c2;
182    G4double eTotMomentum = sqrt(tkin*(tkin + electron_mass_c2*2));
183    G4ThreeVector eDir = G4RandomDirection();
184    lv1 = lv0;
185    lv2.set(eTotMomentum*eDir.x(),eTotMomentum*eDir.y(),
186            eTotMomentum*eDir.z(),eEnergy);
187    G4ThreeVector bst = lv2.boostVector();
188    lv1.boost(-bst);
189
190    // In the rest frame of an electron
191    // The scattered gamma energy is sampled according to Klein - Nishina formula.
192    // The random number techniques of Butcher & Messel are used
193    // (Nuc Phys 20(1960),15).
194 
195    G4double gamEnergy0 = lv1.e();
196    G4double E0_m = gamEnergy0 / electron_mass_c2 ;
197
198    G4ThreeVector gamDirection0 = (lv1.vect()).unit();
199
200    //
201    // sample the energy rate of the scattered gamma
202    //
203
204    G4double epsilon, epsilonsq, onecost, sint2, greject ;
205
206    G4double epsilon0   = 1./(1. + 2.*E0_m);
207    G4double epsilon0sq = epsilon0*epsilon0;
208    G4double alpha1     = - log(epsilon0);
209    G4double alpha2     = 0.5*(1.- epsilon0sq);
210
211    do {
212      if ( alpha1/(alpha1+alpha2) > G4UniformRand() ) {
213        epsilon   = exp(-alpha1*G4UniformRand());   // epsilon0**r
214        epsilonsq = epsilon*epsilon; 
215
216      } else {
217        epsilonsq = epsilon0sq + (1.- epsilon0sq)*G4UniformRand();
218        epsilon   = sqrt(epsilonsq);
219      };
220
221      onecost = (1.- epsilon)/(epsilon*E0_m);
222      sint2   = onecost*(2.-onecost);
223      greject = 1. - epsilon*sint2/(1.+ epsilonsq);
224
225    } while (greject < G4UniformRand());
226 
227    //
228    // scattered gamma angles. ( Z - axis along the parent gamma)
229    //
230
231    G4double cosTeta = 1. - onecost; 
232    G4double sinTeta = sqrt (sint2);
233    G4double Phi     = twopi * G4UniformRand();
234    G4double dirx = sinTeta*cos(Phi), diry = sinTeta*sin(Phi), dirz = cosTeta;
235
236    //
237    // update G4VParticleChange for the scattered gamma
238    //
239   
240    G4ThreeVector gamDirection1 ( dirx,diry,dirz );
241    gamDirection1.rotateUz(gamDirection0);
242    gamEnergy1 = epsilon*gamEnergy0;
243
244    // before scattering
245    lv2.set(0.0,0.0,0.0,electron_mass_c2);
246    lv2 += lv1;
247 
248    // after scattering
249    lv1.set(gamEnergy1*gamDirection1.x(),gamEnergy1*gamDirection1.y(),
250            gamEnergy1*gamDirection1.z(),gamEnergy1);
251    lv2 -= lv1;
252    lv2.boost(bst);
253    lv1.boost(bst);
254    eKinEnergy = lv2.e() - electron_mass_c2 - bindingEnergy;
255  } while ( eKinEnergy < 0.0 );
256
257  // gamma kinematics
258  gamEnergy1 = lv1.e();
259  if(gamEnergy1 > lowestGammaEnergy) {
260    fParticleChange->SetProposedKineticEnergy(gamEnergy1);
261    fParticleChange->ProposeMomentumDirection((lv1.vect()).unit());
262  } else { 
263    fParticleChange->ProposeTrackStatus(fStopAndKill);
264    fParticleChange->ProposeLocalEnergyDeposit(gamEnergy1);
265  }
266
267  //
268  // kinematic of the scattered electron
269  //
270  if(eKinEnergy > DBL_MIN) {
271    G4ThreeVector eDirection = (lv2.vect()).unit();
272    G4DynamicParticle* dp = new G4DynamicParticle(theElectron,eDirection,eKinEnergy);
273    fvect->push_back(dp);
274  }
275  // sample deexcitation
276  //
277  if(fAtomDeexcitation) {
278    G4AtomicShellEnumerator as = G4AtomicShellEnumerator(i);
279    const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as);   
280    fAtomDeexcitation->GenerateParticles(fvect, shell, Z, couple->GetIndex());
281  }
282}
283
284//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
285
Note: See TracBrowser for help on using the repository browser.