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

Last change on this file since 1250 was 1228, checked in by garnier, 15 years ago

update geant4.9.3 tag

File size: 23.4 KB
Line 
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: G4DNABornIonisationModel.cc,v 1.14 2009/11/12 03:08:58 sincerti Exp $
27// GEANT4 tag $Name: geant4-09-03 $
28//
29
30#include "G4DNABornIonisationModel.hh"
31//#include "G4DynamicMolecule.hh"
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 
51  if( verboseLevel>0 ) 
52  { 
53    G4cout << "Born ionisation model is constructed " << G4endl;
54  }
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
107    lowEnergyLimit[electron] = 11. * eV; 
108    highEnergyLimit[electron] = 1. * MeV;
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;
161    highEnergyLimit[proton] = 100. * MeV;
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
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 
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  // Test if water material
242
243  flagMaterialIsWater= false;
244  densityWater = 0;
245
246  const G4ProductionCutsTable* theCoupleTable = G4ProductionCutsTable::GetProductionCutsTable();
247
248  if(theCoupleTable) 
249  {
250    G4int numOfCouples = theCoupleTable->GetTableSize();
251 
252    if(numOfCouples>0) 
253    {
254          for (G4int i=0; i<numOfCouples; i++) 
255          {
256            const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(i);
257            const G4Material* material = couple->GetMaterial();
258
259            if (material->GetName() == "G4_WATER") 
260            {
261              G4double density = material->GetAtomicNumDensityVector()[1];
262              flagMaterialIsWater = true; 
263              densityWater = density; 
264             
265              if (verboseLevel > 3) 
266              G4cout << "****** Water material is found with density(cm^-3)=" << density/(cm*cm*cm) << G4endl;
267            }
268 
269          }
270
271    } // if(numOfCouples>0)
272
273  } // if (theCoupleTable)
274
275}
276
277//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
278
279G4double G4DNABornIonisationModel::CrossSectionPerVolume(const G4Material*,
280                                           const G4ParticleDefinition* particleDefinition,
281                                           G4double ekin,
282                                           G4double,
283                                           G4double)
284{
285  if (verboseLevel > 3)
286    G4cout << "Calling CrossSectionPerVolume() of G4DNABornIonisationModel" << G4endl;
287
288  if (
289      particleDefinition != G4Proton::ProtonDefinition()
290      &&
291      particleDefinition != G4Electron::ElectronDefinition()
292     )
293           
294    return 0;
295 
296  // Calculate total cross section for model
297
298  G4double lowLim = 0;
299  G4double highLim = 0;
300  G4double sigma=0;
301
302  if (flagMaterialIsWater)
303  {
304    const G4String& particleName = particleDefinition->GetParticleName();
305 
306    std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
307    pos1 = lowEnergyLimit.find(particleName);
308    if (pos1 != lowEnergyLimit.end())
309    {
310      lowLim = pos1->second;
311    }
312 
313    std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
314    pos2 = highEnergyLimit.find(particleName);
315    if (pos2 != highEnergyLimit.end())
316    {
317      highLim = pos2->second;
318    }
319
320    if (ekin >= lowLim && ekin < highLim)
321    {
322      std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
323      pos = tableData.find(particleName);
324       
325      if (pos != tableData.end())
326      {
327        G4DNACrossSectionDataSet* table = pos->second;
328        if (table != 0)
329        {
330          sigma = table->FindValue(ekin);
331        }
332      }
333      else
334      {
335        G4Exception("G4DNABornIonisationModel::CrossSectionPerVolume: attempting to calculate cross section for wrong particle");
336      }
337    }
338
339    if (verboseLevel > 3)
340    {
341      G4cout << "---> Kinetic energy(eV)=" << ekin/eV << G4endl;
342      G4cout << " - Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
343      G4cout << " - Cross section per water molecule (cm^-1)=" << sigma*densityWater/(1./cm) << G4endl;
344    } 
345 
346  } // if (waterMaterial)
347 
348 return sigma*densityWater;               
349
350}
351
352//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
353
354void G4DNABornIonisationModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
355                                              const G4MaterialCutsCouple* /*couple*/,
356                                              const G4DynamicParticle* particle,
357                                              G4double,
358                                              G4double)
359{
360
361  if (verboseLevel > 3)
362    G4cout << "Calling SampleSecondaries() of G4DNABornIonisationModel" << G4endl;
363
364  G4double lowLim = 0;
365  G4double highLim = 0;
366
367  G4double k = particle->GetKineticEnergy();
368
369  const G4String& particleName = particle->GetDefinition()->GetParticleName();
370
371  std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
372  pos1 = lowEnergyLimit.find(particleName);
373
374  if (pos1 != lowEnergyLimit.end())
375  {
376    lowLim = pos1->second;
377  }
378
379  std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
380  pos2 = highEnergyLimit.find(particleName);
381
382  if (pos2 != highEnergyLimit.end())
383  {
384    highLim = pos2->second;
385  }
386
387  if (k >= lowLim && k < highLim)
388  {
389    G4ParticleMomentum primaryDirection = particle->GetMomentumDirection();
390    G4double particleMass = particle->GetDefinition()->GetPDGMass();
391    G4double totalEnergy = k + particleMass;
392    G4double pSquare = k * (totalEnergy + particleMass);
393    G4double totalMomentum = std::sqrt(pSquare);
394
395    G4int ionizationShell = RandomSelect(k,particleName);
396 
397    G4double secondaryKinetic = RandomizeEjectedElectronEnergy(particle->GetDefinition(),k,ionizationShell);
398 
399    G4double bindingEnergy = waterStructure.IonisationEnergy(ionizationShell);
400
401    G4double cosTheta = 0.;
402    G4double phi = 0.; 
403    RandomizeEjectedElectronDirection(particle->GetDefinition(), k,secondaryKinetic, cosTheta, phi);
404
405    G4double sinTheta = std::sqrt(1.-cosTheta*cosTheta);
406    G4double dirX = sinTheta*std::cos(phi);
407    G4double dirY = sinTheta*std::sin(phi);
408    G4double dirZ = cosTheta;
409    G4ThreeVector deltaDirection(dirX,dirY,dirZ);
410    deltaDirection.rotateUz(primaryDirection);
411
412    G4double deltaTotalMomentum = std::sqrt(secondaryKinetic*(secondaryKinetic + 2.*electron_mass_c2 ));
413
414    G4double finalPx = totalMomentum*primaryDirection.x() - deltaTotalMomentum*deltaDirection.x();
415    G4double finalPy = totalMomentum*primaryDirection.y() - deltaTotalMomentum*deltaDirection.y();
416    G4double finalPz = totalMomentum*primaryDirection.z() - deltaTotalMomentum*deltaDirection.z();
417    G4double finalMomentum = std::sqrt(finalPx*finalPx + finalPy*finalPy + finalPz*finalPz);
418    finalPx /= finalMomentum;
419    finalPy /= finalMomentum;
420    finalPz /= finalMomentum;
421   
422    G4ThreeVector direction;
423    direction.set(finalPx,finalPy,finalPz);
424   
425    fParticleChangeForGamma->ProposeMomentumDirection(direction.unit()) ;
426    fParticleChangeForGamma->SetProposedKineticEnergy(k-bindingEnergy-secondaryKinetic);
427    fParticleChangeForGamma->ProposeLocalEnergyDeposit(bindingEnergy);
428
429    G4DynamicParticle* dp = new G4DynamicParticle (G4Electron::Electron(),deltaDirection,secondaryKinetic) ;
430    fvect->push_back(dp);
431    /*
432    // creating neutral water molechule...
433
434    G4DNAGenericMoleculeManager *instance;
435    instance = G4DNAGenericMoleculeManager::Instance();
436    G4ParticleDefinition* waterDef = NULL;
437    G4Molecule* water = instance->GetMolecule("H2O");
438    waterDef = (G4ParticleDefinition*)water;
439
440    direction.set(0.,0.,0.);
441
442    //G4DynamicParticle* dynamicWater = new G4DynamicParticle(waterDef, direction, bindingEnergy);
443        G4DynamicMolecule* dynamicWater = new G4DynamicMolecule(water, direction, bindingEnergy);
444
445
446        //dynamicWater->RemoveElectron(ionizationShell, 1);
447
448    G4DynamicMolecule* dynamicWater2 = new G4DynamicMolecule(water, direction, bindingEnergy);
449    G4DynamicMolecule* dynamicWater3 = new G4DynamicMolecule(water, direction, bindingEnergy);
450
451    fvect->push_back(dynamicWater);
452    fvect->push_back(dynamicWater2);
453    fvect->push_back(dynamicWater3);
454    */
455  }
456
457}
458
459//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
460
461G4double G4DNABornIonisationModel::RandomizeEjectedElectronEnergy(G4ParticleDefinition* particleDefinition, 
462G4double k, G4int shell)
463{
464  if (particleDefinition == G4Electron::ElectronDefinition()) 
465  {
466    G4double maximumEnergyTransfer=0.;
467    if ((k+waterStructure.IonisationEnergy(shell))/2. > k) maximumEnergyTransfer=k;
468    else maximumEnergyTransfer = (k+waterStructure.IonisationEnergy(shell))/2.;
469   
470// SI : original method
471/*
472    G4double crossSectionMaximum = 0.;
473    for(G4double value=waterStructure.IonisationEnergy(shell); value<=maximumEnergyTransfer; value+=0.1*eV)
474    {
475      G4double differentialCrossSection = DifferentialCrossSection(particleDefinition, k/eV, value/eV, shell);
476      if(differentialCrossSection >= crossSectionMaximum) crossSectionMaximum = differentialCrossSection;
477    }
478*/
479
480 
481// SI : alternative method
482
483    G4double crossSectionMaximum = 0.;
484
485    G4double minEnergy = waterStructure.IonisationEnergy(shell);
486    G4double maxEnergy = maximumEnergyTransfer;
487    G4int nEnergySteps = 50;
488
489    G4double value(minEnergy);
490    G4double stpEnergy(std::pow(maxEnergy/value, 1./static_cast<G4double>(nEnergySteps-1)));
491    G4int step(nEnergySteps);
492    while (step>0)
493    {
494      step--;
495      G4double differentialCrossSection = DifferentialCrossSection(particleDefinition, k/eV, value/eV, shell);
496      if(differentialCrossSection >= crossSectionMaximum) crossSectionMaximum = differentialCrossSection;
497      value*=stpEnergy;
498    }
499//
500 
501    G4double secondaryElectronKineticEnergy=0.;
502    do 
503    {
504      secondaryElectronKineticEnergy = G4UniformRand() * (maximumEnergyTransfer-waterStructure.IonisationEnergy(shell));
505    } while(G4UniformRand()*crossSectionMaximum >
506      DifferentialCrossSection(particleDefinition, k/eV,(secondaryElectronKineticEnergy+waterStructure.IonisationEnergy(shell))/eV,shell));
507
508    return secondaryElectronKineticEnergy;
509 
510  }
511 
512  if (particleDefinition == G4Proton::ProtonDefinition()) 
513  {
514    G4double maximumKineticEnergyTransfer = 4.* (electron_mass_c2 / proton_mass_c2) * k - (waterStructure.IonisationEnergy(shell));
515
516    G4double crossSectionMaximum = 0.;
517    for (G4double value = waterStructure.IonisationEnergy(shell); 
518         value<=4.*waterStructure.IonisationEnergy(shell) ; 
519         value+=0.1*eV)
520    {
521      G4double differentialCrossSection = DifferentialCrossSection(particleDefinition, k/eV, value/eV, shell);
522      if (differentialCrossSection >= crossSectionMaximum) crossSectionMaximum = differentialCrossSection;
523    }
524
525    G4double secondaryElectronKineticEnergy = 0.;
526    do
527    {
528      secondaryElectronKineticEnergy = G4UniformRand() * maximumKineticEnergyTransfer;
529    } while(G4UniformRand()*crossSectionMaximum >= 
530              DifferentialCrossSection(particleDefinition, k/eV,(secondaryElectronKineticEnergy+waterStructure.IonisationEnergy(shell))/eV,shell));
531
532    return secondaryElectronKineticEnergy;
533  }
534
535  return 0;
536}
537
538//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
539
540void G4DNABornIonisationModel::RandomizeEjectedElectronDirection(G4ParticleDefinition* particleDefinition, 
541                                                                   G4double k, 
542                                                                   G4double secKinetic, 
543                                                                   G4double & cosTheta, 
544                                                                   G4double & phi )
545{
546  if (particleDefinition == G4Electron::ElectronDefinition()) 
547  {
548    phi = twopi * G4UniformRand();
549    if (secKinetic < 50.*eV) cosTheta = (2.*G4UniformRand())-1.;
550    else if (secKinetic <= 200.*eV)     
551    {
552      if (G4UniformRand() <= 0.1) cosTheta = (2.*G4UniformRand())-1.;
553      else cosTheta = G4UniformRand()*(std::sqrt(2.)/2);
554    }
555    else       
556    {
557      G4double sin2O = (1.-secKinetic/k) / (1.+secKinetic/(2.*electron_mass_c2));
558      cosTheta = std::sqrt(1.-sin2O);
559    }
560  }
561 
562  if (particleDefinition == G4Proton::ProtonDefinition()) 
563  {
564    G4double maxSecKinetic = 4.* (electron_mass_c2 / proton_mass_c2) * k;
565    phi = twopi * G4UniformRand();
566    cosTheta = std::sqrt(secKinetic / maxSecKinetic);
567  }                     
568}
569
570//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
571
572double G4DNABornIonisationModel::DifferentialCrossSection(G4ParticleDefinition * particleDefinition, 
573                                                            G4double k, 
574                                                            G4double energyTransfer, 
575                                                            G4int ionizationLevelIndex)
576{
577  G4double sigma = 0.;
578
579  if (energyTransfer >= waterStructure.IonisationEnergy(ionizationLevelIndex))
580  {
581    G4double valueT1 = 0;
582    G4double valueT2 = 0;
583    G4double valueE21 = 0;
584    G4double valueE22 = 0;
585    G4double valueE12 = 0;
586    G4double valueE11 = 0;
587
588    G4double xs11 = 0;   
589    G4double xs12 = 0; 
590    G4double xs21 = 0; 
591    G4double xs22 = 0; 
592
593    if (particleDefinition == G4Electron::ElectronDefinition()) 
594    {
595      // k should be in eV and energy transfer eV also
596
597      std::vector<double>::iterator t2 = std::upper_bound(eTdummyVec.begin(),eTdummyVec.end(), k);
598
599      std::vector<double>::iterator t1 = t2-1;
600
601      // SI : the following condition avoids situations where energyTransfer >last vector element
602      if (energyTransfer <= eVecm[(*t1)].back() && energyTransfer <= eVecm[(*t2)].back() )
603      {
604        std::vector<double>::iterator e12 = std::upper_bound(eVecm[(*t1)].begin(),eVecm[(*t1)].end(), energyTransfer);
605        std::vector<double>::iterator e11 = e12-1;
606
607        std::vector<double>::iterator e22 = std::upper_bound(eVecm[(*t2)].begin(),eVecm[(*t2)].end(), energyTransfer);
608        std::vector<double>::iterator e21 = e22-1;
609
610        valueT1  =*t1;
611        valueT2  =*t2;
612        valueE21 =*e21;
613        valueE22 =*e22;
614        valueE12 =*e12;
615        valueE11 =*e11;
616
617        xs11 = eDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE11];
618        xs12 = eDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE12];
619        xs21 = eDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE21];
620        xs22 = eDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE22];
621      }
622
623    }
624 
625   if (particleDefinition == G4Proton::ProtonDefinition()) 
626   {
627      // k should be in eV and energy transfer eV also
628      std::vector<double>::iterator t2 = std::upper_bound(pTdummyVec.begin(),pTdummyVec.end(), k);
629      std::vector<double>::iterator t1 = t2-1;
630     
631        std::vector<double>::iterator e12 = std::upper_bound(pVecm[(*t1)].begin(),pVecm[(*t1)].end(), energyTransfer);
632        std::vector<double>::iterator e11 = e12-1;
633
634        std::vector<double>::iterator e22 = std::upper_bound(pVecm[(*t2)].begin(),pVecm[(*t2)].end(), energyTransfer);
635        std::vector<double>::iterator e21 = e22-1;
636 
637        valueT1  =*t1;
638        valueT2  =*t2;
639        valueE21 =*e21;
640        valueE22 =*e22;
641        valueE12 =*e12;
642        valueE11 =*e11;
643
644        xs11 = pDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE11];
645        xs12 = pDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE12];
646        xs21 = pDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE21];
647        xs22 = pDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE22];
648
649   }
650 
651   G4double xsProduct = xs11 * xs12 * xs21 * xs22;
652   if (xsProduct != 0.)
653   {
654     sigma = QuadInterpolator(     valueE11, valueE12, 
655                                   valueE21, valueE22, 
656                                   xs11, xs12, 
657                                   xs21, xs22, 
658                                   valueT1, valueT2, 
659                                   k, energyTransfer);
660   }
661 
662 }
663 
664  return sigma;
665}
666
667//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
668
669G4double G4DNABornIonisationModel::LogLogInterpolate(G4double e1, 
670                                                       G4double e2, 
671                                                       G4double e, 
672                                                       G4double xs1, 
673                                                       G4double xs2)
674{
675  G4double a = (std::log10(xs2)-std::log10(xs1)) / (std::log10(e2)-std::log10(e1));
676  G4double b = std::log10(xs2) - a*std::log10(e2);
677  G4double sigma = a*std::log10(e) + b;
678  G4double value = (std::pow(10.,sigma));
679  return value;
680}
681
682//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
683
684G4double G4DNABornIonisationModel::QuadInterpolator(G4double e11, G4double e12, 
685                                                      G4double e21, G4double e22, 
686                                                      G4double xs11, G4double xs12, 
687                                                      G4double xs21, G4double xs22, 
688                                                      G4double t1, G4double t2, 
689                                                      G4double t, G4double e)
690{
691  G4double interpolatedvalue1 = LogLogInterpolate(e11, e12, e, xs11, xs12);
692  G4double interpolatedvalue2 = LogLogInterpolate(e21, e22, e, xs21, xs22);
693  G4double value = LogLogInterpolate(t1, t2, t, interpolatedvalue1, interpolatedvalue2);
694  return value;
695}
696
697//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
698
699G4int G4DNABornIonisationModel::RandomSelect(G4double k, const G4String& particle )
700{   
701  G4int level = 0;
702
703  std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
704  pos = tableData.find(particle);
705
706  if (pos != tableData.end())
707  {
708    G4DNACrossSectionDataSet* table = pos->second;
709
710    if (table != 0)
711    {
712      G4double* valuesBuffer = new G4double[table->NumberOfComponents()];
713      const size_t n(table->NumberOfComponents());
714      size_t i(n);
715      G4double value = 0.;
716           
717      while (i>0)
718      { 
719        i--;
720        valuesBuffer[i] = table->GetComponent(i)->FindValue(k);
721        value += valuesBuffer[i];
722      }
723           
724      value *= G4UniformRand();
725   
726      i = n;
727           
728      while (i > 0)
729      {
730        i--;
731
732        if (valuesBuffer[i] > value)
733        {
734          delete[] valuesBuffer;
735          return i;
736        }
737        value -= valuesBuffer[i];
738      }
739           
740      if (valuesBuffer) delete[] valuesBuffer;
741   
742    }
743  }
744  else
745  {
746    G4Exception("G4DNABornIonisationModel::RandomSelect attempting to calculate cross section for wrong particle");
747  }
748     
749  return level;
750}
751
Note: See TracBrowser for help on using the repository browser.