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