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

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

update

  • Property svn:executable set to *
File size: 10.3 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.2 2008/07/14 20:47:34 sincerti Exp $
27// GEANT4 tag $Name: geant4-09-02 $
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  G4int z = 10;
118
119  G4double numerator;
120  numerator = std::pow(z * aj[excitationLevel], omegaj[excitationLevel]) * 
121    std::pow(tCorrected - waterExcitation.ExcitationEnergy(excitationLevel), nu);
122
123  G4double power;
124  power = omegaj[excitationLevel] + nu;
125
126  G4double denominator;
127  denominator = std::pow(jj[excitationLevel], power) + std::pow(tCorrected, power);
128
129  G4double zEff = particleDefinition->GetPDGCharge() / eplus + particleDefinition->GetLeptonNumber();
130
131  zEff -= ( sCoefficient[0][particleTypeIndex] * S_1s(k, waterExcitation.ExcitationEnergy(excitationLevel), slaterEffectiveCharge[0][particleTypeIndex], 1.) +
132            sCoefficient[1][particleTypeIndex] * S_2s(k, waterExcitation.ExcitationEnergy(excitationLevel), slaterEffectiveCharge[1][particleTypeIndex], 2.) +
133            sCoefficient[2][particleTypeIndex] * S_2p(k, waterExcitation.ExcitationEnergy(excitationLevel), slaterEffectiveCharge[2][particleTypeIndex], 2.) );
134
135  G4double cross = sigma0 * zEff * zEff * numerator / denominator;
136
137  return cross;
138}
139
140//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
141
142G4int G4CrossSectionExcitationMillerGreenPartial::RandomSelect(G4double k, 
143                                                               const G4ParticleDefinition* particle)
144{
145  G4int i = nLevels;
146  G4double value = 0.;
147  std::deque<double> values;
148 
149  // ---- MGP ---- The following algorithm is wrong: it works if the cross section
150  // is a monotone increasing function.
151  // The algorithm should be corrected by building the cumulative function
152  // of the cross section and comparing a random number in the range 0-1 against
153  // the cumulative value at each bin
154
155  G4DNAGenericIonsManager *instance;
156  instance = G4DNAGenericIonsManager::Instance();
157
158  // ELECTRON CORRECTION
159 
160  if ( particle == instance->GetIon("alpha++"))
161  {  while (i > 0)
162     {
163        i--;
164        G4double partial = CrossSection(k,i,particle);
165        values.push_front(partial);
166        value += partial;
167     }
168
169     value *= G4UniformRand();
170   
171     i = nLevels;
172
173     while (i > 0)
174     {
175        i--;
176        if (values[i] > value) return i;
177        value -= values[i];
178     }
179  } 
180
181  // add ONE or TWO electron-water excitation for alpha+ and helium
182   
183  if ( particle == instance->GetIon("alpha+") 
184       ||
185       particle == instance->GetIon("helium")
186     ) 
187  {
188    while (i>0)
189    {
190          i--;
191         
192          G4CrossSectionExcitationEmfietzoglouPartial* excitationXS = 
193            new G4CrossSectionExcitationEmfietzoglouPartial();
194         
195          G4double sigmaExcitation=0;
196          if (k*0.511/3728 > 7.4*eV && k*0.511/3728 < 10*keV) sigmaExcitation = excitationXS->CrossSection(k*0.511/3728,i);
197 
198          G4double partial = CrossSection(k,i,particle);
199          if (particle == instance->GetIon("alpha+")) partial = CrossSection(k,i,particle) + sigmaExcitation;
200          if (particle == instance->GetIon("helium")) partial = CrossSection(k,i,particle) + 2*sigmaExcitation;
201          values.push_front(partial);
202          value += partial;
203          delete excitationXS;
204    }
205 
206    value*=G4UniformRand();
207 
208    i=5;
209    while (i>0)
210    {
211          i--;
212   
213          if (values[i]>value) return i;
214 
215          value-=values[i];
216    }
217  }   
218  //   
219
220  return 0;
221}
222
223//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
224
225G4double G4CrossSectionExcitationMillerGreenPartial::Sum(G4double k, const G4ParticleDefinition* particle)
226{
227  G4double totalCrossSection = 0.;
228
229  for (G4int i=0; i<nLevels; i++)
230  {
231    totalCrossSection += CrossSection(k,i,particle);
232  }
233  return totalCrossSection;
234}
235
236//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
237
238G4double G4CrossSectionExcitationMillerGreenPartial::S_1s(G4double t,
239                                                          G4double energyTransferred,
240                                                          G4double slaterEffectiveCharge,
241                                                          G4double shellNumber)
242{
243  // 1 - e^(-2r) * ( 1 + 2 r + 2 r^2)
244  // Dingfelder, in Chattanooga 2005 proceedings, formula (7)
245 
246  G4double r = R(t, energyTransferred, slaterEffectiveCharge, shellNumber);
247  G4double value = 1. - std::exp(-2 * r) * ( ( 2. * r + 2. ) * r + 1. );
248 
249  return value;
250}
251
252
253//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
254
255G4double G4CrossSectionExcitationMillerGreenPartial::S_2s(G4double t,
256                                                          G4double energyTransferred,
257                                                          G4double slaterEffectiveCharge,
258                                                          G4double shellNumber)
259{
260  // 1 - e^(-2 r) * ( 1 + 2 r + 2 r^2 + 2 r^4)
261  // Dingfelder, in Chattanooga 2005 proceedings, formula (8)
262
263  G4double r = R(t, energyTransferred, slaterEffectiveCharge, shellNumber);
264  G4double value =  1. - std::exp(-2 * r) * (((2. * r * r + 2.) * r + 2.) * r + 1.);
265
266  return value;
267 
268}
269
270//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
271
272G4double G4CrossSectionExcitationMillerGreenPartial::S_2p(G4double t,
273                                                          G4double energyTransferred,
274                                                          G4double slaterEffectiveCharge,
275                                                          G4double shellNumber)
276{
277  // 1 - e^(-2 r) * ( 1 + 2 r + 2 r^2 + 4/3 r^3 + 2/3 r^4)
278  // Dingfelder, in Chattanooga 2005 proceedings, formula (9)
279
280  G4double r = R(t, energyTransferred, slaterEffectiveCharge, shellNumber);
281  G4double value =  1. - std::exp(-2 * r) * (((( 2./3. * r + 4./3.) * r + 2.) * r + 2.) * r  + 1.);
282
283  return value;
284}
285
286//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
287
288G4double G4CrossSectionExcitationMillerGreenPartial::R(G4double t,
289                                                       G4double energyTransferred,
290                                                       G4double slaterEffectiveCharge,
291                                                       G4double shellNumber) 
292{
293  // tElectron = m_electron / m_alpha * t
294  // Hardcoded in Riccardo's implementation; to be corrected
295  // Dingfelder, in Chattanooga 2005 proceedings, p 4
296
297  G4double tElectron = 0.511/3728. * t;
298  G4double value = 2. * tElectron * slaterEffectiveCharge / (energyTransferred * shellNumber);
299 
300  return value;
301}
302
Note: See TracBrowser for help on using the repository browser.