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