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