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

Last change on this file since 1197 was 1196, checked in by garnier, 16 years ago

update CVS release candidate geant4.9.3.01

  • Property svn:executable set to *
File size: 9.5 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.7 2009/06/11 15:47:08 mantero Exp $
27// GEANT4 tag $Name: geant4-09-03-cand-01 $
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 G4cout << G4endl;
158 G4cout << "*******************************************************************************" << G4endl;
159 G4cout << "*******************************************************************************" << G4endl;
160 G4cout << " The class G4CrossSectionIonisationRudd is NOT SUPPORTED ANYMORE. " << G4endl;
161 G4cout << " It will be REMOVED with the next major release of Geant4. " << G4endl;
162 G4cout << " Please consult: https://twiki.cern.ch/twiki/bin/view/Geant4/LoweProcesses" << G4endl;
163 G4cout << "*******************************************************************************" << G4endl;
164 G4cout << "*******************************************************************************" << G4endl;
165 G4cout << G4endl;
166
167}
168
169//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
170
171G4CrossSectionIonisationRudd::~G4CrossSectionIonisationRudd()
172{
173 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
174 for (pos = tableData.begin(); pos != tableData.end(); ++pos)
175 {
176 G4DNACrossSectionDataSet* table = pos->second;
177 delete table;
178 }
179}
180
181//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
182
183G4double G4CrossSectionIonisationRudd::CrossSection(const G4Track& track )
184{
185 const G4DynamicParticle* particle = track.GetDynamicParticle();
186 G4double k = particle->GetKineticEnergy();
187
188 G4DNAGenericIonsManager *instance;
189 instance = G4DNAGenericIonsManager::Instance();
190
191 if (
192 track.GetDefinition() != G4Proton::ProtonDefinition()
193 &&
194 track.GetDefinition() != instance->GetIon("hydrogen")
195 &&
196 track.GetDefinition() != instance->GetIon("alpha++")
197 &&
198 track.GetDefinition() != instance->GetIon("alpha+")
199 &&
200 track.GetDefinition() != instance->GetIon("helium")
201 )
202
203 G4Exception("G4CrossSectionIonisationRudd: attempting to calculate cross section for wrong particle");
204
205 G4double sigma = 0.;
206
207 G4double lowLim = lowEnergyLimitDefault;
208 G4double highLim = highEnergyLimitDefault;
209
210 const G4String& particleName = particle->GetDefinition()->GetParticleName();
211
212 std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
213 pos1 = lowEnergyLimit.find(particleName);
214
215 if (pos1 != lowEnergyLimit.end())
216 {
217 lowLim = pos1->second;
218 }
219
220 std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
221 pos2 = highEnergyLimit.find(particleName);
222
223 if (pos2 != highEnergyLimit.end())
224 {
225 highLim = pos2->second;
226 }
227
228 if (k >= lowLim && k <= highLim)
229 {
230 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
231 pos = tableData.find(particleName);
232
233 if (pos != tableData.end())
234 {
235 G4DNACrossSectionDataSet* table = pos->second;
236 if (table != 0)
237 {
238 sigma = table->FindValue(k);
239
240 // BEGIN ELECTRON CORRECTION
241 // add ONE or TWO electron-water excitation for alpha+ and helium
242
243 if ( particle->GetDefinition() == instance->GetIon("alpha+")
244 ||
245 particle->GetDefinition() == instance->GetIon("helium")
246 )
247 {
248
249 G4DNACrossSectionDataSet* electronDataset = new G4DNACrossSectionDataSet
250 (new G4LogLogInterpolation, eV, (1./3.343e22)*m*m);
251
252 electronDataset->LoadData("dna/sigma_ionisation_e_born");
253
254 G4double kElectron = k * 0.511/3728;
255
256 if ( particle->GetDefinition() == instance->GetIon("alpha+") )
257 {
258 G4double tmp1 = table->FindValue(k) + electronDataset->FindValue(kElectron);
259 delete electronDataset;
260 return tmp1;
261 }
262
263 if ( particle->GetDefinition() == instance->GetIon("helium") )
264 {
265 G4double tmp2 = table->FindValue(k) + 2. * electronDataset->FindValue(kElectron);
266 delete electronDataset;
267 return tmp2;
268 }
269 }
270
271 // END ELECTRON CORRECTION
272 }
273 }
274 else
275 {
276 G4Exception("G4CrossSectionIonisationRudd: attempting to calculate cross section for wrong particle");
277 }
278 }
279
280 return sigma;
281}
282
Note: See TracBrowser for help on using the repository browser.