source: trunk/source/processes/electromagnetic/standard/src/G4IonFluctuations.cc @ 1340

Last change on this file since 1340 was 1340, checked in by garnier, 14 years ago

update ti head

File size: 14.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: G4IonFluctuations.cc,v 1.27 2010/10/25 19:13:23 vnivanch Exp $
27// GEANT4 tag $Name: emstand-V09-03-24 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class file
32//
33//
34// File name:     G4IonFluctuation
35//
36// Author:        Vladimir Ivanchenko
37//
38// Creation date: 03.01.2002
39//
40// Modifications:
41//
42// 28-12-02 add method Dispersion (V.Ivanchenko)
43// 07-02-03 change signature (V.Ivanchenko)
44// 13-02-03 Add name (V.Ivanchenko)
45// 23-05-03 Add control on parthalogical cases (V.Ivanchenko)
46// 16-10-03 Changed interface to Initialisation (V.Ivanchenko)
47// 27-09-07 Use FermiEnergy from material, add cut dependence (V.Ivanchenko)
48// 01-02-08 Add protection for small energies and optimise the code (V.Ivanchenko)
49// 01-06-08 Added initialisation of effective charge prestep (V.Ivanchenko)
50//
51// Class Description:
52//
53// -------------------------------------------------------------------
54//
55
56//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
57//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
58
59#include "G4IonFluctuations.hh"
60#include "Randomize.hh"
61#include "G4Poisson.hh"
62#include "G4Material.hh"
63#include "G4DynamicParticle.hh"
64
65//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
66
67using namespace std;
68
69G4IonFluctuations::G4IonFluctuations(const G4String& nam)
70  : G4VEmFluctuationModel(nam),
71    particle(0),
72    particleMass(proton_mass_c2),
73    charge(1.0),
74    chargeSquare(1.0),
75    effChargeSquare(1.0),
76    parameter(10.0*CLHEP::MeV/CLHEP::proton_mass_c2),
77    minNumberInteractionsBohr(0.0),
78    theBohrBeta2(50.0*keV/CLHEP::proton_mass_c2),
79    minFraction(0.2),
80    xmin(0.2),
81    minLoss(0.001*eV)
82{
83  kineticEnergy = 0.0;
84  beta2 = 0.0;
85}
86
87//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
88
89G4IonFluctuations::~G4IonFluctuations()
90{}
91
92//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
93
94void G4IonFluctuations::InitialiseMe(const G4ParticleDefinition* part)
95{
96  particle       = part;
97  particleMass   = part->GetPDGMass();
98  charge         = part->GetPDGCharge()/eplus;
99  chargeSquare   = charge*charge;
100  effChargeSquare= chargeSquare;
101  uniFluct.InitialiseMe(part);
102}
103
104//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
105
106G4double G4IonFluctuations::SampleFluctuations(const G4Material* material,
107                                               const G4DynamicParticle* dp,
108                                               G4double& tmax,
109                                               G4double& length,
110                                               G4double& meanLoss)
111{
112  //  G4cout << "### meanLoss= " << meanLoss << G4endl;
113  if(meanLoss <= minLoss) return meanLoss;
114
115  //G4cout << "G4IonFluctuations::SampleFluctuations E(MeV)= " << dp->GetKineticEnergy()
116  //     << "  Elim(MeV)= " << parameter*charge*particleMass << G4endl;
117
118  // Vavilov fluctuations
119  if(dp->GetKineticEnergy() > parameter*charge*particleMass) {
120    return uniFluct.SampleFluctuations(material,dp,tmax,length,meanLoss);
121  }
122
123  G4double siga = Dispersion(material,dp,tmax,length);
124  G4double loss = meanLoss;
125 
126  G4double navr = minNumberInteractionsBohr;
127  navr = meanLoss*meanLoss/siga;
128  //G4cout << "### siga= " << sqrt(siga) << "  navr= " << navr << G4endl;
129
130  // Gaussian fluctuation
131  if (navr >= minNumberInteractionsBohr) {
132
133    // Increase fluctuations for big fractional energy loss
134    //G4cout << "siga= " << siga << G4endl;
135    if ( meanLoss > minFraction*kineticEnergy ) {
136      G4double gam = (kineticEnergy - meanLoss)/particleMass + 1.0;
137      G4double b2  = 1.0 - 1.0/(gam*gam);
138      if(b2 < xmin*beta2) b2 = xmin*beta2;
139      G4double x   = b2/beta2;
140      G4double x3  = 1.0/(x*x*x);
141      siga *= 0.25*(1.0 + x)*(x3 + (1.0/b2 - 0.5)/(1.0/beta2 - 0.5) );
142    }
143    //       G4cout << "siga= " << siga << G4endl;
144    siga = sqrt(siga);
145    G4double lossmax = meanLoss+meanLoss;
146
147    if(siga > 5.0*meanLoss) {
148      loss = lossmax*G4UniformRand();
149    } else {
150      do {
151        loss = G4RandGauss::shoot(meanLoss,siga);
152      } while (0.0 > loss || loss > lossmax);
153    }
154  // Poisson fluctuations
155  } else {
156
157    G4double n    = (G4double)(G4Poisson(navr));
158    loss = meanLoss*n/navr;
159  }
160
161  //G4cout << "meanLoss= " << meanLoss << " loss= " << loss << G4endl;
162  return loss;
163}
164
165//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
166
167G4double G4IonFluctuations::Dispersion(const G4Material* material,
168                                       const G4DynamicParticle* dp,
169                                       G4double& tmax,
170                                       G4double& length)
171{
172  kineticEnergy = dp->GetKineticEnergy();
173  G4double etot = kineticEnergy + particleMass;
174  beta2 = kineticEnergy*(kineticEnergy + 2.*particleMass)/(etot*etot);
175
176  G4double electronDensity = material->GetElectronDensity();
177
178  /*
179  G4cout << "e= " <<  kineticEnergy << " m= " << particleMass
180         << " tmax= " << tmax << " l= " << length
181         << " q^2= " << effChargeSquare << " beta2=" << beta2<< G4endl;
182  */
183  G4double siga = (1. - beta2*0.5)*tmax*length*electronDensity*
184    twopi_mc2_rcl2*chargeSquare/beta2;
185
186  // Low velocity - additional ion charge fluctuations according to
187  // Q.Yang et al., NIM B61(1991)149-155.
188  //G4cout << "sigE= " << sqrt(siga) << " charge= " << charge <<G4endl;
189
190  G4double Z = electronDensity/material->GetTotNbOfAtomsPerVolume();
191 
192  G4double fac = Factor(material, Z);
193
194  // heavy ion correction
195//  G4double f1 = 1.065e-4*chargeSquare;
196//  if(beta2 > theBohrBeta2)  f1/= beta2;
197//  else                      f1/= theBohrBeta2;
198//  if(f1 > 2.5) f1 = 2.5;
199//  fac *= (1.0 + f1);
200
201  // taking into account the cut
202  G4double fac_cut = 1.0 + (fac - 1.0)*2.0*electron_mass_c2*beta2/(tmax*(1.0 - beta2));
203  if(fac_cut > 0.01 && fac > 0.01) {
204    siga *= fac_cut;
205  }
206
207  //G4cout << "siga(keV)= " << sqrt(siga)/keV << " fac= " << fac
208  //     << "  f1= " << f1 << G4endl;
209
210  return siga;
211}
212
213//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
214
215G4double G4IonFluctuations::Factor(const G4Material* material, G4double Z)
216{
217  // The aproximation of energy loss fluctuations
218  // Q.Yang et al., NIM B61(1991)149-155.
219
220  // Reduced energy in MeV/AMU
221  G4double energy = kineticEnergy *amu_c2/(particleMass*MeV) ;
222
223  // simple approximation for higher beta2
224  G4double s1 = RelativisticFactor(material, Z);
225
226  // tabulation for lower beta2
227  if( beta2 < 3.0*theBohrBeta2*Z ) {
228
229    static G4double a[96][4] = {
230 {-0.3291, -0.8312,  0.2460, -1.0220},
231 {-0.5615, -0.5898,  0.5205, -0.7258},
232 {-0.5280, -0.4981,  0.5519, -0.5865},
233 {-0.5125, -0.4625,  0.5660, -0.5190},
234 {-0.5127, -0.8595,  0.5626, -0.8721},
235 {-0.5174, -1.1930,  0.5565, -1.1980},
236 {-0.5179, -1.1850,  0.5560, -1.2070},
237 {-0.5209, -0.9355,  0.5590, -1.0250},
238 {-0.5255, -0.7766,  0.5720, -0.9412},
239
240 {-0.5776, -0.6665,  0.6598, -0.8484},
241 {-0.6013, -0.6045,  0.7321, -0.7671},
242 {-0.5781, -0.5518,  0.7605, -0.6919},
243 {-0.5587, -0.4981,  0.7835, -0.6195},
244 {-0.5466, -0.4656,  0.7978, -0.5771},
245 {-0.5406, -0.4690,  0.8031, -0.5718},
246 {-0.5391, -0.5061,  0.8024, -0.5974},
247 {-0.5380, -0.6483,  0.7962, -0.6970},
248 {-0.5355, -0.7722,  0.7962, -0.7839},
249 {-0.5329, -0.7720,  0.7988, -0.7846},
250
251 {-0.5335, -0.7671,  0.7984, -0.7933},
252 {-0.5324, -0.7612,  0.7998, -0.8031},
253 {-0.5305, -0.7300,  0.8031, -0.7990},
254 {-0.5307, -0.7178,  0.8049, -0.8216},
255 {-0.5248, -0.6621,  0.8165, -0.7919},
256 {-0.5180, -0.6502,  0.8266, -0.7986},
257 {-0.5084, -0.6408,  0.8396, -0.8048},
258 {-0.4967, -0.6331,  0.8549, -0.8093},
259 {-0.4861, -0.6508,  0.8712, -0.8432},
260 {-0.4700, -0.6186,  0.8961, -0.8132},
261
262 {-0.4545, -0.5720,  0.9227, -0.7710},
263 {-0.4404, -0.5226,  0.9481, -0.7254},
264 {-0.4288, -0.4778,  0.9701, -0.6850},
265 {-0.4199, -0.4425,  0.9874, -0.6539},
266 {-0.4131, -0.4188,  0.9998, -0.6332},
267 {-0.4089, -0.4057,  1.0070, -0.6218},
268 {-0.4039, -0.3913,  1.0150, -0.6107},
269 {-0.3987, -0.3698,  1.0240, -0.5938},
270 {-0.3977, -0.3608,  1.0260, -0.5852},
271 {-0.3972, -0.3600,  1.0260, -0.5842},
272
273 {-0.3985, -0.3803,  1.0200, -0.6013},
274 {-0.3985, -0.3979,  1.0150, -0.6168},
275 {-0.3968, -0.3990,  1.0160, -0.6195},
276 {-0.3971, -0.4432,  1.0050, -0.6591},
277 {-0.3944, -0.4665,  1.0010, -0.6825},
278 {-0.3924, -0.5109,  0.9921, -0.7235},
279 {-0.3882, -0.5158,  0.9947, -0.7343},
280 {-0.3838, -0.5125,  0.9999, -0.7370},
281 {-0.3786, -0.4976,  1.0090, -0.7310},
282 {-0.3741, -0.4738,  1.0200, -0.7155},
283
284 {-0.3969, -0.4496,  1.0320, -0.6982},
285 {-0.3663, -0.4297,  1.0430, -0.6828},
286 {-0.3630, -0.4120,  1.0530, -0.6689},
287 {-0.3597, -0.3964,  1.0620, -0.6564},
288 {-0.3555, -0.3809,  1.0720, -0.6454},
289 {-0.3525, -0.3607,  1.0820, -0.6289},
290 {-0.3505, -0.3465,  1.0900, -0.6171},
291 {-0.3397, -0.3570,  1.1020, -0.6384},
292 {-0.3314, -0.3552,  1.1130, -0.6441},
293 {-0.3235, -0.3531,  1.1230, -0.6498},
294
295 {-0.3150, -0.3483,  1.1360, -0.6539},
296 {-0.3060, -0.3441,  1.1490, -0.6593},
297 {-0.2968, -0.3396,  1.1630, -0.6649},
298 {-0.2935, -0.3225,  1.1760, -0.6527},
299 {-0.2797, -0.3262,  1.1940, -0.6722},
300 {-0.2704, -0.3202,  1.2100, -0.6770},
301 {-0.2815, -0.3227,  1.2480, -0.6775},
302 {-0.2880, -0.3245,  1.2810, -0.6801},
303 {-0.3034, -0.3263,  1.3270, -0.6778},
304 {-0.2936, -0.3215,  1.3430, -0.6835},
305
306 {-0.3282, -0.3200,  1.3980, -0.6650},
307 {-0.3260, -0.3070,  1.4090, -0.6552},
308 {-0.3511, -0.3074,  1.4470, -0.6442},
309 {-0.3501, -0.3064,  1.4500, -0.6442},
310 {-0.3490, -0.3027,  1.4550, -0.6418},
311 {-0.3487, -0.3048,  1.4570, -0.6447},
312 {-0.3478, -0.3074,  1.4600, -0.6483},
313 {-0.3501, -0.3283,  1.4540, -0.6669},
314 {-0.3494, -0.3373,  1.4550, -0.6765},
315 {-0.3485, -0.3373,  1.4570, -0.6774},
316
317 {-0.3462, -0.3300,  1.4630, -0.6728},
318 {-0.3462, -0.3225,  1.4690, -0.6662},
319 {-0.3453, -0.3094,  1.4790, -0.6553},
320 {-0.3844, -0.3134,  1.5240, -0.6412},
321 {-0.3848, -0.3018,  1.5310, -0.6303},
322 {-0.3862, -0.2955,  1.5360, -0.6237},
323 {-0.4262, -0.2991,  1.5860, -0.6115},
324 {-0.4278, -0.2910,  1.5900, -0.6029},
325 {-0.4303, -0.2817,  1.5940, -0.5927},
326 {-0.4315, -0.2719,  1.6010, -0.5829},
327
328 {-0.4359, -0.2914,  1.6050, -0.6010},
329 {-0.4365, -0.2982,  1.6080, -0.6080},
330 {-0.4253, -0.3037,  1.6120, -0.6150},
331 {-0.4335, -0.3245,  1.6160, -0.6377},
332 {-0.4307, -0.3292,  1.6210, -0.6447},
333 {-0.4284, -0.3204,  1.6290, -0.6380},
334 {-0.4227, -0.3217,  1.6360, -0.6438}
335    } ;
336
337    G4int iz = G4int(Z) - 2;
338    if( 0 > iz )      iz = 0;
339    else if(95 < iz ) iz = 95;
340
341    G4double ss = 1.0 + a[iz][0]*pow(energy,a[iz][1])+
342      + a[iz][2]*pow(energy,a[iz][3]);
343 
344    // protection for the validity range for low beta
345    G4double slim = 0.001;
346    if(ss < slim) s1 = 1.0/slim;
347    // for high value of beta
348    else if(s1*ss < 1.0) s1 = 1.0/ss;
349  }
350
351  G4int i = 0 ;
352  G4double factor = 1.0 ;
353
354  // The index of set of parameters i = 0 for protons(hadrons) in gases
355  //                                    1 for protons(hadrons) in solids
356  //                                    2 for ions in atomic gases
357  //                                    3 for ions in molecular gases
358  //                                    4 for ions in solids
359  static G4double b[5][4] = {
360  {0.1014,  0.3700,  0.9642,  3.987},
361  {0.1955,  0.6941,  2.522,   1.040},
362  {0.05058, 0.08975, 0.1419, 10.80},
363  {0.05009, 0.08660, 0.2751,  3.787},
364  {0.01273, 0.03458, 0.3951,  3.812}
365  } ;
366
367  // protons (hadrons)
368  if(1.5 > charge) {
369    if( kStateGas != material->GetState() ) i = 1 ;
370
371  // ions
372  } else {
373
374    factor = charge * pow(charge/Z, 0.33333333);
375
376    if( kStateGas == material->GetState() ) {
377      energy /= (charge * sqrt(charge)) ;
378
379      if(1 == (material->GetNumberOfElements())) {
380        i = 2 ;
381      } else {
382        i = 3 ;
383      }
384
385    } else {
386      energy /= (charge * sqrt(charge*Z)) ;
387      i = 4 ;
388    }
389  }
390
391  G4double x = b[i][2];
392  G4double y = energy * b[i][3];
393  if(y <= 0.2) x *= (y*(1.0 - 0.5*y));
394  else         x *= (1.0 - exp(-y));
395
396  y = energy - b[i][1];
397
398  G4double s2 = factor * x * b[i][0] / (y*y + x*x);
399  /* 
400  G4cout << "s1= " << s1 << " s2= " << s2 << " q^2= " << effChargeSquare
401         << " e= " << energy << G4endl;
402  */
403  return s1*effChargeSquare/chargeSquare + s2;
404}
405
406//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
407
408G4double G4IonFluctuations::RelativisticFactor(const G4Material* mat, 
409                                               G4double Z)
410{
411  G4double eF = mat->GetIonisation()->GetFermiEnergy();
412  G4double I  = mat->GetIonisation()->GetMeanExcitationEnergy();
413
414  // H.Geissel et al. NIM B, 195 (2002) 3.
415  G4double bF2= 2.0*eF/electron_mass_c2;
416  G4double f  = 0.4*(1.0 - beta2)/((1.0 - 0.5*beta2)*Z);
417  if(beta2 > bF2) f *= log(2.0*electron_mass_c2*beta2/I)*bF2/beta2;
418  else            f *= log(4.0*eF/I);
419
420  //  G4cout << "f= " << f << " beta2= " << beta2
421  //     << " bf2= " << bF2 << " q^2= " << chargeSquare << G4endl;
422
423  return 1.0 + f;
424}
425
426//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
427
428void G4IonFluctuations::SetParticleAndCharge(const G4ParticleDefinition* part,
429                                             G4double q2)
430{
431  if(part != particle) {
432    particle       = part;
433    particleMass   = part->GetPDGMass();
434    charge         = part->GetPDGCharge()/eplus;
435    chargeSquare   = charge*charge;
436  }
437  effChargeSquare  = q2;
438  uniFluct.SetParticleAndCharge(part, q2);
439}
440
441//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
Note: See TracBrowser for help on using the repository browser.