CMS 3D CMS Logo

LandauFP420.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 // -------------------------------------------------------------------
24 //
25 // GEANT4 Class file
26 //
27 //
28 // File name: G4UniversalFluctuation
29 //
30 // Author: Vladimir Ivanchenko
31 //
32 // Creation date: 03.01.2002
33 //
34 // Modifications:
35 //
36 // 28-12-02 add method Dispersion (V.Ivanchenko)
37 // 07-02-03 change signature (V.Ivanchenko)
38 // 13-02-03 Add name (V.Ivanchenko)
39 // Modified for standalone use in ORCA. d.k. 6/04
40 //
41 // -------------------------------------------------------------------
42 
44 #include "CLHEP/Units/GlobalSystemOfUnits.h"
45 #include "CLHEP/Units/GlobalPhysicalConstants.h"
46 #include "CLHEP/Random/RandGaussQ.h"
47 #include "CLHEP/Random/RandPoisson.h"
48 #include "CLHEP/Random/RandFlat.h"
49 
50 #include <cstdio>
51 #include <gsl/gsl_fit.h>
52 #include <cmath>
53 #include<vector>
54 
55 //#include "Randomize.hh"
56 //#include "G4Poisson.hh"
57 //#include "G4Step.hh"
58 //#include "G4Material.hh"
59 //#include "G4DynamicParticle.hh"
60 //#include "G4ParticleDefinition.hh"
61 using namespace std;
62 
63 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
64 // The constructor setups various constants pluc eloss parameters
65 // for silicon.
67  :minNumberInteractionsBohr(10.0),
68  theBohrBeta2(50.0*keV/proton_mass_c2),
69  minLoss(0.000001*eV),
70  problim(0.01),
71  alim(10.),
72  nmaxCont1(4),
73  nmaxCont2(16)
74 {
75  sumalim = -log(problim);
76 
77  chargeSquare = 1.; //Assume all particles have charge 1
78  // Taken from Geant4 printout, HARDWIRED for Silicon.
79  ipotFluct = 0.0001736; //material->GetIonisation()->GetMeanExcitationEnergy();
80  electronDensity = 6.797E+20; // material->GetElectronDensity();
81  f1Fluct = 0.8571; // material->GetIonisation()->GetF1fluct();
82  f2Fluct = 0.1429; //material->GetIonisation()->GetF2fluct();
83  e1Fluct = 0.000116;// material->GetIonisation()->GetEnergy1fluct();
84  e2Fluct = 0.00196; //material->GetIonisation()->GetEnergy2fluct();
85  e1LogFluct = -9.063; //material->GetIonisation()->GetLogEnergy1fluct();
86  e2LogFluct = -6.235; //material->GetIonisation()->GetLogEnergy2fluct();
87  rateFluct = 0.4; //material->GetIonisation()->GetRateionexcfluct();
88  ipotLogFluct = -8.659; //material->GetIonisation()->GetLogMeanExcEnergy();
89  e0 = 1.E-5; //material->GetIonisation()->GetEnergy0fluct();
90 
91 }
92 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
94 {}
95 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
96 // The main dedx fluctuation routine.
97 // Arguments: momentum in MeV/c, mass in MeV, delta ray cut (tmax) in
98 // MeV, silicon thickness in mm, mean eloss in MeV.
99 double LandauFP420::SampleFluctuations(const double momentum,
100  const double mass,
101  double& tmax,
102  const double length,
103  const double meanLoss)
104 {
105 // calculate actual loss from the mean loss
106 // The model used to get the fluctuation is essentially the same
107 // as in Glandz in Geant3.
108 
109  // shortcut for very very small loss
110  if(meanLoss < minLoss) return meanLoss;
111 
112  //if(dp->GetDefinition() != particle) {
113  particleMass = mass; // dp->GetMass();
114  //G4double q = dp->GetCharge();
115  //chargeSquare = q*q;
116  //}
117 
118  //double gam = (dp->GetKineticEnergy())/particleMass + 1.0;
119  //double gam2 = gam*gam;
120  double gam2 = (momentum*momentum)/(particleMass*particleMass) + 1.0;
121  double beta2 = 1.0 - 1.0/gam2;
122 
123  // Validity range for delta electron cross section
124  double loss, siga;
125  // Gaussian fluctuation
126  if(meanLoss >= minNumberInteractionsBohr*tmax ||
128  {
129  siga = (1.0/beta2 - 0.5) * twopi_mc2_rcl2 * tmax * length
131  siga = sqrt(siga);
132  do {
133  //loss = G4RandGauss::shoot(meanLoss,siga);
134  loss = CLHEP::RandGaussQ::shoot(meanLoss,siga);
135  } while (loss < 0. || loss > 2.*meanLoss);
136 
137  return loss;
138  }
139 
140  // Non Gaussian fluctuation
141  double suma,w1,w2,C,lossc,w;
142  double a1,a2,a3;
143  int p1,p2,p3;
144  int nb;
145  double corrfac, na,alfa,rfac,namean,sa,alfa1,ea,sea;
146  double dp3;
147 
148  w1 = tmax/ipotFluct;
149  w2 = log(2.*electron_mass_c2*(gam2 - 1.0));
150 
151  C = meanLoss*(1.-rateFluct)/(w2-ipotLogFluct-beta2);
152 
153  a1 = C*f1Fluct*(w2-e1LogFluct-beta2)/e1Fluct;
154  a2 = C*f2Fluct*(w2-e2LogFluct-beta2)/e2Fluct;
155  a3 = rateFluct*meanLoss*(tmax-ipotFluct)/(ipotFluct*tmax*log(w1));
156  if(a1 < 0.) a1 = 0.;
157  if(a2 < 0.) a2 = 0.;
158  if(a3 < 0.) a3 = 0.;
159 
160  suma = a1+a2+a3;
161 
162  loss = 0. ;
163 
164  if(suma < sumalim) // very small Step
165  {
166  //e0 = material->GetIonisation()->GetEnergy0fluct();//Hardwired in const
167  if(tmax == ipotFluct)
168  {
169  a3 = meanLoss/e0;
170 
171  if(a3>alim)
172  {
173  siga=sqrt(a3) ;
174  //p3 = G4std::max(0,int(G4RandGauss::shoot(a3,siga)+0.5));
175  p3 = std::max(0,int(CLHEP::RandGaussQ::shoot(a3,siga)+0.5));
176  }
177  else
178  p3 = CLHEP::RandPoisson::shoot(a3);
179  //p3 = G4Poisson(a3);
180 
181  loss = p3*e0 ;
182 
183  if(p3 > 0)
184  //loss += (1.-2.*G4UniformRand())*e0 ;
185  loss += (1.-2.*CLHEP::RandFlat::shoot())*e0 ;
186 
187  }
188  else
189  {
190  tmax = tmax-ipotFluct+e0 ;
191  a3 = meanLoss*(tmax-e0)/(tmax*e0*log(tmax/e0));
192 
193  if(a3>alim)
194  {
195  siga=sqrt(a3) ;
196  //p3 = G4std::max(0,int(G4RandGauss::shoot(a3,siga)+0.5));
197  p3 = std::max(0,int(CLHEP::RandGaussQ::shoot(a3,siga)+0.5));
198  }
199  else
200  p3 = CLHEP::RandPoisson::shoot(a3);
201  //p3 = G4Poisson(a3);
202 
203  if(p3 > 0)
204  {
205  w = (tmax-e0)/tmax ;
206  if(p3 > nmaxCont2)
207  {
208  dp3 = float(p3) ;
209  corrfac = dp3/float(nmaxCont2) ;
210  p3 = nmaxCont2 ;
211  }
212  else
213  corrfac = 1. ;
214 
215  //for(int i=0; i<p3; i++) loss += 1./(1.-w*G4UniformRand()) ;
216  for(int i=0; i<p3; i++) loss += 1./(1.-w*CLHEP::RandFlat::shoot()) ;
217  loss *= e0*corrfac ;
218  }
219  }
220  }
221 
222  else // not so small Step
223  {
224  // excitation type 1
225  if(a1>alim)
226  {
227  siga=sqrt(a1) ;
228  //p1 = std::max(0,int(G4RandGauss::shoot(a1,siga)+0.5));
229  p1 = std::max(0,int(CLHEP::RandGaussQ::shoot(a1,siga)+0.5));
230  }
231  else
232  p1 = CLHEP::RandPoisson::shoot(a1);
233  //p1 = G4Poisson(a1);
234 
235  // excitation type 2
236  if(a2>alim)
237  {
238  siga=sqrt(a2) ;
239  //p2 = std::max(0,int(G4RandGauss::shoot(a2,siga)+0.5));
240  p2 = std::max(0,int(CLHEP::RandGaussQ::shoot(a2,siga)+0.5));
241  }
242  else
243  p2 = CLHEP::RandPoisson::shoot(a2);
244  //p2 = G4Poisson(a2);
245 
246  loss = p1*e1Fluct+p2*e2Fluct;
247 
248  // smearing to avoid unphysical peaks
249  if(p2 > 0)
250  //loss += (1.-2.*G4UniformRand())*e2Fluct;
251  loss += (1.-2.*CLHEP::RandFlat::shoot())*e2Fluct;
252  else if (loss>0.)
253  loss += (1.-2.*CLHEP::RandFlat::shoot())*e1Fluct;
254 
255  // ionisation .......................................
256  if(a3 > 0.)
257  {
258  if(a3>alim)
259  {
260  siga=sqrt(a3) ;
261  p3 = std::max(0,int(CLHEP::RandGaussQ::shoot(a3,siga)+0.5));
262  }
263  else
264  p3 = CLHEP::RandPoisson::shoot(a3);
265 
266  lossc = 0.;
267  if(p3 > 0)
268  {
269  na = 0.;
270  alfa = 1.;
271  if (p3 > nmaxCont2)
272  {
273  dp3 = float(p3);
274  rfac = dp3/(float(nmaxCont2)+dp3);
275  namean = float(p3)*rfac;
276  sa = float(nmaxCont1)*rfac;
277  na = CLHEP::RandGaussQ::shoot(namean,sa);
278  if (na > 0.)
279  {
280  alfa = w1*float(nmaxCont2+p3)/
281  (w1*float(nmaxCont2)+float(p3));
282  alfa1 = alfa*log(alfa)/(alfa-1.);
283  ea = na*ipotFluct*alfa1;
284  sea = ipotFluct*sqrt(na*(alfa-alfa1*alfa1));
285  lossc += CLHEP::RandGaussQ::shoot(ea,sea);
286  }
287  }
288 
289  nb = int(float(p3)-na);
290  if (nb > 0)
291  {
292  w2 = alfa*ipotFluct;
293  w = (tmax-w2)/tmax;
294  for (int k=0; k<nb; k++) lossc += w2/(1.-w*CLHEP::RandFlat::shoot());
295  }
296  }
297  loss += lossc;
298  }
299  }
300 
301  return loss;
302 }
303 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
double e1Fluct
Definition: LandauFP420.h:94
common ppss p3p6s2 common epss epspn46 common const1 w2
Definition: inclppp.h:1
const double w
Definition: UKUtility.cc:23
double e2LogFluct
Definition: LandauFP420.h:98
double e1LogFluct
Definition: LandauFP420.h:97
double chargeSquare
Definition: LandauFP420.h:85
double rateFluct
Definition: LandauFP420.h:96
double minLoss
Definition: LandauFP420.h:104
double alim
Definition: LandauFP420.h:107
double sumalim
Definition: LandauFP420.h:106
double electronDensity
Definition: LandauFP420.h:89
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 minNumberInteractionsBohr
Definition: LandauFP420.h:102
double p1[4]
Definition: TauolaWrapper.h:89
double ipotFluct
Definition: LandauFP420.h:88
double particleMass
Definition: LandauFP420.h:84
double e2Fluct
Definition: LandauFP420.h:95
double SampleFluctuations(const double momentum, const double mass, double &tmax, const double length, const double meanLoss)
Definition: LandauFP420.cc:99
double f1Fluct
Definition: LandauFP420.h:92
double f2Fluct
Definition: LandauFP420.h:93
double problim
Definition: LandauFP420.h:105
double ipotLogFluct
Definition: LandauFP420.h:99
double p3[4]
Definition: TauolaWrapper.h:91