source: trunk/source/processes/electromagnetic/lowenergy/src/G4CrossSectionExcitationMillerGreenPartial.cc@ 1199

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

update CVS release candidate geant4.9.3.01

  • Property svn:executable set to *
File size: 10.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: G4CrossSectionExcitationMillerGreenPartial.cc,v 1.4 2009/06/10 13:32:36 mantero Exp $
27// GEANT4 tag $Name: geant4-09-03-cand-01 $
28
29#include "G4CrossSectionExcitationMillerGreenPartial.hh"
30
31//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
32
33G4CrossSectionExcitationMillerGreenPartial::G4CrossSectionExcitationMillerGreenPartial()
34{
35 nLevels = waterExcitation.NumberOfLevels();
36
37 //PROTON
38 kineticEnergyCorrection[0] = 1.;
39 slaterEffectiveCharge[0][0] = 0.;
40 slaterEffectiveCharge[1][0] = 0.;
41 slaterEffectiveCharge[2][0] = 0.;
42 sCoefficient[0][0] = 0.;
43 sCoefficient[1][0] = 0.;
44 sCoefficient[2][0] = 0.;
45
46 //ALPHA++
47 kineticEnergyCorrection[1] = 0.9382723/3.727417;
48 slaterEffectiveCharge[0][1]=0.;
49 slaterEffectiveCharge[1][1]=0.;
50 slaterEffectiveCharge[2][1]=0.;
51 sCoefficient[0][1]=0.;
52 sCoefficient[1][1]=0.;
53 sCoefficient[2][1]=0.;
54
55 // ALPHA+
56 kineticEnergyCorrection[2] = 0.9382723/3.727417;
57 slaterEffectiveCharge[0][2]=2.0;
58 slaterEffectiveCharge[1][2]=1.15;
59 slaterEffectiveCharge[2][2]=1.15;
60 sCoefficient[0][2]=0.7;
61 sCoefficient[1][2]=0.15;
62 sCoefficient[2][2]=0.15;
63
64 // HELIUM
65 kineticEnergyCorrection[3] = 0.9382723/3.727417;
66 slaterEffectiveCharge[0][3]=1.7;
67 slaterEffectiveCharge[1][3]=1.15;
68 slaterEffectiveCharge[2][3]=1.15;
69 sCoefficient[0][3]=0.5;
70 sCoefficient[1][3]=0.25;
71 sCoefficient[2][3]=0.25;
72
73}
74
75//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
76
77G4CrossSectionExcitationMillerGreenPartial::~G4CrossSectionExcitationMillerGreenPartial()
78{}
79
80//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
81
82G4double G4CrossSectionExcitationMillerGreenPartial::CrossSection(G4double k, G4int excitationLevel,
83 const G4ParticleDefinition* particleDefinition)
84{
85 // ( ( z * aj ) ^ omegaj ) * ( t - ej ) ^ nu
86 // sigma(t) = zEff^2 * sigma0 * --------------------------------------------
87 // jj ^ ( omegaj + nu ) + t ^ ( omegaj + nu )
88 //
89 // where t is the kinetic energy corrected by Helium mass over proton mass for Helium ions
90 //
91 // zEff is:
92 // 1 for protons
93 // 2 for alpha++
94 // and 2 - c1 S_1s - c2 S_2s - c3 S_2p for alpha+ and He
95 //
96 // Dingfelder et al., RPC 59, 255-275, 2000 from Miller and Green (1973)
97 // Formula (34) and Table 2
98
99 const G4double sigma0(1.E+8 * barn);
100 const G4double nu(1.);
101 const G4double aj[]={876.*eV, 2084.* eV, 1373.*eV, 692.*eV, 900.*eV};
102 const G4double jj[]={19820.*eV, 23490.*eV, 27770.*eV, 30830.*eV, 33080.*eV};
103 const G4double omegaj[]={0.85, 0.88, 0.88, 0.78, 0.78};
104
105 G4int particleTypeIndex = 0;
106 G4DNAGenericIonsManager* instance;
107 instance = G4DNAGenericIonsManager::Instance();
108
109 if (particleDefinition == G4Proton::ProtonDefinition()) particleTypeIndex=0;
110 if (particleDefinition == instance->GetIon("alpha++")) particleTypeIndex=1;
111 if (particleDefinition == instance->GetIon("alpha+")) particleTypeIndex=2;
112 if (particleDefinition == instance->GetIon("helium")) particleTypeIndex=3;
113
114 G4double tCorrected;
115 tCorrected = k * kineticEnergyCorrection[particleTypeIndex];
116
117 // SI - added protection
118 if (tCorrected < waterExcitation.ExcitationEnergy(excitationLevel)) return 0;
119 //
120
121 G4int z = 10;
122
123 G4double numerator;
124 numerator = std::pow(z * aj[excitationLevel], omegaj[excitationLevel]) *
125 std::pow(tCorrected - waterExcitation.ExcitationEnergy(excitationLevel), nu);
126
127 G4double power;
128 power = omegaj[excitationLevel] + nu;
129
130 G4double denominator;
131 denominator = std::pow(jj[excitationLevel], power) + std::pow(tCorrected, power);
132
133 G4double zEff = particleDefinition->GetPDGCharge() / eplus + particleDefinition->GetLeptonNumber();
134
135 zEff -= ( sCoefficient[0][particleTypeIndex] * S_1s(k, waterExcitation.ExcitationEnergy(excitationLevel), slaterEffectiveCharge[0][particleTypeIndex], 1.) +
136 sCoefficient[1][particleTypeIndex] * S_2s(k, waterExcitation.ExcitationEnergy(excitationLevel), slaterEffectiveCharge[1][particleTypeIndex], 2.) +
137 sCoefficient[2][particleTypeIndex] * S_2p(k, waterExcitation.ExcitationEnergy(excitationLevel), slaterEffectiveCharge[2][particleTypeIndex], 2.) );
138
139 G4double cross = sigma0 * zEff * zEff * numerator / denominator;
140
141 return cross;
142}
143
144//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
145
146G4int G4CrossSectionExcitationMillerGreenPartial::RandomSelect(G4double k,
147 const G4ParticleDefinition* particle)
148{
149 G4int i = nLevels;
150 G4double value = 0.;
151 std::deque<double> values;
152
153 // ---- MGP ---- The following algorithm is wrong: it works if the cross section
154 // is a monotone increasing function.
155 // The algorithm should be corrected by building the cumulative function
156 // of the cross section and comparing a random number in the range 0-1 against
157 // the cumulative value at each bin
158
159 G4DNAGenericIonsManager *instance;
160 instance = G4DNAGenericIonsManager::Instance();
161
162 // ELECTRON CORRECTION
163
164 if ( particle == instance->GetIon("alpha++")||
165 particle == G4Proton::ProtonDefinition())
166
167 { while (i > 0)
168 {
169 i--;
170 G4double partial = CrossSection(k,i,particle);
171 values.push_front(partial);
172 value += partial;
173 }
174
175 value *= G4UniformRand();
176
177 i = nLevels;
178
179 while (i > 0)
180 {
181 i--;
182 if (values[i] > value) return i;
183 value -= values[i];
184 }
185 }
186
187 // add ONE or TWO electron-water excitation for alpha+ and helium
188
189 if ( particle == instance->GetIon("alpha+")
190 ||
191 particle == instance->GetIon("helium")
192 )
193 {
194 while (i>0)
195 {
196 i--;
197
198 G4CrossSectionExcitationEmfietzoglouPartial* excitationXS =
199 new G4CrossSectionExcitationEmfietzoglouPartial();
200
201 G4double sigmaExcitation=0;
202 if (k*0.511/3728 > 7.4*eV && k*0.511/3728 < 10*keV) sigmaExcitation = excitationXS->CrossSection(k*0.511/3728,i);
203
204 G4double partial = CrossSection(k,i,particle);
205 if (particle == instance->GetIon("alpha+")) partial = CrossSection(k,i,particle) + sigmaExcitation;
206 if (particle == instance->GetIon("helium")) partial = CrossSection(k,i,particle) + 2*sigmaExcitation;
207 values.push_front(partial);
208 value += partial;
209 delete excitationXS;
210 }
211
212 value*=G4UniformRand();
213
214 i=5;
215 while (i>0)
216 {
217 i--;
218
219 if (values[i]>value) return i;
220
221 value-=values[i];
222 }
223 }
224 //
225
226 return 0;
227}
228
229//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
230
231G4double G4CrossSectionExcitationMillerGreenPartial::Sum(G4double k, const G4ParticleDefinition* particle)
232{
233 G4double totalCrossSection = 0.;
234
235 for (G4int i=0; i<nLevels; i++)
236 {
237 totalCrossSection += CrossSection(k,i,particle);
238 }
239 return totalCrossSection;
240}
241
242//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
243
244G4double G4CrossSectionExcitationMillerGreenPartial::S_1s(G4double t,
245 G4double energyTransferred,
246 G4double slaterEffectiveCharge,
247 G4double shellNumber)
248{
249 // 1 - e^(-2r) * ( 1 + 2 r + 2 r^2)
250 // Dingfelder, in Chattanooga 2005 proceedings, formula (7)
251
252 G4double r = R(t, energyTransferred, slaterEffectiveCharge, shellNumber);
253 G4double value = 1. - std::exp(-2 * r) * ( ( 2. * r + 2. ) * r + 1. );
254
255 return value;
256}
257
258
259//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
260
261G4double G4CrossSectionExcitationMillerGreenPartial::S_2s(G4double t,
262 G4double energyTransferred,
263 G4double slaterEffectiveCharge,
264 G4double shellNumber)
265{
266 // 1 - e^(-2 r) * ( 1 + 2 r + 2 r^2 + 2 r^4)
267 // Dingfelder, in Chattanooga 2005 proceedings, formula (8)
268
269 G4double r = R(t, energyTransferred, slaterEffectiveCharge, shellNumber);
270 G4double value = 1. - std::exp(-2 * r) * (((2. * r * r + 2.) * r + 2.) * r + 1.);
271
272 return value;
273
274}
275
276//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
277
278G4double G4CrossSectionExcitationMillerGreenPartial::S_2p(G4double t,
279 G4double energyTransferred,
280 G4double slaterEffectiveCharge,
281 G4double shellNumber)
282{
283 // 1 - e^(-2 r) * ( 1 + 2 r + 2 r^2 + 4/3 r^3 + 2/3 r^4)
284 // Dingfelder, in Chattanooga 2005 proceedings, formula (9)
285
286 G4double r = R(t, energyTransferred, slaterEffectiveCharge, shellNumber);
287 G4double value = 1. - std::exp(-2 * r) * (((( 2./3. * r + 4./3.) * r + 2.) * r + 2.) * r + 1.);
288
289 return value;
290}
291
292//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
293
294G4double G4CrossSectionExcitationMillerGreenPartial::R(G4double t,
295 G4double energyTransferred,
296 G4double slaterEffectiveCharge,
297 G4double shellNumber)
298{
299 // tElectron = m_electron / m_alpha * t
300 // Hardcoded in Riccardo's implementation; to be corrected
301 // Dingfelder, in Chattanooga 2005 proceedings, p 4
302
303 G4double tElectron = 0.511/3728. * t;
304 G4double value = 2. * tElectron * slaterEffectiveCharge / (energyTransferred * shellNumber);
305
306 return value;
307}
308
Note: See TracBrowser for help on using the repository browser.