source: trunk/source/processes/electromagnetic/lowenergy/src/G4FinalStateIonisationBorn.cc @ 924

Last change on this file since 924 was 819, checked in by garnier, 16 years ago

import all except CVS

  • Property svn:executable set to *
File size: 16.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//
27// $Id: G4FinalStateIonisationBorn.cc,v 1.9 2007/11/26 17:27:09 pia Exp $
28// GEANT4 tag $Name:  $
29//
30// Contact Author: Sebastien Incerti (incerti@cenbg.in2p3.fr)
31//                 Maria Grazia Pia  (Maria.Grazia.Pia@cern.ch)
32//
33// Reference: TNS Geant4-DNA paper
34// Reference for implementation model: NIM. 155, pp. 145-156, 1978
35
36// History:
37// -----------
38// Date         Name              Modification
39// 28 Apr 2007  M.G. Pia          Created in compliance with design described in TNS paper
40//    Nov 2007  S. Incerti        Implementation
41// 26 Nov 2007  MGP               Cleaned up std::
42//
43// -------------------------------------------------------------------
44
45// Class description:
46// Reference: TNS Geant4-DNA paper
47// S. Chauvie et al., Geant4 physics processes for microdosimetry simulation:
48// design foundation and implementation of the first set of models,
49// IEEE Trans. Nucl. Sci., vol. 54, no. 6, Dec. 2007.
50// Further documentation available from http://www.ge.infn.it/geant4/dna
51
52// -------------------------------------------------------------------
53
54
55#include "G4FinalStateIonisationBorn.hh"
56#include "G4Track.hh"
57#include "G4Step.hh"
58#include "G4DynamicParticle.hh"
59#include "Randomize.hh"
60
61#include "G4ParticleTypes.hh"
62#include "G4ParticleDefinition.hh"
63#include "G4Electron.hh"
64#include "G4Proton.hh"
65#include "G4SystemOfUnits.hh"
66#include "G4ParticleMomentum.hh"
67
68
69G4FinalStateIonisationBorn::G4FinalStateIonisationBorn()
70{
71
72  name = "IonisationBorn";
73
74  // NEW
75  // Factor to scale microscopic/macroscopic cross section data in water
76 
77  G4double scaleFactor = (1.e-22 / 3.343) * m*m;
78
79  // Energy limits
80  G4ParticleDefinition* electronDef = G4Electron::ElectronDefinition();
81  G4ParticleDefinition* protonDef = G4Proton::ProtonDefinition();
82
83  G4String electron;
84  G4String proton;
85
86  // Default energy limits (defined for protection against anomalous behaviour only)
87  lowEnergyLimitDefault = 25 * eV;
88  highEnergyLimitDefault = 10 * MeV;
89
90  char *path = getenv("G4LEDATA");
91 
92  if (!path)
93    G4Exception("G4DNACrossSectionDataSet::FullFileName - G4LEDATA environment variable not set");
94
95  // Data members for electrons
96
97  if (electronDef != 0)
98    {
99      electron = electronDef->GetParticleName();
100      lowEnergyLimit[electron] = 25. * eV;
101      highEnergyLimit[electron] = 30. * keV;
102
103      std::ostringstream eFullFileName;
104      eFullFileName << path << "/dna/sigmadiff_ionisation_e_born.dat";
105      std::ifstream eDiffCrossSection(eFullFileName.str().c_str());
106      // eDiffCrossSection(eFullFileName.str().c_str());
107      if (!eDiffCrossSection)
108        { 
109          // G4cout << "ERROR OPENING DATA FILE IN ELECTRON BORN IONIZATION !!! " << G4endl;
110          G4Exception("G4FinalStateIonisationBorn::ERROR OPENING electron DATA FILE");
111          while(1); // ---- MGP ---- What is this?
112        }
113     
114      eTdummyVec.push_back(0.);
115      while(!eDiffCrossSection.eof())
116        {
117          double tDummy;
118          double eDummy;
119          eDiffCrossSection>>tDummy>>eDummy;
120          if (tDummy != eTdummyVec.back()) eTdummyVec.push_back(tDummy);
121          for (int j=0; j<5; j++)
122            {
123              eDiffCrossSection>>eDiffCrossSectionData[j][tDummy][eDummy];
124              eDiffCrossSectionData[j][tDummy][eDummy]*=scaleFactor;
125              eVecm[tDummy].push_back(eDummy);
126            }
127        }
128
129    }
130  else
131    {
132      G4Exception("G4FinalStateIonisationBorn Constructor: electron is not defined");
133    }
134
135  // Data members for protons
136
137  if (protonDef != 0)
138    {
139      proton = protonDef->GetParticleName();
140      lowEnergyLimit[proton] = 500. * keV;
141      highEnergyLimit[proton] = 10. * MeV;
142
143      std::ostringstream pFullFileName;
144      pFullFileName << path << "/dna/sigmadiff_ionisation_p_born.dat";
145      std::ifstream pDiffCrossSection(pFullFileName.str().c_str());
146      //     pDiffCrossSection(pFullFileName.str().c_str());
147      if (!pDiffCrossSection)
148        { 
149          // G4cout<<"ERROR OPENING DATA FILE IN PROTON BORN IONIZATION !!! "<<G4endl;
150          G4Exception("G4FinalStateIonisationBorn::ERROR OPENING proton DATA FILE");
151          while(1); // ---- MGP ---- What is this?
152        }
153     
154      pTdummyVec.push_back(0.);
155      while(!pDiffCrossSection.eof())
156        {
157          double tDummy;
158          double eDummy;
159          pDiffCrossSection>>tDummy>>eDummy;
160          if (tDummy != pTdummyVec.back()) pTdummyVec.push_back(tDummy);
161          for (int j=0; j<5; j++)
162            {
163              pDiffCrossSection>>pDiffCrossSectionData[j][tDummy][eDummy];
164              pDiffCrossSectionData[j][tDummy][eDummy]*=scaleFactor;
165              //G4cout << "j=" << j << " Tdum=" << tDummy << " Edum=" << eDummy << " pDiff=" << pDiffCrossSectionData[j][tDummy][eDummy] << G4endl;
166              pVecm[tDummy].push_back(eDummy);
167            }
168        }
169    }
170  else
171    {
172      G4Exception("G4FinalStateIonisationBorn Constructor: proton is not defined");
173    }
174}
175
176
177G4FinalStateIonisationBorn::~G4FinalStateIonisationBorn()
178{
179  eVecm.clear();
180  pVecm.clear();
181}
182
183
184const G4FinalStateProduct& G4FinalStateIonisationBorn::GenerateFinalState(const G4Track& track, const G4Step& /* step */)
185{
186  // Clear previous secondaries, energy deposit and particle kill status
187  product.Clear();
188
189  const G4DynamicParticle* particle = track.GetDynamicParticle();
190
191  G4double lowLim = lowEnergyLimitDefault;
192  G4double highLim = highEnergyLimitDefault;
193
194  G4double k = particle->GetKineticEnergy();
195
196  const G4String& particleName = particle->GetDefinition()->GetParticleName();
197
198  // Retrieve energy limits for the current particle type
199
200  std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
201  pos1 = lowEnergyLimit.find(particleName);
202
203  // Lower limit
204  if (pos1 != lowEnergyLimit.end())
205    {
206      lowLim = pos1->second;
207    }
208
209  // Upper limit
210  std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
211  pos2 = highEnergyLimit.find(particleName);
212
213  if (pos2 != highEnergyLimit.end())
214    {
215      highLim = pos2->second;
216    }
217
218  // Verify that the current track is within the energy limits of validity of the cross section model
219
220  if (k >= lowLim && k <= highLim)
221    {
222      // Kinetic energy of primary particle
223
224      G4ParticleMomentum primaryDirection = particle->GetMomentumDirection();
225      G4double particleMass = particle->GetDefinition()->GetPDGMass();
226      G4double totalEnergy = k + particleMass;
227      G4double pSquare = k * (totalEnergy + particleMass);
228      G4double totalMomentum = std::sqrt(pSquare);
229
230      const G4String& particleName = particle->GetDefinition()->GetParticleName();
231 
232      G4int ionizationShell = cross.RandomSelect(k,particleName);
233 
234      G4double secondaryKinetic = RandomizeEjectedElectronEnergy(particle->GetDefinition(),k,ionizationShell);
235 
236      G4double bindingEnergy = waterStructure.IonisationEnergy(ionizationShell);
237
238      G4double cosTheta = 0.;
239      G4double phi = 0.; 
240      RandomizeEjectedElectronDirection(track.GetDefinition(), k,secondaryKinetic, cosTheta, phi);
241
242      G4double sinTheta = std::sqrt(1.-cosTheta*cosTheta);
243      G4double dirX = sinTheta*std::cos(phi);
244      G4double dirY = sinTheta*std::sin(phi);
245      G4double dirZ = cosTheta;
246      G4ThreeVector deltaDirection(dirX,dirY,dirZ);
247      deltaDirection.rotateUz(primaryDirection);
248
249      G4double deltaTotalMomentum = std::sqrt(secondaryKinetic*(secondaryKinetic + 2.*electron_mass_c2 ));
250
251      //Primary Particle Direction
252      G4double finalPx = totalMomentum*primaryDirection.x() - deltaTotalMomentum*deltaDirection.x();
253      G4double finalPy = totalMomentum*primaryDirection.y() - deltaTotalMomentum*deltaDirection.y();
254      G4double finalPz = totalMomentum*primaryDirection.z() - deltaTotalMomentum*deltaDirection.z();
255      G4double finalMomentum = std::sqrt(finalPx*finalPx + finalPy*finalPy + finalPz*finalPz);
256      finalPx /= finalMomentum;
257      finalPy /= finalMomentum;
258      finalPz /= finalMomentum;
259
260      product.ModifyPrimaryParticle(finalPx,finalPy,finalPz,k-bindingEnergy-secondaryKinetic);
261      product.AddEnergyDeposit(bindingEnergy);
262
263      G4DynamicParticle* aElectron = new G4DynamicParticle(G4Electron::Electron(),deltaDirection,secondaryKinetic);
264      product.AddSecondary(aElectron);
265    }
266
267  return product;
268}
269
270
271G4double G4FinalStateIonisationBorn::RandomizeEjectedElectronEnergy(G4ParticleDefinition* particleDefinition, 
272                                                                    G4double k, 
273                                                                    G4int shell)
274{
275
276  if (particleDefinition == G4Electron::ElectronDefinition()) 
277    {
278
279      G4double maximumEnergyTransfer=0.;
280      if ((k+waterStructure.IonisationEnergy(shell))/2. > k) maximumEnergyTransfer=k;
281      else maximumEnergyTransfer = (k+waterStructure.IonisationEnergy(shell))/2.;
282   
283      G4double crossSectionMaximum = 0.;
284      for(G4double value=waterStructure.IonisationEnergy(shell); value<=maximumEnergyTransfer; value+=0.1*eV)
285        {
286          G4double differentialCrossSection = DifferentialCrossSection(particleDefinition, k/eV, value/eV, shell);
287          if(differentialCrossSection >= crossSectionMaximum) crossSectionMaximum = differentialCrossSection;
288        }
289 
290      G4double secondaryElectronKineticEnergy=0.;
291      do 
292        {
293          secondaryElectronKineticEnergy = G4UniformRand() * (maximumEnergyTransfer-waterStructure.IonisationEnergy(shell));
294        } while(G4UniformRand()*crossSectionMaximum >
295              DifferentialCrossSection(particleDefinition, k/eV,(secondaryElectronKineticEnergy+waterStructure.IonisationEnergy(shell))/eV,shell));
296
297      return secondaryElectronKineticEnergy;
298 
299    }
300 
301  if (particleDefinition == G4Proton::ProtonDefinition()) 
302    {
303      G4double maximumKineticEnergyTransfer = 4.* (electron_mass_c2 / proton_mass_c2) * k - (waterStructure.IonisationEnergy(shell));
304
305      G4double crossSectionMaximum = 0.;
306      for (G4double value = waterStructure.IonisationEnergy(shell); 
307           value<=4.*waterStructure.IonisationEnergy(shell) ; 
308           value+=0.1*eV)
309        {
310          G4double differentialCrossSection = DifferentialCrossSection(particleDefinition, k/eV, value/eV, shell);
311          if (differentialCrossSection >= crossSectionMaximum) crossSectionMaximum = differentialCrossSection;
312        }
313
314      G4double secondaryElectronKineticEnergy = 0.;
315      do
316        {
317          secondaryElectronKineticEnergy = G4UniformRand() * maximumKineticEnergyTransfer;
318        } while(G4UniformRand()*crossSectionMaximum >= 
319              DifferentialCrossSection(particleDefinition, k/eV,(secondaryElectronKineticEnergy+waterStructure.IonisationEnergy(shell))/eV,shell));
320
321      return secondaryElectronKineticEnergy;
322    }
323
324  return 0;
325}
326
327
328void G4FinalStateIonisationBorn::RandomizeEjectedElectronDirection(G4ParticleDefinition* particleDefinition, 
329                                                                   G4double k, 
330                                                                   G4double secKinetic, 
331                                                                   G4double & cosTheta, 
332                                                                   G4double & phi )
333{
334  if (particleDefinition == G4Electron::ElectronDefinition()) 
335    {
336
337      phi = twopi * G4UniformRand();
338      if (secKinetic < 50.*eV) cosTheta = (2.*G4UniformRand())-1.;
339      else if (secKinetic <= 200.*eV)   
340        {
341          if (G4UniformRand() <= 0.1) cosTheta = (2.*G4UniformRand())-1.;
342          else cosTheta = G4UniformRand()*(std::sqrt(2.)/2);
343        }
344      else     
345        {
346          G4double sin2O = (1.-secKinetic/k) / (1.+secKinetic/(2.*electron_mass_c2));
347          cosTheta = std::sqrt(1.-sin2O);
348        }
349    }
350 
351  if (particleDefinition == G4Proton::ProtonDefinition()) 
352    {
353      G4double maxSecKinetic = 4.* (electron_mass_c2 / proton_mass_c2) * k;
354      phi = twopi * G4UniformRand();
355      cosTheta = std::sqrt(secKinetic / maxSecKinetic);
356    }                   
357}
358
359
360double G4FinalStateIonisationBorn::DifferentialCrossSection(G4ParticleDefinition * particleDefinition, 
361                                                            G4double k, 
362                                                            G4double energyTransfer, 
363                                                            G4int ionizationLevelIndex)
364{
365  G4double sigma = 0.;
366
367  if (energyTransfer >= waterStructure.IonisationEnergy(ionizationLevelIndex))
368    {
369      G4double valueT1 = 0;
370      G4double valueT2 = 0;
371      G4double valueE21 = 0;
372      G4double valueE22 = 0;
373      G4double valueE12 = 0;
374      G4double valueE11 = 0;
375
376      G4double xs11  =  0;   
377      G4double xs12 = 0; 
378      G4double xs21 = 0; 
379      G4double xs22 = 0; 
380
381 
382      if (particleDefinition == G4Electron::ElectronDefinition()) 
383        {
384          // k should be in eV and energy transfer eV also
385          std::vector<double>::iterator t2 = std::upper_bound(eTdummyVec.begin(),eTdummyVec.end(), k);
386          std::vector<double>::iterator t1 = t2-1;
387          std::vector<double>::iterator e12 = std::upper_bound(eVecm[(*t1)].begin(),eVecm[(*t1)].end(), energyTransfer);
388          std::vector<double>::iterator e11 = e12-1;
389
390          std::vector<double>::iterator e22 = std::upper_bound(eVecm[(*t2)].begin(),eVecm[(*t2)].end(), energyTransfer);
391          std::vector<double>::iterator e21 = e22-1;
392
393          valueT1  =*t1;
394          valueT2  =*t2;
395          valueE21 =*e21;
396          valueE22 =*e22;
397          valueE12 =*e12;
398          valueE11 =*e11;
399
400          xs11 = eDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE11];
401          xs12 = eDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE12];
402          xs21 = eDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE21];
403          xs22 = eDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE22];
404
405        }
406 
407      if (particleDefinition == G4Proton::ProtonDefinition()) 
408        {
409          // k should be in eV and energy transfer eV also
410          std::vector<double>::iterator t2 = std::upper_bound(pTdummyVec.begin(),pTdummyVec.end(), k);
411          std::vector<double>::iterator t1 = t2-1;
412          std::vector<double>::iterator e12 = std::upper_bound(pVecm[(*t1)].begin(),pVecm[(*t1)].end(), energyTransfer);
413          std::vector<double>::iterator e11 = e12-1;
414
415          std::vector<double>::iterator e22 = std::upper_bound(pVecm[(*t2)].begin(),pVecm[(*t2)].end(), energyTransfer);
416          std::vector<double>::iterator e21 = e22-1;
417 
418          valueT1  =*t1;
419          valueT2  =*t2;
420          valueE21 =*e21;
421          valueE22 =*e22;
422          valueE12 =*e12;
423          valueE11 =*e11;
424
425          xs11 = pDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE11];
426          xs12 = pDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE12];
427          xs21 = pDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE21];
428          xs22 = pDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE22];
429        }
430 
431      G4double xsProduct = xs11 * xs12 * xs21 * xs22;
432      // if (xs11==0 || xs12==0 ||xs21==0 ||xs22==0) return (0.);
433      if (xsProduct != 0.)
434        {
435          sigma = QuadInterpolator(valueE11, valueE12, 
436                                   valueE21, valueE22, 
437                                   xs11, xs12, 
438                                   xs21, xs22, 
439                                   valueT1, valueT2, 
440                                   k, energyTransfer);
441        }
442    }
443  return sigma;
444}
445
446
447G4double G4FinalStateIonisationBorn::LogLogInterpolate(G4double e1, 
448                                                       G4double e2, 
449                                                       G4double e, 
450                                                       G4double xs1, 
451                                                       G4double xs2)
452{
453  G4double a = (std::log10(xs2)-std::log10(xs1)) / (std::log10(e2)-std::log10(e1));
454  G4double b = std::log10(xs2) - a*std::log10(e2);
455  G4double sigma = a*std::log10(e) + b;
456  G4double value = (std::pow(10.,sigma));
457  return value;
458}
459
460
461G4double G4FinalStateIonisationBorn::QuadInterpolator(G4double e11, G4double e12, 
462                                                      G4double e21, G4double e22, 
463                                                      G4double xs11, G4double xs12, 
464                                                      G4double xs21, G4double xs22, 
465                                                      G4double t1, G4double t2, 
466                                                      G4double t, G4double e)
467{
468  G4double interpolatedvalue1 = LogLogInterpolate(e11, e12, e, xs11, xs12);
469  G4double interpolatedvalue2 = LogLogInterpolate(e21, e22, e, xs21, xs22);
470  G4double value = LogLogInterpolate(t1, t2, t, interpolatedvalue1, interpolatedvalue2);
471  return value;
472}
473
474
Note: See TracBrowser for help on using the repository browser.