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: G4DNAEmfietzoglouExcitationModel.cc,v 1.5 2009/02/16 12:46:58 sincerti Exp $ |
---|
27 | // GEANT4 tag $Name: geant4-09-03-beta-cand-01 $ |
---|
28 | // |
---|
29 | |
---|
30 | #include "G4DNAEmfietzoglouExcitationModel.hh" |
---|
31 | |
---|
32 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
33 | |
---|
34 | using namespace std; |
---|
35 | |
---|
36 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
37 | |
---|
38 | G4DNAEmfietzoglouExcitationModel::G4DNAEmfietzoglouExcitationModel(const G4ParticleDefinition*, |
---|
39 | const G4String& nam) |
---|
40 | :G4VEmModel(nam),isInitialised(false) |
---|
41 | { |
---|
42 | |
---|
43 | lowEnergyLimit = 8.23 * eV; |
---|
44 | highEnergyLimit = 10 * MeV; |
---|
45 | SetLowEnergyLimit(lowEnergyLimit); |
---|
46 | SetHighEnergyLimit(highEnergyLimit); |
---|
47 | |
---|
48 | verboseLevel= 0; |
---|
49 | // Verbosity scale: |
---|
50 | // 0 = nothing |
---|
51 | // 1 = warning for energy non-conservation |
---|
52 | // 2 = details of energy budget |
---|
53 | // 3 = calculation of cross sections, file openings, sampling of atoms |
---|
54 | // 4 = entering in methods |
---|
55 | |
---|
56 | if (verboseLevel > 3) |
---|
57 | |
---|
58 | G4cout << "Emfietzoglou Excitation model is constructed " << G4endl |
---|
59 | << "Energy range: " |
---|
60 | << lowEnergyLimit / eV << " eV - " |
---|
61 | << highEnergyLimit / MeV << " MeV" |
---|
62 | << G4endl; |
---|
63 | |
---|
64 | } |
---|
65 | |
---|
66 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
67 | |
---|
68 | G4DNAEmfietzoglouExcitationModel::~G4DNAEmfietzoglouExcitationModel() |
---|
69 | {} |
---|
70 | |
---|
71 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
72 | |
---|
73 | void G4DNAEmfietzoglouExcitationModel::Initialise(const G4ParticleDefinition* /*particle*/, |
---|
74 | const G4DataVector& /*cuts*/) |
---|
75 | { |
---|
76 | |
---|
77 | if (verboseLevel > 3) |
---|
78 | G4cout << "Calling G4DNAEmfietzoglouExcitationModel::Initialise()" << G4endl; |
---|
79 | |
---|
80 | // Energy limits |
---|
81 | |
---|
82 | if (LowEnergyLimit() < lowEnergyLimit) |
---|
83 | { |
---|
84 | G4cout << "G4DNAEmfietzoglouExcitationModel: low energy limit increased from " << |
---|
85 | LowEnergyLimit()/eV << " eV to " << lowEnergyLimit/eV << " eV" << G4endl; |
---|
86 | SetLowEnergyLimit(lowEnergyLimit); |
---|
87 | } |
---|
88 | |
---|
89 | if (HighEnergyLimit() > highEnergyLimit) |
---|
90 | { |
---|
91 | G4cout << "G4DNAEmfietzoglouExcitationModel: high energy limit decreased from " << |
---|
92 | HighEnergyLimit()/MeV << " MeV to " << highEnergyLimit/MeV << " MeV" << G4endl; |
---|
93 | SetHighEnergyLimit(highEnergyLimit); |
---|
94 | } |
---|
95 | |
---|
96 | // |
---|
97 | |
---|
98 | nLevels = waterExcitation.NumberOfLevels(); |
---|
99 | |
---|
100 | // |
---|
101 | |
---|
102 | G4cout << "Emfietzoglou Excitation model is initialized " << G4endl |
---|
103 | << "Energy range: " |
---|
104 | << LowEnergyLimit() / eV << " eV - " |
---|
105 | << HighEnergyLimit() / MeV << " MeV" |
---|
106 | << G4endl; |
---|
107 | |
---|
108 | if(!isInitialised) |
---|
109 | { |
---|
110 | isInitialised = true; |
---|
111 | |
---|
112 | if(pParticleChange) |
---|
113 | fParticleChangeForGamma = reinterpret_cast<G4ParticleChangeForGamma*>(pParticleChange); |
---|
114 | else |
---|
115 | fParticleChangeForGamma = new G4ParticleChangeForGamma(); |
---|
116 | } |
---|
117 | |
---|
118 | // InitialiseElementSelectors(particle,cuts); |
---|
119 | |
---|
120 | // Test if water material |
---|
121 | |
---|
122 | flagMaterialIsWater= false; |
---|
123 | densityWater = 0; |
---|
124 | |
---|
125 | const G4ProductionCutsTable* theCoupleTable = G4ProductionCutsTable::GetProductionCutsTable(); |
---|
126 | |
---|
127 | if(theCoupleTable) |
---|
128 | { |
---|
129 | G4int numOfCouples = theCoupleTable->GetTableSize(); |
---|
130 | |
---|
131 | if(numOfCouples>0) |
---|
132 | { |
---|
133 | for (G4int i=0; i<numOfCouples; i++) |
---|
134 | { |
---|
135 | const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(i); |
---|
136 | const G4Material* material = couple->GetMaterial(); |
---|
137 | |
---|
138 | size_t j = material->GetNumberOfElements(); |
---|
139 | while (j>0) |
---|
140 | { |
---|
141 | j--; |
---|
142 | const G4Element* element(material->GetElement(j)); |
---|
143 | if (element->GetZ() == 8.) |
---|
144 | { |
---|
145 | G4double density = material->GetAtomicNumDensityVector()[j]; |
---|
146 | if (density > 0.) |
---|
147 | { |
---|
148 | flagMaterialIsWater = true; |
---|
149 | densityWater = density; |
---|
150 | |
---|
151 | if (verboseLevel > 3) |
---|
152 | G4cout << "Water material is found with density(cm^-3)=" << density/(cm*cm*cm) << G4endl; |
---|
153 | } |
---|
154 | } |
---|
155 | } |
---|
156 | |
---|
157 | } |
---|
158 | } // if(numOfCouples>0) |
---|
159 | |
---|
160 | } // if (theCoupleTable) |
---|
161 | |
---|
162 | } |
---|
163 | |
---|
164 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
165 | |
---|
166 | G4double G4DNAEmfietzoglouExcitationModel::CrossSectionPerVolume(const G4Material*, |
---|
167 | const G4ParticleDefinition* particleDefinition, |
---|
168 | G4double ekin, |
---|
169 | G4double, |
---|
170 | G4double) |
---|
171 | { |
---|
172 | if (verboseLevel > 3) |
---|
173 | G4cout << "Calling CrossSectionPerVolume() of G4DNAEmfietzoglouExcitationModel" << G4endl; |
---|
174 | |
---|
175 | // Calculate total cross section for model |
---|
176 | |
---|
177 | G4double sigma=0; |
---|
178 | |
---|
179 | if (flagMaterialIsWater) |
---|
180 | { |
---|
181 | |
---|
182 | if (particleDefinition == G4Electron::ElectronDefinition()) |
---|
183 | { |
---|
184 | if (ekin >= lowEnergyLimit && ekin < highEnergyLimit) |
---|
185 | { |
---|
186 | sigma = Sum(ekin); |
---|
187 | } |
---|
188 | } |
---|
189 | |
---|
190 | if (verboseLevel > 3) |
---|
191 | { |
---|
192 | G4cout << "---> Kinetic energy(eV)=" << ekin/eV << G4endl; |
---|
193 | G4cout << " - Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl; |
---|
194 | G4cout << " - Cross section per water molecule (cm^-1)=" << sigma*densityWater/(1./cm) << G4endl; |
---|
195 | } |
---|
196 | |
---|
197 | } // if (flagMaterialIsWater) |
---|
198 | |
---|
199 | return sigma*densityWater; |
---|
200 | } |
---|
201 | |
---|
202 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
203 | |
---|
204 | void G4DNAEmfietzoglouExcitationModel::SampleSecondaries(std::vector<G4DynamicParticle*>* /*fvect*/, |
---|
205 | const G4MaterialCutsCouple* /*couple*/, |
---|
206 | const G4DynamicParticle* aDynamicElectron, |
---|
207 | G4double, |
---|
208 | G4double) |
---|
209 | { |
---|
210 | |
---|
211 | if (verboseLevel > 3) |
---|
212 | G4cout << "Calling SampleSecondaries() of G4DNAEmfietzoglouExcitationModel" << G4endl; |
---|
213 | |
---|
214 | G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy(); |
---|
215 | |
---|
216 | G4int level = RandomSelect(electronEnergy0); |
---|
217 | |
---|
218 | G4double excitationEnergy = waterExcitation.ExcitationEnergy(level); |
---|
219 | G4double newEnergy = electronEnergy0 - excitationEnergy; |
---|
220 | |
---|
221 | if (electronEnergy0 < highEnergyLimit) |
---|
222 | { |
---|
223 | if (newEnergy >= lowEnergyLimit) |
---|
224 | { |
---|
225 | fParticleChangeForGamma->ProposeMomentumDirection(aDynamicElectron->GetMomentumDirection()); |
---|
226 | fParticleChangeForGamma->SetProposedKineticEnergy(newEnergy); |
---|
227 | fParticleChangeForGamma->ProposeLocalEnergyDeposit(excitationEnergy); |
---|
228 | } |
---|
229 | |
---|
230 | else |
---|
231 | { |
---|
232 | fParticleChangeForGamma->ProposeTrackStatus(fStopAndKill); |
---|
233 | fParticleChangeForGamma->ProposeLocalEnergyDeposit(electronEnergy0); |
---|
234 | } |
---|
235 | } |
---|
236 | } |
---|
237 | |
---|
238 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
239 | |
---|
240 | G4double G4DNAEmfietzoglouExcitationModel::PartialCrossSection(G4double t, G4int level) |
---|
241 | { |
---|
242 | // Aj T |
---|
243 | // Sigma(T) = ------------- (Bj / T) ln(Cj ---) [1 - Bj / T]^Pj |
---|
244 | // 2 pi alpha0 R |
---|
245 | // |
---|
246 | // Sigma is the macroscopic cross section = N sigma, where N = number of target particles per unit volume |
---|
247 | // and sigma is the microscopic cross section |
---|
248 | // T is the incoming electron kinetic energy |
---|
249 | // alpha0 is the Bohr Radius (Bohr_radius) |
---|
250 | // Aj, Bj, Cj & Pj are parameters that can be found in Emfietzoglou's papers |
---|
251 | // |
---|
252 | // From Phys. Med. Biol. 48 (2003) 2355-2371, D.Emfietzoglou, |
---|
253 | // Monte Carlo Simulation of the energy loss of low energy electrons in liquid Water |
---|
254 | // |
---|
255 | // Scaling for macroscopic cross section: number of water moleculs per unit volume |
---|
256 | // const G4double sigma0 = (10. / 3.343e22) * cm2; |
---|
257 | |
---|
258 | const G4double density = 3.34192e+19 * mm3; |
---|
259 | |
---|
260 | const G4double aj[]={0.0205, 0.0209, 0.0130, 0.0026, 0.0025}; |
---|
261 | const G4double cj[]={4.9801, 3.3850, 2.8095, 1.9242, 3.4624}; |
---|
262 | const G4double pj[]={0.4757, 0.3483, 0.4443, 0.3429, 0.4379}; |
---|
263 | const G4double r = 13.6 * eV; |
---|
264 | |
---|
265 | G4double sigma = 0.; |
---|
266 | |
---|
267 | G4double exc = waterExcitation.ExcitationEnergy(level); |
---|
268 | |
---|
269 | if (t >= exc) |
---|
270 | { |
---|
271 | G4double excitationSigma = ( aj[level] / (2.*pi*Bohr_radius)) |
---|
272 | * (exc / t) |
---|
273 | * std::log(cj[level]*(t/r)) |
---|
274 | * std::pow((1.- (exc/t)), pj[level]); |
---|
275 | sigma = excitationSigma / density; |
---|
276 | } |
---|
277 | return sigma; |
---|
278 | } |
---|
279 | |
---|
280 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
281 | |
---|
282 | G4int G4DNAEmfietzoglouExcitationModel::RandomSelect(G4double k) |
---|
283 | { |
---|
284 | G4int i = nLevels; |
---|
285 | G4double value = 0.; |
---|
286 | std::deque<double> values; |
---|
287 | |
---|
288 | while (i > 0) |
---|
289 | { |
---|
290 | i--; |
---|
291 | G4double partial = PartialCrossSection(k,i); |
---|
292 | values.push_front(partial); |
---|
293 | value += partial; |
---|
294 | } |
---|
295 | |
---|
296 | value *= G4UniformRand(); |
---|
297 | |
---|
298 | i = nLevels; |
---|
299 | |
---|
300 | while (i > 0) |
---|
301 | { |
---|
302 | i--; |
---|
303 | if (values[i] > value) return i; |
---|
304 | value -= values[i]; |
---|
305 | } |
---|
306 | |
---|
307 | return 0; |
---|
308 | } |
---|
309 | |
---|
310 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
311 | |
---|
312 | G4double G4DNAEmfietzoglouExcitationModel::Sum(G4double k) |
---|
313 | { |
---|
314 | G4double totalCrossSection = 0.; |
---|
315 | |
---|
316 | for (G4int i=0; i<nLevels; i++) |
---|
317 | { |
---|
318 | totalCrossSection += PartialCrossSection(k,i); |
---|
319 | } |
---|
320 | return totalCrossSection; |
---|
321 | } |
---|
322 | |
---|