source: trunk/source/processes/electromagnetic/lowenergy/src/G4CrossSectionElasticScreenedRutherfordLE.cc@ 1005

Last change on this file since 1005 was 991, checked in by garnier, 17 years ago

update

File size: 4.4 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: G4CrossSectionElasticScreenedRutherfordLE.cc,v 1.2 2008/07/14 20:47:34 sincerti Exp $
27// GEANT4 tag $Name: geant4-09-02 $
28
29#include "G4CrossSectionElasticScreenedRutherfordLE.hh"
30
31//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
32
33G4CrossSectionElasticScreenedRutherfordLE::G4CrossSectionElasticScreenedRutherfordLE()
34{
35 lowEnergyLimit = 0 * eV;
36 highEnergyLimit = 200 * eV;
37}
38
39//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
40
41G4CrossSectionElasticScreenedRutherfordLE::~G4CrossSectionElasticScreenedRutherfordLE()
42{}
43
44//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
45
46G4double G4CrossSectionElasticScreenedRutherfordLE::CrossSection(const G4Track& track)
47{
48 const G4DynamicParticle* particle = track.GetDynamicParticle();
49 G4double k = particle->GetKineticEnergy();
50
51 G4double screenedCrossSection = 0.;
52
53 if (k > lowEnergyLimit && k <= highEnergyLimit)
54 {
55 G4double z = 10.;
56 G4double n = ScreeningFactor(k,z);
57 G4double crossSection = RutherfordCrossSection(k, z);
58 screenedCrossSection = pi * crossSection / (n * (n + 1.));
59 }
60
61 return screenedCrossSection;
62}
63
64//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
65
66G4double G4CrossSectionElasticScreenedRutherfordLE::RutherfordCrossSection(G4double k, G4double z)
67{
68 //
69 // e^4 / K + m_e c^2 \^2
70 // sigma_Ruth(K) = Z (Z+1) -------------------- | --------------------- |
71 // (4 pi epsilon_0)^2 \ K * (K + 2 m_e c^2) /
72 //
73 // Where K is the electron non-relativistic kinetic energy
74 //
75 // NIM 155, pp. 145-156, 1978
76
77 G4double length =(e_squared * (k + electron_mass_c2)) / (4 * pi *epsilon0 * k * ( k + 2 * electron_mass_c2));
78 G4double cross = z * ( z + 1) * length * length;
79
80 return cross;
81}
82
83//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
84
85G4double G4CrossSectionElasticScreenedRutherfordLE::ScreeningFactor(G4double k, G4double z)
86{
87 //
88 // alpha_1 + beta_1 ln(K/eV) constK Z^(2/3)
89 // n(T) = -------------------------- -----------------
90 // K/(m_e c^2) 2 + K/(m_e c^2)
91 //
92 // Where K is the electron non-relativistic kinetic energy
93 //
94 // n(T) > 0 for T < ~ 400 MeV
95 //
96 // NIM 155, pp. 145-156, 1978
97 // Formulae (2) and (5)
98
99 const G4double alpha_1(1.64);
100 const G4double beta_1(-0.0825);
101 const G4double constK(1.7E-5);
102
103 G4double numerator = (alpha_1 + beta_1 * std::log(k/eV)) * constK * std::pow(z, 2./3.);
104
105 k /= electron_mass_c2;
106
107 G4double denominator = k * (2 + k);
108
109 G4double value = 0.;
110 if (denominator > 0.) value = numerator / denominator;
111
112 return value;
113}
Note: See TracBrowser for help on using the repository browser.