Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4UniversalFluctuation.cc
Go to the documentation of this file.
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: G4UniversalFluctuation.cc 79188 2014-02-20 09:22:48Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4UniversalFluctuation
34 //
35 // Author: Laszlo Urban
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 // 16-10-03 Changed interface to Initialisation (V.Ivanchenko)
45 // 07-11-03 Fix problem of rounding of double in G4UniversalFluctuations
46 // 06-02-04 Add control on big sigma > 2*meanLoss (V.Ivanchenko)
47 // 26-04-04 Comment out the case of very small step (V.Ivanchenko)
48 // 07-02-05 define problim = 5.e-3 (mma)
49 // 03-05-05 conditions of Gaussian fluctuation changed (bugfix)
50 // + smearing for very small loss (L.Urban)
51 // 03-10-05 energy dependent rate -> cut dependence of the
52 // distribution is much weaker (L.Urban)
53 // 17-10-05 correction for very small loss (L.Urban)
54 // 20-03-07 'GLANDZ' part rewritten completely, no 'very small loss'
55 // regime any more (L.Urban)
56 // 03-04-07 correction to get better width of eloss distr.(L.Urban)
57 // 13-07-07 add protection for very small step or low-density material (VI)
58 // 19-03-09 new width correction (does not depend on previous steps) (L.Urban)
59 // 20-03-09 modification in the width correction (L.Urban)
60 // 14-06-10 fixed tail distribution - do not use uniform function (L.Urban)
61 // 08-08-10 width correction algorithm has bee modified -->
62 // better results for thin targets (L.Urban)
63 // 06-02-11 correction for very small losses (L.Urban)
64 //
65 
66 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
67 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
68 
70 #include "G4PhysicalConstants.hh"
71 #include "G4SystemOfUnits.hh"
72 #include "Randomize.hh"
73 #include "G4Poisson.hh"
74 #include "G4Step.hh"
75 #include "G4Material.hh"
76 #include "G4MaterialCutsCouple.hh"
77 #include "G4DynamicParticle.hh"
78 #include "G4ParticleDefinition.hh"
79 #include "G4Log.hh"
80 #include "G4Exp.hh"
81 
82 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
83 
84 using namespace std;
85 
88  particle(0),
89  minNumberInteractionsBohr(10.0),
90  theBohrBeta2(50.0*keV/proton_mass_c2),
91  minLoss(10.*eV),
92  nmaxCont(16.),
93  rate(0.55),
94  fw(4.)
95 {
96  lastMaterial = 0;
97 
98  particleMass = chargeSquare = ipotFluct = electronDensity = f1Fluct = f2Fluct
99  = e1Fluct = e2Fluct = e1LogFluct = e2LogFluct = ipotLogFluct = e0 = esmall
100  = e1 = e2 = 0;
101  m_Inv_particleMass = m_massrate = DBL_MAX;
102 }
103 
104 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
105 
107 {}
108 
109 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
110 
112 {
113  particle = part;
114  particleMass = part->GetPDGMass();
115  G4double q = part->GetPDGCharge()/eplus;
116 
117  // Derived quantities
118  m_Inv_particleMass = 1.0 / particleMass;
119  m_massrate = electron_mass_c2 * m_Inv_particleMass ;
120  chargeSquare = q*q;
121 }
122 
123 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
124 
125 G4double
127  const G4DynamicParticle* dp,
128  G4double tmax,
129  G4double length,
130  G4double averageLoss)
131 {
132  // Calculate actual loss from the mean loss.
133  // The model used to get the fluctuations is essentially the same
134  // as in Glandz in Geant3 (Cern program library W5013, phys332).
135  // L. Urban et al. NIM A362, p.416 (1995) and Geant4 Physics Reference Manual
136 
137  // shortcut for very small loss or from a step nearly equal to the range
138  // (out of validity of the model)
139  //
140  G4double meanLoss = averageLoss;
141  G4double tkin = dp->GetKineticEnergy();
142  //G4cout<< "Emean= "<< meanLoss<< " tmax= "<< tmax<< " L= "<<length<<G4endl;
143  if (meanLoss < minLoss) { return meanLoss; }
144 
145  if(dp->GetDefinition() != particle) { InitialiseMe(dp->GetDefinition()); }
146 
147  G4double tau = tkin * m_Inv_particleMass;
148  G4double gam = tau + 1.0;
149  G4double gam2 = gam*gam;
150  G4double beta2 = tau*(tau + 2.0)/gam2;
151 
152  G4double loss(0.), siga(0.);
153 
154  const G4Material* material = couple->GetMaterial();
155 
156  // Gaussian regime
157  // for heavy particles only and conditions
158  // for Gauusian fluct. has been changed
159  //
160  if ((particleMass > electron_mass_c2) &&
161  (meanLoss >= minNumberInteractionsBohr*tmax))
162  {
163  G4double tmaxkine = 2.*electron_mass_c2*beta2*gam2/
164  (1.+m_massrate*(2.*gam+m_massrate)) ;
165  if (tmaxkine <= 2.*tmax)
166  {
167  electronDensity = material->GetElectronDensity();
168  siga = sqrt((1.0/beta2 - 0.5) * twopi_mc2_rcl2 * tmax * length
169  * electronDensity * chargeSquare);
170 
171  G4double sn = meanLoss/siga;
172 
173  // thick target case
174  if (sn >= 2.0) {
175 
176  G4double twomeanLoss = meanLoss + meanLoss;
177  do {
178  loss = G4RandGauss::shoot(meanLoss,siga);
179  } while (0.0 > loss || twomeanLoss < loss);
180 
181  // Gamma distribution
182  } else {
183 
184  G4double neff = sn*sn;
185  loss = meanLoss*G4RandGamma::shoot(neff,1.0)/neff;
186  }
187  //G4cout << "Gauss: " << loss << G4endl;
188  return loss;
189  }
190  }
191 
192  // Glandz regime : initialisation
193  //
194  if (material != lastMaterial) {
195  f1Fluct = material->GetIonisation()->GetF1fluct();
196  f2Fluct = material->GetIonisation()->GetF2fluct();
197  e1Fluct = material->GetIonisation()->GetEnergy1fluct();
198  e2Fluct = material->GetIonisation()->GetEnergy2fluct();
199  e1LogFluct = material->GetIonisation()->GetLogEnergy1fluct();
200  e2LogFluct = material->GetIonisation()->GetLogEnergy2fluct();
201  ipotFluct = material->GetIonisation()->GetMeanExcitationEnergy();
202  ipotLogFluct = material->GetIonisation()->GetLogMeanExcEnergy();
203  e0 = material->GetIonisation()->GetEnergy0fluct();
204  esmall = 0.5*sqrt(e0*ipotFluct);
205  lastMaterial = material;
206  }
207 
208  // very small step or low-density material
209  if(tmax <= e0) { return meanLoss; }
210 
211  G4double losstot = 0.;
212  G4int nstep = 1;
213  if(meanLoss < 25.*ipotFluct)
214  {
215  if(G4UniformRand()*ipotFluct< 0.04*meanLoss)
216  { nstep = 1; }
217  else
218  {
219  nstep = 2;
220  meanLoss *= 0.5;
221  }
222  }
223 
224  for (G4int istep=0; istep < nstep; ++istep) {
225 
226  loss = 0.;
227 
228  G4double a1 = 0. , a2 = 0., a3 = 0. ;
229 
230  if(tmax > ipotFluct) {
231  G4double w2 = G4Log(2.*electron_mass_c2*beta2*gam2)-beta2;
232 
233  if(w2 > ipotLogFluct) {
234  G4double C = meanLoss*(1.-rate)/(w2-ipotLogFluct);
235  a1 = C*f1Fluct*(w2-e1LogFluct)/e1Fluct;
236  if(w2 > e2LogFluct) {
237  a2 = C*f2Fluct*(w2-e2LogFluct)/e2Fluct;
238  }
239  if(a1 < nmaxCont) {
240  //small energy loss
241  G4double sa1 = sqrt(a1);
242  if(G4UniformRand() < G4Exp(-sa1))
243  {
244  e1 = esmall;
245  a1 = meanLoss*(1.-rate)/e1;
246  a2 = 0.;
247  e2 = e2Fluct;
248  }
249  else
250  {
251  a1 = sa1 ;
252  e1 = sa1*e1Fluct;
253  e2 = e2Fluct;
254  }
255 
256  } else {
257  //not small energy loss
258  //correction to get better fwhm value
259  a1 /= fw;
260  e1 = fw*e1Fluct;
261  e2 = e2Fluct;
262  }
263  }
264  }
265 
266  G4double w1 = tmax/e0;
267  if(tmax > e0) {
268  a3 = rate*meanLoss*(tmax-e0)/(e0*tmax*G4Log(w1));
269  if(a1+a2 <= 0.) {
270  a3 /= rate;
271  }
272  }
273  //'nearly' Gaussian fluctuation if a1>nmaxCont&&a2>nmaxCont&&a3>nmaxCont
274  G4double emean = 0.;
275  G4double sig2e = 0., sige = 0.;
276  G4double p1 = 0., p2 = 0., p3 = 0.;
277 
278  // excitation of type 1
279  if(a1 > nmaxCont)
280  {
281  emean += a1*e1;
282  sig2e += a1*e1*e1;
283  }
284  else if(a1 > 0.)
285  {
286  p1 = G4double(G4Poisson(a1));
287  loss += p1*e1;
288  if(p1 > 0.) {
289  loss += (1.-2.*G4UniformRand())*e1;
290  }
291  }
292 
293  // excitation of type 2
294  if(a2 > nmaxCont)
295  {
296  emean += a2*e2;
297  sig2e += a2*e2*e2;
298  }
299  else if(a2 > 0.)
300  {
301  p2 = G4double(G4Poisson(a2));
302  loss += p2*e2;
303  if(p2 > 0.)
304  loss += (1.-2.*G4UniformRand())*e2;
305  }
306  if(emean > 0.)
307  {
308  sige = sqrt(sig2e);
309  loss += max(0.,G4RandGauss::shoot(emean,sige));
310  }
311 
312  // ionisation
313  G4double lossc = 0.;
314  if(a3 > 0.) {
315  emean = 0.;
316  sig2e = 0.;
317  sige = 0.;
318  p3 = a3;
319  G4double alfa = 1.;
320  if(a3 > nmaxCont)
321  {
322  alfa = w1*(nmaxCont+a3)/(w1*nmaxCont+a3);
323  G4double alfa1 = alfa*G4Log(alfa)/(alfa-1.);
324  G4double namean = a3*w1*(alfa-1.)/((w1-1.)*alfa);
325  emean += namean*e0*alfa1;
326  sig2e += e0*e0*namean*(alfa-alfa1*alfa1);
327  p3 = a3-namean;
328  }
329 
330  G4double w2 = alfa*e0;
331  G4double w = (tmax-w2)/tmax;
332  G4int nb = G4Poisson(p3);
333  if(nb > 0) {
334  for (G4int k=0; k<nb; k++) { lossc += w2/(1.-w*G4UniformRand()); }
335  }
336 
337  if(emean > 0.)
338  {
339  sige = sqrt(sig2e);
340  lossc += max(0.,G4RandGauss::shoot(emean,sige));
341  }
342  }
343 
344  loss += lossc;
345 
346  losstot += loss;
347  }
348  //G4cout << "Vavilov: " << losstot << " Nstep= " << nstep << G4endl;
349 
350  return losstot;
351 
352 }
353 
354 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
355 
356 
358  const G4Material* material,
359  const G4DynamicParticle* dp,
360  G4double tmax,
361  G4double length)
362 {
363  if(dp->GetDefinition() != particle) { InitialiseMe(dp->GetDefinition()); }
364 
365  electronDensity = material->GetElectronDensity();
366 
367  G4double gam = (dp->GetKineticEnergy())*m_Inv_particleMass + 1.0;
368  G4double beta2 = 1.0 - 1.0/(gam*gam);
369 
370  G4double siga = (1.0/beta2 - 0.5) * twopi_mc2_rcl2 * tmax * length
371  * electronDensity * chargeSquare;
372 
373  return siga;
374 }
375 
376 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
377 
378 void
380  G4double q2)
381 {
382  if(part != particle) {
383  particle = part;
384  particleMass = part->GetPDGMass();
385 
386  // Derived quantities
387  if( particleMass != 0.0 ){
388  m_Inv_particleMass = 1.0 / particleMass;
389  m_massrate = electron_mass_c2 * m_Inv_particleMass ;
390  }else{
391  m_Inv_particleMass = DBL_MAX;
392  m_massrate = DBL_MAX;
393  }
394  }
395  chargeSquare = q2;
396 }
397 
398 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
ThreeVector shoot(const G4int Ap, const G4int Af)
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:224
G4long G4Poisson(G4double mean)
Definition: G4Poisson.hh:51
G4double GetKineticEnergy() const
G4double GetEnergy2fluct() const
G4double GetLogEnergy2fluct() const
G4UniversalFluctuation(const G4String &nam="UniFluc")
G4ParticleDefinition * GetDefinition() const
G4double GetLogMeanExcEnergy() const
int G4int
Definition: G4Types.hh:78
virtual G4double Dispersion(const G4Material *, const G4DynamicParticle *, G4double, G4double)
string material
Definition: eplot.py:19
G4double GetEnergy0fluct() const
#define G4UniformRand()
Definition: Randomize.hh:87
virtual void SetParticleAndCharge(const G4ParticleDefinition *, G4double q2)
G4double GetElectronDensity() const
Definition: G4Material.hh:215
float proton_mass_c2
Definition: hepunit.py:275
float electron_mass_c2
Definition: hepunit.py:274
virtual void InitialiseMe(const G4ParticleDefinition *)
virtual G4double SampleFluctuations(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double, G4double, G4double)
G4double G4Log(G4double x)
Definition: G4Log.hh:227
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
G4double GetLogEnergy1fluct() const
G4double GetPDGMass() const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4double GetMeanExcitationEnergy() const
G4double GetF2fluct() const
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
#define DBL_MAX
Definition: templates.hh:83
G4double GetF1fluct() const
const G4Material * GetMaterial() const
G4double GetEnergy1fluct() const