CMS 3D CMS Logo

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