source: trunk/source/processes/electromagnetic/lowenergy/src/G4CrossSectionElasticScreenedRutherford.cc@ 968

Last change on this file since 968 was 961, checked in by garnier, 17 years ago

update processes

File size: 5.8 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//
27// $Id: G4CrossSectionElasticScreenedRutherford.cc,v 1.1 2007/10/12 23:11:41 pia Exp $
28// GEANT4 tag $Name: geant4-09-02-ref-02 $
29//
30// Contact Author: Maria Grazia Pia (Maria.Grazia.Pia@cern.ch)
31//
32// Reference: TNS Geant4-DNA paper
33// Reference for implementation model: NIM. 155, pp. 145-156, 1978
34
35// History:
36// -----------
37// Date Name Modification
38// 28 Apr 2007 M.G. Pia Created in compliance with design described in TNS paper
39//
40// -------------------------------------------------------------------
41
42// Class description:
43// Geant4-DNA Cross total cross section for electron elastic scattering in water
44// Reference: TNS Geant4-DNA paper
45// S. Chauvie et al., Geant4 physics processes for microdosimetry simulation:
46// design foundation and implementation of the first set of models,
47// IEEE Trans. Nucl. Sci., vol. 54, no. 6, Dec. 2007.
48// Further documentation available from http://www.ge.infn.it/geant4/dna
49
50// -------------------------------------------------------------------
51
52
53#include "G4CrossSectionElasticScreenedRutherford.hh"
54#include "G4Track.hh"
55#include "G4DynamicParticle.hh"
56#include "G4ParticleDefinition.hh"
57#include "G4Electron.hh"
58
59
60G4CrossSectionElasticScreenedRutherford::G4CrossSectionElasticScreenedRutherford()
61{
62
63 name = "CrossSectionElasticScreenedRutherford";
64 lowEnergyLimit = 7. * eV;
65 highEnergyLimit = 10 * MeV;
66
67// if (verboseLevel > 0)
68// {
69// G4cout << name << " is created " << G4endl
70// << "Energy range: "
71// << lowEnergyLimit / keV << " keV - "
72// << highEnergyLimit / GeV << " GeV"
73// << G4endl;
74// }
75}
76
77
78G4CrossSectionElasticScreenedRutherford::~G4CrossSectionElasticScreenedRutherford()
79{ }
80
81
82G4double G4CrossSectionElasticScreenedRutherford::CrossSection(const G4Track& track)
83{
84 const G4DynamicParticle* particle = track.GetDynamicParticle();
85 G4double k = particle->GetKineticEnergy();
86
87 // Cross section = 0 outside the energy validity limits set in the constructor
88 // ---- MGP ---- Better handling of these limits to be set in a following design iteration
89
90 G4double screenedCrossSection = 0.;
91
92 if (k > lowEnergyLimit && k < highEnergyLimit)
93 {
94 // G4Material* material = track.GetMaterial();
95
96 // Assume that the material is water; proper algorithm to calculate z correctly for any material to be inserted here
97 // For H20 Z = 10 (total number of electrons)
98 G4double z = 10.;
99
100 G4double n = ScreeningFactor(k,z);
101 G4double crossSection = RutherfordCrossSection(k, z);
102 screenedCrossSection = pi * crossSection / (n * (n + 1.));
103 }
104
105 return screenedCrossSection;
106}
107
108G4double G4CrossSectionElasticScreenedRutherford::RutherfordCrossSection(G4double k, G4double z)
109{
110 //
111 // e^4 / K + m_e c^2 \^2
112 // sigma_Ruth(K) = Z (Z+1) -------------------- | --------------------- |
113 // (4 pi epsilon_0)^2 \ K * (K + 2 m_e c^2) /
114 //
115 // Where K is the electron non-relativistic kinetic energy
116 //
117 // NIM 155, pp. 145-156, 1978
118
119 G4double length =(e_squared * (k + electron_mass_c2)) / (4 * pi *epsilon0 * k * ( k + 2 * electron_mass_c2));
120 G4double cross = z * ( z + 1) * length * length;
121
122 return cross;
123}
124
125//G4bool G4CrossSectionElasticScreenedRutherford::IsApplicable(const G4ParticleDefinition& particle)
126//{
127// return ( &particle == G4Electron::Electron() );
128//}
129
130
131G4double G4CrossSectionElasticScreenedRutherford::ScreeningFactor(G4double k, G4double z)
132{
133 //
134 // alpha_1 + beta_1 ln(K/eV) constK Z^(2/3)
135 // n(T) = -------------------------- -----------------
136 // K/(m_e c^2) 2 + K/(m_e c^2)
137 //
138 // Where K is the electron non-relativistic kinetic energy
139 //
140 // n(T) > 0 for T < ~ 400 MeV
141 //
142 // NIM 155, pp. 145-156, 1978
143 // Formulae (2) and (5)
144
145 const G4double alpha_1(1.64);
146 const G4double beta_1(-0.0825);
147 const G4double constK(1.7E-5);
148
149 G4double numerator = (alpha_1 + beta_1 * std::log(k/eV)) * constK * std::pow(z, 2./3.);
150
151 k /= electron_mass_c2;
152
153 G4double denominator = k * (2 + k);
154
155 G4double value = 0.;
156 if (denominator > 0.) value = numerator / denominator;
157
158 return value;
159
160}
161
Note: See TracBrowser for help on using the repository browser.