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 //#include "G4UniversalFluctuation.hh"
64 //#include "Randomize.hh"
65 //#include "G4Poisson.hh"
66 //#include "G4Step.hh"
67 //#include "G4Material.hh"
68 //#include "G4DynamicParticle.hh"
69 //#include "G4ParticleDefinition.hh"
70 
71 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
72 
73 using namespace std;
74 
76  :minNumberInteractionsBohr(10.0),
77  theBohrBeta2(50.0*keV/proton_mass_c2),
78  minLoss(10.*eV),
79  problim(5.e-3),
80  alim(10.),
81  nmaxCont1(4.),
82  nmaxCont2(16.)
83 {
84  sumalim = -log(problim);
85  //lastMaterial = 0;
86 
87  // Add these definitions d.k.
88  chargeSquare = 1.; //Assume all particles have charge 1
89  // Taken from Geant4 printout, HARDWIRED for Silicon.
90  ipotFluct = 0.0001736; //material->GetIonisation()->GetMeanExcitationEnergy();
91  electronDensity = 6.797E+20; // material->GetElectronDensity();
92  f1Fluct = 0.8571; // material->GetIonisation()->GetF1fluct();
93  f2Fluct = 0.1429; //material->GetIonisation()->GetF2fluct();
94  e1Fluct = 0.000116;// material->GetIonisation()->GetEnergy1fluct();
95  e2Fluct = 0.00196; //material->GetIonisation()->GetEnergy2fluct();
96  e1LogFluct = -9.063; //material->GetIonisation()->GetLogEnergy1fluct();
97  e2LogFluct = -6.235; //material->GetIonisation()->GetLogEnergy2fluct();
98  rateFluct = 0.4; //material->GetIonisation()->GetRateionexcfluct();
99  ipotLogFluct = -8.659; //material->GetIonisation()->GetLogMeanExcEnergy();
100  e0 = 1.E-5; //material->GetIonisation()->GetEnergy0fluct();
101 
102  //cout << " init new fluct +++++++++++++++++++++++++++++++++++++++++"<<endl;
103 }
104 
105 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
106 // The main dedx fluctuation routine.
107 // Arguments: momentum in MeV/c, mass in MeV, delta ray cut (tmax) in
108 // MeV, silicon thickness in mm, mean eloss in MeV.
109 
111 {
112 }
113 
114 
115 double SiG4UniversalFluctuation::SampleFluctuations(const double momentum,
116  const double mass,
117  double& tmax,
118  const double length,
119  const double meanLoss,
120  CLHEP::HepRandomEngine* engine)
121 {
122 // Calculate actual loss from the mean loss.
123 // The model used to get the fluctuations is essentially the same
124 // as in Glandz in Geant3 (Cern program library W5013, phys332).
125 // L. Urban et al. NIM A362, p.416 (1995) and Geant4 Physics Reference Manual
126 
127  // shortcut for very very small loss (out of validity of the model)
128  //
129  if (meanLoss < minLoss) return meanLoss;
130 
131  //if(!particle) InitialiseMe(dp->GetDefinition());
132  //G4double tau = dp->GetKineticEnergy()/particleMass;
133  //G4double gam = tau + 1.0;
134  //G4double gam2 = gam*gam;
135  //G4double beta2 = tau*(tau + 2.0)/gam2;
136 
137  particleMass = mass; // dp->GetMass();
138  double gam2 = (momentum*momentum)/(particleMass*particleMass) + 1.0;
139  double beta2 = 1.0 - 1.0/gam2;
140  double gam = sqrt(gam2);
141 
142  double loss(0.), siga(0.);
143 
144  // Gaussian regime
145  // for heavy particles only and conditions
146  // for Gauusian fluct. has been changed
147  //
148  if ((particleMass > electron_mass_c2) &&
149  (meanLoss >= minNumberInteractionsBohr*tmax))
150  {
151  double massrate = electron_mass_c2/particleMass ;
152  double tmaxkine = 2.*electron_mass_c2*beta2*gam2/
153  (1.+massrate*(2.*gam+massrate)) ;
154  if (tmaxkine <= 2.*tmax)
155  {
156  //electronDensity = material->GetElectronDensity();
157  siga = (1.0/beta2 - 0.5) * twopi_mc2_rcl2 * tmax * length
159  siga = sqrt(siga);
160  double twomeanLoss = meanLoss + meanLoss;
161  if (twomeanLoss < siga) {
162  double x;
163  do {
164  loss = twomeanLoss*CLHEP::RandFlat::shoot(engine);
165  x = (loss - meanLoss)/siga;
166  } while (1.0 - 0.5*x*x < CLHEP::RandFlat::shoot(engine));
167  } else {
168  do {
169  loss = CLHEP::RandGaussQ::shoot(engine, meanLoss, siga);
170  } while (loss < 0. || loss > twomeanLoss);
171  }
172  return loss;
173  }
174  }
175 
176  // Glandz regime : initialisation
177  //
178 // if (material != lastMaterial) {
179 // f1Fluct = material->GetIonisation()->GetF1fluct();
180 // f2Fluct = material->GetIonisation()->GetF2fluct();
181 // e1Fluct = material->GetIonisation()->GetEnergy1fluct();
182 // e2Fluct = material->GetIonisation()->GetEnergy2fluct();
183 // e1LogFluct = material->GetIonisation()->GetLogEnergy1fluct();
184 // e2LogFluct = material->GetIonisation()->GetLogEnergy2fluct();
185 // rateFluct = material->GetIonisation()->GetRateionexcfluct();
186 // ipotFluct = material->GetIonisation()->GetMeanExcitationEnergy();
187 // ipotLogFluct = material->GetIonisation()->GetLogMeanExcEnergy();
188 // lastMaterial = material;
189 // }
190 
191  double a1 = 0. , a2 = 0., a3 = 0. ;
192  double p1,p2,p3;
193  double rate = rateFluct ;
194 
195  double w1 = tmax/ipotFluct;
196  double w2 = vdt::fast_log(2.*electron_mass_c2*beta2*gam2)-beta2;
197 
198  if(w2 > ipotLogFluct)
199  {
200  double C = meanLoss*(1.-rateFluct)/(w2-ipotLogFluct);
201  a1 = C*f1Fluct*(w2-e1LogFluct)/e1Fluct;
202  a2 = C*f2Fluct*(w2-e2LogFluct)/e2Fluct;
203  if(a2 < 0.)
204  {
205  a1 = 0. ;
206  a2 = 0. ;
207  rate = 1. ;
208  }
209  }
210  else
211  {
212  rate = 1. ;
213  }
214 
215  // added
216  if(tmax > ipotFluct) {
217  a3 = rate*meanLoss*(tmax-ipotFluct)/(ipotFluct*tmax*vdt::fast_log(w1));
218  }
219  double suma = a1+a2+a3;
220 
221  // Glandz regime
222  //
223  if (suma > sumalim)
224  {
225  p1 = 0., p2 = 0 ;
226  if((a1+a2) > 0.)
227  {
228  // excitation type 1
229  if (a1>alim) {
230  siga=sqrt(a1) ;
231  p1 = max(0., CLHEP::RandGaussQ::shoot(engine, a1, siga) + 0.5);
232  } else {
233  p1 = double(CLHEP::RandPoissonQ::shoot(engine,a1));
234  }
235 
236  // excitation type 2
237  if (a2>alim) {
238  siga=sqrt(a2) ;
239  p2 = max(0., CLHEP::RandGaussQ::shoot(engine, a2, siga) + 0.5);
240  } else {
241  p2 = double(CLHEP::RandPoissonQ::shoot(engine,a2));
242  }
243 
244  loss = p1*e1Fluct+p2*e2Fluct;
245 
246  // smearing to avoid unphysical peaks
247  if (p2 > 0.)
248  loss += (1.-2.*CLHEP::RandFlat::shoot(engine))*e2Fluct;
249  else if (loss>0.)
250  loss += (1.-2.*CLHEP::RandFlat::shoot(engine))*e1Fluct;
251  if (loss < 0.) loss = 0.0;
252  }
253 
254  // ionisation
255  if (a3 > 0.) {
256  if (a3>alim) {
257  siga=sqrt(a3) ;
258  p3 = max(0., CLHEP::RandGaussQ::shoot(engine, a3, siga) + 0.5);
259  } else {
260  p3 = double(CLHEP::RandPoissonQ::shoot(engine,a3));
261  }
262  double lossc = 0.;
263  if (p3 > 0) {
264  double na = 0.;
265  double alfa = 1.;
266  if (p3 > nmaxCont2) {
267  double rfac = p3/(nmaxCont2+p3);
268  double namean = p3*rfac;
269  double sa = nmaxCont1*rfac;
270  na = CLHEP::RandGaussQ::shoot(engine, namean, sa);
271  if (na > 0.) {
272  alfa = w1*(nmaxCont2+p3)/(w1*nmaxCont2+p3);
273  double alfa1 = alfa*vdt::fast_log(alfa)/(alfa-1.);
274  double ea = na*ipotFluct*alfa1;
275  double sea = ipotFluct*sqrt(na*(alfa-alfa1*alfa1));
276  lossc += CLHEP::RandGaussQ::shoot(engine, ea, sea);
277  }
278  }
279 
280  if (p3 > na) {
281  w2 = alfa*ipotFluct;
282  double w = (tmax-w2)/tmax;
283  int nb = int(p3-na);
284  for (int k=0; k<nb; k++) lossc += w2/(1.-w*CLHEP::RandFlat::shoot(engine));
285  }
286  }
287  loss += lossc;
288  }
289  return loss;
290  }
291 
292  // suma < sumalim; very small energy loss;
293  //
294  //double e0 = material->GetIonisation()->GetEnergy0fluct();
295 
296  a3 = meanLoss*(tmax-e0)/(tmax*e0*vdt::fast_log(tmax/e0));
297  if (a3 > alim)
298  {
299  siga=sqrt(a3);
300  p3 = max(0., CLHEP::RandGaussQ::shoot(engine, a3, siga) + 0.5);
301  } else {
302  p3 = double(CLHEP::RandPoissonQ::shoot(engine,a3));
303  }
304  if (p3 > 0.) {
305  double w = (tmax-e0)/tmax;
306  double corrfac = 1.;
307  if (p3 > nmaxCont2) {
308  corrfac = p3/nmaxCont2;
309  p3 = nmaxCont2;
310  }
311  int ip3 = (int)p3;
312  for (int i=0; i<ip3; i++) loss += 1./(1.-w*CLHEP::RandFlat::shoot(engine));
313  loss *= e0*corrfac;
314  // smearing for losses near to e0
315  if(p3 <= 2.)
316  loss += e0*(1.-2.*CLHEP::RandFlat::shoot(engine)) ;
317  }
318 
319  return loss;
320 }
321 
322 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
323 // G4double SiG4UniversalFluctuation::Dispersion(
324 // const G4Material* material,
325 // const G4DynamicParticle* dp,
326 // G4double& tmax,
327 // G4double& length)
328 // {
329 // if(!particle) InitialiseMe(dp->GetDefinition());
330 
331 // electronDensity = material->GetElectronDensity();
332 
333 // G4double gam = (dp->GetKineticEnergy())/particleMass + 1.0;
334 // G4double beta2 = 1.0 - 1.0/(gam*gam);
335 
336 // G4double siga = (1.0/beta2 - 0.5) * twopi_mc2_rcl2 * tmax * length
337 // * electronDensity * chargeSquare;
338 
339 // return siga;
340 // }
341 
342 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
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