source: trunk/source/processes/electromagnetic/lowenergy/src/G4CrossSectionIonisationBornPartial.cc@ 968

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

update processes

File size: 6.4 KB
RevLine 
[819]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//
[961]26// $Id: G4CrossSectionIonisationBornPartial.cc,v 1.5 2009/01/20 07:40:53 sincerti Exp $
27// GEANT4 tag $Name: geant4-09-02-ref-02 $
[819]28
29#include "G4CrossSectionIonisationBornPartial.hh"
30
[961]31//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
[819]32
33G4CrossSectionIonisationBornPartial::G4CrossSectionIonisationBornPartial()
34{
[961]35 lowEnergyLimitDefault = 12.61 * eV;
[819]36 highEnergyLimitDefault = 30 * keV;
37
38 G4String fileElectron("dna/sigma_ionisation_e_born");
39 G4String fileProton("dna/sigma_ionisation_p_born");
40
41 G4ParticleDefinition* electronDef = G4Electron::ElectronDefinition();
42 G4ParticleDefinition* protonDef = G4Proton::ProtonDefinition();
43
44 G4String electron;
45 G4String proton;
46
47 G4double scaleFactor = (1.e-22 / 3.343) * m*m;
48
49 if (electronDef != 0)
[961]50 {
51 electron = electronDef->GetParticleName();
52 tableFile[electron] = fileElectron;
[819]53
[961]54 lowEnergyLimit[electron] = 12.61 * eV;
55 highEnergyLimit[electron] = 30. * keV;
[819]56
[961]57 G4DNACrossSectionDataSet* tableE = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, eV,scaleFactor );
58 tableE->LoadData(fileElectron);
[819]59
[961]60 tableData[electron] = tableE;
61 }
[819]62 else
[961]63 {
64 G4Exception("G4CrossSectionIonisationBornPartial Constructor: electron is not defined");
65 }
[819]66
67 if (protonDef != 0)
[961]68 {
69 proton = protonDef->GetParticleName();
70 tableFile[proton] = fileProton;
[819]71
[961]72 lowEnergyLimit[proton] = 500. * keV;
73 highEnergyLimit[proton] = 10. * MeV;
[819]74
[961]75 G4DNACrossSectionDataSet* tableP = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, eV,scaleFactor );
76 tableP->LoadData(fileProton);
[819]77
[961]78 tableData[proton] = tableP;
79 }
[819]80 else
[961]81 {
82 G4Exception("G4CrossSectionIonisationBornPartial Constructor: proton is not defined");
83 }
[819]84}
85
[961]86//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
[819]87
88G4CrossSectionIonisationBornPartial::~G4CrossSectionIonisationBornPartial()
89{
90 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
91 for (pos = tableData.begin(); pos != tableData.end(); ++pos)
[961]92 {
93 G4DNACrossSectionDataSet* table = pos->second;
94 delete table;
95 }
[819]96}
97
[961]98//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
[819]99
100G4int G4CrossSectionIonisationBornPartial::RandomSelect(G4double k, const G4String& particle )
101{
102 G4int level = 0;
103
[961]104 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
105 pos = tableData.find(particle);
[819]106
[961]107 if (pos != tableData.end())
108 {
109 G4DNACrossSectionDataSet* table = pos->second;
[819]110
[961]111 if (table != 0)
112 {
113 G4double* valuesBuffer = new G4double[table->NumberOfComponents()];
114 const size_t n(table->NumberOfComponents());
115 size_t i(n);
116 G4double value = 0.;
[819]117
[961]118 while (i>0)
119 {
120 i--;
121 valuesBuffer[i] = table->GetComponent(i)->FindValue(k);
122 value += valuesBuffer[i];
123 }
[819]124
[961]125 value *= G4UniformRand();
126
127 i = n;
[819]128
[961]129 while (i > 0)
130 {
131 i--;
[819]132
[961]133 if (valuesBuffer[i] > value)
134 {
135 delete[] valuesBuffer;
136 return i;
137 }
138 value -= valuesBuffer[i];
139 }
[819]140
[961]141 if (valuesBuffer) delete[] valuesBuffer;
142
143 }
144 }
145 else
146 {
147 G4Exception("G4CrossSectionIonisationBornPartial: attempting to calculate cross section for wrong particle");
148 }
[819]149
150 return level;
151}
152
[961]153//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
[819]154
155G4double G4CrossSectionIonisationBornPartial::CrossSection(const G4Track& track )
156{
157 G4double sigma = 0.;
158
159 const G4DynamicParticle* particle = track.GetDynamicParticle();
160 G4double k = particle->GetKineticEnergy();
161
162 G4double lowLim = lowEnergyLimitDefault;
163 G4double highLim = highEnergyLimitDefault;
164
165 const G4String& particleName = particle->GetDefinition()->GetParticleName();
166
[961]167 std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
168 pos1 = lowEnergyLimit.find(particleName);
[819]169
[961]170 if (pos1 != lowEnergyLimit.end())
171 {
172 lowLim = pos1->second;
173 }
[819]174
[961]175 std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
176 pos2 = highEnergyLimit.find(particleName);
[819]177
[961]178 if (pos2 != highEnergyLimit.end())
179 {
180 highLim = pos2->second;
181 }
[819]182
[961]183 if (k > lowLim && k < highLim)
184 {
185 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
186 pos = tableData.find(particleName);
[819]187
[961]188 if (pos != tableData.end())
189 {
190 G4DNACrossSectionDataSet* table = pos->second;
191 if (table != 0)
192 {
193 sigma = table->FindValue(k);
194 }
195 }
196 else
197 {
198 G4Exception("G4CrossSectionIonisationBornPartial: attempting to calculate cross section for wrong particle");
199 }
200 }
[819]201
[961]202 return sigma;
[819]203}
204
[961]205//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
[819]206
207G4double G4CrossSectionIonisationBornPartial::Sum(G4double /* energy */, const G4String& /* particle */)
208{
209 return 0;
210}
Note: See TracBrowser for help on using the repository browser.