CMS 3D CMS Logo

SiG4UniversalFluctuation.cc
Go to the documentation of this file.
1 
2 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
3 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
4 
6 #include "CLHEP/Units/GlobalSystemOfUnits.h"
7 #include "CLHEP/Units/GlobalPhysicalConstants.h"
8 #include "CLHEP/Random/RandGaussQ.h"
9 #include "CLHEP/Random/RandPoissonQ.h"
10 #include "CLHEP/Random/RandFlat.h"
11 #include <math.h>
12 #include "vdt/log.h"
13 
14 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
15 
16 using namespace std;
17 
19  :minNumberInteractionsBohr(10.0),
20  theBohrBeta2(50.0*keV/proton_mass_c2),
21  minLoss(10.*eV),
22  problim(5.e-3),
23  alim(10.),
24  nmaxCont1(4.),
25  nmaxCont2(16.)
26 {
27  sumalim = -log(problim);
28 
29  // Add these definitions d.k.
30  chargeSquare = 1.; //Assume all particles have charge 1
31  // Taken from Geant4 printout, HARDWIRED for Silicon.
32  ipotFluct = 0.0001736; //material->GetIonisation()->GetMeanExcitationEnergy();
33  electronDensity = 6.797E+20; // material->GetElectronDensity();
34  f1Fluct = 0.8571; // material->GetIonisation()->GetF1fluct();
35  f2Fluct = 0.1429; //material->GetIonisation()->GetF2fluct();
36  e1Fluct = 0.000116;// material->GetIonisation()->GetEnergy1fluct();
37  e2Fluct = 0.00196; //material->GetIonisation()->GetEnergy2fluct();
38  e1LogFluct = -9.063; //material->GetIonisation()->GetLogEnergy1fluct();
39  e2LogFluct = -6.235; //material->GetIonisation()->GetLogEnergy2fluct();
40  rateFluct = 0.4; //material->GetIonisation()->GetRateionexcfluct();
41  ipotLogFluct = -8.659; //material->GetIonisation()->GetLogMeanExcEnergy();
42  e0 = 1.E-5; //material->GetIonisation()->GetEnergy0fluct();
43 
44  //cout << " init new fluct +++++++++++++++++++++++++++++++++++++++++"<<endl;
45 }
46 
47 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
48 // The main dedx fluctuation routine.
49 // Arguments: momentum in MeV/c, mass in MeV, delta ray cut (tmax) in
50 // MeV, silicon thickness in mm, mean eloss in MeV.
51 
53 {
54 }
55 
56 double SiG4UniversalFluctuation::SampleFluctuations(const double momentum,
57  const double mass,
58  double& tmax,
59  const double length,
60  const double meanLoss,
61  CLHEP::HepRandomEngine* engine)
62 {
63 // Calculate actual loss from the mean loss.
64 // The model used to get the fluctuations is essentially the same
65 // as in Glandz in Geant3 (Cern program library W5013, phys332).
66 // L. Urban et al. NIM A362, p.416 (1995) and Geant4 Physics Reference Manual
67 
68  // shortcut for very very small loss (out of validity of the model)
69  //
70  if (meanLoss < minLoss) return meanLoss;
71 
72  particleMass = mass;
73  double gam2 = (momentum*momentum)/(particleMass*particleMass) + 1.0;
74  double beta2 = 1.0 - 1.0/gam2;
75  double gam = sqrt(gam2);
76 
77  double loss(0.), siga(0.);
78 
79  // Gaussian regime
80  // for heavy particles only and conditions
81  // for Gauusian fluct. has been changed
82  //
83  if ((particleMass > electron_mass_c2) &&
84  (meanLoss >= minNumberInteractionsBohr*tmax))
85  {
86  double massrate = electron_mass_c2/particleMass ;
87  double tmaxkine = 2.*electron_mass_c2*beta2*gam2/
88  (1.+massrate*(2.*gam+massrate)) ;
89  if (tmaxkine <= 2.*tmax)
90  {
91  siga = (1.0/beta2 - 0.5) * twopi_mc2_rcl2 * tmax * length
93  siga = sqrt(siga);
94  double twomeanLoss = meanLoss + meanLoss;
95  if (twomeanLoss < siga) {
96  double x;
97  do {
98  loss = twomeanLoss*CLHEP::RandFlat::shoot(engine);
99  x = (loss - meanLoss)/siga;
100  } while (1.0 - 0.5*x*x < CLHEP::RandFlat::shoot(engine));
101  } else {
102  do {
103  loss = CLHEP::RandGaussQ::shoot(engine, meanLoss, siga);
104  } while (loss < 0. || loss > twomeanLoss);
105  }
106  return loss;
107  }
108  }
109 
110  double a1 = 0., a2 = 0., a3 = 0.;
111  double p3;
112  double rate = rateFluct ;
113 
114  double w1 = tmax/ipotFluct;
115  double w2 = vdt::fast_log(2.*electron_mass_c2*beta2*gam2)-beta2;
116 
117  if(w2 > ipotLogFluct)
118  {
119  double C = meanLoss*(1.-rateFluct)/(w2-ipotLogFluct);
120  a1 = C*f1Fluct*(w2-e1LogFluct)/e1Fluct;
121  a2 = C*f2Fluct*(w2-e2LogFluct)/e2Fluct;
122  if(a2 < 0.)
123  {
124  a1 = 0. ;
125  a2 = 0. ;
126  rate = 1. ;
127  }
128  }
129  else
130  {
131  rate = 1. ;
132  }
133 
134  // added
135  if(tmax > ipotFluct) {
136  a3 = rate*meanLoss*(tmax-ipotFluct)/(ipotFluct*tmax*vdt::fast_log(w1));
137  }
138  double suma = a1+a2+a3;
139 
140  // Glandz regime
141  //
142  if (suma > sumalim)
143  {
144  if((a1+a2) > 0.)
145  {
146  double p1, p2;
147  // excitation type 1
148  if (a1>alim) {
149  siga=sqrt(a1) ;
150  p1 = max(0., CLHEP::RandGaussQ::shoot(engine, a1, siga) + 0.5);
151  } else {
152  p1 = double(CLHEP::RandPoissonQ::shoot(engine,a1));
153  }
154 
155  // excitation type 2
156  if (a2>alim) {
157  siga=sqrt(a2) ;
158  p2 = max(0., CLHEP::RandGaussQ::shoot(engine, a2, siga) + 0.5);
159  } else {
160  p2 = double(CLHEP::RandPoissonQ::shoot(engine,a2));
161  }
162 
163  loss = p1*e1Fluct+p2*e2Fluct;
164 
165  // smearing to avoid unphysical peaks
166  if (p2 > 0.)
167  loss += (1.-2.*CLHEP::RandFlat::shoot(engine))*e2Fluct;
168  else if (loss>0.)
169  loss += (1.-2.*CLHEP::RandFlat::shoot(engine))*e1Fluct;
170  if (loss < 0.) loss = 0.0;
171  }
172 
173  // ionisation
174  if (a3 > 0.) {
175  if (a3>alim) {
176  siga=sqrt(a3) ;
177  p3 = max(0., CLHEP::RandGaussQ::shoot(engine, a3, siga) + 0.5);
178  } else {
179  p3 = double(CLHEP::RandPoissonQ::shoot(engine,a3));
180  }
181  double lossc = 0.;
182  if (p3 > 0) {
183  double na = 0.;
184  double alfa = 1.;
185  if (p3 > nmaxCont2) {
186  double rfac = p3/(nmaxCont2+p3);
187  double namean = p3*rfac;
188  double sa = nmaxCont1*rfac;
189  na = CLHEP::RandGaussQ::shoot(engine, namean, sa);
190  if (na > 0.) {
191  alfa = w1*(nmaxCont2+p3)/(w1*nmaxCont2+p3);
192  double alfa1 = alfa*vdt::fast_log(alfa)/(alfa-1.);
193  double ea = na*ipotFluct*alfa1;
194  double sea = ipotFluct*sqrt(na*(alfa-alfa1*alfa1));
195  lossc += CLHEP::RandGaussQ::shoot(engine, ea, sea);
196  }
197  }
198 
199  if (p3 > na) {
200  w2 = alfa*ipotFluct;
201  double w = (tmax-w2)/tmax;
202  int nb = int(p3-na);
203  for (int k=0; k<nb; k++) lossc += w2/(1.-w*CLHEP::RandFlat::shoot(engine));
204  }
205  }
206  loss += lossc;
207  }
208  return loss;
209  }
210 
211  // suma < sumalim; very small energy loss;
212  a3 = meanLoss*(tmax-e0)/(tmax*e0*vdt::fast_log(tmax/e0));
213  if (a3 > alim)
214  {
215  siga=sqrt(a3);
216  p3 = max(0., CLHEP::RandGaussQ::shoot(engine, a3, siga) + 0.5);
217  } else {
218  p3 = double(CLHEP::RandPoissonQ::shoot(engine,a3));
219  }
220  if (p3 > 0.) {
221  double w = (tmax-e0)/tmax;
222  double corrfac = 1.;
223  if (p3 > nmaxCont2) {
224  corrfac = p3/nmaxCont2;
225  p3 = nmaxCont2;
226  }
227  int ip3 = (int)p3;
228  for (int i=0; i<ip3; i++) loss += 1./(1.-w*CLHEP::RandFlat::shoot(engine));
229  loss *= e0*corrfac;
230  // smearing for losses near to e0
231  if(p3 <= 2.)
232  loss += e0*(1.-2.*CLHEP::RandFlat::shoot(engine)) ;
233  }
234 
235  return loss;
236 }
237 
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:18
double p2[4]
Definition: TauolaWrapper.h:90
int k[5][pyjets_maxn]
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