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 // $Id: SiG4UniversalFluctuation.cc,v 1.7 2009/05/25 07:38:46 fabiocos Exp $
24 // GEANT4 tag $Name: CMSSW_5_3_11_patch3 $
25 //
26 // -------------------------------------------------------------------
27 //
28 // GEANT4 Class file
29 //
30 //
31 // File name: G4UniversalFluctuation
32 //
33 // Author: Vladimir Ivanchenko
34 //
35 // Creation date: 03.01.2002
36 //
37 // Modifications:
38 //
39 // 28-12-02 add method Dispersion (V.Ivanchenko)
40 // 07-02-03 change signature (V.Ivanchenko)
41 // 13-02-03 Add name (V.Ivanchenko)
42 // 16-10-03 Changed interface to Initialisation (V.Ivanchenko)
43 // 07-11-03 Fix problem of rounding of double in G4UniversalFluctuations
44 // 06-02-04 Add control on big sigma > 2*meanLoss (V.Ivanchenko)
45 // 26-04-04 Comment out the case of very small step (V.Ivanchenko)
46 // 07-02-05 define problim = 5.e-3 (mma)
47 // 03-05-05 conditions of Gaussian fluctuation changed (bugfix)
48 // + smearing for very small loss (L.Urban)
49 //
50 // Modified for standalone use in CMSSW. danek k. 2/06
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/RandPoisson.h"
60 #include "CLHEP/Random/RandFlat.h"
61 #include <math.h>
62 //#include "G4UniversalFluctuation.hh"
63 //#include "Randomize.hh"
64 //#include "G4Poisson.hh"
65 //#include "G4Step.hh"
66 //#include "G4Material.hh"
67 //#include "G4DynamicParticle.hh"
68 //#include "G4ParticleDefinition.hh"
69 
70 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
71 
72 using namespace std;
73 
75  :rndEngine(eng),
76  gaussQDistribution(0),
77  poissonDistribution(0),
78  flatDistribution(0),
79  minNumberInteractionsBohr(10.0),
80  theBohrBeta2(50.0*keV/proton_mass_c2),
81  minLoss(10.*eV),
82  problim(5.e-3),
83  alim(10.),
84  nmaxCont1(4.),
85  nmaxCont2(16.)
86 {
87  sumalim = -log(problim);
88  //lastMaterial = 0;
89 
90  // Add these definitions d.k.
91  chargeSquare = 1.; //Assume all particles have charge 1
92  // Taken from Geant4 printout, HARDWIRED for Silicon.
93  ipotFluct = 0.0001736; //material->GetIonisation()->GetMeanExcitationEnergy();
94  electronDensity = 6.797E+20; // material->GetElectronDensity();
95  f1Fluct = 0.8571; // material->GetIonisation()->GetF1fluct();
96  f2Fluct = 0.1429; //material->GetIonisation()->GetF2fluct();
97  e1Fluct = 0.000116;// material->GetIonisation()->GetEnergy1fluct();
98  e2Fluct = 0.00196; //material->GetIonisation()->GetEnergy2fluct();
99  e1LogFluct = -9.063; //material->GetIonisation()->GetLogEnergy1fluct();
100  e2LogFluct = -6.235; //material->GetIonisation()->GetLogEnergy2fluct();
101  rateFluct = 0.4; //material->GetIonisation()->GetRateionexcfluct();
102  ipotLogFluct = -8.659; //material->GetIonisation()->GetLogMeanExcEnergy();
103  e0 = 1.E-5; //material->GetIonisation()->GetEnergy0fluct();
104 
105  gaussQDistribution = new CLHEP::RandGaussQ(rndEngine);
106  poissonDistribution = new CLHEP::RandPoisson(rndEngine);
107  flatDistribution = new CLHEP::RandFlat(rndEngine);
108 
109  //cout << " init new fluct +++++++++++++++++++++++++++++++++++++++++"<<endl;
110 }
111 
112 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
113 // The main dedx fluctuation routine.
114 // Arguments: momentum in MeV/c, mass in MeV, delta ray cut (tmax) in
115 // MeV, silicon thickness in mm, mean eloss in MeV.
116 
118 {
119  delete gaussQDistribution;
120  delete poissonDistribution;
121  delete flatDistribution;
122 
123 }
124 
125 
126 double SiG4UniversalFluctuation::SampleFluctuations(const double momentum,
127  const double mass,
128  double& tmax,
129  const double length,
130  const double meanLoss)
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 very small loss (out of validity of the model)
138  //
139  if (meanLoss < minLoss) return meanLoss;
140 
141  //if(!particle) InitialiseMe(dp->GetDefinition());
142  //G4double tau = dp->GetKineticEnergy()/particleMass;
143  //G4double gam = tau + 1.0;
144  //G4double gam2 = gam*gam;
145  //G4double beta2 = tau*(tau + 2.0)/gam2;
146 
147  particleMass = mass; // dp->GetMass();
148  double gam2 = (momentum*momentum)/(particleMass*particleMass) + 1.0;
149  double beta2 = 1.0 - 1.0/gam2;
150  double gam = sqrt(gam2);
151 
152  double loss(0.), siga(0.);
153 
154  // Gaussian regime
155  // for heavy particles only and conditions
156  // for Gauusian fluct. has been changed
157  //
158  if ((particleMass > electron_mass_c2) &&
159  (meanLoss >= minNumberInteractionsBohr*tmax))
160  {
161  double massrate = electron_mass_c2/particleMass ;
162  double tmaxkine = 2.*electron_mass_c2*beta2*gam2/
163  (1.+massrate*(2.*gam+massrate)) ;
164  if (tmaxkine <= 2.*tmax)
165  {
166  //electronDensity = material->GetElectronDensity();
167  siga = (1.0/beta2 - 0.5) * twopi_mc2_rcl2 * tmax * length
169  siga = sqrt(siga);
170  double twomeanLoss = meanLoss + meanLoss;
171  if (twomeanLoss < siga) {
172  double x;
173  do {
174  loss = twomeanLoss*flatDistribution->fire();
175  x = (loss - meanLoss)/siga;
176  } while (1.0 - 0.5*x*x < flatDistribution->fire());
177  } else {
178  do {
179  loss = gaussQDistribution->fire(meanLoss,siga);
180  } while (loss < 0. || loss > twomeanLoss);
181  }
182  return loss;
183  }
184  }
185 
186  // Glandz regime : initialisation
187  //
188 // if (material != lastMaterial) {
189 // f1Fluct = material->GetIonisation()->GetF1fluct();
190 // f2Fluct = material->GetIonisation()->GetF2fluct();
191 // e1Fluct = material->GetIonisation()->GetEnergy1fluct();
192 // e2Fluct = material->GetIonisation()->GetEnergy2fluct();
193 // e1LogFluct = material->GetIonisation()->GetLogEnergy1fluct();
194 // e2LogFluct = material->GetIonisation()->GetLogEnergy2fluct();
195 // rateFluct = material->GetIonisation()->GetRateionexcfluct();
196 // ipotFluct = material->GetIonisation()->GetMeanExcitationEnergy();
197 // ipotLogFluct = material->GetIonisation()->GetLogMeanExcEnergy();
198 // lastMaterial = material;
199 // }
200 
201  double a1 = 0. , a2 = 0., a3 = 0. ;
202  double p1,p2,p3;
203  double rate = rateFluct ;
204 
205  double w1 = tmax/ipotFluct;
206  double w2 = log(2.*electron_mass_c2*beta2*gam2)-beta2;
207 
208  if(w2 > ipotLogFluct)
209  {
210  double C = meanLoss*(1.-rateFluct)/(w2-ipotLogFluct);
211  a1 = C*f1Fluct*(w2-e1LogFluct)/e1Fluct;
212  a2 = C*f2Fluct*(w2-e2LogFluct)/e2Fluct;
213  if(a2 < 0.)
214  {
215  a1 = 0. ;
216  a2 = 0. ;
217  rate = 1. ;
218  }
219  }
220  else
221  {
222  rate = 1. ;
223  }
224 
225  a3 = rate*meanLoss*(tmax-ipotFluct)/(ipotFluct*tmax*log(w1));
226 
227  double suma = a1+a2+a3;
228 
229  // Glandz regime
230  //
231  if (suma > sumalim)
232  {
233  p1 = 0., p2 = 0 ;
234  if((a1+a2) > 0.)
235  {
236  // excitation type 1
237  if (a1>alim) {
238  siga=sqrt(a1) ;
239  p1 = max(0.,gaussQDistribution->fire(a1,siga)+0.5);
240  } else {
241  p1 = double(poissonDistribution->fire(a1));
242  }
243 
244  // excitation type 2
245  if (a2>alim) {
246  siga=sqrt(a2) ;
247  p2 = max(0.,gaussQDistribution->fire(a2,siga)+0.5);
248  } else {
249  p2 = double(poissonDistribution->fire(a2));
250  }
251 
252  loss = p1*e1Fluct+p2*e2Fluct;
253 
254  // smearing to avoid unphysical peaks
255  if (p2 > 0.)
256  loss += (1.-2.*flatDistribution->fire())*e2Fluct;
257  else if (loss>0.)
258  loss += (1.-2.*flatDistribution->fire())*e1Fluct;
259  if (loss < 0.) loss = 0.0;
260  }
261 
262  // ionisation
263  if (a3 > 0.) {
264  if (a3>alim) {
265  siga=sqrt(a3) ;
266  p3 = max(0.,gaussQDistribution->fire(a3,siga)+0.5);
267  } else {
268  p3 = double(poissonDistribution->fire(a3));
269  }
270  double lossc = 0.;
271  if (p3 > 0) {
272  double na = 0.;
273  double alfa = 1.;
274  if (p3 > nmaxCont2) {
275  double rfac = p3/(nmaxCont2+p3);
276  double namean = p3*rfac;
277  double sa = nmaxCont1*rfac;
278  na = gaussQDistribution->fire(namean,sa);
279  if (na > 0.) {
280  alfa = w1*(nmaxCont2+p3)/(w1*nmaxCont2+p3);
281  double alfa1 = alfa*log(alfa)/(alfa-1.);
282  double ea = na*ipotFluct*alfa1;
283  double sea = ipotFluct*sqrt(na*(alfa-alfa1*alfa1));
284  lossc += gaussQDistribution->fire(ea,sea);
285  }
286  }
287 
288  if (p3 > na) {
289  w2 = alfa*ipotFluct;
290  double w = (tmax-w2)/tmax;
291  int nb = int(p3-na);
292  for (int k=0; k<nb; k++) lossc += w2/(1.-w*flatDistribution->fire());
293  }
294  }
295  loss += lossc;
296  }
297  return loss;
298  }
299 
300  // suma < sumalim; very small energy loss;
301  //
302  //double e0 = material->GetIonisation()->GetEnergy0fluct();
303 
304  a3 = meanLoss*(tmax-e0)/(tmax*e0*log(tmax/e0));
305  if (a3 > alim)
306  {
307  siga=sqrt(a3);
308  p3 = max(0.,gaussQDistribution->fire(a3,siga)+0.5);
309  } else {
310  p3 = double(poissonDistribution->fire(a3));
311  }
312  if (p3 > 0.) {
313  double w = (tmax-e0)/tmax;
314  double corrfac = 1.;
315  if (p3 > nmaxCont2) {
316  corrfac = p3/nmaxCont2;
317  p3 = nmaxCont2;
318  }
319  int ip3 = (int)p3;
320  for (int i=0; i<ip3; i++) loss += 1./(1.-w*flatDistribution->fire());
321  loss *= e0*corrfac;
322  // smearing for losses near to e0
323  if(p3 <= 2.)
324  loss += e0*(1.-2.*flatDistribution->fire()) ;
325  }
326 
327  return loss;
328 }
329 
330 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
331 // G4double SiG4UniversalFluctuation::Dispersion(
332 // const G4Material* material,
333 // const G4DynamicParticle* dp,
334 // G4double& tmax,
335 // G4double& length)
336 // {
337 // if(!particle) InitialiseMe(dp->GetDefinition());
338 
339 // electronDensity = material->GetElectronDensity();
340 
341 // G4double gam = (dp->GetKineticEnergy())/particleMass + 1.0;
342 // G4double beta2 = 1.0 - 1.0/(gam*gam);
343 
344 // G4double siga = (1.0/beta2 - 0.5) * twopi_mc2_rcl2 * tmax * length
345 // * electronDensity * chargeSquare;
346 
347 // return siga;
348 // }
349 
350 //....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:46
double p2[4]
Definition: TauolaWrapper.h:90
list rate
Definition: scaleCards.py:77
int k[5][pyjets_maxn]
static const double tmax[3]
CLHEP::RandGaussQ * gaussQDistribution
double p1[4]
Definition: TauolaWrapper.h:89
tuple mass
Definition: scaleCards.py:27
Definition: DDAxes.h:10
CLHEP::RandPoisson * poissonDistribution
CLHEP::HepRandomEngine & rndEngine
T w() const
double p3[4]
Definition: TauolaWrapper.h:91