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: G4CrossSectionExcitationEmfietzoglouPartial.cc,v 1.3 2008/07/14 20:47:34 sincerti Exp $ |
---|
27 | // GEANT4 tag $Name: geant4-09-02-ref-02 $ |
---|
28 | |
---|
29 | #include "G4CrossSectionExcitationEmfietzoglouPartial.hh" |
---|
30 | |
---|
31 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
32 | |
---|
33 | G4CrossSectionExcitationEmfietzoglouPartial::G4CrossSectionExcitationEmfietzoglouPartial() |
---|
34 | { |
---|
35 | nLevels = waterExcitation.NumberOfLevels(); |
---|
36 | } |
---|
37 | |
---|
38 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
39 | |
---|
40 | G4CrossSectionExcitationEmfietzoglouPartial::~G4CrossSectionExcitationEmfietzoglouPartial() |
---|
41 | {} |
---|
42 | |
---|
43 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
44 | |
---|
45 | G4double G4CrossSectionExcitationEmfietzoglouPartial::CrossSection(G4double t, G4int level) |
---|
46 | { |
---|
47 | // Aj T |
---|
48 | // Sigma(T) = ------------- (Bj / T) ln(Cj ---) [1 - Bj / T]^Pj |
---|
49 | // 2 pi alpha0 R |
---|
50 | // |
---|
51 | // Sigma is the macroscopic cross section = N sigma, where N = number of target particles per unit volume |
---|
52 | // and sigma is the microscopic cross section |
---|
53 | // T is the incoming electron kinetic energy |
---|
54 | // alpha0 is the Bohr Radius (Bohr_radius) |
---|
55 | // Aj, Bj, Cj & Pj are parameters that can be found in Emfietzoglou's papers |
---|
56 | // |
---|
57 | // From Phys. Med. Biol. 48 (2003) 2355-2371, D.Emfietzoglou, |
---|
58 | // Monte Carlo Simulation of the energy loss of low energy electrons in liquid Water |
---|
59 | // |
---|
60 | // Scaling for macroscopic cross section: number of water moleculs per unit volume |
---|
61 | // const G4double sigma0 = (10. / 3.343e22) * cm2; |
---|
62 | |
---|
63 | const G4double density = 3.34192e+19 * mm3; |
---|
64 | |
---|
65 | const G4double aj[]={0.0205, 0.0209, 0.0130, 0.0026, 0.0025}; |
---|
66 | const G4double cj[]={4.9801, 3.3850, 2.8095, 1.9242, 3.4624}; |
---|
67 | const G4double pj[]={0.4757, 0.3483, 0.4443, 0.3429, 0.4379}; |
---|
68 | const G4double r = 13.6 * eV; |
---|
69 | |
---|
70 | G4double sigma = 0.; |
---|
71 | |
---|
72 | G4double exc = waterExcitation.ExcitationEnergy(level); |
---|
73 | |
---|
74 | if (t >= exc) |
---|
75 | { |
---|
76 | G4double excitationSigma = ( aj[level] / (2.*pi*Bohr_radius)) |
---|
77 | * (exc / t) |
---|
78 | * std::log(cj[level]*(t/r)) |
---|
79 | * std::pow((1.- (exc/t)), pj[level]); |
---|
80 | sigma = excitationSigma / density; |
---|
81 | } |
---|
82 | return sigma; |
---|
83 | } |
---|
84 | |
---|
85 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
86 | |
---|
87 | G4int G4CrossSectionExcitationEmfietzoglouPartial::RandomSelect(G4double k) |
---|
88 | { |
---|
89 | G4int i = nLevels; |
---|
90 | G4double value = 0.; |
---|
91 | std::deque<double> values; |
---|
92 | |
---|
93 | // ---- MGP ---- The following algorithm is wrong: it works if the cross section |
---|
94 | // is a monotone increasing function. |
---|
95 | // The algorithm should be corrected by building the cumulative function |
---|
96 | // of the cross section and comparing a random number in the range 0-1 against |
---|
97 | // the cumulative value at each bin |
---|
98 | |
---|
99 | while (i > 0) |
---|
100 | { |
---|
101 | i--; |
---|
102 | G4double partial = CrossSection(k,i); |
---|
103 | values.push_front(partial); |
---|
104 | value += partial; |
---|
105 | } |
---|
106 | |
---|
107 | value *= G4UniformRand(); |
---|
108 | |
---|
109 | i = nLevels; |
---|
110 | |
---|
111 | while (i > 0) |
---|
112 | { |
---|
113 | i--; |
---|
114 | if (values[i] > value) return i; |
---|
115 | value -= values[i]; |
---|
116 | } |
---|
117 | |
---|
118 | return 0; |
---|
119 | } |
---|
120 | |
---|
121 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
122 | |
---|
123 | G4double G4CrossSectionExcitationEmfietzoglouPartial::Sum(G4double k) |
---|
124 | { |
---|
125 | G4double totalCrossSection = 0.; |
---|
126 | |
---|
127 | for (G4int i=0; i<nLevels; i++) |
---|
128 | { |
---|
129 | totalCrossSection += CrossSection(k,i); |
---|
130 | } |
---|
131 | return totalCrossSection; |
---|
132 | } |
---|