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

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

geant4 tag 9.4

  • Property svn:executable set to *
File size: 15.8 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: G4FinalStateIonisationBorn.cc,v 1.19 2009/06/11 15:47:08 mantero Exp $
27// GEANT4 tag $Name: geant4-09-04-ref-00 $
28
29#include "G4FinalStateIonisationBorn.hh"
30
31//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
32
33G4FinalStateIonisationBorn::G4FinalStateIonisationBorn()
34{
35  G4double scaleFactor = (1.e-22 / 3.343) * m*m;
36
37  G4ParticleDefinition* electronDef = G4Electron::ElectronDefinition();
38  G4ParticleDefinition* protonDef = G4Proton::ProtonDefinition();
39
40  G4String electron;
41  G4String proton;
42
43  lowEnergyLimitDefault = 12.61 * eV; // SI: i/o 25 eV
44  highEnergyLimitDefault = 10 * MeV;
45
46  char *path = getenv("G4LEDATA");
47 
48  if (!path)
49    G4Exception("G4DNACrossSectionDataSet::FullFileName - G4LEDATA environment variable not set");
50
51  if (electronDef != 0)
52  {
53    electron = electronDef->GetParticleName();
54    lowEnergyLimit[electron] = 12.61 * eV; // SI: i/o 25 eV
55    highEnergyLimit[electron] = 30. * keV;
56
57    std::ostringstream eFullFileName;
58    eFullFileName << path << "/dna/sigmadiff_ionisation_e_born.dat";
59    std::ifstream eDiffCrossSection(eFullFileName.str().c_str());
60    if (!eDiffCrossSection)
61    { 
62      G4Exception("G4FinalStateIonisationBorn::ERROR OPENING electron DATA FILE");
63    }
64     
65    eTdummyVec.push_back(0.);
66    while(!eDiffCrossSection.eof())
67    {
68      double tDummy;
69      double eDummy;
70      eDiffCrossSection>>tDummy>>eDummy;
71      if (tDummy != eTdummyVec.back()) eTdummyVec.push_back(tDummy);
72      for (int j=0; j<5; j++)
73      {
74        eDiffCrossSection>>eDiffCrossSectionData[j][tDummy][eDummy];
75
76        // SI - only if eof is not reached !
77        if (!eDiffCrossSection.eof()) eDiffCrossSectionData[j][tDummy][eDummy]*=scaleFactor;
78
79        eVecm[tDummy].push_back(eDummy);
80
81      }
82    }
83
84  }
85  else
86  {
87    G4Exception("G4FinalStateIonisationBorn Constructor: electron is not defined");
88  }
89
90  if (protonDef != 0)
91  {
92    proton = protonDef->GetParticleName();
93    lowEnergyLimit[proton] = 500. * keV;
94    highEnergyLimit[proton] = 10. * MeV;
95
96    std::ostringstream pFullFileName;
97    pFullFileName << path << "/dna/sigmadiff_ionisation_p_born.dat";
98    std::ifstream pDiffCrossSection(pFullFileName.str().c_str());
99    if (!pDiffCrossSection)
100    { 
101      G4Exception("G4FinalStateIonisationBorn::ERROR OPENING proton DATA FILE");
102    }
103     
104    pTdummyVec.push_back(0.);
105    while(!pDiffCrossSection.eof())
106    {
107      double tDummy;
108      double eDummy;
109      pDiffCrossSection>>tDummy>>eDummy;
110      if (tDummy != pTdummyVec.back()) pTdummyVec.push_back(tDummy);
111      for (int j=0; j<5; j++)
112      {
113        pDiffCrossSection>>pDiffCrossSectionData[j][tDummy][eDummy];
114
115        // SI - only if eof is not reached !
116        if (!pDiffCrossSection.eof()) pDiffCrossSectionData[j][tDummy][eDummy]*=scaleFactor;
117
118        pVecm[tDummy].push_back(eDummy);
119      }
120    }
121  }
122  else
123  {
124    G4Exception("G4FinalStateIonisationBorn Constructor: proton is not defined");
125  }
126
127   G4cout << G4endl;
128   G4cout << "*******************************************************************************" << G4endl;
129   G4cout << "*******************************************************************************" << G4endl;
130   G4cout << "   The class G4FinalStateIonisationBorn is NOT SUPPORTED ANYMORE. " << G4endl;
131   G4cout << "   It will be REMOVED with the next major release of Geant4. " << G4endl;
132   G4cout << "   Please consult: https://twiki.cern.ch/twiki/bin/view/Geant4/LoweProcesses" << G4endl;
133   G4cout << "*******************************************************************************" << G4endl;
134   G4cout << "*******************************************************************************" << G4endl;
135   G4cout << G4endl;
136}
137
138//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
139
140G4FinalStateIonisationBorn::~G4FinalStateIonisationBorn()
141{
142  eVecm.clear();
143  pVecm.clear();
144}
145
146//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
147
148const G4FinalStateProduct& G4FinalStateIonisationBorn::GenerateFinalState(const G4Track& track, const G4Step& /* step */)
149{
150  product.Clear();
151
152  const G4DynamicParticle* particle = track.GetDynamicParticle();
153
154  G4double lowLim = lowEnergyLimitDefault;
155  G4double highLim = highEnergyLimitDefault;
156
157  G4double k = particle->GetKineticEnergy();
158
159  const G4String& particleName = particle->GetDefinition()->GetParticleName();
160
161  std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
162  pos1 = lowEnergyLimit.find(particleName);
163
164  if (pos1 != lowEnergyLimit.end())
165  {
166    lowLim = pos1->second;
167  }
168
169  std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
170  pos2 = highEnergyLimit.find(particleName);
171
172  if (pos2 != highEnergyLimit.end())
173  {
174    highLim = pos2->second;
175  }
176
177  if (k >= lowLim && k <= highLim)
178  {
179    G4ParticleMomentum primaryDirection = particle->GetMomentumDirection();
180    G4double particleMass = particle->GetDefinition()->GetPDGMass();
181    G4double totalEnergy = k + particleMass;
182    G4double pSquare = k * (totalEnergy + particleMass);
183    G4double totalMomentum = std::sqrt(pSquare);
184
185    const G4String& particleName = particle->GetDefinition()->GetParticleName();
186 
187    G4int ionizationShell = cross.RandomSelect(k,particleName);
188 
189    G4double secondaryKinetic = RandomizeEjectedElectronEnergy(particle->GetDefinition(),k,ionizationShell);
190 
191    G4double bindingEnergy = waterStructure.IonisationEnergy(ionizationShell);
192
193    G4double cosTheta = 0.;
194    G4double phi = 0.; 
195    RandomizeEjectedElectronDirection(track.GetDefinition(), k,secondaryKinetic, cosTheta, phi);
196
197    G4double sinTheta = std::sqrt(1.-cosTheta*cosTheta);
198    G4double dirX = sinTheta*std::cos(phi);
199    G4double dirY = sinTheta*std::sin(phi);
200    G4double dirZ = cosTheta;
201    G4ThreeVector deltaDirection(dirX,dirY,dirZ);
202    deltaDirection.rotateUz(primaryDirection);
203
204    G4double deltaTotalMomentum = std::sqrt(secondaryKinetic*(secondaryKinetic + 2.*electron_mass_c2 ));
205
206    G4double finalPx = totalMomentum*primaryDirection.x() - deltaTotalMomentum*deltaDirection.x();
207    G4double finalPy = totalMomentum*primaryDirection.y() - deltaTotalMomentum*deltaDirection.y();
208    G4double finalPz = totalMomentum*primaryDirection.z() - deltaTotalMomentum*deltaDirection.z();
209    G4double finalMomentum = std::sqrt(finalPx*finalPx + finalPy*finalPy + finalPz*finalPz);
210    finalPx /= finalMomentum;
211    finalPy /= finalMomentum;
212    finalPz /= finalMomentum;
213
214    product.ModifyPrimaryParticle(finalPx,finalPy,finalPz,k-bindingEnergy-secondaryKinetic);
215    product.AddEnergyDeposit(bindingEnergy);
216
217    G4DynamicParticle* aElectron = new G4DynamicParticle(G4Electron::Electron(),deltaDirection,secondaryKinetic);
218    product.AddSecondary(aElectron);
219  }
220
221  return product;
222}
223
224//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
225
226G4double G4FinalStateIonisationBorn::RandomizeEjectedElectronEnergy(G4ParticleDefinition* particleDefinition, 
227G4double k, G4int shell)
228{
229  if (particleDefinition == G4Electron::ElectronDefinition()) 
230  {
231    G4double maximumEnergyTransfer=0.;
232    if ((k+waterStructure.IonisationEnergy(shell))/2. > k) maximumEnergyTransfer=k;
233    else maximumEnergyTransfer = (k+waterStructure.IonisationEnergy(shell))/2.;
234   
235    G4double crossSectionMaximum = 0.;
236    for(G4double value=waterStructure.IonisationEnergy(shell); value<=maximumEnergyTransfer; value+=0.1*eV)
237    {
238      G4double differentialCrossSection = DifferentialCrossSection(particleDefinition, k/eV, value/eV, shell);
239      if(differentialCrossSection >= crossSectionMaximum) crossSectionMaximum = differentialCrossSection;
240    }
241 
242    G4double secondaryElectronKineticEnergy=0.;
243    do 
244    {
245      secondaryElectronKineticEnergy = G4UniformRand() * (maximumEnergyTransfer-waterStructure.IonisationEnergy(shell));
246    } while(G4UniformRand()*crossSectionMaximum >
247      DifferentialCrossSection(particleDefinition, k/eV,(secondaryElectronKineticEnergy+waterStructure.IonisationEnergy(shell))/eV,shell));
248
249    return secondaryElectronKineticEnergy;
250 
251  }
252 
253  if (particleDefinition == G4Proton::ProtonDefinition()) 
254  {
255    G4double maximumKineticEnergyTransfer = 4.* (electron_mass_c2 / proton_mass_c2) * k - (waterStructure.IonisationEnergy(shell));
256
257    G4double crossSectionMaximum = 0.;
258    for (G4double value = waterStructure.IonisationEnergy(shell); 
259         value<=4.*waterStructure.IonisationEnergy(shell) ; 
260         value+=0.1*eV)
261    {
262      G4double differentialCrossSection = DifferentialCrossSection(particleDefinition, k/eV, value/eV, shell);
263      if (differentialCrossSection >= crossSectionMaximum) crossSectionMaximum = differentialCrossSection;
264    }
265
266    G4double secondaryElectronKineticEnergy = 0.;
267    do
268    {
269      secondaryElectronKineticEnergy = G4UniformRand() * maximumKineticEnergyTransfer;
270    } while(G4UniformRand()*crossSectionMaximum >= 
271              DifferentialCrossSection(particleDefinition, k/eV,(secondaryElectronKineticEnergy+waterStructure.IonisationEnergy(shell))/eV,shell));
272
273    return secondaryElectronKineticEnergy;
274  }
275
276  return 0;
277}
278
279//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
280
281void G4FinalStateIonisationBorn::RandomizeEjectedElectronDirection(G4ParticleDefinition* particleDefinition, 
282                                                                   G4double k, 
283                                                                   G4double secKinetic, 
284                                                                   G4double & cosTheta, 
285                                                                   G4double & phi )
286{
287  if (particleDefinition == G4Electron::ElectronDefinition()) 
288  {
289    phi = twopi * G4UniformRand();
290    if (secKinetic < 50.*eV) cosTheta = (2.*G4UniformRand())-1.;
291    else if (secKinetic <= 200.*eV)     
292    {
293      if (G4UniformRand() <= 0.1) cosTheta = (2.*G4UniformRand())-1.;
294      else cosTheta = G4UniformRand()*(std::sqrt(2.)/2);
295    }
296    else       
297    {
298      G4double sin2O = (1.-secKinetic/k) / (1.+secKinetic/(2.*electron_mass_c2));
299      cosTheta = std::sqrt(1.-sin2O);
300    }
301  }
302 
303  if (particleDefinition == G4Proton::ProtonDefinition()) 
304  {
305    G4double maxSecKinetic = 4.* (electron_mass_c2 / proton_mass_c2) * k;
306    phi = twopi * G4UniformRand();
307    cosTheta = std::sqrt(secKinetic / maxSecKinetic);
308  }                     
309}
310
311//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
312
313double G4FinalStateIonisationBorn::DifferentialCrossSection(G4ParticleDefinition * particleDefinition, 
314                                                            G4double k, 
315                                                            G4double energyTransfer, 
316                                                            G4int ionizationLevelIndex)
317{
318  G4double sigma = 0.;
319
320  if (energyTransfer >= waterStructure.IonisationEnergy(ionizationLevelIndex))
321  {
322    G4double valueT1 = 0;
323    G4double valueT2 = 0;
324    G4double valueE21 = 0;
325    G4double valueE22 = 0;
326    G4double valueE12 = 0;
327    G4double valueE11 = 0;
328
329    G4double xs11 = 0;   
330    G4double xs12 = 0; 
331    G4double xs21 = 0; 
332    G4double xs22 = 0; 
333
334    if (particleDefinition == G4Electron::ElectronDefinition()) 
335    {
336      // k should be in eV and energy transfer eV also
337
338      std::vector<double>::iterator t2 = std::upper_bound(eTdummyVec.begin(),eTdummyVec.end(), k);
339
340      std::vector<double>::iterator t1 = t2-1;
341
342      // SI : the following condition avoids situations where energyTransfer >last vector element
343      if (energyTransfer <= eVecm[(*t1)].back())
344      {
345        std::vector<double>::iterator e12 = std::upper_bound(eVecm[(*t1)].begin(),eVecm[(*t1)].end(), energyTransfer);
346        std::vector<double>::iterator e11 = e12-1;
347
348        std::vector<double>::iterator e22 = std::upper_bound(eVecm[(*t2)].begin(),eVecm[(*t2)].end(), energyTransfer);
349        std::vector<double>::iterator e21 = e22-1;
350
351        valueT1  =*t1;
352        valueT2  =*t2;
353        valueE21 =*e21;
354        valueE22 =*e22;
355        valueE12 =*e12;
356        valueE11 =*e11;
357
358        xs11 = eDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE11];
359        xs12 = eDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE12];
360        xs21 = eDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE21];
361        xs22 = eDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE22];
362      }
363
364    }
365 
366   if (particleDefinition == G4Proton::ProtonDefinition()) 
367   {
368      // k should be in eV and energy transfer eV also
369      std::vector<double>::iterator t2 = std::upper_bound(pTdummyVec.begin(),pTdummyVec.end(), k);
370      std::vector<double>::iterator t1 = t2-1;
371     
372        std::vector<double>::iterator e12 = std::upper_bound(pVecm[(*t1)].begin(),pVecm[(*t1)].end(), energyTransfer);
373        std::vector<double>::iterator e11 = e12-1;
374
375        std::vector<double>::iterator e22 = std::upper_bound(pVecm[(*t2)].begin(),pVecm[(*t2)].end(), energyTransfer);
376        std::vector<double>::iterator e21 = e22-1;
377 
378        valueT1  =*t1;
379        valueT2  =*t2;
380        valueE21 =*e21;
381        valueE22 =*e22;
382        valueE12 =*e12;
383        valueE11 =*e11;
384
385        xs11 = pDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE11];
386        xs12 = pDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE12];
387        xs21 = pDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE21];
388        xs22 = pDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE22];
389
390   }
391 
392   G4double xsProduct = xs11 * xs12 * xs21 * xs22;
393   if (xsProduct != 0.)
394   {
395     sigma = QuadInterpolator(     valueE11, valueE12, 
396                                   valueE21, valueE22, 
397                                   xs11, xs12, 
398                                   xs21, xs22, 
399                                   valueT1, valueT2, 
400                                   k, energyTransfer);
401   }
402 
403 }
404 
405  return sigma;
406}
407
408//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
409
410G4double G4FinalStateIonisationBorn::LogLogInterpolate(G4double e1, 
411                                                       G4double e2, 
412                                                       G4double e, 
413                                                       G4double xs1, 
414                                                       G4double xs2)
415{
416  G4double a = (std::log10(xs2)-std::log10(xs1)) / (std::log10(e2)-std::log10(e1));
417  G4double b = std::log10(xs2) - a*std::log10(e2);
418  G4double sigma = a*std::log10(e) + b;
419  G4double value = (std::pow(10.,sigma));
420  return value;
421}
422
423//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
424
425G4double G4FinalStateIonisationBorn::QuadInterpolator(G4double e11, G4double e12, 
426                                                      G4double e21, G4double e22, 
427                                                      G4double xs11, G4double xs12, 
428                                                      G4double xs21, G4double xs22, 
429                                                      G4double t1, G4double t2, 
430                                                      G4double t, G4double e)
431{
432  G4double interpolatedvalue1 = LogLogInterpolate(e11, e12, e, xs11, xs12);
433  G4double interpolatedvalue2 = LogLogInterpolate(e21, e22, e, xs21, xs22);
434  G4double value = LogLogInterpolate(t1, t2, t, interpolatedvalue1, interpolatedvalue2);
435  return value;
436}
437
438
Note: See TracBrowser for help on using the repository browser.