CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SiG4UniversalFluctuation.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * DISCLAIMER *
4 // * *
5 // * The following disclaimer summarizes all the specific disclaimers *
6 // * of contributors to this software. The specific disclaimers,which *
7 // * govern, are listed with their locations in: *
8 // * http://cern.ch/geant4/license *
9 // * *
10 // * Neither the authors of this software system, nor their employing *
11 // * institutes,nor the agencies providing financial support for this *
12 // * work make any representation or warranty, express or implied, *
13 // * regarding this software system or assume any liability for its *
14 // * use. *
15 // * *
16 // * This code implementation is the intellectual property of the *
17 // * GEANT4 collaboration. *
18 // * By copying, distributing or modifying the Program (or any work *
19 // * based on the Program) you indicate your acceptance of this *
20 // * statement, and all its terms. *
21 // ********************************************************************
22 //
23 // GEANT4 tag $Name: $
24 //
25 // -------------------------------------------------------------------
26 //
27 // GEANT4 Class file
28 //
29 //
30 // File name: G4UniversalFluctuation
31 //
32 // Author: Vladimir Ivanchenko
33 //
34 // Creation date: 03.01.2002
35 //
36 // Modifications:
37 //
38 // 28-12-02 add method Dispersion (V.Ivanchenko)
39 // 07-02-03 change signature (V.Ivanchenko)
40 // 13-02-03 Add name (V.Ivanchenko)
41 // 16-10-03 Changed interface to Initialisation (V.Ivanchenko)
42 // 07-11-03 Fix problem of rounding of double in G4UniversalFluctuations
43 // 06-02-04 Add control on big sigma > 2*meanLoss (V.Ivanchenko)
44 // 26-04-04 Comment out the case of very small step (V.Ivanchenko)
45 // 07-02-05 define problim = 5.e-3 (mma)
46 // 03-05-05 conditions of Gaussian fluctuation changed (bugfix)
47 // + smearing for very small loss (L.Urban)
48 //
49 // Modified for standalone use in CMSSW. danek k. 2/06
50 // 25-04-13 Used vdt::log, added check a3>0 (V.Ivanchenko & D. Nikolopoulos)
51 
52 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
53 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
54 
56 #include "CLHEP/Units/GlobalSystemOfUnits.h"
57 #include "CLHEP/Units/GlobalPhysicalConstants.h"
58 #include "CLHEP/Random/RandGaussQ.h"
59 #include "CLHEP/Random/RandPoissonQ.h"
60 #include "CLHEP/Random/RandFlat.h"
61 #include <math.h>
62 #include "vdt/log.h"
63 
64 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
65 
66 using namespace std;
67 
69  :minNumberInteractionsBohr(10.0),
70  theBohrBeta2(50.0*keV/proton_mass_c2),
71  minLoss(10.*eV),
72  problim(5.e-3),
73  alim(10.),
74  nmaxCont1(4.),
75  nmaxCont2(16.)
76 {
77  sumalim = -log(problim);
78  //lastMaterial = 0;
79 
80  // Add these definitions d.k.
81  chargeSquare = 1.; //Assume all particles have charge 1
82  // Taken from Geant4 printout, HARDWIRED for Silicon.
83  ipotFluct = 0.0001736; //material->GetIonisation()->GetMeanExcitationEnergy();
84  electronDensity = 6.797E+20; // material->GetElectronDensity();
85  f1Fluct = 0.8571; // material->GetIonisation()->GetF1fluct();
86  f2Fluct = 0.1429; //material->GetIonisation()->GetF2fluct();
87  e1Fluct = 0.000116;// material->GetIonisation()->GetEnergy1fluct();
88  e2Fluct = 0.00196; //material->GetIonisation()->GetEnergy2fluct();
89  e1LogFluct = -9.063; //material->GetIonisation()->GetLogEnergy1fluct();
90  e2LogFluct = -6.235; //material->GetIonisation()->GetLogEnergy2fluct();
91  rateFluct = 0.4; //material->GetIonisation()->GetRateionexcfluct();
92  ipotLogFluct = -8.659; //material->GetIonisation()->GetLogMeanExcEnergy();
93  e0 = 1.E-5; //material->GetIonisation()->GetEnergy0fluct();
94 
95  //cout << " init new fluct +++++++++++++++++++++++++++++++++++++++++"<<endl;
96 }
97 
98 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
99 // The main dedx fluctuation routine.
100 // Arguments: momentum in MeV/c, mass in MeV, delta ray cut (tmax) in
101 // MeV, silicon thickness in mm, mean eloss in MeV.
102 
104 {
105 }
106 
107 
108 double SiG4UniversalFluctuation::SampleFluctuations(const double momentum,
109  const double mass,
110  double& tmax,
111  const double length,
112  const double meanLoss,
113  CLHEP::HepRandomEngine* engine)
114 {
115 // Calculate actual loss from the mean loss.
116 // The model used to get the fluctuations is essentially the same
117 // as in Glandz in Geant3 (Cern program library W5013, phys332).
118 // L. Urban et al. NIM A362, p.416 (1995) and Geant4 Physics Reference Manual
119 
120  // shortcut for very very small loss (out of validity of the model)
121  //
122  if (meanLoss < minLoss) return meanLoss;
123 
124  particleMass = mass;
125  double gam2 = (momentum*momentum)/(particleMass*particleMass) + 1.0;
126  double beta2 = 1.0 - 1.0/gam2;
127  double gam = sqrt(gam2);
128 
129  double loss(0.), siga(0.);
130 
131  // Gaussian regime
132  // for heavy particles only and conditions
133  // for Gauusian fluct. has been changed
134  //
135  if ((particleMass > electron_mass_c2) &&
136  (meanLoss >= minNumberInteractionsBohr*tmax))
137  {
138  double massrate = electron_mass_c2/particleMass ;
139  double tmaxkine = 2.*electron_mass_c2*beta2*gam2/
140  (1.+massrate*(2.*gam+massrate)) ;
141  if (tmaxkine <= 2.*tmax)
142  {
143  //electronDensity = material->GetElectronDensity();
144  siga = (1.0/beta2 - 0.5) * twopi_mc2_rcl2 * tmax * length
146  siga = sqrt(siga);
147  double twomeanLoss = meanLoss + meanLoss;
148  if (twomeanLoss < siga) {
149  double x;
150  do {
151  loss = twomeanLoss*CLHEP::RandFlat::shoot(engine);
152  x = (loss - meanLoss)/siga;
153  } while (1.0 - 0.5*x*x < CLHEP::RandFlat::shoot(engine));
154  } else {
155  do {
156  loss = CLHEP::RandGaussQ::shoot(engine, meanLoss, siga);
157  } while (loss < 0. || loss > twomeanLoss);
158  }
159  return loss;
160  }
161  }
162 
163  double a1 = 0., a2 = 0., a3 = 0.;
164  double p3;
165  double rate = rateFluct ;
166 
167  double w1 = tmax/ipotFluct;
168  double w2 = vdt::fast_log(2.*electron_mass_c2*beta2*gam2)-beta2;
169 
170  if(w2 > ipotLogFluct)
171  {
172  double C = meanLoss*(1.-rateFluct)/(w2-ipotLogFluct);
173  a1 = C*f1Fluct*(w2-e1LogFluct)/e1Fluct;
174  a2 = C*f2Fluct*(w2-e2LogFluct)/e2Fluct;
175  if(a2 < 0.)
176  {
177  a1 = 0. ;
178  a2 = 0. ;
179  rate = 1. ;
180  }
181  }
182  else
183  {
184  rate = 1. ;
185  }
186 
187  // added
188  if(tmax > ipotFluct) {
189  a3 = rate*meanLoss*(tmax-ipotFluct)/(ipotFluct*tmax*vdt::fast_log(w1));
190  }
191  double suma = a1+a2+a3;
192 
193  // Glandz regime
194  //
195  if (suma > sumalim)
196  {
197  if((a1+a2) > 0.)
198  {
199  double p1, p2;
200  // excitation type 1
201  if (a1>alim) {
202  siga=sqrt(a1) ;
203  p1 = max(0., CLHEP::RandGaussQ::shoot(engine, a1, siga) + 0.5);
204  } else {
205  p1 = double(CLHEP::RandPoissonQ::shoot(engine,a1));
206  }
207 
208  // excitation type 2
209  if (a2>alim) {
210  siga=sqrt(a2) ;
211  p2 = max(0., CLHEP::RandGaussQ::shoot(engine, a2, siga) + 0.5);
212  } else {
213  p2 = double(CLHEP::RandPoissonQ::shoot(engine,a2));
214  }
215 
216  loss = p1*e1Fluct+p2*e2Fluct;
217 
218  // smearing to avoid unphysical peaks
219  if (p2 > 0.)
220  loss += (1.-2.*CLHEP::RandFlat::shoot(engine))*e2Fluct;
221  else if (loss>0.)
222  loss += (1.-2.*CLHEP::RandFlat::shoot(engine))*e1Fluct;
223  if (loss < 0.) loss = 0.0;
224  }
225 
226  // ionisation
227  if (a3 > 0.) {
228  if (a3>alim) {
229  siga=sqrt(a3) ;
230  p3 = max(0., CLHEP::RandGaussQ::shoot(engine, a3, siga) + 0.5);
231  } else {
232  p3 = double(CLHEP::RandPoissonQ::shoot(engine,a3));
233  }
234  double lossc = 0.;
235  if (p3 > 0) {
236  double na = 0.;
237  double alfa = 1.;
238  if (p3 > nmaxCont2) {
239  double rfac = p3/(nmaxCont2+p3);
240  double namean = p3*rfac;
241  double sa = nmaxCont1*rfac;
242  na = CLHEP::RandGaussQ::shoot(engine, namean, sa);
243  if (na > 0.) {
244  alfa = w1*(nmaxCont2+p3)/(w1*nmaxCont2+p3);
245  double alfa1 = alfa*vdt::fast_log(alfa)/(alfa-1.);
246  double ea = na*ipotFluct*alfa1;
247  double sea = ipotFluct*sqrt(na*(alfa-alfa1*alfa1));
248  lossc += CLHEP::RandGaussQ::shoot(engine, ea, sea);
249  }
250  }
251 
252  if (p3 > na) {
253  w2 = alfa*ipotFluct;
254  double w = (tmax-w2)/tmax;
255  int nb = int(p3-na);
256  for (int k=0; k<nb; k++) lossc += w2/(1.-w*CLHEP::RandFlat::shoot(engine));
257  }
258  }
259  loss += lossc;
260  }
261  return loss;
262  }
263 
264  // suma < sumalim; very small energy loss;
265  a3 = meanLoss*(tmax-e0)/(tmax*e0*vdt::fast_log(tmax/e0));
266  if (a3 > alim)
267  {
268  siga=sqrt(a3);
269  p3 = max(0., CLHEP::RandGaussQ::shoot(engine, a3, siga) + 0.5);
270  } else {
271  p3 = double(CLHEP::RandPoissonQ::shoot(engine,a3));
272  }
273  if (p3 > 0.) {
274  double w = (tmax-e0)/tmax;
275  double corrfac = 1.;
276  if (p3 > nmaxCont2) {
277  corrfac = p3/nmaxCont2;
278  p3 = nmaxCont2;
279  }
280  int ip3 = (int)p3;
281  for (int i=0; i<ip3; i++) loss += 1./(1.-w*CLHEP::RandFlat::shoot(engine));
282  loss *= e0*corrfac;
283  // smearing for losses near to e0
284  if(p3 <= 2.)
285  loss += e0*(1.-2.*CLHEP::RandFlat::shoot(engine)) ;
286  }
287 
288  return loss;
289 }
290 
int i
Definition: DBlmapReader.cc:9
common ppss p3p6s2 common epss epspn46 common const1 w2
Definition: inclppp.h:1
const double w
Definition: UKUtility.cc:23
T sqrt(T t)
Definition: SSEVec.h:48
double p2[4]
Definition: TauolaWrapper.h:90
static const double tmax[3]
double SampleFluctuations(const double momentum, const double mass, double &tmax, const double length, const double meanLoss, CLHEP::HepRandomEngine *)
double rate(double x)
Definition: Constants.cc:3
double p1[4]
Definition: TauolaWrapper.h:89
double p3[4]
Definition: TauolaWrapper.h:91
tuple log
Definition: cmsBatch.py:341