CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
GflashHadronShowerProfile.cc
Go to the documentation of this file.
2 
6 
7 #include <CLHEP/GenericFunctions/IncompleteGamma.hh>
8 #include <CLHEP/GenericFunctions/LogGamma.hh>
9 #include <CLHEP/Units/PhysicalConstants.h>
10 #include <CLHEP/Units/SystemOfUnits.h>
11 #include <CLHEP/Random/Randomize.h>
12 #include <CLHEP/Random/RandGaussQ.h>
13 #include <CLHEP/Random/RandPoissonQ.h>
14 
15 #include <math.h>
16 
18 {
19  theBField = parSet.getParameter<double>("bField");
20  theGflashHcalOuter = parSet.getParameter<bool>("GflashHcalOuter");
21 
22  theShowino = new GflashShowino();
24 }
25 
27 {
28  if(theShowino) delete theShowino;
29 }
30 
31 void GflashHadronShowerProfile::initialize(int showerType, double energy, double globalTime,double charge,
32  Gflash3Vector &position,Gflash3Vector &momentum) {
33 
34  //initialize GflashShowino for this track
35  theShowino->initialize(showerType, energy, globalTime, charge,
36  position, momentum, theBField);
37 
38 }
39 
41 {
42  // The skeleton of this method is based on the fortran code gfshow.F originally written
43  // by S. Peters and G. Grindhammer (also see NIM A290 (1990) 469-488), but longitudinal
44  // parameterizations of hadron showers are significantly modified for the CMS calorimeter
45 
46  // unit convention: energy in [GeV] and length in [cm]
47  // intrinsic properties of hadronic showers (lateral shower profile)
48 
49  // The step size of showino along the helix trajectory in cm unit
50  double showerDepthR50 = 0.0;
51  bool firstHcalHit = true;
52 
53  //initial valuses that will be changed as the shower developes
54  double stepLengthLeft = theShowino->getStepLengthToOut();
55 
56  double deltaStep = 0.0;
57  double showerDepth = 0.0;
58 
60 
61  theGflashHitList.clear();
62 
63  GflashHit aHit;
64 
65  while(stepLengthLeft > 0.0) {
66 
67  // update shower depth and stepLengthLeft
68  if ( stepLengthLeft < Gflash::divisionStep ) {
69  deltaStep = stepLengthLeft;
70  stepLengthLeft = 0.0;
71  }
72  else {
73  deltaStep = Gflash::divisionStep;
74  stepLengthLeft -= deltaStep;
75  }
76 
77  showerDepth += deltaStep;
78  showerDepthR50 += deltaStep;
79 
80  //update GflashShowino
81  theShowino->updateShowino(deltaStep);
82 
83  //evaluate energy in this deltaStep along the longitudinal shower profile
84  double heightProfile = 0.;
85  double deltaEnergy = 0.;
86 
88 
89  //skip if Showino is outside envelopes
90  if(whichCalor == Gflash::kNULL ) continue;
91 
92  heightProfile = longitudinalProfile();
93 
94  //skip if the delta energy for this step will be very small
95  if(heightProfile < 1.00e-08 ) continue;
96 
97  //get energy deposition for this step
98  deltaEnergy = heightProfile*Gflash::divisionStep*energyScale[whichCalor];
99  theShowino->addEnergyDeposited(deltaEnergy);
100 
101  //apply sampling fluctuation if showino is inside the sampling calorimeter
102  double hadronicFraction = 1.0;
103  double fluctuatedEnergy = deltaEnergy;
104 
105  int nSpotsInStep = std::max(1,static_cast<int>(getNumberOfSpots(whichCalor)*(deltaEnergy/energyScale[whichCalor]) ));
106  double sampleSpotEnergy = hadronicFraction*fluctuatedEnergy/nSpotsInStep;
107 
108  // Lateral shape and fluctuations
109 
110  if((whichCalor==Gflash::kHB || whichCalor==Gflash::kHE) && firstHcalHit) {
111  firstHcalHit = false;
112  //reset the showerDepth used in the lateral parameterization inside Hcal
113  showerDepthR50 = Gflash::divisionStep;
114  }
115 
116  //evaluate the fluctuated median of the lateral distribution, R50
117  double R50 = medianLateralArm(showerDepthR50,whichCalor);
118 
119  double hitEnergy = sampleSpotEnergy*CLHEP::GeV;
120  double hitTime = theShowino->getGlobalTime()*CLHEP::nanosecond;
121 
123 
124  for (int ispot = 0 ; ispot < nSpotsInStep ; ispot++) {
125 
126  // Compute global position of generated spots with taking into account magnetic field
127  // Divide deltaStep into nSpotsInStep and give a spot a global position
128  double incrementPath = theShowino->getPathLength()
129  + (deltaStep/nSpotsInStep)*(ispot+0.5 - 0.5*nSpotsInStep);
130 
131  // trajectoryPoint give a spot an imaginary point along the shower development
132  GflashTrajectoryPoint trajectoryPoint;
133  theShowino->getHelix()->getGflashTrajectoryPoint(trajectoryPoint,incrementPath);
134 
135  Gflash3Vector hitPosition = locateHitPosition(trajectoryPoint,R50);
136 
137  hitCalor = Gflash::getCalorimeterNumber(hitPosition);
138 
139  if( hitCalor == Gflash::kNULL) continue;
140 
141  hitPosition *= CLHEP::cm;
142 
143  aHit.setTime(hitTime);
144  aHit.setEnergy(hitEnergy);
145  aHit.setPosition(hitPosition);
146  theGflashHitList.push_back(aHit);
147 
148  } // end of for spot iteration
149  } // end of while for longitudinal integration
150 
151  //HO parameterization
152 
154  fabs(theShowino->getPositionAtShower().pseudoRapidity()) < Gflash::EtaMax[Gflash::kHO] &&
156 
157  //non zero ho fraction to simulate based on geant4
158  double nonzeroProb = 0.7*fTanh(theShowino->getEnergy(),Gflash::ho_nonzero);
159  double r0 = CLHEP::HepUniformRand();
160  double leftoverE = theShowino->getEnergy() - theShowino->getEnergyDeposited();
161 
162  //@@@ nonzeroProb is not random - need further correlation for non-zero HO energy
163 
164  if( r0 < nonzeroProb && leftoverE > 0.0) {
165 
166  //starting path Length and stepLengthLeft
168  stepLengthLeft = theShowino->getHelix()->getPathLengthAtRhoEquals(Gflash::Rmax[Gflash::kHO]+10) - pathLength;
169  showerDepth = pathLength - theShowino->getPathLengthAtShower();
170 
171  theShowino->setPathLength(pathLength);
172 
175 
176  while(stepLengthLeft > 0.0) {
177 
178  // update shower depth and stepLengthLeft
179  if ( stepLengthLeft < Gflash::divisionStep ) {
180  deltaStep = stepLengthLeft;
181  stepLengthLeft = 0.0;
182  }
183  else {
184  deltaStep = Gflash::divisionStep;
185  stepLengthLeft -= deltaStep;
186  }
187 
188  showerDepth += deltaStep;
189 
190  //update GflashShowino
191  theShowino->updateShowino(deltaStep);
192 
193  //evaluate energy in this deltaStep along the longitudinal shower profile
194  double heightProfile = 0.;
195  double deltaEnergy = 0.;
196 
197  double hoScale = leftoverE*(pathLengthx-pathLengthy)/(pathLengthx- theShowino->getPathLengthAtShower());
198  double refDepth = theShowino->getPathLength()
200 
201  if( refDepth > 0) {
202  heightProfile = hoProfile(theShowino->getPathLength(),refDepth);
203  deltaEnergy = heightProfile*Gflash::divisionStep*hoScale;
204  }
205 
206  int nSpotsInStep = std::max(50,static_cast<int>((160.+40* CLHEP::RandGaussQ::shoot())*std::log(theShowino->getEnergy())+50.));
207 
208  double hoFraction = 1.00;
209  double poissonProb = CLHEP::RandPoissonQ::shoot(1.0);
210 
211  double fluctuatedEnergy = deltaEnergy*poissonProb;
212  double sampleSpotEnergy = hoFraction*fluctuatedEnergy/nSpotsInStep;
213 
214  // Lateral shape and fluctuations
215 
216  //evaluate the fluctuated median of the lateral distribution, R50
217  double R50 = medianLateralArm(showerDepth,Gflash::kHB);
218 
219  double hitEnergy = sampleSpotEnergy*CLHEP::GeV;
220  double hitTime = theShowino->getGlobalTime()*CLHEP::nanosecond;
221 
222  for (int ispot = 0 ; ispot < nSpotsInStep ; ispot++) {
223 
224  double incrementPath = theShowino->getPathLength()
225  + (deltaStep/nSpotsInStep)*(ispot+0.5 - 0.5*nSpotsInStep);
226 
227  // trajectoryPoint give a spot an imaginary point along the shower development
228  GflashTrajectoryPoint trajectoryPoint;
229  theShowino->getHelix()->getGflashTrajectoryPoint(trajectoryPoint,incrementPath);
230 
231  Gflash3Vector hitPosition = locateHitPosition(trajectoryPoint,R50);
232  hitPosition *= CLHEP::cm;
233 
234  if(std::fabs(hitPosition.getZ()/CLHEP::cm) > Gflash::Zmax[Gflash::kHO] ) continue;
235 
236  aHit.setTime(hitTime);
237  aHit.setEnergy(hitEnergy);
238  aHit.setPosition(hitPosition);
239  theGflashHitList.push_back(aHit);
240 
241  } // end of for HO spot iteration
242  } // end of while for HO longitudinal integration
243  }
244  }
245 
246  // delete theGflashNavigator;
247 
248 }
249 
251 {
252  edm::LogInfo("SimGeneralGFlash") << "GflashHadronShowerProfile::loadParameters() "
253  << "should be implimented for each particle type";
254 }
255 
257 {
258  double lateralArm = 0.0;
259  if(kCalor != Gflash::kNULL) {
260 
261  double showerDepthR50X = std::min(showerDepthR50/22.4, Gflash::maxShowerDepthforR50);
262  double R50 = lateralPar[kCalor][0] + std::max(0.0,lateralPar[kCalor][1]) * showerDepthR50X;
263  double varinanceR50 = std::pow((lateralPar[kCalor][2] + lateralPar[kCalor][3] * showerDepthR50X) * R50, 2);
264 
265  // Simulation of lognormal distribution
266 
267  if(R50>0) {
268  double sigmaSq = std::log(varinanceR50/(R50*R50)+1.0);
269  double sigmaR50 = std::sqrt(sigmaSq);
270  double meanR50 = std::log(R50) - (sigmaSq/2.);
271 
272  lateralArm = std::exp(meanR50 + sigmaR50* CLHEP::RandGaussQ::shoot());
273  }
274  }
275  return lateralArm;
276 }
277 
279 {
280  // Smearing in r according to f(r)= 2.*r*R50**2/(r**2+R50**2)**2
281  double rnunif = CLHEP::HepUniformRand();
282  double rxPDF = std::sqrt(rnunif/(1.-rnunif));
283  double rShower = lateralArm*rxPDF;
284 
285  //rShower within maxLateralArmforR50
286  rShower = std::min(Gflash::maxLateralArmforR50,rShower);
287 
288  // Uniform smearing in phi
289  double azimuthalAngle = CLHEP::twopi*CLHEP::HepUniformRand();
290 
291  // actual spot position by adding a radial vector to a trajectoryPoint
293  rShower*std::cos(azimuthalAngle)*point.getOrthogonalUnitVector() +
294  rShower*std::sin(azimuthalAngle)*point.getCrossUnitVector();
295 
296  //@@@debugging histograms
297  if(theHisto->getStoreFlag()) {
298  theHisto->rshower->Fill(rShower);
299  theHisto->lateralx->Fill(rShower*std::cos(azimuthalAngle));
300  theHisto->lateraly->Fill(rShower*std::sin(azimuthalAngle));
301  }
302  return position;
303 }
304 
305 //double GflashHadronShowerProfile::longitudinalProfile(double showerDepth, double pathLength) {
307 
308  double heightProfile = 0.0;
309 
311  int showerType = theShowino->getShowerType();
312  double showerDepth = theShowino->getDepth();
313  double transDepth = theShowino->getStepLengthToHcal();
314 
315  // Energy in a delta step (dz) = (energy to deposite)*[Gamma(z+dz)-Gamma(z)]*dz
316  // where the incomplete Gamma function gives an intergrate probability of the longitudinal
317  // shower up to the shower depth (z).
318  // Instead, we use approximated energy; energy in dz = (energy to deposite)*gamma(z)*dz
319  // where gamma is the Gamma-distributed probability function
320 
322 
323  if(showerType == 0 || showerType == 4 ) {
325  if(shiftDepth > 0 ) {
326  heightProfile = twoGammaProfile(longEcal,showerDepth-shiftDepth,whichCalor);
327  }
328  else {
329  heightProfile = 0.;
330  // std::cout << "negative shiftDepth for showerType 0 " << shiftDepth << std::endl;
331  }
332  }
333  else if(showerType == 1 || showerType == 5 ) {
334  if(whichCalor == Gflash::kESPM || whichCalor == Gflash::kENCA ) {
335  heightProfile = twoGammaProfile(longEcal,showerDepth,whichCalor);
336  }
337  else if(whichCalor == Gflash::kHB || whichCalor == Gflash::kHE) {
338  heightProfile = twoGammaProfile(longHcal,showerDepth-transDepth,whichCalor);
339  }
340  else heightProfile = 0.;
341  }
342  else if (showerType == 2 || showerType == 6 ) {
343  //two gammas between crystal and Hcal
344  if((showerDepth - transDepth) > 0.0) {
345  heightProfile = twoGammaProfile(longHcal,showerDepth-transDepth,Gflash::kHB);
346  }
347  else heightProfile = 0.;
348  }
349  else if (showerType == 3 || showerType == 7 ) {
350  //two gammas inside Hcal
351  heightProfile = twoGammaProfile(longHcal,showerDepth,Gflash::kHB);
352  }
353 
354  return heightProfile;
355 }
356 
357 double GflashHadronShowerProfile::hoProfile(double pathLength, double refDepth) {
358 
359  double heightProfile = 0;
360 
361  GflashTrajectoryPoint tempPoint;
362  theShowino->getHelix()->getGflashTrajectoryPoint(tempPoint,pathLength);
363 
364  double dint = 1.4*Gflash::intLength[Gflash::kHO]*std::sin(tempPoint.getPosition().getTheta());
365  heightProfile = std::exp(-1.0*refDepth/dint);
366 
367  return heightProfile;
368 }
369 
370 void GflashHadronShowerProfile::getFluctuationVector(double *lowTriangle, double *correlationVector) {
371 
372  const int dim = Gflash::NPar;
373 
374  double **xr = new double *[dim];
375  double **xrho = new double *[dim];
376 
377  for(int j=0;j<dim;j++) {
378  xr[j] = new double [dim];
379  xrho[j] = new double [dim];
380  }
381 
382  for(int i = 0; i < dim; i++) {
383  for(int j = 0; j < i+1 ; j++) {
384  if(j==i) xrho[i][j] = 1.0;
385  else {
386  xrho[i][j] = lowTriangle[i*(i-1)/2 + j];
387  xrho[j][i] = xrho[i][j];
388  }
389  }
390  }
391 
392  doCholeskyReduction(xrho,xr,dim);
393 
394  for(int i = 0 ; i < dim ; i++) {
395  for (int j = 0 ; j < i+1 ; j++){
396  correlationVector[i*(i+1)/2 + j] = xr[i][j];
397  }
398  }
399 
400  for(int j=0;j<dim;j++) delete [] xr[j];
401  delete [] xr;
402  for(int j=0;j<dim;j++) delete [] xrho[j];
403  delete [] xrho;
404 }
405 
406 void GflashHadronShowerProfile::doCholeskyReduction(double **vv, double **cc, const int ndim) {
407 
408  double sumCjkSquare;
409  double vjjLess;
410  double sumCikjk;
411 
412  cc[0][0] = std::sqrt(vv[0][0]);
413 
414  for(int j=1 ; j < ndim ; j++) {
415  cc[j][0] = vv[j][0]/cc[0][0];
416  }
417 
418  for(int j=1 ; j < ndim ; j++) {
419 
420  sumCjkSquare = 0.0;
421  for (int k=0 ; k < j ; k++) sumCjkSquare += cc[j][k]*cc[j][k];
422 
423  vjjLess = vv[j][j] - sumCjkSquare;
424 
425  //check for the case that vjjLess is negative
426  cc[j][j] = std::sqrt(std::fabs(vjjLess));
427 
428  for (int i=j+1 ; i < ndim ; i++) {
429  sumCikjk = 0.;
430  for(int k=0 ; k < j ; k++) sumCikjk += cc[i][k]*cc[j][k];
431  cc[i][j] = (vv[i][j] - sumCikjk)/cc[j][j];
432  }
433  }
434 }
435 
437  //generator number of spots: energy dependent Gamma distribution of Nspots based on Geant4
438  //replacing old parameterization of H1,
439  //int numberOfSpots = std::max( 50, static_cast<int>(80.*std::log(einc)+50.));
440 
441  double einc = theShowino->getEnergy();
442  int showerType = theShowino->getShowerType();
443 
444  int numberOfSpots = 0;
445  double nmean = 0.0;
446  double nsigma = 0.0;
447 
448  if(showerType == 0 || showerType == 1 || showerType == 4 || showerType == 5 ) {
449  if(kCalor == Gflash::kESPM || kCalor == Gflash::kENCA) {
450  nmean = 10000 + 5000*log(einc);
451  nsigma = 1000;
452  }
453  if(kCalor == Gflash::kHB || kCalor == Gflash::kHE) {
454  nmean = 5000 + 2500*log(einc);
455  nsigma = 500;
456  }
457  }
458  else if (showerType == 2 || showerType == 3 || showerType == 6 || showerType == 7 ) {
459  if(kCalor == Gflash::kHB || kCalor == Gflash::kHE) {
460  nmean = 5000 + 2500*log(einc);
461  nsigma = 500;
462  }
463  else {
464  nmean = 10000;
465  nsigma = 1000;
466  }
467  }
468  //@@@need correlation and individual fluctuation on alphaNspots and betaNspots here:
469  //evaluating covariance should be straight forward since the distribution is 'one' Gamma
470 
471  numberOfSpots = std::max(500,static_cast<int> (nmean+nsigma* CLHEP::RandGaussQ::shoot()));
472 
473  //until we optimize the reduction scale in the number of Nspots
474 
475  if( kCalor == Gflash::kESPM || kCalor == Gflash::kENCA) {
476  numberOfSpots = static_cast<int>(numberOfSpots/100);
477  }
478  else {
479  numberOfSpots = static_cast<int>(numberOfSpots/3.0);
480  }
481 
482  return numberOfSpots;
483 }
484 
485 double GflashHadronShowerProfile::fTanh(double einc, const double *par) {
486  double func = 0.0;
487  if(einc>0.0) func = par[0]+par[1]*std::tanh(par[2]*(std::log(einc)-par[3])) + par[4]*std::log(einc);
488  return func;
489 }
490 
491 double GflashHadronShowerProfile::fLnE1(double einc, const double *par) {
492  double func = 0.0;
493  if(einc>0.0) func = par[0]+par[1]*std::log(einc);
494  return func;
495 }
496 
497 double GflashHadronShowerProfile::depthScale(double ssp, double ssp0, double length) {
498  double func = 0.0;
499  if(length>0.0) func = std::pow((ssp-ssp0)/length,2.0);
500  return func;
501 }
502 
503 double GflashHadronShowerProfile::gammaProfile(double alpha, double beta, double showerDepth, double lengthUnit) {
504  double gamma = 0.0;
505  // if(alpha > 0 && beta > 0 && lengthUnit > 0) {
506  if(showerDepth>0.0) {
507  Genfun::LogGamma lgam;
508  double x = showerDepth*(beta/lengthUnit);
509  gamma = (beta/lengthUnit)*std::pow(x,alpha-1.0)*std::exp(-x)/std::exp(lgam(alpha));
510  }
511  return gamma;
512 }
513 
514 double GflashHadronShowerProfile::twoGammaProfile(double *longPar, double depth, Gflash::CalorimeterNumber kIndex) {
515  double twoGamma = 0.0;
516 
517  longPar[0] = std::min(1.0,longPar[0]);
518  longPar[0] = std::max(0.0,longPar[0]);
519 
520  if(longPar[3] > 4.0 || longPar[4] > 4.0) {
521  double rfactor = 2.0/std::max(longPar[3],longPar[4]);
522  longPar[3] = rfactor*(longPar[3]+1.0);
523  longPar[4] *= rfactor;
524  }
525 
526  twoGamma = longPar[0]*gammaProfile(exp(longPar[1]),exp(longPar[2]),depth,Gflash::radLength[kIndex])
527  +(1-longPar[0])*gammaProfile(exp(longPar[3]),exp(longPar[4]),depth,Gflash::intLength[kIndex]);
528  return twoGamma;
529 }
Gflash3Vector locateHitPosition(GflashTrajectoryPoint &point, double lateralArm)
const double beta
double fTanh(double einc, const double *par)
T getParameter(std::string const &) const
static GflashHistogram * instance()
int i
Definition: DBlmapReader.cc:9
void setEnergy(const double energy)
Definition: GflashHit.h:19
float alpha
Definition: AMPTWrapper.h:95
const double divisionStep
double lateralPar[Gflash::kNumberCalorimeter][Gflash::Nrpar]
const double Zmax[kNumberCalorimeter]
Gflash3Vector getCrossUnitVector()
double twoGammaProfile(double *par, double depth, Gflash::CalorimeterNumber kIndex)
Gflash3Vector & getPosition()
Definition: GflashShowino.h:34
GflashTrajectory * getHelix()
Definition: GflashShowino.h:30
double medianLateralArm(double depth, Gflash::CalorimeterNumber kCalor)
double getEnergy()
Definition: GflashShowino.h:24
int getShowerType()
Definition: GflashShowino.h:23
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
double getPathLength()
Definition: GflashShowino.h:33
void setPosition(const Gflash3Vector &pos)
Definition: GflashHit.h:20
const double EtaMax[kNumberCalorimeter]
double charge(const std::vector< uint8_t > &Ampls)
static int position[TOTALCHAMBERS][3]
Definition: ReadPGInfo.cc:509
void doCholeskyReduction(double **cc, double **vv, const int ndim)
void addEnergyDeposited(double energy)
Definition: GflashShowino.h:41
double getStepLengthToHcal()
Definition: GflashShowino.h:28
void setTime(const double time)
Definition: GflashHit.h:18
const double ho_nonzero[5]
CalorimeterNumber getCalorimeterNumber(const Gflash3Vector &position)
const T & max(const T &a, const T &b)
double getPathLengthAtShower()
Definition: GflashShowino.h:26
double fLnE1(double einc, const double *par)
T sqrt(T t)
Definition: SSEVec.h:48
double getEnergyDeposited()
Definition: GflashShowino.h:35
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
const double MinEnergyCutOffForHO
int j
Definition: DBlmapReader.cc:9
void setPathLength(double pathLength)
Definition: GflashShowino.h:39
double gammaProfile(double alpha, double beta, double depth, double lengthUnit)
Gflash3Vector & getPositionAtShower()
Definition: GflashShowino.h:27
const int NPar
Gflash3Vector getOrthogonalUnitVector()
int k[5][pyjets_maxn]
const double intLength[kNumberCalorimeter]
void getGflashTrajectoryPoint(GflashTrajectoryPoint &point, double s) const
std::vector< GflashHit > theGflashHitList
const double maxShowerDepthforR50
double getDepth()
Definition: GflashShowino.h:36
const double radLength[kNumberCalorimeter]
GflashHadronShowerProfile(const edm::ParameterSet &parSet)
CLHEP::Hep3Vector Gflash3Vector
Definition: Gflash3Vector.h:6
double getGlobalTime()
Definition: GflashShowino.h:32
double getPathLengthOnEcal()
Definition: GflashShowino.h:25
Gflash3Vector & getPosition()
const double maxLateralArmforR50
void getFluctuationVector(double *lowTriangle, double *correlationVector)
double depthScale(double ssp, double ssp0, double length)
const double Rmax[kNumberCalorimeter]
Definition: DDAxes.h:10
int getNumberOfSpots(Gflash::CalorimeterNumber kCalor)
double hoProfile(double pathLength, double refDepth)
double getStepLengthToOut()
Definition: GflashShowino.h:29
const double Rmin[kNumberCalorimeter]
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
void initialize(int showerType, double energy, double globalTime, double charge, Gflash3Vector &position, Gflash3Vector &momentum)
double getPathLengthAtRhoEquals(double rho) const
void initialize(int showerType, double energy, double globalTime, double charge, Gflash3Vector &position, Gflash3Vector &momentum, double magneticField)
double energyScale[Gflash::kNumberCalorimeter]
void updateShowino(double deltaStep)