source: trunk/source/processes/electromagnetic/lowenergy/src/G4eCrossSectionScreenedRutherford.cc @ 1228

Last change on this file since 1228 was 1228, checked in by garnier, 14 years ago

update geant4.9.3 tag

File size: 5.5 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: G4eCrossSectionScreenedRutherford.cc,v 1.3 2007/10/12 12:27:19 pia Exp $
28// GEANT4 tag $Name: geant4-09-03 $
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// Further documentation available from http://www.ge.infn.it/geant4/dna
45
46// -------------------------------------------------------------------
47
48
49#include "G4eCrossSectionScreenedRutherford.hh"
50#include "G4Track.hh"
51#include "G4DynamicParticle.hh"
52#include "G4ParticleDefinition.hh"
53#include "G4Electron.hh"
54
55
56G4eCrossSectionScreenedRutherford::G4eCrossSectionScreenedRutherford()
57{
58
59  name = "eCrossSectionScreenedRutherford";
60  lowEnergyLimit = 7. * eV;
61  highEnergyLimit = 10 * MeV;
62
63//  if (verboseLevel > 0)
64//  {
65//    G4cout << name << " is created " << G4endl
66//     << "Energy range: "
67//     << lowEnergyLimit / keV << " keV - "
68//     << highEnergyLimit / GeV << " GeV"
69//     << G4endl;
70//  }
71}
72
73
74G4eCrossSectionScreenedRutherford::~G4eCrossSectionScreenedRutherford()
75{ }
76 
77
78G4double G4eCrossSectionScreenedRutherford::CrossSection(const G4Track& track)
79{
80  const G4DynamicParticle* particle = track.GetDynamicParticle();
81  G4double k = particle->GetKineticEnergy();
82
83  // Cross section = 0 outside the energy validity limits set in the constructor
84  // ---- MGP ---- Better handling of these limits to be set in a following design iteration
85
86  G4double screenedCrossSection = 0.;
87
88  if (k > lowEnergyLimit && k < highEnergyLimit)
89    {     
90      // G4Material* material = track.GetMaterial();
91
92      // Assume that the material is water; proper algorithm to calculate z correctly for any material to be inserted here
93      // For H20 Z = 10 (total number of electrons)
94      G4double z = 10.;
95     
96      G4double n = ScreeningFactor(k,z);
97      G4double crossSection = RutherfordCrossSection(k, z);
98      screenedCrossSection = pi *  crossSection / (n * (n + 1.));
99    }   
100
101  return screenedCrossSection;
102}
103 
104G4double G4eCrossSectionScreenedRutherford::RutherfordCrossSection(G4double k, G4double z)
105{
106  //   
107  //                               e^4         /      K + m_e c^2      \^2
108  // sigma_Ruth(K) = Z (Z+1) -------------------- | --------------------- |
109  //                          (4 pi epsilon_0)^2  \  K * (K + 2 m_e c^2)  /
110  //
111  // Where K is the electron non-relativistic kinetic energy
112  //
113  // NIM 155, pp. 145-156, 1978
114 
115  G4double length =(e_squared * (k + electron_mass_c2)) / (4 * pi *epsilon0 * k * ( k + 2 * electron_mass_c2));
116  G4double cross = z * ( z + 1) * length * length;
117 
118  return cross;
119}
120
121//G4bool G4eCrossSectionScreenedRutherford::IsApplicable(const G4ParticleDefinition& particle)
122//{
123//  return ( &particle == G4Electron::Electron() );
124//}
125
126
127G4double G4eCrossSectionScreenedRutherford::ScreeningFactor(G4double k, G4double z)
128{
129  //
130  //         alpha_1 + beta_1 ln(K/eV)   constK Z^(2/3)
131  // n(T) = -------------------------- -----------------
132  //              K/(m_e c^2)            2 + K/(m_e c^2)
133  //
134  // Where K is the electron non-relativistic kinetic energy
135  //
136  // n(T) > 0 for T < ~ 400 MeV
137  //
138  // NIM 155, pp. 145-156, 1978
139  // Formulae (2) and (5)
140
141  const G4double alpha_1(1.64);
142  const G4double beta_1(-0.0825);
143  const G4double constK(1.7E-5);
144
145  G4double numerator = (alpha_1 + beta_1 * std::log(k/eV)) * constK * std::pow(z, 2./3.);
146
147  k /= electron_mass_c2;
148
149  G4double denominator = k * (2 + k);
150
151  G4double value = 0.;
152  if (denominator > 0.) value = numerator / denominator;
153
154  return value;
155
156}
157
Note: See TracBrowser for help on using the repository browser.