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

Last change on this file since 1316 was 1315, checked in by garnier, 15 years ago

update geant4-09-04-beta-cand-01 interfaces-V09-03-09 vis-V09-03-08

File size: 21.7 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//
[1315]26// $Id: G4DNABornIonisationModel.cc,v 1.16 2010/03/26 18:10:43 sincerti Exp $
27// GEANT4 tag $Name: geant4-09-04-beta-cand-01 $
[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
378 G4double deltaTotalMomentum = std::sqrt(secondaryKinetic*(secondaryKinetic + 2.*electron_mass_c2 ));
379
380 G4double finalPx = totalMomentum*primaryDirection.x() - deltaTotalMomentum*deltaDirection.x();
381 G4double finalPy = totalMomentum*primaryDirection.y() - deltaTotalMomentum*deltaDirection.y();
382 G4double finalPz = totalMomentum*primaryDirection.z() - deltaTotalMomentum*deltaDirection.z();
383 G4double finalMomentum = std::sqrt(finalPx*finalPx + finalPy*finalPy + finalPz*finalPz);
384 finalPx /= finalMomentum;
385 finalPy /= finalMomentum;
386 finalPz /= finalMomentum;
387
388 G4ThreeVector direction;
389 direction.set(finalPx,finalPy,finalPz);
390
391 fParticleChangeForGamma->ProposeMomentumDirection(direction.unit()) ;
392 fParticleChangeForGamma->SetProposedKineticEnergy(k-bindingEnergy-secondaryKinetic);
393 fParticleChangeForGamma->ProposeLocalEnergyDeposit(bindingEnergy);
394
395 G4DynamicParticle* dp = new G4DynamicParticle (G4Electron::Electron(),deltaDirection,secondaryKinetic) ;
396 fvect->push_back(dp);
[1192]397
[1058]398 }
399
400}
401
402//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
403
404G4double G4DNABornIonisationModel::RandomizeEjectedElectronEnergy(G4ParticleDefinition* particleDefinition,
405G4double k, G4int shell)
406{
407 if (particleDefinition == G4Electron::ElectronDefinition())
408 {
409 G4double maximumEnergyTransfer=0.;
410 if ((k+waterStructure.IonisationEnergy(shell))/2. > k) maximumEnergyTransfer=k;
411 else maximumEnergyTransfer = (k+waterStructure.IonisationEnergy(shell))/2.;
412
[1192]413// SI : original method
414/*
[1058]415 G4double crossSectionMaximum = 0.;
416 for(G4double value=waterStructure.IonisationEnergy(shell); value<=maximumEnergyTransfer; value+=0.1*eV)
417 {
418 G4double differentialCrossSection = DifferentialCrossSection(particleDefinition, k/eV, value/eV, shell);
419 if(differentialCrossSection >= crossSectionMaximum) crossSectionMaximum = differentialCrossSection;
420 }
[1192]421*/
422
[1058]423
[1192]424// SI : alternative method
425
426 G4double crossSectionMaximum = 0.;
427
428 G4double minEnergy = waterStructure.IonisationEnergy(shell);
429 G4double maxEnergy = maximumEnergyTransfer;
430 G4int nEnergySteps = 50;
431
432 G4double value(minEnergy);
433 G4double stpEnergy(std::pow(maxEnergy/value, 1./static_cast<G4double>(nEnergySteps-1)));
434 G4int step(nEnergySteps);
435 while (step>0)
436 {
437 step--;
438 G4double differentialCrossSection = DifferentialCrossSection(particleDefinition, k/eV, value/eV, shell);
439 if(differentialCrossSection >= crossSectionMaximum) crossSectionMaximum = differentialCrossSection;
440 value*=stpEnergy;
441 }
442//
443
[1058]444 G4double secondaryElectronKineticEnergy=0.;
445 do
446 {
447 secondaryElectronKineticEnergy = G4UniformRand() * (maximumEnergyTransfer-waterStructure.IonisationEnergy(shell));
448 } while(G4UniformRand()*crossSectionMaximum >
449 DifferentialCrossSection(particleDefinition, k/eV,(secondaryElectronKineticEnergy+waterStructure.IonisationEnergy(shell))/eV,shell));
450
451 return secondaryElectronKineticEnergy;
452
453 }
454
455 if (particleDefinition == G4Proton::ProtonDefinition())
456 {
[1315]457 G4double maximumKineticEnergyTransfer = 4.* (electron_mass_c2 / proton_mass_c2) * k;
[1058]458
459 G4double crossSectionMaximum = 0.;
460 for (G4double value = waterStructure.IonisationEnergy(shell);
461 value<=4.*waterStructure.IonisationEnergy(shell) ;
462 value+=0.1*eV)
463 {
464 G4double differentialCrossSection = DifferentialCrossSection(particleDefinition, k/eV, value/eV, shell);
465 if (differentialCrossSection >= crossSectionMaximum) crossSectionMaximum = differentialCrossSection;
466 }
467
468 G4double secondaryElectronKineticEnergy = 0.;
469 do
470 {
471 secondaryElectronKineticEnergy = G4UniformRand() * maximumKineticEnergyTransfer;
472 } while(G4UniformRand()*crossSectionMaximum >=
473 DifferentialCrossSection(particleDefinition, k/eV,(secondaryElectronKineticEnergy+waterStructure.IonisationEnergy(shell))/eV,shell));
474
475 return secondaryElectronKineticEnergy;
476 }
477
478 return 0;
479}
480
481//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
482
483void G4DNABornIonisationModel::RandomizeEjectedElectronDirection(G4ParticleDefinition* particleDefinition,
484 G4double k,
485 G4double secKinetic,
486 G4double & cosTheta,
487 G4double & phi )
488{
489 if (particleDefinition == G4Electron::ElectronDefinition())
490 {
491 phi = twopi * G4UniformRand();
492 if (secKinetic < 50.*eV) cosTheta = (2.*G4UniformRand())-1.;
493 else if (secKinetic <= 200.*eV)
494 {
495 if (G4UniformRand() <= 0.1) cosTheta = (2.*G4UniformRand())-1.;
496 else cosTheta = G4UniformRand()*(std::sqrt(2.)/2);
497 }
498 else
499 {
500 G4double sin2O = (1.-secKinetic/k) / (1.+secKinetic/(2.*electron_mass_c2));
501 cosTheta = std::sqrt(1.-sin2O);
502 }
503 }
504
505 if (particleDefinition == G4Proton::ProtonDefinition())
506 {
507 G4double maxSecKinetic = 4.* (electron_mass_c2 / proton_mass_c2) * k;
508 phi = twopi * G4UniformRand();
509 cosTheta = std::sqrt(secKinetic / maxSecKinetic);
510 }
511}
512
513//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
514
515double G4DNABornIonisationModel::DifferentialCrossSection(G4ParticleDefinition * particleDefinition,
516 G4double k,
517 G4double energyTransfer,
518 G4int ionizationLevelIndex)
519{
520 G4double sigma = 0.;
521
522 if (energyTransfer >= waterStructure.IonisationEnergy(ionizationLevelIndex))
523 {
524 G4double valueT1 = 0;
525 G4double valueT2 = 0;
526 G4double valueE21 = 0;
527 G4double valueE22 = 0;
528 G4double valueE12 = 0;
529 G4double valueE11 = 0;
530
531 G4double xs11 = 0;
532 G4double xs12 = 0;
533 G4double xs21 = 0;
534 G4double xs22 = 0;
535
536 if (particleDefinition == G4Electron::ElectronDefinition())
537 {
538 // k should be in eV and energy transfer eV also
539
540 std::vector<double>::iterator t2 = std::upper_bound(eTdummyVec.begin(),eTdummyVec.end(), k);
541
542 std::vector<double>::iterator t1 = t2-1;
543
544 // SI : the following condition avoids situations where energyTransfer >last vector element
[1192]545 if (energyTransfer <= eVecm[(*t1)].back() && energyTransfer <= eVecm[(*t2)].back() )
[1058]546 {
547 std::vector<double>::iterator e12 = std::upper_bound(eVecm[(*t1)].begin(),eVecm[(*t1)].end(), energyTransfer);
548 std::vector<double>::iterator e11 = e12-1;
549
550 std::vector<double>::iterator e22 = std::upper_bound(eVecm[(*t2)].begin(),eVecm[(*t2)].end(), energyTransfer);
551 std::vector<double>::iterator e21 = e22-1;
552
553 valueT1 =*t1;
554 valueT2 =*t2;
555 valueE21 =*e21;
556 valueE22 =*e22;
557 valueE12 =*e12;
558 valueE11 =*e11;
559
560 xs11 = eDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE11];
561 xs12 = eDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE12];
562 xs21 = eDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE21];
563 xs22 = eDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE22];
564 }
565
566 }
567
568 if (particleDefinition == G4Proton::ProtonDefinition())
569 {
570 // k should be in eV and energy transfer eV also
571 std::vector<double>::iterator t2 = std::upper_bound(pTdummyVec.begin(),pTdummyVec.end(), k);
572 std::vector<double>::iterator t1 = t2-1;
573
574 std::vector<double>::iterator e12 = std::upper_bound(pVecm[(*t1)].begin(),pVecm[(*t1)].end(), energyTransfer);
575 std::vector<double>::iterator e11 = e12-1;
576
577 std::vector<double>::iterator e22 = std::upper_bound(pVecm[(*t2)].begin(),pVecm[(*t2)].end(), energyTransfer);
578 std::vector<double>::iterator e21 = e22-1;
579
580 valueT1 =*t1;
581 valueT2 =*t2;
582 valueE21 =*e21;
583 valueE22 =*e22;
584 valueE12 =*e12;
585 valueE11 =*e11;
586
587 xs11 = pDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE11];
588 xs12 = pDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE12];
589 xs21 = pDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE21];
590 xs22 = pDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE22];
591
592 }
593
594 G4double xsProduct = xs11 * xs12 * xs21 * xs22;
595 if (xsProduct != 0.)
596 {
597 sigma = QuadInterpolator( valueE11, valueE12,
598 valueE21, valueE22,
599 xs11, xs12,
600 xs21, xs22,
601 valueT1, valueT2,
602 k, energyTransfer);
603 }
604
605 }
606
607 return sigma;
608}
609
610//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
611
612G4double G4DNABornIonisationModel::LogLogInterpolate(G4double e1,
613 G4double e2,
614 G4double e,
615 G4double xs1,
616 G4double xs2)
617{
618 G4double a = (std::log10(xs2)-std::log10(xs1)) / (std::log10(e2)-std::log10(e1));
619 G4double b = std::log10(xs2) - a*std::log10(e2);
620 G4double sigma = a*std::log10(e) + b;
621 G4double value = (std::pow(10.,sigma));
622 return value;
623}
624
625//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
626
627G4double G4DNABornIonisationModel::QuadInterpolator(G4double e11, G4double e12,
628 G4double e21, G4double e22,
629 G4double xs11, G4double xs12,
630 G4double xs21, G4double xs22,
631 G4double t1, G4double t2,
632 G4double t, G4double e)
633{
634 G4double interpolatedvalue1 = LogLogInterpolate(e11, e12, e, xs11, xs12);
635 G4double interpolatedvalue2 = LogLogInterpolate(e21, e22, e, xs21, xs22);
636 G4double value = LogLogInterpolate(t1, t2, t, interpolatedvalue1, interpolatedvalue2);
637 return value;
638}
639
640//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
641
642G4int G4DNABornIonisationModel::RandomSelect(G4double k, const G4String& particle )
643{
644 G4int level = 0;
645
646 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
647 pos = tableData.find(particle);
648
649 if (pos != tableData.end())
650 {
651 G4DNACrossSectionDataSet* table = pos->second;
652
653 if (table != 0)
654 {
655 G4double* valuesBuffer = new G4double[table->NumberOfComponents()];
656 const size_t n(table->NumberOfComponents());
657 size_t i(n);
658 G4double value = 0.;
659
660 while (i>0)
661 {
662 i--;
663 valuesBuffer[i] = table->GetComponent(i)->FindValue(k);
664 value += valuesBuffer[i];
665 }
666
667 value *= G4UniformRand();
668
669 i = n;
670
671 while (i > 0)
672 {
673 i--;
674
675 if (valuesBuffer[i] > value)
676 {
677 delete[] valuesBuffer;
678 return i;
679 }
680 value -= valuesBuffer[i];
681 }
682
683 if (valuesBuffer) delete[] valuesBuffer;
684
685 }
686 }
687 else
688 {
689 G4Exception("G4DNABornIonisationModel::RandomSelect attempting to calculate cross section for wrong particle");
690 }
691
692 return level;
693}
694
Note: See TracBrowser for help on using the repository browser.