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: G4CrossSectionElasticScreenedRutherfordHE.cc,v 1.5 2009/06/11 15:47:08 mantero Exp $ |
---|
27 | // GEANT4 tag $Name: geant4-09-04-beta-cand-01 $ |
---|
28 | |
---|
29 | #include "G4CrossSectionElasticScreenedRutherfordHE.hh" |
---|
30 | |
---|
31 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
32 | |
---|
33 | G4CrossSectionElasticScreenedRutherfordHE::G4CrossSectionElasticScreenedRutherfordHE() |
---|
34 | { |
---|
35 | lowEnergyLimit = 200. * eV; |
---|
36 | highEnergyLimit = 10. * MeV; |
---|
37 | |
---|
38 | G4cout << G4endl; |
---|
39 | G4cout << "*******************************************************************************" << G4endl; |
---|
40 | G4cout << "*******************************************************************************" << G4endl; |
---|
41 | G4cout << " The class G4CrossSectionElasticScreenedRutherfordHE is NOT SUPPORTED ANYMORE. " << G4endl; |
---|
42 | G4cout << " It will be REMOVED with the next major release of Geant4. " << G4endl; |
---|
43 | G4cout << " Please consult: https://twiki.cern.ch/twiki/bin/view/Geant4/LoweProcesses" << G4endl; |
---|
44 | G4cout << "*******************************************************************************" << G4endl; |
---|
45 | G4cout << "*******************************************************************************" << G4endl; |
---|
46 | G4cout << G4endl; |
---|
47 | } |
---|
48 | |
---|
49 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
50 | |
---|
51 | G4CrossSectionElasticScreenedRutherfordHE::~G4CrossSectionElasticScreenedRutherfordHE() |
---|
52 | {} |
---|
53 | |
---|
54 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
55 | |
---|
56 | G4double G4CrossSectionElasticScreenedRutherfordHE::CrossSection(const G4Track& track) |
---|
57 | { |
---|
58 | const G4DynamicParticle* particle = track.GetDynamicParticle(); |
---|
59 | G4double k = particle->GetKineticEnergy(); |
---|
60 | |
---|
61 | G4double screenedCrossSection = 0.; |
---|
62 | |
---|
63 | if (k > lowEnergyLimit && k <= highEnergyLimit) |
---|
64 | { |
---|
65 | G4double z = 10.; |
---|
66 | G4double n = ScreeningFactor(k,z); |
---|
67 | G4double crossSection = RutherfordCrossSection(k, z); |
---|
68 | screenedCrossSection = pi * crossSection / (n * (n + 1.)); |
---|
69 | } |
---|
70 | |
---|
71 | return screenedCrossSection; |
---|
72 | } |
---|
73 | |
---|
74 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
75 | |
---|
76 | G4double G4CrossSectionElasticScreenedRutherfordHE::RutherfordCrossSection(G4double k, G4double z) |
---|
77 | { |
---|
78 | // |
---|
79 | // e^4 / K + m_e c^2 \^2 |
---|
80 | // sigma_Ruth(K) = Z (Z+1) -------------------- | --------------------- | |
---|
81 | // (4 pi epsilon_0)^2 \ K * (K + 2 m_e c^2) / |
---|
82 | // |
---|
83 | // Where K is the electron non-relativistic kinetic energy |
---|
84 | // |
---|
85 | // NIM 155, pp. 145-156, 1978 |
---|
86 | |
---|
87 | G4double length =(e_squared * (k + electron_mass_c2)) / (4 * pi *epsilon0 * k * ( k + 2 * electron_mass_c2)); |
---|
88 | G4double cross = z * ( z + 1) * length * length; |
---|
89 | |
---|
90 | return cross; |
---|
91 | } |
---|
92 | |
---|
93 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
94 | |
---|
95 | G4double G4CrossSectionElasticScreenedRutherfordHE::ScreeningFactor(G4double k, G4double z) |
---|
96 | { |
---|
97 | // |
---|
98 | // alpha_1 + beta_1 ln(K/eV) constK Z^(2/3) |
---|
99 | // n(T) = -------------------------- ----------------- |
---|
100 | // K/(m_e c^2) 2 + K/(m_e c^2) |
---|
101 | // |
---|
102 | // Where K is the electron non-relativistic kinetic energy |
---|
103 | // |
---|
104 | // n(T) > 0 for T < ~ 400 MeV |
---|
105 | // |
---|
106 | // NIM 155, pp. 145-156, 1978 |
---|
107 | // Formulae (2) and (5) |
---|
108 | |
---|
109 | const G4double alpha_1(1.64); |
---|
110 | const G4double beta_1(-0.0825); |
---|
111 | const G4double constK(1.7E-5); |
---|
112 | |
---|
113 | G4double numerator = (alpha_1 + beta_1 * std::log(k/eV)) * constK * std::pow(z, 2./3.); |
---|
114 | |
---|
115 | k /= electron_mass_c2; |
---|
116 | |
---|
117 | G4double denominator = k * (2 + k); |
---|
118 | |
---|
119 | G4double value = 0.; |
---|
120 | if (denominator > 0.) value = numerator / denominator; |
---|
121 | |
---|
122 | return value; |
---|
123 | } |
---|
124 | |
---|