source: trunk/source/processes/electromagnetic/lowenergy/src/G4CrossSectionIonisationRudd.cc @ 961

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

update processes

  • Property svn:executable set to *
File size: 8.7 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: G4CrossSectionIonisationRudd.cc,v 1.4 2008/07/14 20:47:34 sincerti Exp $
27// GEANT4 tag $Name: geant4-09-02-ref-02 $
28
29#include "G4CrossSectionIonisationRudd.hh"
30
31//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
32
33G4CrossSectionIonisationRudd::G4CrossSectionIonisationRudd()
34{
35  lowEnergyLimitDefault = 0 * eV;
36  highEnergyLimitDefault = 100 * MeV;
37
38  G4String fileProton("dna/sigma_ionisation_p_rudd");
39  G4String fileHydrogen("dna/sigma_ionisation_h_rudd");
40  G4String fileAlphaPlusPlus("dna/sigma_ionisation_alphaplusplus_rudd");
41  G4String fileAlphaPlus("dna/sigma_ionisation_alphaplus_rudd");
42  G4String fileHelium("dna/sigma_ionisation_he_rudd");
43
44  G4DNAGenericIonsManager *instance;
45  instance = G4DNAGenericIonsManager::Instance();
46  G4ParticleDefinition* protonDef = G4Proton::ProtonDefinition();
47  G4ParticleDefinition* hydrogenDef = instance->GetIon("hydrogen");
48  G4ParticleDefinition* alphaPlusPlusDef = instance->GetIon("alpha++");
49  G4ParticleDefinition* alphaPlusDef = instance->GetIon("alpha+");
50  G4ParticleDefinition* heliumDef = instance->GetIon("helium");
51
52  G4String proton;
53  G4String hydrogen;
54  G4String alphaPlusPlus;
55  G4String alphaPlus;
56  G4String helium;
57
58  G4double scaleFactor = 1 * m*m;
59
60  if (protonDef != 0)
61  {
62    proton = protonDef->GetParticleName();
63    tableFile[proton] = fileProton;
64
65    lowEnergyLimit[proton] = 0. * eV;
66    highEnergyLimit[proton] = 500. * keV;
67
68    G4DNACrossSectionDataSet* tableProton = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, 
69                                                                           eV,
70                                                                           scaleFactor );
71    tableProton->LoadData(fileProton);
72    tableData[proton] = tableProton;
73  }
74  else
75  {
76    G4Exception("G4CrossSectionIonisationRudd Constructor: proton is not defined");
77  }
78
79  if (hydrogenDef != 0)
80  {
81    hydrogen = hydrogenDef->GetParticleName();
82    tableFile[hydrogen] = fileHydrogen;
83
84    lowEnergyLimit[hydrogen] = 0. * eV;
85    highEnergyLimit[hydrogen] = 100. * MeV;
86
87    G4DNACrossSectionDataSet* tableHydrogen = new G4DNACrossSectionDataSet(new G4LogLogInterpolation,
88                                                                             eV,
89                                                                             scaleFactor );
90    tableHydrogen->LoadData(fileHydrogen);
91     
92    tableData[hydrogen] = tableHydrogen;
93  }
94  else
95  {
96    G4Exception("G4CrossSectionIonisationRudd Constructor: hydrogen is not defined");
97  }
98
99  if (alphaPlusPlusDef != 0)
100  {
101    alphaPlusPlus = alphaPlusPlusDef->GetParticleName();
102    tableFile[alphaPlusPlus] = fileAlphaPlusPlus;
103
104    lowEnergyLimit[alphaPlusPlus] = 0. * eV;
105    highEnergyLimit[alphaPlusPlus] = 10. * MeV;
106
107    G4DNACrossSectionDataSet* tableAlphaPlusPlus = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, 
108                                                                                  eV,
109                                                                                  scaleFactor );
110    tableAlphaPlusPlus->LoadData(fileAlphaPlusPlus);
111     
112    tableData[alphaPlusPlus] = tableAlphaPlusPlus;
113  }
114  else
115  {
116    G4Exception("G4CrossSectionIonisationRudd Constructor: alphaPlusPlus is not defined");
117  }
118
119  if (alphaPlusDef != 0)
120  {
121    alphaPlus = alphaPlusDef->GetParticleName();
122    tableFile[alphaPlus] = fileAlphaPlus;
123
124    lowEnergyLimit[alphaPlus] = 0. * eV;
125    highEnergyLimit[alphaPlus] = 10. * MeV;
126
127    G4DNACrossSectionDataSet* tableAlphaPlus = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, 
128                                                                              eV,
129                                                                              scaleFactor );
130    tableAlphaPlus->LoadData(fileAlphaPlus);
131    tableData[alphaPlus] = tableAlphaPlus;
132  }
133  else
134  {
135    G4Exception("G4CrossSectionIonisationRudd Constructor: alphaPlus is not defined");
136  }
137
138  if (heliumDef != 0)
139  {
140    helium = heliumDef->GetParticleName();
141    tableFile[helium] = fileHelium;
142
143    lowEnergyLimit[helium] = 0. * eV;
144    highEnergyLimit[helium] = 10. * MeV;
145
146    G4DNACrossSectionDataSet* tableHelium = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, 
147                                                                           eV,
148                                                                           scaleFactor );
149    tableHelium->LoadData(fileHelium);
150    tableData[helium] = tableHelium;
151    }
152  else
153  {
154    G4Exception("G4CrossSectionIonisationRudd Constructor: helium is not defined");
155  }
156}
157
158//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
159
160G4CrossSectionIonisationRudd::~G4CrossSectionIonisationRudd()
161{
162  std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
163  for (pos = tableData.begin(); pos != tableData.end(); ++pos)
164  {
165    G4DNACrossSectionDataSet* table = pos->second;
166    delete table;
167  }
168}
169
170//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
171
172G4double G4CrossSectionIonisationRudd::CrossSection(const G4Track& track )
173{
174  const G4DynamicParticle* particle = track.GetDynamicParticle();
175  G4double k = particle->GetKineticEnergy();
176 
177  G4DNAGenericIonsManager *instance;
178  instance = G4DNAGenericIonsManager::Instance();
179
180  if (
181      track.GetDefinition() != G4Proton::ProtonDefinition()
182      &&
183      track.GetDefinition() != instance->GetIon("hydrogen")
184      &&
185      track.GetDefinition() != instance->GetIon("alpha++")
186      &&
187      track.GetDefinition() != instance->GetIon("alpha+")
188      &&
189      track.GetDefinition() != instance->GetIon("helium")
190     )
191           
192    G4Exception("G4CrossSectionIonisationRudd: attempting to calculate cross section for wrong particle");
193
194  G4double sigma = 0.;
195
196  G4double lowLim = lowEnergyLimitDefault;
197  G4double highLim = highEnergyLimitDefault;
198
199  const G4String& particleName = particle->GetDefinition()->GetParticleName();
200
201  std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
202  pos1 = lowEnergyLimit.find(particleName);
203
204  if (pos1 != lowEnergyLimit.end())
205  {
206    lowLim = pos1->second;
207  }
208
209  std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
210  pos2 = highEnergyLimit.find(particleName);
211
212  if (pos2 != highEnergyLimit.end())
213  {
214    highLim = pos2->second;
215  }
216
217  if (k >= lowLim && k <= highLim)
218  {
219    std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
220    pos = tableData.find(particleName);
221       
222    if (pos != tableData.end())
223    {
224       G4DNACrossSectionDataSet* table = pos->second;
225       if (table != 0)
226       {
227              sigma = table->FindValue(k);
228
229              // BEGIN ELECTRON CORRECTION
230              // add ONE or TWO electron-water excitation for alpha+ and helium
231   
232              if ( particle->GetDefinition() == instance->GetIon("alpha+") 
233                   ||
234                   particle->GetDefinition() == instance->GetIon("helium")
235                   ) 
236              {
237     
238                  G4DNACrossSectionDataSet* electronDataset = new G4DNACrossSectionDataSet
239                    (new G4LogLogInterpolation, eV, (1./3.343e22)*m*m);
240       
241                  electronDataset->LoadData("dna/sigma_ionisation_e_born");
242
243                  G4double kElectron = k * 0.511/3728;
244
245                  if ( particle->GetDefinition() == instance->GetIon("alpha+") ) 
246                  {
247                      G4double tmp1 = table->FindValue(k) + electronDataset->FindValue(kElectron);
248                      delete electronDataset;
249                      return tmp1;
250                  }
251
252                  if ( particle->GetDefinition() == instance->GetIon("helium") ) 
253                  {
254                      G4double tmp2 = table->FindValue(k) +  2. * electronDataset->FindValue(kElectron);
255                      delete electronDataset;
256                      return tmp2;
257                  }
258              }     
259
260              // END ELECTRON CORRECTION
261       }
262    }
263    else
264    {
265        G4Exception("G4CrossSectionIonisationRudd: attempting to calculate cross section for wrong particle");
266    }
267  }
268
269  return sigma;
270}
271
Note: See TracBrowser for help on using the repository browser.