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