source: trunk/source/processes/electromagnetic/lowenergy/src/G4DNABornIonisationModel.cc

Last change on this file was 1347, checked in by garnier, 15 years ago

geant4 tag 9.4

File size: 22.1 KB
RevLine 
[1192]1//
[1058]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//
[1340]26// $Id: G4DNABornIonisationModel.cc,v 1.18 2010/11/03 12:22:36 sincerti Exp $
[1347]27// GEANT4 tag $Name: geant4-09-04-ref-00 $
[1058]28//
29
30#include "G4DNABornIonisationModel.hh"
[1192]31//#include "G4DynamicMolecule.hh"
[1058]32
33//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
34
35using namespace std;
36
37//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
38
39G4DNABornIonisationModel::G4DNABornIonisationModel(const G4ParticleDefinition*,
40 const G4String& nam)
41:G4VEmModel(nam),isInitialised(false)
42{
43 verboseLevel= 0;
44 // Verbosity scale:
45 // 0 = nothing
46 // 1 = warning for energy non-conservation
47 // 2 = details of energy budget
48 // 3 = calculation of cross sections, file openings, sampling of atoms
49 // 4 = entering in methods
50
[1192]51 if( verboseLevel>0 )
52 {
53 G4cout << "Born ionisation model is constructed " << G4endl;
54 }
[1058]55}
56
57//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
58
59G4DNABornIonisationModel::~G4DNABornIonisationModel()
60{
61 // Cross section
62
63 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
64 for (pos = tableData.begin(); pos != tableData.end(); ++pos)
65 {
66 G4DNACrossSectionDataSet* table = pos->second;
67 delete table;
68 }
69
70 // Final state
71
72 eVecm.clear();
73 pVecm.clear();
74
75}
76
77//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
78
79void G4DNABornIonisationModel::Initialise(const G4ParticleDefinition* particle,
80 const G4DataVector& /*cuts*/)
81{
82
83 if (verboseLevel > 3)
84 G4cout << "Calling G4DNABornIonisationModel::Initialise()" << G4endl;
85
86 // Energy limits
87
88 G4String fileElectron("dna/sigma_ionisation_e_born");
89 G4String fileProton("dna/sigma_ionisation_p_born");
90
91 G4ParticleDefinition* electronDef = G4Electron::ElectronDefinition();
92 G4ParticleDefinition* protonDef = G4Proton::ProtonDefinition();
93
94 G4String electron;
95 G4String proton;
96
97 G4double scaleFactor = (1.e-22 / 3.343) * m*m;
98
99 char *path = getenv("G4LEDATA");
100
101 if (electronDef != 0)
102 {
103 electron = electronDef->GetParticleName();
104
105 tableFile[electron] = fileElectron;
106
[1192]107 lowEnergyLimit[electron] = 11. * eV;
108 highEnergyLimit[electron] = 1. * MeV;
[1058]109
110 // Cross section
111
112 G4DNACrossSectionDataSet* tableE = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, eV,scaleFactor );
113 tableE->LoadData(fileElectron);
114
115 tableData[electron] = tableE;
116
117 // Final state
118
119 std::ostringstream eFullFileName;
120 eFullFileName << path << "/dna/sigmadiff_ionisation_e_born.dat";
121 std::ifstream eDiffCrossSection(eFullFileName.str().c_str());
122
123 if (!eDiffCrossSection)
124 {
125 G4Exception("G4DNABornIonisationModel::ERROR OPENING electron DATA FILE");
126 }
127
128 eTdummyVec.push_back(0.);
129 while(!eDiffCrossSection.eof())
130 {
131 double tDummy;
132 double eDummy;
133 eDiffCrossSection>>tDummy>>eDummy;
134 if (tDummy != eTdummyVec.back()) eTdummyVec.push_back(tDummy);
135 for (int j=0; j<5; j++)
136 {
137 eDiffCrossSection>>eDiffCrossSectionData[j][tDummy][eDummy];
138
139 // SI - only if eof is not reached !
140 if (!eDiffCrossSection.eof()) eDiffCrossSectionData[j][tDummy][eDummy]*=scaleFactor;
141
142 eVecm[tDummy].push_back(eDummy);
143
144 }
145 }
146
147 //
148 }
149 else
150 {
151 G4Exception("G4DNABornIonisationModel::Initialise(): electron is not defined");
152 }
153
154 if (protonDef != 0)
155 {
156 proton = protonDef->GetParticleName();
157
158 tableFile[proton] = fileProton;
159
160 lowEnergyLimit[proton] = 500. * keV;
[1192]161 highEnergyLimit[proton] = 100. * MeV;
[1058]162
163 // Cross section
164
165 G4DNACrossSectionDataSet* tableP = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, eV,scaleFactor );
166 tableP->LoadData(fileProton);
167
168 tableData[proton] = tableP;
169
170 // Final state
171
172 std::ostringstream pFullFileName;
173 pFullFileName << path << "/dna/sigmadiff_ionisation_p_born.dat";
174 std::ifstream pDiffCrossSection(pFullFileName.str().c_str());
175
176 if (!pDiffCrossSection)
177 {
178 G4Exception("G4DNABornIonisationModel::ERROR OPENING proton DATA FILE");
179 }
180
181 pTdummyVec.push_back(0.);
182 while(!pDiffCrossSection.eof())
183 {
184 double tDummy;
185 double eDummy;
186 pDiffCrossSection>>tDummy>>eDummy;
187 if (tDummy != pTdummyVec.back()) pTdummyVec.push_back(tDummy);
188 for (int j=0; j<5; j++)
189 {
190 pDiffCrossSection>>pDiffCrossSectionData[j][tDummy][eDummy];
191
192 // SI - only if eof is not reached !
193 if (!pDiffCrossSection.eof()) pDiffCrossSectionData[j][tDummy][eDummy]*=scaleFactor;
194
195 pVecm[tDummy].push_back(eDummy);
196 }
197 }
198
199 }
200 else
201 {
202 G4Exception("G4DNABornIonisationModel::Initialise(): proton is not defined");
203 }
204
205 if (particle==electronDef)
206 {
207 SetLowEnergyLimit(lowEnergyLimit[electron]);
208 SetHighEnergyLimit(highEnergyLimit[electron]);
209 }
210
211 if (particle==protonDef)
212 {
213 SetLowEnergyLimit(lowEnergyLimit[proton]);
214 SetHighEnergyLimit(highEnergyLimit[proton]);
215 }
216
[1192]217 if( verboseLevel>0 )
218 {
219 G4cout << "Born ionisation model is initialized " << G4endl
220 << "Energy range: "
221 << LowEnergyLimit() / eV << " eV - "
222 << HighEnergyLimit() / keV << " keV for "
223 << particle->GetParticleName()
224 << G4endl;
225 }
226
[1058]227 //
228
229 if(!isInitialised)
230 {
231 isInitialised = true;
232
233 if(pParticleChange)
234 fParticleChangeForGamma = reinterpret_cast<G4ParticleChangeForGamma*>(pParticleChange);
235 else
236 fParticleChangeForGamma = new G4ParticleChangeForGamma();
237 }
238
239 // InitialiseElementSelectors(particle,cuts);
240
241}
242
243//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
244
[1315]245G4double G4DNABornIonisationModel::CrossSectionPerVolume(const G4Material* material,
[1058]246 const G4ParticleDefinition* particleDefinition,
247 G4double ekin,
248 G4double,
249 G4double)
250{
251 if (verboseLevel > 3)
252 G4cout << "Calling CrossSectionPerVolume() of G4DNABornIonisationModel" << G4endl;
253
254 if (
255 particleDefinition != G4Proton::ProtonDefinition()
256 &&
257 particleDefinition != G4Electron::ElectronDefinition()
258 )
259
260 return 0;
261
262 // Calculate total cross section for model
263
264 G4double lowLim = 0;
265 G4double highLim = 0;
266 G4double sigma=0;
267
[1315]268 if (material->GetName() == "G4_WATER")
[1058]269 {
270 const G4String& particleName = particleDefinition->GetParticleName();
271
272 std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
273 pos1 = lowEnergyLimit.find(particleName);
274 if (pos1 != lowEnergyLimit.end())
275 {
276 lowLim = pos1->second;
277 }
278
279 std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
280 pos2 = highEnergyLimit.find(particleName);
281 if (pos2 != highEnergyLimit.end())
282 {
283 highLim = pos2->second;
284 }
285
286 if (ekin >= lowLim && ekin < highLim)
287 {
288 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
289 pos = tableData.find(particleName);
290
291 if (pos != tableData.end())
292 {
293 G4DNACrossSectionDataSet* table = pos->second;
294 if (table != 0)
295 {
296 sigma = table->FindValue(ekin);
297 }
298 }
299 else
300 {
301 G4Exception("G4DNABornIonisationModel::CrossSectionPerVolume: attempting to calculate cross section for wrong particle");
302 }
303 }
304
305 if (verboseLevel > 3)
306 {
307 G4cout << "---> Kinetic energy(eV)=" << ekin/eV << G4endl;
308 G4cout << " - Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
[1315]309 G4cout << " - Cross section per water molecule (cm^-1)=" << sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
[1058]310 }
311
312 } // if (waterMaterial)
313
[1315]314 return sigma*material->GetAtomicNumDensityVector()[1];
[1058]315
316}
317
318//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
319
320void G4DNABornIonisationModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
321 const G4MaterialCutsCouple* /*couple*/,
322 const G4DynamicParticle* particle,
323 G4double,
324 G4double)
325{
326
327 if (verboseLevel > 3)
328 G4cout << "Calling SampleSecondaries() of G4DNABornIonisationModel" << G4endl;
329
330 G4double lowLim = 0;
331 G4double highLim = 0;
332
333 G4double k = particle->GetKineticEnergy();
334
335 const G4String& particleName = particle->GetDefinition()->GetParticleName();
336
337 std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
338 pos1 = lowEnergyLimit.find(particleName);
339
340 if (pos1 != lowEnergyLimit.end())
341 {
342 lowLim = pos1->second;
343 }
344
345 std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
346 pos2 = highEnergyLimit.find(particleName);
347
348 if (pos2 != highEnergyLimit.end())
349 {
350 highLim = pos2->second;
351 }
352
353 if (k >= lowLim && k < highLim)
354 {
355 G4ParticleMomentum primaryDirection = particle->GetMomentumDirection();
356 G4double particleMass = particle->GetDefinition()->GetPDGMass();
357 G4double totalEnergy = k + particleMass;
358 G4double pSquare = k * (totalEnergy + particleMass);
359 G4double totalMomentum = std::sqrt(pSquare);
360
361 G4int ionizationShell = RandomSelect(k,particleName);
362
363 G4double secondaryKinetic = RandomizeEjectedElectronEnergy(particle->GetDefinition(),k,ionizationShell);
364
365 G4double bindingEnergy = waterStructure.IonisationEnergy(ionizationShell);
366
367 G4double cosTheta = 0.;
368 G4double phi = 0.;
369 RandomizeEjectedElectronDirection(particle->GetDefinition(), k,secondaryKinetic, cosTheta, phi);
370
371 G4double sinTheta = std::sqrt(1.-cosTheta*cosTheta);
372 G4double dirX = sinTheta*std::cos(phi);
373 G4double dirY = sinTheta*std::sin(phi);
374 G4double dirZ = cosTheta;
375 G4ThreeVector deltaDirection(dirX,dirY,dirZ);
376 deltaDirection.rotateUz(primaryDirection);
377
[1340]378 if (particle->GetDefinition() == G4Electron::ElectronDefinition())
379 {
380 G4double deltaTotalMomentum = std::sqrt(secondaryKinetic*(secondaryKinetic + 2.*electron_mass_c2 ));
[1058]381
[1340]382 G4double finalPx = totalMomentum*primaryDirection.x() - deltaTotalMomentum*deltaDirection.x();
383 G4double finalPy = totalMomentum*primaryDirection.y() - deltaTotalMomentum*deltaDirection.y();
384 G4double finalPz = totalMomentum*primaryDirection.z() - deltaTotalMomentum*deltaDirection.z();
385 G4double finalMomentum = std::sqrt(finalPx*finalPx + finalPy*finalPy + finalPz*finalPz);
386 finalPx /= finalMomentum;
387 finalPy /= finalMomentum;
388 finalPz /= finalMomentum;
[1058]389
[1340]390 G4ThreeVector direction;
391 direction.set(finalPx,finalPy,finalPz);
[1058]392
[1340]393 fParticleChangeForGamma->ProposeMomentumDirection(direction.unit()) ;
394 }
395
396 else fParticleChangeForGamma->ProposeMomentumDirection(primaryDirection) ;
397
[1058]398 fParticleChangeForGamma->SetProposedKineticEnergy(k-bindingEnergy-secondaryKinetic);
399 fParticleChangeForGamma->ProposeLocalEnergyDeposit(bindingEnergy);
400
401 G4DynamicParticle* dp = new G4DynamicParticle (G4Electron::Electron(),deltaDirection,secondaryKinetic) ;
402 fvect->push_back(dp);
[1192]403
[1058]404 }
405
406}
407
408//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
409
410G4double G4DNABornIonisationModel::RandomizeEjectedElectronEnergy(G4ParticleDefinition* particleDefinition,
411G4double k, G4int shell)
412{
413 if (particleDefinition == G4Electron::ElectronDefinition())
414 {
415 G4double maximumEnergyTransfer=0.;
416 if ((k+waterStructure.IonisationEnergy(shell))/2. > k) maximumEnergyTransfer=k;
417 else maximumEnergyTransfer = (k+waterStructure.IonisationEnergy(shell))/2.;
418
[1192]419// SI : original method
420/*
[1058]421 G4double crossSectionMaximum = 0.;
422 for(G4double value=waterStructure.IonisationEnergy(shell); value<=maximumEnergyTransfer; value+=0.1*eV)
423 {
424 G4double differentialCrossSection = DifferentialCrossSection(particleDefinition, k/eV, value/eV, shell);
425 if(differentialCrossSection >= crossSectionMaximum) crossSectionMaximum = differentialCrossSection;
426 }
[1192]427*/
428
[1058]429
[1192]430// SI : alternative method
431
432 G4double crossSectionMaximum = 0.;
433
434 G4double minEnergy = waterStructure.IonisationEnergy(shell);
435 G4double maxEnergy = maximumEnergyTransfer;
436 G4int nEnergySteps = 50;
437
438 G4double value(minEnergy);
439 G4double stpEnergy(std::pow(maxEnergy/value, 1./static_cast<G4double>(nEnergySteps-1)));
440 G4int step(nEnergySteps);
441 while (step>0)
442 {
443 step--;
444 G4double differentialCrossSection = DifferentialCrossSection(particleDefinition, k/eV, value/eV, shell);
445 if(differentialCrossSection >= crossSectionMaximum) crossSectionMaximum = differentialCrossSection;
446 value*=stpEnergy;
447 }
448//
449
[1058]450 G4double secondaryElectronKineticEnergy=0.;
451 do
452 {
453 secondaryElectronKineticEnergy = G4UniformRand() * (maximumEnergyTransfer-waterStructure.IonisationEnergy(shell));
454 } while(G4UniformRand()*crossSectionMaximum >
455 DifferentialCrossSection(particleDefinition, k/eV,(secondaryElectronKineticEnergy+waterStructure.IonisationEnergy(shell))/eV,shell));
456
457 return secondaryElectronKineticEnergy;
458
459 }
460
461 if (particleDefinition == G4Proton::ProtonDefinition())
462 {
[1315]463 G4double maximumKineticEnergyTransfer = 4.* (electron_mass_c2 / proton_mass_c2) * k;
[1058]464
465 G4double crossSectionMaximum = 0.;
466 for (G4double value = waterStructure.IonisationEnergy(shell);
467 value<=4.*waterStructure.IonisationEnergy(shell) ;
468 value+=0.1*eV)
469 {
470 G4double differentialCrossSection = DifferentialCrossSection(particleDefinition, k/eV, value/eV, shell);
471 if (differentialCrossSection >= crossSectionMaximum) crossSectionMaximum = differentialCrossSection;
472 }
473
474 G4double secondaryElectronKineticEnergy = 0.;
475 do
476 {
477 secondaryElectronKineticEnergy = G4UniformRand() * maximumKineticEnergyTransfer;
478 } while(G4UniformRand()*crossSectionMaximum >=
479 DifferentialCrossSection(particleDefinition, k/eV,(secondaryElectronKineticEnergy+waterStructure.IonisationEnergy(shell))/eV,shell));
480
481 return secondaryElectronKineticEnergy;
482 }
483
484 return 0;
485}
486
487//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
488
489void G4DNABornIonisationModel::RandomizeEjectedElectronDirection(G4ParticleDefinition* particleDefinition,
490 G4double k,
491 G4double secKinetic,
492 G4double & cosTheta,
493 G4double & phi )
494{
495 if (particleDefinition == G4Electron::ElectronDefinition())
496 {
497 phi = twopi * G4UniformRand();
498 if (secKinetic < 50.*eV) cosTheta = (2.*G4UniformRand())-1.;
499 else if (secKinetic <= 200.*eV)
500 {
501 if (G4UniformRand() <= 0.1) cosTheta = (2.*G4UniformRand())-1.;
502 else cosTheta = G4UniformRand()*(std::sqrt(2.)/2);
503 }
504 else
505 {
506 G4double sin2O = (1.-secKinetic/k) / (1.+secKinetic/(2.*electron_mass_c2));
507 cosTheta = std::sqrt(1.-sin2O);
508 }
509 }
510
511 if (particleDefinition == G4Proton::ProtonDefinition())
512 {
513 G4double maxSecKinetic = 4.* (electron_mass_c2 / proton_mass_c2) * k;
514 phi = twopi * G4UniformRand();
[1340]515
516 // cosTheta = std::sqrt(secKinetic / maxSecKinetic);
517
518 // Restriction below 100 eV from Emfietzoglou (2000)
519
520 if (secKinetic>100*eV) cosTheta = std::sqrt(secKinetic / maxSecKinetic);
521 else cosTheta = (2.*G4UniformRand())-1.;
522
[1058]523 }
524}
525
526//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
527
528double G4DNABornIonisationModel::DifferentialCrossSection(G4ParticleDefinition * particleDefinition,
529 G4double k,
530 G4double energyTransfer,
531 G4int ionizationLevelIndex)
532{
533 G4double sigma = 0.;
534
535 if (energyTransfer >= waterStructure.IonisationEnergy(ionizationLevelIndex))
536 {
537 G4double valueT1 = 0;
538 G4double valueT2 = 0;
539 G4double valueE21 = 0;
540 G4double valueE22 = 0;
541 G4double valueE12 = 0;
542 G4double valueE11 = 0;
543
544 G4double xs11 = 0;
545 G4double xs12 = 0;
546 G4double xs21 = 0;
547 G4double xs22 = 0;
548
549 if (particleDefinition == G4Electron::ElectronDefinition())
550 {
551 // k should be in eV and energy transfer eV also
552
553 std::vector<double>::iterator t2 = std::upper_bound(eTdummyVec.begin(),eTdummyVec.end(), k);
554
555 std::vector<double>::iterator t1 = t2-1;
556
557 // SI : the following condition avoids situations where energyTransfer >last vector element
[1192]558 if (energyTransfer <= eVecm[(*t1)].back() && energyTransfer <= eVecm[(*t2)].back() )
[1058]559 {
560 std::vector<double>::iterator e12 = std::upper_bound(eVecm[(*t1)].begin(),eVecm[(*t1)].end(), energyTransfer);
561 std::vector<double>::iterator e11 = e12-1;
562
563 std::vector<double>::iterator e22 = std::upper_bound(eVecm[(*t2)].begin(),eVecm[(*t2)].end(), energyTransfer);
564 std::vector<double>::iterator e21 = e22-1;
565
566 valueT1 =*t1;
567 valueT2 =*t2;
568 valueE21 =*e21;
569 valueE22 =*e22;
570 valueE12 =*e12;
571 valueE11 =*e11;
572
573 xs11 = eDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE11];
574 xs12 = eDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE12];
575 xs21 = eDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE21];
576 xs22 = eDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE22];
577 }
578
579 }
580
581 if (particleDefinition == G4Proton::ProtonDefinition())
582 {
583 // k should be in eV and energy transfer eV also
584 std::vector<double>::iterator t2 = std::upper_bound(pTdummyVec.begin(),pTdummyVec.end(), k);
585 std::vector<double>::iterator t1 = t2-1;
586
587 std::vector<double>::iterator e12 = std::upper_bound(pVecm[(*t1)].begin(),pVecm[(*t1)].end(), energyTransfer);
588 std::vector<double>::iterator e11 = e12-1;
589
590 std::vector<double>::iterator e22 = std::upper_bound(pVecm[(*t2)].begin(),pVecm[(*t2)].end(), energyTransfer);
591 std::vector<double>::iterator e21 = e22-1;
592
593 valueT1 =*t1;
594 valueT2 =*t2;
595 valueE21 =*e21;
596 valueE22 =*e22;
597 valueE12 =*e12;
598 valueE11 =*e11;
599
600 xs11 = pDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE11];
601 xs12 = pDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE12];
602 xs21 = pDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE21];
603 xs22 = pDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE22];
604
605 }
606
607 G4double xsProduct = xs11 * xs12 * xs21 * xs22;
608 if (xsProduct != 0.)
609 {
610 sigma = QuadInterpolator( valueE11, valueE12,
611 valueE21, valueE22,
612 xs11, xs12,
613 xs21, xs22,
614 valueT1, valueT2,
615 k, energyTransfer);
616 }
617
618 }
619
620 return sigma;
621}
622
623//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
624
625G4double G4DNABornIonisationModel::LogLogInterpolate(G4double e1,
626 G4double e2,
627 G4double e,
628 G4double xs1,
629 G4double xs2)
630{
631 G4double a = (std::log10(xs2)-std::log10(xs1)) / (std::log10(e2)-std::log10(e1));
632 G4double b = std::log10(xs2) - a*std::log10(e2);
633 G4double sigma = a*std::log10(e) + b;
634 G4double value = (std::pow(10.,sigma));
635 return value;
636}
637
638//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
639
640G4double G4DNABornIonisationModel::QuadInterpolator(G4double e11, G4double e12,
641 G4double e21, G4double e22,
642 G4double xs11, G4double xs12,
643 G4double xs21, G4double xs22,
644 G4double t1, G4double t2,
645 G4double t, G4double e)
646{
647 G4double interpolatedvalue1 = LogLogInterpolate(e11, e12, e, xs11, xs12);
648 G4double interpolatedvalue2 = LogLogInterpolate(e21, e22, e, xs21, xs22);
649 G4double value = LogLogInterpolate(t1, t2, t, interpolatedvalue1, interpolatedvalue2);
650 return value;
651}
652
653//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
654
655G4int G4DNABornIonisationModel::RandomSelect(G4double k, const G4String& particle )
656{
657 G4int level = 0;
658
659 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
660 pos = tableData.find(particle);
661
662 if (pos != tableData.end())
663 {
664 G4DNACrossSectionDataSet* table = pos->second;
665
666 if (table != 0)
667 {
668 G4double* valuesBuffer = new G4double[table->NumberOfComponents()];
669 const size_t n(table->NumberOfComponents());
670 size_t i(n);
671 G4double value = 0.;
672
673 while (i>0)
674 {
675 i--;
676 valuesBuffer[i] = table->GetComponent(i)->FindValue(k);
677 value += valuesBuffer[i];
678 }
679
680 value *= G4UniformRand();
681
682 i = n;
683
684 while (i > 0)
685 {
686 i--;
687
688 if (valuesBuffer[i] > value)
689 {
690 delete[] valuesBuffer;
691 return i;
692 }
693 value -= valuesBuffer[i];
694 }
695
696 if (valuesBuffer) delete[] valuesBuffer;
697
698 }
699 }
700 else
701 {
702 G4Exception("G4DNABornIonisationModel::RandomSelect attempting to calculate cross section for wrong particle");
703 }
704
705 return level;
706}
707
Note: See TracBrowser for help on using the repository browser.