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: G4DNAChampionElasticModel.cc,v 1.15 2010/10/17 11:28:51 sincerti Exp $ |
---|
27 | // GEANT4 tag $Name: emlowen-V09-03-54 $ |
---|
28 | // |
---|
29 | |
---|
30 | #include "G4DNAChampionElasticModel.hh" |
---|
31 | |
---|
32 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
33 | |
---|
34 | using namespace std; |
---|
35 | |
---|
36 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
37 | |
---|
38 | G4DNAChampionElasticModel::G4DNAChampionElasticModel(const G4ParticleDefinition*, |
---|
39 | const G4String& nam) |
---|
40 | :G4VEmModel(nam),isInitialised(false) |
---|
41 | { |
---|
42 | |
---|
43 | killBelowEnergy = 0.025*eV; // Minimum e- energy for energy loss by excitation |
---|
44 | lowEnergyLimit = 0 * eV; |
---|
45 | lowEnergyLimitOfModel = 0.025 * eV; |
---|
46 | highEnergyLimit = 1. * MeV; |
---|
47 | SetLowEnergyLimit(lowEnergyLimit); |
---|
48 | SetHighEnergyLimit(highEnergyLimit); |
---|
49 | |
---|
50 | verboseLevel= 0; |
---|
51 | // Verbosity scale: |
---|
52 | // 0 = nothing |
---|
53 | // 1 = warning for energy non-conservation |
---|
54 | // 2 = details of energy budget |
---|
55 | // 3 = calculation of cross sections, file openings, sampling of atoms |
---|
56 | // 4 = entering in methods |
---|
57 | |
---|
58 | if( verboseLevel>0 ) |
---|
59 | { |
---|
60 | G4cout << "Champion Elastic model is constructed " << G4endl |
---|
61 | << "Energy range: " |
---|
62 | << lowEnergyLimit / eV << " eV - " |
---|
63 | << highEnergyLimit / MeV << " MeV" |
---|
64 | << G4endl; |
---|
65 | } |
---|
66 | } |
---|
67 | |
---|
68 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
69 | |
---|
70 | G4DNAChampionElasticModel::~G4DNAChampionElasticModel() |
---|
71 | { |
---|
72 | // For total cross section |
---|
73 | |
---|
74 | std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos; |
---|
75 | for (pos = tableData.begin(); pos != tableData.end(); ++pos) |
---|
76 | { |
---|
77 | G4DNACrossSectionDataSet* table = pos->second; |
---|
78 | delete table; |
---|
79 | } |
---|
80 | |
---|
81 | // For final state |
---|
82 | |
---|
83 | eVecm.clear(); |
---|
84 | |
---|
85 | } |
---|
86 | |
---|
87 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
88 | |
---|
89 | void G4DNAChampionElasticModel::Initialise(const G4ParticleDefinition* /*particle*/, |
---|
90 | const G4DataVector& /*cuts*/) |
---|
91 | { |
---|
92 | |
---|
93 | if (verboseLevel > 3) |
---|
94 | G4cout << "Calling G4DNAChampionElasticModel::Initialise()" << G4endl; |
---|
95 | |
---|
96 | // Energy limits |
---|
97 | |
---|
98 | if (LowEnergyLimit() < lowEnergyLimit) |
---|
99 | { |
---|
100 | G4cout << "G4DNAChampionElasticModel: low energy limit increased from " << |
---|
101 | LowEnergyLimit()/eV << " eV to " << lowEnergyLimit/eV << " eV" << G4endl; |
---|
102 | SetLowEnergyLimit(lowEnergyLimit); |
---|
103 | } |
---|
104 | |
---|
105 | if (HighEnergyLimit() > highEnergyLimit) |
---|
106 | { |
---|
107 | G4cout << "G4DNAChampionElasticModel: high energy limit decreased from " << |
---|
108 | HighEnergyLimit()/MeV << " MeV to " << highEnergyLimit/MeV << " MeV" << G4endl; |
---|
109 | SetHighEnergyLimit(highEnergyLimit); |
---|
110 | } |
---|
111 | |
---|
112 | // Reading of data files |
---|
113 | |
---|
114 | G4double scaleFactor = 1e-16*cm*cm; |
---|
115 | |
---|
116 | G4String fileElectron("dna/sigma_elastic_e_champion"); |
---|
117 | |
---|
118 | G4ParticleDefinition* electronDef = G4Electron::ElectronDefinition(); |
---|
119 | G4String electron; |
---|
120 | |
---|
121 | if (electronDef != 0) |
---|
122 | { |
---|
123 | // For total cross section |
---|
124 | |
---|
125 | electron = electronDef->GetParticleName(); |
---|
126 | |
---|
127 | tableFile[electron] = fileElectron; |
---|
128 | |
---|
129 | G4DNACrossSectionDataSet* tableE = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, eV,scaleFactor ); |
---|
130 | tableE->LoadData(fileElectron); |
---|
131 | tableData[electron] = tableE; |
---|
132 | |
---|
133 | // For final state |
---|
134 | |
---|
135 | char *path = getenv("G4LEDATA"); |
---|
136 | |
---|
137 | if (!path) |
---|
138 | G4Exception("G4FinalStateElasticChampion::Initialise: G4LEDATA environment variable not set"); |
---|
139 | |
---|
140 | std::ostringstream eFullFileName; |
---|
141 | eFullFileName << path << "/dna/sigmadiff_cumulatedshort_elastic_e_champion.dat"; |
---|
142 | std::ifstream eDiffCrossSection(eFullFileName.str().c_str()); |
---|
143 | |
---|
144 | if (!eDiffCrossSection) G4Exception("G4DNAChampionElasticModel::Initialise: error opening electron DATA FILE"); |
---|
145 | |
---|
146 | eTdummyVec.push_back(0.); |
---|
147 | |
---|
148 | while(!eDiffCrossSection.eof()) |
---|
149 | { |
---|
150 | double tDummy; |
---|
151 | double eDummy; |
---|
152 | eDiffCrossSection>>tDummy>>eDummy; |
---|
153 | |
---|
154 | // SI : mandatory eVecm initialization |
---|
155 | |
---|
156 | if (tDummy != eTdummyVec.back()) |
---|
157 | { |
---|
158 | eTdummyVec.push_back(tDummy); |
---|
159 | eVecm[tDummy].push_back(0.); |
---|
160 | } |
---|
161 | |
---|
162 | eDiffCrossSection>>eDiffCrossSectionData[tDummy][eDummy]; |
---|
163 | |
---|
164 | if (eDummy != eVecm[tDummy].back()) eVecm[tDummy].push_back(eDummy); |
---|
165 | |
---|
166 | } |
---|
167 | |
---|
168 | // End final state |
---|
169 | |
---|
170 | } |
---|
171 | else G4Exception("G4DNAChampionElasticModel::Initialise: electron is not defined"); |
---|
172 | |
---|
173 | if (verboseLevel > 2) |
---|
174 | G4cout << "Loaded cross section files for Champion Elastic model" << G4endl; |
---|
175 | |
---|
176 | if( verboseLevel>0 ) |
---|
177 | { |
---|
178 | G4cout << "Champion Elastic model is initialized " << G4endl |
---|
179 | << "Energy range: " |
---|
180 | << LowEnergyLimit() / eV << " eV - " |
---|
181 | << HighEnergyLimit() / MeV << " MeV" |
---|
182 | << G4endl; |
---|
183 | } |
---|
184 | |
---|
185 | if(!isInitialised) |
---|
186 | { |
---|
187 | isInitialised = true; |
---|
188 | |
---|
189 | if(pParticleChange) |
---|
190 | fParticleChangeForGamma = reinterpret_cast<G4ParticleChangeForGamma*>(pParticleChange); |
---|
191 | else |
---|
192 | fParticleChangeForGamma = new G4ParticleChangeForGamma(); |
---|
193 | } |
---|
194 | |
---|
195 | // InitialiseElementSelectors(particle,cuts); |
---|
196 | |
---|
197 | } |
---|
198 | |
---|
199 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
200 | |
---|
201 | G4double G4DNAChampionElasticModel::CrossSectionPerVolume(const G4Material* material, |
---|
202 | const G4ParticleDefinition* p, |
---|
203 | G4double ekin, |
---|
204 | G4double, |
---|
205 | G4double) |
---|
206 | { |
---|
207 | if (verboseLevel > 3) |
---|
208 | G4cout << "Calling CrossSectionPerVolume() of G4DNAChampionElasticModel" << G4endl; |
---|
209 | |
---|
210 | // Calculate total cross section for model |
---|
211 | |
---|
212 | G4double sigma=0; |
---|
213 | |
---|
214 | if (material->GetName() == "G4_WATER") |
---|
215 | { |
---|
216 | const G4String& particleName = p->GetParticleName(); |
---|
217 | |
---|
218 | if (ekin < highEnergyLimit) |
---|
219 | { |
---|
220 | //SI : XS must not be zero otherwise sampling of secondaries method ignored |
---|
221 | if (ekin < lowEnergyLimitOfModel) ekin = lowEnergyLimitOfModel; |
---|
222 | // |
---|
223 | |
---|
224 | std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos; |
---|
225 | pos = tableData.find(particleName); |
---|
226 | |
---|
227 | if (pos != tableData.end()) |
---|
228 | { |
---|
229 | G4DNACrossSectionDataSet* table = pos->second; |
---|
230 | if (table != 0) |
---|
231 | { |
---|
232 | sigma = table->FindValue(ekin); |
---|
233 | } |
---|
234 | } |
---|
235 | else |
---|
236 | { |
---|
237 | G4Exception("G4DNAChampionElasticModel::ComputeCrossSectionPerVolume: attempting to calculate cross section for wrong particle"); |
---|
238 | } |
---|
239 | } |
---|
240 | |
---|
241 | if (verboseLevel > 3) |
---|
242 | { |
---|
243 | G4cout << "---> Kinetic energy(eV)=" << ekin/eV << G4endl; |
---|
244 | G4cout << " - Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl; |
---|
245 | G4cout << " - Cross section per water molecule (cm^-1)=" << sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl; |
---|
246 | } |
---|
247 | |
---|
248 | } |
---|
249 | |
---|
250 | return sigma*material->GetAtomicNumDensityVector()[1]; |
---|
251 | } |
---|
252 | |
---|
253 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
254 | |
---|
255 | void G4DNAChampionElasticModel::SampleSecondaries(std::vector<G4DynamicParticle*>* /*fvect*/, |
---|
256 | const G4MaterialCutsCouple* /*couple*/, |
---|
257 | const G4DynamicParticle* aDynamicElectron, |
---|
258 | G4double, |
---|
259 | G4double) |
---|
260 | { |
---|
261 | |
---|
262 | if (verboseLevel > 3) |
---|
263 | G4cout << "Calling SampleSecondaries() of G4DNAChampionElasticModel" << G4endl; |
---|
264 | |
---|
265 | G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy(); |
---|
266 | |
---|
267 | if (electronEnergy0 < killBelowEnergy) |
---|
268 | { |
---|
269 | fParticleChangeForGamma->ProposeTrackStatus(fStopAndKill); |
---|
270 | fParticleChangeForGamma->ProposeLocalEnergyDeposit(electronEnergy0); |
---|
271 | return ; |
---|
272 | } |
---|
273 | |
---|
274 | if (electronEnergy0>= killBelowEnergy && electronEnergy0 < highEnergyLimit) |
---|
275 | { |
---|
276 | |
---|
277 | G4double cosTheta = RandomizeCosTheta(electronEnergy0); |
---|
278 | |
---|
279 | G4double phi = 2. * pi * G4UniformRand(); |
---|
280 | |
---|
281 | G4ThreeVector zVers = aDynamicElectron->GetMomentumDirection(); |
---|
282 | G4ThreeVector xVers = zVers.orthogonal(); |
---|
283 | G4ThreeVector yVers = zVers.cross(xVers); |
---|
284 | |
---|
285 | G4double xDir = std::sqrt(1. - cosTheta*cosTheta); |
---|
286 | G4double yDir = xDir; |
---|
287 | xDir *= std::cos(phi); |
---|
288 | yDir *= std::sin(phi); |
---|
289 | |
---|
290 | G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers)); |
---|
291 | |
---|
292 | fParticleChangeForGamma->ProposeMomentumDirection(zPrimeVers.unit()) ; |
---|
293 | |
---|
294 | fParticleChangeForGamma->SetProposedKineticEnergy(electronEnergy0); |
---|
295 | } |
---|
296 | |
---|
297 | } |
---|
298 | |
---|
299 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
300 | |
---|
301 | G4double G4DNAChampionElasticModel::Theta |
---|
302 | (G4ParticleDefinition * particleDefinition, G4double k, G4double integrDiff) |
---|
303 | { |
---|
304 | G4double theta = 0.; |
---|
305 | G4double valueT1 = 0; |
---|
306 | G4double valueT2 = 0; |
---|
307 | G4double valueE21 = 0; |
---|
308 | G4double valueE22 = 0; |
---|
309 | G4double valueE12 = 0; |
---|
310 | G4double valueE11 = 0; |
---|
311 | G4double xs11 = 0; |
---|
312 | G4double xs12 = 0; |
---|
313 | G4double xs21 = 0; |
---|
314 | G4double xs22 = 0; |
---|
315 | |
---|
316 | if (particleDefinition == G4Electron::ElectronDefinition()) |
---|
317 | { |
---|
318 | std::vector<double>::iterator t2 = std::upper_bound(eTdummyVec.begin(),eTdummyVec.end(), k); |
---|
319 | std::vector<double>::iterator t1 = t2-1; |
---|
320 | |
---|
321 | std::vector<double>::iterator e12 = std::upper_bound(eVecm[(*t1)].begin(),eVecm[(*t1)].end(), integrDiff); |
---|
322 | std::vector<double>::iterator e11 = e12-1; |
---|
323 | |
---|
324 | std::vector<double>::iterator e22 = std::upper_bound(eVecm[(*t2)].begin(),eVecm[(*t2)].end(), integrDiff); |
---|
325 | std::vector<double>::iterator e21 = e22-1; |
---|
326 | |
---|
327 | valueT1 =*t1; |
---|
328 | valueT2 =*t2; |
---|
329 | valueE21 =*e21; |
---|
330 | valueE22 =*e22; |
---|
331 | valueE12 =*e12; |
---|
332 | valueE11 =*e11; |
---|
333 | |
---|
334 | xs11 = eDiffCrossSectionData[valueT1][valueE11]; |
---|
335 | xs12 = eDiffCrossSectionData[valueT1][valueE12]; |
---|
336 | xs21 = eDiffCrossSectionData[valueT2][valueE21]; |
---|
337 | xs22 = eDiffCrossSectionData[valueT2][valueE22]; |
---|
338 | } |
---|
339 | |
---|
340 | if (xs11==0 && xs12==0 && xs21==0 && xs22==0) return (0.); |
---|
341 | |
---|
342 | theta = QuadInterpolator ( valueE11, valueE12, |
---|
343 | valueE21, valueE22, |
---|
344 | xs11, xs12, |
---|
345 | xs21, xs22, |
---|
346 | valueT1, valueT2, |
---|
347 | k, integrDiff ); |
---|
348 | |
---|
349 | return theta; |
---|
350 | } |
---|
351 | |
---|
352 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
353 | |
---|
354 | G4double G4DNAChampionElasticModel::LinLogInterpolate(G4double e1, |
---|
355 | G4double e2, |
---|
356 | G4double e, |
---|
357 | G4double xs1, |
---|
358 | G4double xs2) |
---|
359 | { |
---|
360 | G4double d1 = std::log(xs1); |
---|
361 | G4double d2 = std::log(xs2); |
---|
362 | G4double value = std::exp(d1 + (d2 - d1)*(e - e1)/ (e2 - e1)); |
---|
363 | return value; |
---|
364 | } |
---|
365 | |
---|
366 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
367 | |
---|
368 | G4double G4DNAChampionElasticModel::LinLinInterpolate(G4double e1, |
---|
369 | G4double e2, |
---|
370 | G4double e, |
---|
371 | G4double xs1, |
---|
372 | G4double xs2) |
---|
373 | { |
---|
374 | G4double d1 = xs1; |
---|
375 | G4double d2 = xs2; |
---|
376 | G4double value = (d1 + (d2 - d1)*(e - e1)/ (e2 - e1)); |
---|
377 | return value; |
---|
378 | } |
---|
379 | |
---|
380 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
381 | |
---|
382 | G4double G4DNAChampionElasticModel::LogLogInterpolate(G4double e1, |
---|
383 | G4double e2, |
---|
384 | G4double e, |
---|
385 | G4double xs1, |
---|
386 | G4double xs2) |
---|
387 | { |
---|
388 | G4double a = (std::log10(xs2)-std::log10(xs1)) / (std::log10(e2)-std::log10(e1)); |
---|
389 | G4double b = std::log10(xs2) - a*std::log10(e2); |
---|
390 | G4double sigma = a*std::log10(e) + b; |
---|
391 | G4double value = (std::pow(10.,sigma)); |
---|
392 | return value; |
---|
393 | } |
---|
394 | |
---|
395 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
396 | |
---|
397 | |
---|
398 | G4double G4DNAChampionElasticModel::QuadInterpolator(G4double e11, G4double e12, |
---|
399 | G4double e21, G4double e22, |
---|
400 | G4double xs11, G4double xs12, |
---|
401 | G4double xs21, G4double xs22, |
---|
402 | G4double t1, G4double t2, |
---|
403 | G4double t, G4double e) |
---|
404 | { |
---|
405 | // Log-Log |
---|
406 | /* |
---|
407 | G4double interpolatedvalue1 = LogLogInterpolate(e11, e12, e, xs11, xs12); |
---|
408 | G4double interpolatedvalue2 = LogLogInterpolate(e21, e22, e, xs21, xs22); |
---|
409 | G4double value = LogLogInterpolate(t1, t2, t, interpolatedvalue1, interpolatedvalue2); |
---|
410 | |
---|
411 | |
---|
412 | // Lin-Log |
---|
413 | G4double interpolatedvalue1 = LinLogInterpolate(e11, e12, e, xs11, xs12); |
---|
414 | G4double interpolatedvalue2 = LinLogInterpolate(e21, e22, e, xs21, xs22); |
---|
415 | G4double value = LinLogInterpolate(t1, t2, t, interpolatedvalue1, interpolatedvalue2); |
---|
416 | */ |
---|
417 | |
---|
418 | // Lin-Lin |
---|
419 | G4double interpolatedvalue1 = LinLinInterpolate(e11, e12, e, xs11, xs12); |
---|
420 | G4double interpolatedvalue2 = LinLinInterpolate(e21, e22, e, xs21, xs22); |
---|
421 | G4double value = LinLinInterpolate(t1, t2, t, interpolatedvalue1, interpolatedvalue2); |
---|
422 | |
---|
423 | return value; |
---|
424 | } |
---|
425 | |
---|
426 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
427 | |
---|
428 | G4double G4DNAChampionElasticModel::RandomizeCosTheta(G4double k) |
---|
429 | { |
---|
430 | |
---|
431 | G4double integrdiff=0; |
---|
432 | G4double uniformRand=G4UniformRand(); |
---|
433 | integrdiff = uniformRand; |
---|
434 | |
---|
435 | G4double theta=0.; |
---|
436 | G4double cosTheta=0.; |
---|
437 | theta = Theta(G4Electron::ElectronDefinition(),k/eV,integrdiff); |
---|
438 | |
---|
439 | cosTheta= std::cos(theta*pi/180); |
---|
440 | |
---|
441 | return cosTheta; |
---|
442 | } |
---|