CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
HDShower.cc
Go to the documentation of this file.
1 //updated by Reza Goldouzian
2 //FastSimulation Headers
5 
8 
9 // CMSSW headers
11 
12 #include <cmath>
13 
14 // number attempts for transverse distribution if exit on a spec. condition
15 #define infinity 10000
16 // debugging flag ( 0, 1, 2, 3)
17 #define debug 0
18 
19 using namespace edm;
20 
22  HDShowerParametrization* myParam,
23  EcalHitMaker* myGrid,
24  HcalHitMaker* myHcalHitMaker,
25  int onECAL,
26  double epart,
27  double pmip)
28  : theParam(myParam),
29  theGrid(myGrid),
30  theHcalHitMaker(myHcalHitMaker),
31  onEcal(onECAL),
32  e(epart),
33  // pmip(pmip),
34  random(engine) {
35  // To get an access to constants read in FASTCalorimeter
36  // FASTCalorimeter * myCalorimeter= FASTCalorimeter::instance();
37 
38  // Values taken from FamosGeneric/FamosCalorimeter/src/FASTCalorimeter.cc
39  lossesOpt = myParam->hsParameters()->getHDlossesOpt();
41  nTRsteps = myParam->hsParameters()->getHDnTRsteps();
42  transParam = myParam->hsParameters()->getHDtransParam();
43  eSpotSize = myParam->hsParameters()->getHDeSpotSize();
44  depthStep = myParam->hsParameters()->getHDdepthStep();
47  balanceEH = myParam->hsParameters()->getHDbalanceEH();
49 
50  // Special tr.size fluctuations
51  transParam *= (1. + random->flatShoot());
52 
53  // Special long. fluctuations
54  hcalDepthFactor += 0.05 * (2. * random->flatShoot() - 1.);
55 
56  transFactor = 1.; // normally 1, in HF - might be smaller
57  // to take into account
58  // a narrowness of the HF shower (Cherenkov light)
59 
60  // simple protection ...
61  if (e < 0)
62  e = 0.;
63 
64  // Get the Famos Histos pointer
65  // myHistos = FamosHistos::instance();
66  // std::cout << " Hello FamosShower " << std::endl;
67 
70 
71  double emax = theParam->emax();
72  double emid = theParam->emid();
73  double emin = theParam->emin();
74  double effective = e;
75 
76  if (e < emid) {
77  theParam->setCase(1);
78  // avoid "underflow" below Emin (for parameters calculation only)
79  if (e < emin)
80  effective = emin;
81  } else
82  theParam->setCase(2);
83 
84  // Avoid "overflow" beyond Emax (for parameters)
85  if (effective > 0.5 * emax) {
86  eSpotSize *= 2.5;
87  if (effective > emax) {
88  effective = emax;
89  eSpotSize *= 4.;
90  depthStep *= 2.;
91  if (effective > 1000.)
92  eSpotSize *= 2.;
93  }
94  }
95 
96  if (debug == 2)
97  LogInfo("FastCalorimetry") << " HDShower : " << std::endl
98  << " Energy " << e << std::endl
99  << " lossesOpt " << lossesOpt << std::endl
100  << " nDepthSteps " << nDepthSteps << std::endl
101  << " nTRsteps " << nTRsteps << std::endl
102  << " transParam " << transParam << std::endl
103  << " eSpotSize " << eSpotSize << std::endl
104  << " criticalEnergy " << criticalEnergy << std::endl
105  << " maxTRfactor " << maxTRfactor << std::endl
106  << " balanceEH " << balanceEH << std::endl
107  << "hcalDepthFactor " << hcalDepthFactor << std::endl;
108 
109  double alpEM1 = theParam->alpe1();
110  double alpEM2 = theParam->alpe2();
111 
112  double betEM1 = theParam->bete1();
113  double betEM2 = theParam->bete2();
114 
115  double alpHD1 = theParam->alph1();
116  double alpHD2 = theParam->alph2();
117 
118  double betHD1 = theParam->beth1();
119  double betHD2 = theParam->beth2();
120 
121  double part1 = theParam->part1();
122  double part2 = theParam->part2();
123 
124  aloge = std::log(effective);
125 
126  double edpar = (theParam->e1() + aloge * theParam->e2()) * effective;
127  double aedep = std::log(edpar);
128 
129  if (debug == 2)
130  LogInfo("FastCalorimetry") << " HDShower : " << std::endl
131  << " edpar " << edpar << " aedep " << aedep << std::endl
132  << " alpEM1 " << alpEM1 << std::endl
133  << " alpEM2 " << alpEM2 << std::endl
134  << " betEM1 " << betEM1 << std::endl
135  << " betEM2 " << betEM2 << std::endl
136  << " alpHD1 " << alpHD1 << std::endl
137  << " alpHD2 " << alpHD2 << std::endl
138  << " betHD1 " << betHD1 << std::endl
139  << " betHD2 " << betHD2 << std::endl
140  << " part1 " << part1 << std::endl
141  << " part2 " << part2 << std::endl;
142 
143  // private members to set
144  theR1 = theParam->r1();
145  theR2 = theParam->r2();
146  theR3 = theParam->r3();
147 
148  alpEM = alpEM1 + alpEM2 * aedep;
149  tgamEM = tgamma(alpEM);
150  betEM = betEM1 - betEM2 * aedep;
151  alpHD = alpHD1 + alpHD2 * aedep;
152  tgamHD = tgamma(alpHD);
153  betHD = betHD1 - betHD2 * aedep;
154  part = part1 - part2 * aedep;
155  if (part > 1.)
156  part = 1.; // protection - just in case of
157 
158  if (debug == 2)
159  LogInfo("FastCalorimetry") << " HDShower : " << std::endl
160  << " alpEM " << alpEM << std::endl
161  << " tgamEM " << tgamEM << std::endl
162  << " betEM " << betEM << std::endl
163  << " alpHD " << alpHD << std::endl
164  << " tgamHD " << tgamHD << std::endl
165  << " betHD " << betHD << std::endl
166  << " part " << part << std::endl;
167 
168  if (onECAL) {
171  } else {
172  lambdaEM = 0.;
173  x0EM = 0.;
174  }
177 
178  if (debug == 2)
179  LogInfo("FastCalorimetry") << " HDShower e " << e << std::endl
180  << " x0EM = " << x0EM << std::endl
181  << " x0HD = " << x0HD << std::endl
182  << " lamEM = " << lambdaEM << std::endl
183  << " lamHD = " << lambdaHD << std::endl;
184 
185  // Starting point of the shower
186  // try first with ECAL lambda
187 
188  double sum1 = 0.; // lambda path from the ECAL/HF entrance;
189  double sum2 = 0.; // lambda path from the interaction point;
190  double sum3 = 0.; // x0 path from the interaction point;
191  int nsteps = 0; // full number of longitudinal steps (counter);
192 
193  int nmoresteps; // how many longitudinal steps in addition to
194  // one (if interaction happens there) in ECAL
195 
196  mip = 1; // just to initiate particle as MIP in ECAL
197 
198  if (e < criticalEnergy)
199  nmoresteps = 1;
200  else
201  nmoresteps = nDepthSteps;
202 
203  depthECAL = 0.;
204  depthGAP = 0.;
205  depthGAPx0 = 0.;
206  if (onECAL) {
207  depthECAL = theGrid->ecalTotalL0(); // ECAL depth segment
208  //depthECAL = theGrid->ecalTotalL0() + theGrid->ps1TotalL0() + theGrid->ps2TotalL0() + theGrid->ps2eeTotalL0(); //TEST: include preshower
209  depthGAP = theGrid->ecalHcalGapTotalL0(); // GAP depth segment
210  depthGAPx0 = theGrid->ecalHcalGapTotalX0(); // GAP depth x0
211  }
212 
213  depthHCAL = theGrid->hcalTotalL0(); // HCAL depth segment
215 
216  //---------------------------------------------------------------------------
217  // Depth simulation & various protections, among them
218  // if too deep - get flat random in the allowed region
219  // if no HCAL material behind - force to deposit in ECAL
220  double maxDepth = depthToHCAL + depthHCAL - 1.1 * depthStep;
221  double depthStart = std::log(1. / random->flatShoot()); // starting point lambda unts
222  //std::cout<<"generated depth "<<depthStart<<std::endl;
223  if (pmip == 1) {
224  depthStart = depthToHCAL;
225  } else {
226  depthStart = depthStart * 0.9 * depthECAL / std::log(1. / pmip);
227  }
228  // std::cout<<"modified depth "<<depthStart<<std::endl;
229 
230  if (e < emin) {
231  if (debug)
232  LogInfo("FastCalorimetry") << " FamosHDShower : e <emin -> depthStart = 0" << std::endl;
233  depthStart = 0.;
234  }
235 
236  if (depthStart > maxDepth) {
237  if (debug)
238  LogInfo("FastCalorimetry") << " FamosHDShower : depthStart too big ... = " << depthStart << std::endl;
239  depthStart = maxDepth * random->flatShoot();
240  if (depthStart < 0.)
241  depthStart = 0.;
242  if (debug)
243  LogInfo("FastCalorimetry") << " FamosHDShower : depthStart re-calculated = " << depthStart << std::endl;
244  }
245 
246  if (onECAL && e < emid) {
247  if (depthECAL > depthStep && (depthECAL - depthStart) / depthECAL > 0.2) {
248  depthStart = 0.5 * depthECAL * random->flatShoot();
249  if (debug)
250  LogInfo("FastCalorimetry") << " FamosHDShower : small energy, "
251  << " depthStart reduced to = " << depthStart << std::endl;
252  }
253  }
254 
255  if (depthHCAL < depthStep) {
256  if (debug)
257  LogInfo("FastCalorimetry") << " FamosHDShower : depthHCAL too small ... = " << depthHCAL
258  << " depthStart -> forced to 0 !!!" << std::endl;
259  depthStart = 0.;
260  nmoresteps = 0;
261 
262  if (depthECAL < depthStep) {
263  nsteps = -1;
264  LogInfo("FastCalorimetry") << " FamosHDShower : too small ECAL and HCAL depths - "
265  << " particle is lost !!! " << std::endl;
266  }
267  }
268 
269  if (debug)
270  LogInfo("FastCalorimetry") << " FamosHDShower depths(lam) - " << std::endl
271  << " ECAL = " << depthECAL << std::endl
272  << " GAP = " << depthGAP << std::endl
273  << " HCAL = " << depthHCAL << std::endl
274  << " starting point = " << depthStart << std::endl;
275 
276  if (onEcal) {
277  if (debug)
278  LogInfo("FastCalorimetry") << " FamosHDShower : onECAL" << std::endl;
279  if (depthStart < depthECAL) {
280  if (debug)
281  LogInfo("FastCalorimetry") << " FamosHDShower : depthStart < depthECAL" << std::endl;
282  if (depthECAL > depthStep && (depthECAL - depthStart) / depthECAL > 0.1) {
283  if (debug)
284  LogInfo("FastCalorimetry") << " FamosHDShower : enough space to make ECAL step" << std::endl;
285 
286  // ECAL - one step
287  nsteps++;
288  sum1 += depthECAL; // at the end of step
289  sum2 += depthECAL - depthStart;
290  sum3 += sum2 * lambdaEM / x0EM;
291  lamtotal.push_back(sum1);
292  lamdepth.push_back(sum2);
293  lamcurr.push_back(lambdaEM);
294  lamstep.push_back(depthECAL - depthStart);
295  x0depth.push_back(sum3);
296  x0curr.push_back(x0EM);
297  detector.push_back(1);
298  mip = 0;
299 
300  if (debug)
301  LogInfo("FastCalorimetry") << " FamosHDShower : "
302  << " in ECAL sum1, sum2 " << sum1 << " " << sum2 << std::endl;
303 
304  // Gap - no additional step after ECAL
305  // just move further to HCAL over the gap
306  sum1 += depthGAP;
307  sum2 += depthGAP;
308  sum3 += depthGAPx0;
309  } else { // Just shift starting point to HCAL
310  // cout << " FamosHDShower : not enough space to make ECAL step" << std::endl;
311  if (debug)
312  LogInfo("FastCalorimetry") << " FamosHDShower : goto HCAL" << std::endl;
313 
314  depthStart = depthToHCAL;
315  sum1 += depthStart;
316  }
317  } else { // GAP or HCAL
318  if (depthStart >= depthECAL && depthStart < depthToHCAL) {
319  depthStart = depthToHCAL; // just a shift to HCAL for simplicity
320  }
321  sum1 += depthStart;
322  if (debug)
323  LogInfo("FastCalorimetry") << " FamosHDShower : goto HCAL" << std::endl;
324  }
325  } else { // Forward
326  if (debug)
327  LogInfo("FastCalorimetry") << " FamosHDShower : forward" << std::endl;
328  sum1 += depthStart;
329  transFactor = 0.5; // makes narower tresverse size of shower
330  }
331 
332  for (int i = 0; i < nmoresteps; i++) {
333  sum1 += depthStep;
334  if (sum1 > (depthECAL + depthGAP + depthHCAL))
335  break;
336  sum2 += depthStep;
337  sum3 += sum2 * lambdaHD / x0HD;
338  lamtotal.push_back(sum1);
339  lamdepth.push_back(sum2);
340  lamcurr.push_back(lambdaHD);
341  lamstep.push_back(depthStep);
342  x0depth.push_back(sum3);
343  x0curr.push_back(x0HD);
344  detector.push_back(3);
345  nsteps++;
346  }
347 
348  // Make fractions of energy and transverse radii at each step
349 
350  if (nsteps > 0) {
351  makeSteps(nsteps);
352  }
353 }
354 
355 void HDShower::makeSteps(int nsteps) {
356  double sumes = 0.;
357  double sum = 0.;
358  std::vector<double> temp;
359 
360  if (debug)
361  LogInfo("FastCalorimetry") << " FamosHDShower::makeSteps - "
362  << " nsteps required : " << nsteps << std::endl;
363 
364  int count = 0;
365  for (int i = 0; i < nsteps; i++) {
366  double deplam = lamdepth[i] - 0.5 * lamstep[i];
367  double depx0 = x0depth[i] - 0.5 * lamstep[i] / x0curr[i];
368  double x = betEM * depx0;
369  double y = betHD * deplam;
370 
371  if (debug == 2)
372  LogInfo("FastCalorimetry") << " FamosHDShower::makeSteps "
373  << " - step " << i << " depx0, x = " << depx0 << ", " << x
374  << " deplam, y = " << deplam << ", " << y << std::endl;
375 
376  double est = (part * betEM * gam(x, alpEM) * lamcurr[i] / (x0curr[i] * tgamEM) +
377  (1. - part) * betHD * gam(y, alpHD) / tgamHD) *
378  lamstep[i];
379 
380  // protection ...
381  if (est < 0.) {
382  LogInfo("FastCalorimetry") << "*** FamosHDShower::makeSteps "
383  << " - negative step energy !!!" << std::endl;
384  est = 0.;
385  break;
386  }
387 
388  // for estimates only
389  sum += est;
390  int nPest = (int)(est * e / sum / eSpotSize);
391 
392  if (debug == 2)
393  LogInfo("FastCalorimetry") << " FamosHDShower::makeSteps - nPoints estimate = " << nPest << std::endl;
394 
395  if (nPest <= 1 && count != 0)
396  break;
397 
398  // good step - to proceed
399 
400  temp.push_back(est);
401  sumes += est;
402 
403  rlamStep.push_back(transParam * (theR1 + (theR2 - theR3 * aloge)) * deplam * transFactor);
404  count++;
405  }
406 
407  // fluctuations in ECAL and re-distribution of remaining energy in HCAL
408  if (detector[0] == 1 && count > 1) {
409  double oldECALenergy = temp[0];
410  double oldHCALenergy = sumes - oldECALenergy;
411  double newECALenergy = 2. * sumes;
412  for (int i = 0; newECALenergy > sumes && i < infinity; i++)
413  newECALenergy = 2. * balanceEH * random->flatShoot() * oldECALenergy;
414 
415  if (debug == 2)
416  LogInfo("FastCalorimetry") << "*** FamosHDShower::makeSteps "
417  << " ECAL fraction : old/new - " << oldECALenergy / sumes << "/"
418  << newECALenergy / sumes << std::endl;
419 
420  temp[0] = newECALenergy;
421  double newHCALenergy = sumes - newECALenergy;
422  double newHCALreweight = newHCALenergy / oldHCALenergy;
423 
424  for (int i = 1; i < count; i++) {
425  temp[i] *= newHCALreweight;
426  }
427  }
428 
429  // final re-normalization of the energy fractions
430  for (int i = 0; i < count; i++) {
431  eStep.push_back(temp[i] * e / sumes);
432  nspots.push_back((int)(eStep[i] / eSpotSize) + 1);
433 
434  if (debug)
435  LogInfo("FastCalorimetry") << " step " << i << " det: " << detector[i]
436  << " xO and lamdepth at the end of step = " << x0depth[i] << " " << lamdepth[i]
437  << " Estep func = " << eStep[i] << " Rstep = " << rlamStep[i]
438  << " Nspots = " << nspots[i] << " espot = " << eStep[i] / (double)nspots[i]
439  << std::endl;
440  }
441 
442  // The only step is in ECAL - let's make the size bigger ...
443  if (count == 1 and detector[0] == 1)
444  rlamStep[0] *= 2.;
445 
446  if (debug) {
447  if (eStep[0] > 0.95 * e && detector[0] == 1)
448  LogInfo("FastCalorimetry") << " FamosHDShower::makeSteps - "
449  << "ECAL energy = " << eStep[0] << " out of total = " << e << std::endl;
450  }
451 }
452 
454  // TimeMe theT("FamosHDShower::compute");
455 
456  bool status = false;
457  int numLongit = eStep.size();
458  if (debug)
459  LogInfo("FastCalorimetry") << " FamosHDShower::compute - "
460  << " N_long.steps required : " << numLongit << std::endl;
461 
462  if (numLongit > 0) {
463  status = true;
464  // Prepare the trsanverse probability function
465  std::vector<double> Fhist;
466  std::vector<double> rhist;
467  for (int j = 0; j < nTRsteps + 1; j++) {
468  rhist.push_back(maxTRfactor * j / nTRsteps);
469  Fhist.push_back(transProb(maxTRfactor, 1., rhist[j]));
470  if (debug == 3)
471  LogInfo("FastCalorimetry") << "indexFinder - i, Fhist[i] = " << j << " " << Fhist[j] << std::endl;
472  }
473 
474  // Longitudinal steps
475  for (int i = 0; i < numLongit; i++) {
476  double currentDepthL0 = lamtotal[i] - 0.5 * lamstep[i];
477  // vary the longitudinal profile if needed
478  if (detector[i] != 1)
479  currentDepthL0 *= hcalDepthFactor;
480  if (debug)
481  LogInfo("FastCalorimetry") << " FamosHDShower::compute - detector = " << detector[i]
482  << " currentDepthL0 = " << currentDepthL0 << std::endl;
483 
484  double maxTRsize = maxTRfactor * rlamStep[i]; // in lambda units
485  double rbinsize = maxTRsize / nTRsteps;
486  double espot = eStep[i] / (double)nspots[i]; // re-adjust espot
487 
488  if (espot > 2. || espot < 0.)
489  LogInfo("FastCalorimetry") << " FamosHDShower::compute - unphysical espot = " << espot << std::endl;
490 
491  int ecal = 0;
492  if (detector[i] != 1) {
493  bool setHDdepth = theHcalHitMaker->setDepth(currentDepthL0);
494 
495  if (debug)
496  LogInfo("FastCalorimetry") << " FamosHDShower::compute - status of "
497  << " theHcalHitMaker->setDepth(currentDepthL0) is " << setHDdepth << std::endl;
498 
499  if (!setHDdepth) {
500  currentDepthL0 -= lamstep[i];
501  setHDdepth = theHcalHitMaker->setDepth(currentDepthL0);
502  }
503 
504  if (!setHDdepth)
505  continue;
506 
508 
509  //fill hcal longitudinal distribution histogram
510  } else {
511  ecal = 1;
512  bool status = theGrid->getPads(currentDepthL0);
513 
514  if (debug)
515  LogInfo("FastCalorimetry") << " FamosHDShower::compute - status of Grid = " << status << std::endl;
516 
517  if (!status)
518  continue;
519 
520  int ntry = nspots[i] * 10;
521  if (ntry >= infinity) { // use max allowed in case of too many spots
522  nspots[i] = 0.5 * infinity;
523  espot *= 0.1 * (double)ntry / double(nspots[i]);
524  } else {
525  espot *= 0.1; // fine-grain energy spots in ECAL
526  // to avoid false ECAL clustering
527  nspots[i] = ntry;
528  }
529 
530  theGrid->setSpotEnergy(espot);
531 
532  //fill ecal longitudinal distribution histogram
533  }
534 
535  // Transverse distribition
536  int nok = 0; // counter of OK
537  int count = 0;
538  int inf = infinity;
539  if (lossesOpt)
540  inf = nspots[i]; // if losses are enabled, otherwise
541  // only OK points are counted ...
542  if (nspots[i] > inf)
543  std::cout << " FamosHDShower::compute - at long.step " << i << " too many spots required : " << nspots[i]
544  << " !!! " << std::endl;
545 
546  for (int j = 0; j < inf; j++) {
547  if (nok == nspots[i])
548  break;
549  count++;
550 
551  double prob = random->flatShoot();
552  int index = indexFinder(prob, Fhist);
553  double radius = rlamStep[i] * rhist[index] + random->flatShoot() * rbinsize; // in-bin
554  double phi = 2. * M_PI * random->flatShoot();
555 
556  if (debug == 2)
557  LogInfo("FastCalorimetry") << std::endl
558  << " FamosHDShower::compute "
559  << " r = " << radius << " phi = " << phi << std::endl;
560 
561  bool result;
562  if (ecal) {
563  result = theGrid->addHit(radius, phi, 0);
564 
565  if (debug == 2)
566  LogInfo("FastCalorimetry") << " FamosHDShower::compute - "
567  << " theGrid->addHit result = " << result << std::endl;
568 
569  //fill ecal transverse distribution histogram
570  } else {
571  result = theHcalHitMaker->addHit(radius, phi, 0);
572 
573  if (debug == 2)
574  LogInfo("FastCalorimetry") << " FamosHDShower::compute - "
575  << " theHcalHitMaker->addHit result = " << result << std::endl;
576 
577  //fill hcal transverse distribution histogram
578  }
579 
580  if (result)
581  nok++;
582 
583  } // end of tranverse simulation
584 
585  if (count == infinity) {
586  if (debug)
587  LogInfo("FastCalorimetry") << " FamosHDShower::compute "
588  << " maximum number of"
589  << " transverse points " << count << " is used !!!" << std::endl;
590  }
591 
592  if (debug)
593  LogInfo("FastCalorimetry") << " FamosHDShower::compute "
594  << " long.step No." << i << " Ntry, Nok = " << count << " " << nok << std::endl;
595  } // end of longitudinal steps
596  } // end of no steps
597 
598  return status;
599 }
600 
601 int HDShower::indexFinder(double x, const std::vector<double>& Fhist) {
602  // binary search in the vector of doubles
603  int size = Fhist.size();
604 
605  int curr = size / 2;
606  int step = size / 4;
607  int iter;
608  int prevdir = 0;
609  int actudir = 0;
610 
611  for (iter = 0; iter < size; iter++) {
612  if (curr >= size || curr < 1)
613  LogWarning("FastCalorimetry") << " FamosHDShower::indexFinder - wrong current index = " << curr << " !!!"
614  << std::endl;
615 
616  if ((x <= Fhist[curr]) && (x > Fhist[curr - 1]))
617  break;
618  prevdir = actudir;
619  if (x > Fhist[curr]) {
620  actudir = 1;
621  } else {
622  actudir = -1;
623  }
624  if (prevdir * actudir < 0) {
625  if (step > 1)
626  step /= 2;
627  }
628  curr += actudir * step;
629  if (curr > size)
630  curr = size;
631  else {
632  if (curr < 1) {
633  curr = 1;
634  }
635  }
636 
637  if (debug == 3)
638  LogInfo("FastCalorimetry") << " indexFinder - end of iter." << iter << " curr, F[curr-1], F[curr] = " << curr
639  << " " << Fhist[curr - 1] << " " << Fhist[curr] << std::endl;
640  }
641 
642  if (debug == 3)
643  LogInfo("FastCalorimetry") << " indexFinder x = " << x << " found index = " << curr - 1 << std::endl;
644 
645  return curr - 1;
646 }
HDShower(const RandomEngineAndDistribution *engine, HDShowerParametrization *myParam, EcalHitMaker *myGrid, HcalHitMaker *myHcalHitMaker, int onECAL, double epart, double pmip)
Definition: HDShower.cc:21
double depthStart
Definition: HDShower.h:75
EcalHitMaker * theGrid
Definition: HDShower.h:86
double radLenIncm() const override
Radiation length in cm.
bool addHit(double r, double phi, unsigned layer=0) override
static std::vector< std::string > checklist log
double balanceEH
Definition: HDShower.h:119
double eSpotSize
Definition: HDShower.h:111
double radLenIncm() const override
Radiation length in cm.
double flatShoot(double xmin=0.0, double xmax=1.0) const
double betEM
Definition: HDShower.h:71
double maxTRfactor
Definition: HDShower.h:117
double criticalEnergy
Definition: HDShower.h:115
const ECALProperties * theECALproperties
Definition: HDShower.h:66
int getHDnTRsteps() const
Definition: HSParameters.h:21
double getHDeSpotSize() const
Definition: HSParameters.h:23
double depthGAP
Definition: HDShower.h:127
std::vector< double > lamstep
Definition: HDShower.h:81
int nTRsteps
Definition: HDShower.h:105
double x0EM
Definition: HDShower.h:74
bool compute()
Compute the shower longitudinal and lateral development.
Definition: HDShower.cc:453
list status
Definition: mps_update.py:107
double theR3
Definition: HDShower.h:70
std::vector< double > x0curr
Definition: HDShower.h:80
double getHDhcalDepthFactor() const
Definition: HSParameters.h:28
double interactionLength() const override
Interaction length in cm.
double betHD
Definition: HDShower.h:71
double theR1
Definition: HDShower.h:70
tuple result
Definition: mps_fire.py:311
double theR2
Definition: HDShower.h:70
std::vector< double > lamdepth
Definition: HDShower.h:81
int getHDlossesOpt() const
Definition: HSParameters.h:19
bool addHit(double r, double phi, unsigned layer=0) override
add the hit in the HCAL in local coordinates
Definition: HcalHitMaker.cc:30
std::vector< double > lamcurr
Definition: HDShower.h:81
double hcalDepthFactor
Definition: HDShower.h:121
const HSParameters * hsParameters() const
int getHDnDepthSteps() const
Definition: HSParameters.h:20
double lambdaHD
Definition: HDShower.h:74
double gam(double x, double a) const
Definition: HDShower.h:48
int onEcal
Definition: HDShower.h:92
string inf
Definition: EcalCondDB.py:96
double transProb(double factor, double R, double r)
Definition: HDShower.h:53
double depthECAL
Definition: HDShower.h:127
uint16_t const *__restrict__ x
Definition: gpuClustering.h:39
std::vector< double > x0depth
Definition: HDShower.h:80
HDShowerParametrization * theParam
Definition: HDShower.h:63
int infinity
Definition: HDShower.h:83
int mip
Definition: HDShower.h:95
int nDepthSteps
Definition: HDShower.h:103
void setSpotEnergy(double e) override
Set the spot energy.
Definition: HcalHitMaker.h:26
double depthGAPx0
Definition: HDShower.h:127
double depthHCAL
Definition: HDShower.h:127
double aloge
Definition: HDShower.h:76
double getHDdepthStep() const
Definition: HSParameters.h:24
#define M_PI
double tgamHD
Definition: HDShower.h:71
Log< level::Info, false > LogInfo
double alpEM
Definition: HDShower.h:71
double getHDbalanceEH() const
Definition: HSParameters.h:27
#define debug
Definition: HDRShower.cc:19
double transParam
Definition: HDShower.h:107
double ecalHcalGapTotalX0() const
ECAL-HCAL transition.
Definition: EcalHitMaker.h:70
std::vector< int > nspots
Definition: HDShower.h:78
part
Definition: HCALResponse.h:20
double alpHD
Definition: HDShower.h:71
double depthStep
Definition: HDShower.h:113
HcalHitMaker * theHcalHitMaker
Definition: HDShower.h:89
std::vector< int > detector
Definition: HDShower.h:78
double transFactor
Definition: HDShower.h:109
std::vector< double > lamtotal
Definition: HDShower.h:81
bool getPads(double depth, bool inCm=false)
int lossesOpt
Definition: HDShower.h:101
const HCALProperties * hcalProperties() const
double getHDcriticalEnergy() const
Definition: HSParameters.h:25
const ECALProperties * ecalProperties() const
double getHDmaxTRfactor() const
Definition: HSParameters.h:26
double e
Definition: HDShower.h:98
const RandomEngineAndDistribution * random
Definition: HDShower.h:124
double depthToHCAL
Definition: HDShower.h:127
double getHDtransParam() const
Definition: HSParameters.h:22
const HCALProperties * theHCALproperties
Definition: HDShower.h:67
tuple cout
Definition: gather_cfg.py:144
std::vector< double > eStep
Definition: HDShower.h:79
int indexFinder(double x, const std::vector< double > &Fhist)
Definition: HDShower.cc:601
step
Definition: StallMonitor.cc:94
double lambdaEM
Definition: HDShower.h:74
Log< level::Warning, false > LogWarning
std::vector< double > rlamStep
Definition: HDShower.h:79
double ecalTotalL0() const
in the ECAL
Definition: EcalHitMaker.h:85
double interactionLength() const override
Interaction length in cm: 18.5 for Standard ECAL.
void makeSteps(int nsteps)
Definition: HDShower.cc:355
double tgamEM
Definition: HDShower.h:71
tuple size
Write out results.
bool setDepth(double, bool inCm=false)
set the depth in X0 or Lambda0 units depending on showerType
void setSpotEnergy(double e) override
Definition: EcalHitMaker.h:112
double ecalHcalGapTotalL0() const
ECAL-HCAL transition.
Definition: EcalHitMaker.h:91
double x0HD
Definition: HDShower.h:74
double hcalTotalL0() const
in the HCAL
Definition: EcalHitMaker.h:88