CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HcalSimpleRecAlgo.cc
Go to the documentation of this file.
7 
8 #include <algorithm>
9 #include <cmath>
10 
11 //--- temporary for printouts
12 // #include<iostream>
13 
14 constexpr double MaximumFractionalError = 0.002; // 0.2% error allowed from this source
17 
18 HcalSimpleRecAlgo::HcalSimpleRecAlgo(bool correctForTimeslew, bool correctForPulse, float phaseNS) :
19  correctForTimeslew_(correctForTimeslew),
20  correctForPulse_(correctForPulse),
21  phaseNS_(phaseNS), runnum_(0), setLeakCorrection_(false), puCorrMethod_(0)
22 {
23  pulseCorr_ = std::auto_ptr<HcalPulseContainmentManager>(new HcalPulseContainmentManager(MaximumFractionalError));
24  pedSubFxn_ = std::auto_ptr<PedestalSub>(new PedestalSub());
25  hltOOTpuCorr_ = std::auto_ptr<HcalDeterministicFit>(new HcalDeterministicFit());
26 }
27 
28 
30 {
31  pulseCorr_->beginRun(es);
32 }
33 
34 
36 {
37  pulseCorr_->endRun();
38 }
39 
40 
42 }
43 
44 void HcalSimpleRecAlgo::setRecoParams(bool correctForTimeslew, bool correctForPulse, bool setLeakCorrection, int pileupCleaningID, float phaseNS){
46  correctForPulse_=correctForPulse;
47  phaseNS_=phaseNS;
49  pileupCleaningID_=pileupCleaningID;
50 }
51 
52 void HcalSimpleRecAlgo::setpuCorrParams(bool iPedestalConstraint, bool iTimeConstraint,bool iAddPulseJitter,
53  bool iUnConstrainedFit, bool iApplyTimeSlew,double iTS4Min, double iTS4Max,
54  double iPulseJitter,double iTimeMean,double iTimeSig,double iPedMean,double iPedSig,
55  double iNoise,double iTMin,double iTMax,
56  double its3Chi2,double its4Chi2,double its345Chi2,double iChargeThreshold, int iFitTimes) {
57  if( iPedestalConstraint ) assert ( iPedSig );
58  if( iTimeConstraint ) assert( iTimeSig );
59  psFitOOTpuCorr_->setPUParams(iPedestalConstraint,iTimeConstraint,iAddPulseJitter,iUnConstrainedFit,iApplyTimeSlew,
60  iTS4Min, iTS4Max, iPulseJitter,iTimeMean,iTimeSig,iPedMean,iPedSig,iNoise,iTMin,iTMax,its3Chi2,its4Chi2,its345Chi2,
61  iChargeThreshold,HcalTimeSlew::Medium, iFitTimes);
62 // int shapeNum = HPDShapev3MCNum;
63 // psFitOOTpuCorr_->setPulseShapeTemplate(theHcalPulseShapes_.getShape(shapeNum));
64 }
65 
66 void HcalSimpleRecAlgo::setMeth3Params(int iPedSubMethod, float iPedSubThreshold, int iTimeSlewParsType, std::vector<double> iTimeSlewPars, double irespCorrM3) {
67 
68  pedSubFxn_->init(((PedestalSub::Method)iPedSubMethod), 0, iPedSubThreshold, 0.0);
69  hltOOTpuCorr_->init((HcalTimeSlew::ParaSource)iTimeSlewParsType, HcalTimeSlew::Medium, (HcalDeterministicFit::NegStrategy)2, *pedSubFxn_, iTimeSlewPars,irespCorrM3);
70 }
71 
73  runnum_ = runnum;
74  if( puCorrMethod_ ==2 ){
75  int shapeNum = HPDShapev3MCNum;
76  if( runnum_ > 0 ){
77  shapeNum = HPDShapev3DataNum;
78  }
79  psFitOOTpuCorr_->setPulseShapeTemplate(theHcalPulseShapes_.getShape(shapeNum));
80  }
81 }
82 
84 
86  boost::shared_ptr<AbsOOTPileupCorrection> corr)
87 {
89 }
90 
92  boost::shared_ptr<AbsOOTPileupCorrection> corr)
93 {
95 }
96 
98  boost::shared_ptr<AbsOOTPileupCorrection> corr)
99 {
101 }
102 
103 void HcalSimpleRecAlgo::setBXInfo(const BunchXParameter* info,
104  const unsigned lenInfo)
105 {
107  lenBunchCrossingInfo_ = lenInfo;
108 }
109 
115 static float timeshift_ns_hbheho(float wpksamp);
116 
118 static float timeshift_ns_hf(float wpksamp);
119 
121 static float eCorr(int ieta, int iphi, double ampl, int runnum);
122 
124 static float leakCorr(double energy);
125 
126 
127 namespace HcalSimpleRecAlgoImpl {
128  template<class Digi>
129  inline float recoHFTime(const Digi& digi, const int maxI, const double amp_fC,
130  const bool slewCorrect, double maxA, float t0, float t2)
131  {
132  // Handle negative excursions by moving "zero":
133  float zerocorr=std::min(t0,t2);
134  if (zerocorr<0.f) {
135  t0 -= zerocorr;
136  t2 -= zerocorr;
137  maxA -= zerocorr;
138  }
139 
140  // pair the peak with the larger of the two neighboring time samples
141  float wpksamp=0.f;
142  if (t0>t2) {
143  wpksamp = t0+maxA;
144  if (wpksamp != 0.f) wpksamp = maxA/wpksamp;
145  } else {
146  wpksamp = maxA+t2;
147  if (wpksamp != 0.f) wpksamp = 1.+(t2/wpksamp);
148  }
149 
150  float time = (maxI - digi.presamples())*25.0 + timeshift_ns_hf(wpksamp);
151 
152  if (slewCorrect && amp_fC > 0.0) {
153  // -5.12327 - put in calibs.timecorr()
154  double tslew=exp(0.337681-5.94689e-4*amp_fC)+exp(2.44628-1.34888e-2*amp_fC);
155  time -= (float)tslew;
156  }
157 
158  return time;
159  }
160 
161 
162  template<class Digi>
163  inline void removePileup(const Digi& digi, const HcalCoder& coder,
164  const HcalCalibrations& calibs,
165  const int ifirst, const int n,
166  const bool pulseCorrect,
168  const AbsOOTPileupCorrection* pileupCorrection,
169  const BunchXParameter* bxInfo, const unsigned lenInfo,
170  double* p_maxA, double* p_ampl, double* p_uncorr_ampl,
171  double* p_fc_ampl, int* p_nRead, int* p_maxI,
172  bool* leakCorrApplied, float* p_t0, float* p_t2)
173  {
174  CaloSamples cs;
175  coder.adc2fC(digi,cs);
176  const int nRead = cs.size();
177  const int iStop = std::min(nRead, n + ifirst);
178 
179  // Signal energy will be calculated both with
180  // and without OOT pileup corrections. Try to
181  // arrange the calculations so that we do not
182  // repeat them.
183  double uncorrectedEnergy[CaloSamples::MAXSAMPLES] {}, buf[CaloSamples::MAXSAMPLES] {};
184  double* correctedEnergy = 0;
185  double fc_ampl = 0.0, corr_fc_ampl = 0.0;
186  bool pulseShapeCorrApplied = false, readjustTiming = false;
187  *leakCorrApplied = false;
188 
189  if (pileupCorrection)
190  {
191  correctedEnergy = &buf[0];
192 
193  double correctionInput[CaloSamples::MAXSAMPLES] {};
194  double gains[CaloSamples::MAXSAMPLES] {};
195 
196  for (int i=0; i<nRead; ++i)
197  {
198  const int capid = digi[i].capid();
199  correctionInput[i] = cs[i] - calibs.pedestal(capid);
200  gains[i] = calibs.respcorrgain(capid);
201  }
202 
203  for (int i=ifirst; i<iStop; ++i)
204  fc_ampl += correctionInput[i];
205 
206  const bool useGain = pileupCorrection->inputIsEnergy();
207  for (int i=0; i<nRead; ++i)
208  {
209  uncorrectedEnergy[i] = correctionInput[i]*gains[i];
210  if (useGain)
211  correctionInput[i] = uncorrectedEnergy[i];
212  }
213 
214  pileupCorrection->apply(digi.id(), correctionInput, nRead,
215  bxInfo, lenInfo, ifirst, n,
216  correctedEnergy, CaloSamples::MAXSAMPLES,
217  &pulseShapeCorrApplied, leakCorrApplied,
218  &readjustTiming);
219  if (useGain)
220  {
221  // Gain factors have been already applied.
222  // Divide by them for accumulating corr_fc_ampl.
223  for (int i=ifirst; i<iStop; ++i)
224  if (gains[i])
225  corr_fc_ampl += correctedEnergy[i]/gains[i];
226  }
227  else
228  {
229  for (int i=ifirst; i<iStop; ++i)
230  corr_fc_ampl += correctedEnergy[i];
231  for (int i=0; i<nRead; ++i)
232  correctedEnergy[i] *= gains[i];
233  }
234  }
235  else
236  {
237  correctedEnergy = &uncorrectedEnergy[0];
238 
239  // In this situation, we do not need to process all time slices
240  const int istart = std::max(ifirst - 1, 0);
241  const int iend = std::min(n + ifirst + 1, nRead);
242  for (int i=istart; i<iend; ++i)
243  {
244  const int capid = digi[i].capid();
245  float ta = cs[i] - calibs.pedestal(capid);
246  if (i >= ifirst && i < iStop)
247  fc_ampl += ta;
248  ta *= calibs.respcorrgain(capid);
249  uncorrectedEnergy[i] = ta;
250  }
251  corr_fc_ampl = fc_ampl;
252  }
253 
254  // Uncorrected and corrected energies
255  double ampl = 0.0, corr_ampl = 0.0;
256  for (int i=ifirst; i<iStop; ++i)
257  {
258  ampl += uncorrectedEnergy[i];
259  corr_ampl += correctedEnergy[i];
260  }
261 
262  // Apply phase-based amplitude correction:
263  if (corr && pulseCorrect)
264  {
265  ampl *= corr->getCorrection(fc_ampl);
266  if (pileupCorrection)
267  {
268  if (!pulseShapeCorrApplied)
269  corr_ampl *= corr->getCorrection(corr_fc_ampl);
270  }
271  else
272  corr_ampl = ampl;
273  }
274 
275  // Which energies we want to use for timing?
276  const double *etime = readjustTiming ? &correctedEnergy[0] : &uncorrectedEnergy[0];
277  int maxI = -1; double maxA = -1.e300;
278  for (int i=ifirst; i<iStop; ++i)
279  if (etime[i] > maxA)
280  {
281  maxA = etime[i];
282  maxI = i;
283  }
284 
285  // Fill out the output
286  *p_maxA = maxA;
287  *p_ampl = corr_ampl;
288  *p_uncorr_ampl = ampl;
289  *p_fc_ampl = readjustTiming ? corr_fc_ampl : fc_ampl;
290  *p_nRead = nRead;
291  *p_maxI = maxI;
292 
293  if (maxI <= 0 || maxI >= (nRead-1))
294  {
295  LogDebug("HCAL Pulse") << "HcalSimpleRecAlgoImpl::removePileup :"
296  << " Invalid max amplitude position, "
297  << " max Amplitude: " << maxI
298  << " first: " << ifirst
299  << " last: " << ifirst + n
300  << std::endl;
301  *p_t0 = 0.f;
302  *p_t2 = 0.f;
303  }
304  else
305  {
306  *p_t0 = etime[maxI - 1];
307  *p_t2 = etime[maxI + 1];
308  }
309  }
310 
311 
312  template<class Digi, class RecHit>
313  inline RecHit reco(const Digi& digi, const HcalCoder& coder,
314  const HcalCalibrations& calibs,
315  const int ifirst, const int n, const bool slewCorrect,
316  const bool pulseCorrect, const HcalPulseContainmentCorrection* corr,
317  const HcalTimeSlew::BiasSetting slewFlavor,
318  const int runnum, const bool useLeak,
319  const AbsOOTPileupCorrection* pileupCorrection,
320  const BunchXParameter* bxInfo, const unsigned lenInfo, const int puCorrMethod, const PulseShapeFitOOTPileupCorrection * psFitOOTpuCorr, HcalDeterministicFit * hltOOTpuCorr, PedestalSub * hltPedSub /* whatever don't know what to do with the pointer...*/)// const on end
321  {
322  double fc_ampl =0, ampl =0, uncorr_ampl =0, m3_ampl =0, maxA = -1.e300;
323  int nRead = 0, maxI = -1;
324  bool leakCorrApplied = false;
325  float t0 =0, t2 =0;
326  float time = -9999;
327 
328 // Disable method 1 inside the removePileup function this way!
329 // Some code in removePileup does NOT do pileup correction & to make sure maximum share of code
330  const AbsOOTPileupCorrection * inputAbsOOTpuCorr = ( puCorrMethod == 1 ? pileupCorrection: 0 );
331 
332  removePileup(digi, coder, calibs, ifirst, n,
333  pulseCorrect, corr, inputAbsOOTpuCorr,
334  bxInfo, lenInfo, &maxA, &ampl,
335  &uncorr_ampl, &fc_ampl, &nRead, &maxI,
336  &leakCorrApplied, &t0, &t2);
337 
338  if (maxI > 0 && maxI < (nRead - 1))
339  {
340  // Handle negative excursions by moving "zero":
341  float minA=t0;
342  if (maxA<minA) minA=maxA;
343  if (t2<minA) minA=t2;
344  if (minA<0) { maxA-=minA; t0-=minA; t2-=minA; } // positivizes all samples
345 
346  float wpksamp = (t0 + maxA + t2);
347  if (wpksamp!=0) wpksamp=(maxA + 2.0*t2) / wpksamp;
348  time = (maxI - digi.presamples())*25.0 + timeshift_ns_hbheho(wpksamp);
349 
350  if (slewCorrect) time-=HcalTimeSlew::delay(std::max(1.0,fc_ampl),slewFlavor);
351 
352  time=time-calibs.timecorr(); // time calibration
353  }
354 
355  // Note that uncorr_ampl is always set from outside of method 2!
356  if( puCorrMethod == 2 ){
357  std::vector<double> correctedOutput;
358 
359  CaloSamples cs;
360  coder.adc2fC(digi,cs);
361  std::vector<int> capidvec;
362  for(int ip=0; ip<cs.size(); ip++){
363  const int capid = digi[ip].capid();
364  capidvec.push_back(capid);
365  }
366  psFitOOTpuCorr->apply(cs, capidvec, calibs, correctedOutput);
367  if( correctedOutput.back() == 0 && correctedOutput.size() >1 ){
368  time = correctedOutput[1]; ampl = correctedOutput[0];
369  }
370  }
371 
372  // S. Brandt - Feb 19th : Adding Section for HLT
373  // Run "Method 3" all the time.
374  {
375  std::vector<double> hltCorrOutput;
376 
377  CaloSamples cs;
378  coder.adc2fC(digi,cs);
379  std::vector<int> capidvec;
380  for(int ip=0; ip<cs.size(); ip++){
381  const int capid = digi[ip].capid();
382  capidvec.push_back(capid);
383  }
384  hltOOTpuCorr->apply(cs, capidvec, calibs, digi, hltCorrOutput);
385  if( hltCorrOutput.size() > 1 ){
386  m3_ampl = hltCorrOutput[0];
387  if (puCorrMethod == 3) {
388  time = hltCorrOutput[1]; ampl = hltCorrOutput[0];
389  }
390  }
391  }
392 
393  // Temporary hack to apply energy-dependent corrections to some HB- cells
394  if (runnum > 0) {
395  const HcalDetId& cell = digi.id();
396  if (cell.subdet() == HcalBarrel) {
397  const int ieta = cell.ieta();
398  const int iphi = cell.iphi();
399  uncorr_ampl *= eCorr(ieta, iphi, uncorr_ampl, runnum);
400  ampl *= eCorr(ieta, iphi, ampl, runnum);
401  m3_ampl *= eCorr(ieta, iphi, m3_ampl, runnum);
402  }
403  }
404 
405  // Correction for a leak to pre-sample
406  if(useLeak && !leakCorrApplied) {
407  uncorr_ampl *= leakCorr(uncorr_ampl);
408  if (puCorrMethod < 2)
409  ampl *= leakCorr(ampl);
410  }
411 
412  RecHit rh(digi.id(),ampl,time);
413  setRawEnergy(rh, static_cast<float>(uncorr_ampl));
414  setAuxEnergy(rh, static_cast<float>(m3_ampl));
415  return rh;
416  }
417 
418  template<class Digi, class RecHit>
419  inline RecHit recoHBHE(const Digi& digi, const HcalCoder& coder,
420  const HcalCalibrations& calibs,
421  const int ifirst, const int n, const bool slewCorrect,
422  const bool pulseCorrect, const HcalPulseContainmentCorrection* corr,
423  const HcalTimeSlew::BiasSetting slewFlavor,
424  const int runnum, const bool useLeak,
425  const AbsOOTPileupCorrection* pileupCorrection,
426  const BunchXParameter* bxInfo, const unsigned lenInfo, const int puCorrMethod, const PulseShapeFitOOTPileupCorrection * psFitOOTpuCorr, HcalDeterministicFit * hltOOTpuCorr, PedestalSub * hltPedSub)// const on end
427  {
428  double fc_ampl =0, ampl =0, uncorr_ampl =0, m3_ampl =0, maxA = -1.e300;
429  int nRead = 0, maxI = -1;
430  bool leakCorrApplied = false;
431  float t0 =0, t2 =0;
432  float time = -9999;
433  bool useTriple = false;
434 
435  // Disable method 1 inside the removePileup function this way!
436  // Some code in removePileup does NOT do pileup correction & to make sure maximum share of code
437  const AbsOOTPileupCorrection * inputAbsOOTpuCorr = ( puCorrMethod == 1 ? pileupCorrection: 0 );
438 
439  removePileup(digi, coder, calibs, ifirst, n,
440  pulseCorrect, corr, inputAbsOOTpuCorr,
441  bxInfo, lenInfo, &maxA, &ampl,
442  &uncorr_ampl, &fc_ampl, &nRead, &maxI,
443  &leakCorrApplied, &t0, &t2);
444 
445  if (maxI > 0 && maxI < (nRead - 1))
446  {
447  // Handle negative excursions by moving "zero":
448  float minA=t0;
449  if (maxA<minA) minA=maxA;
450  if (t2<minA) minA=t2;
451  if (minA<0) { maxA-=minA; t0-=minA; t2-=minA; } // positivizes all samples
452 
453  float wpksamp = (t0 + maxA + t2);
454  if (wpksamp!=0) wpksamp=(maxA + 2.0*t2) / wpksamp;
455  time = (maxI - digi.presamples())*25.0 + timeshift_ns_hbheho(wpksamp);
456 
457  if (slewCorrect) time-=HcalTimeSlew::delay(std::max(1.0,fc_ampl),slewFlavor);
458 
459  time=time-calibs.timecorr(); // time calibration
460  }
461 
462  // Note that uncorr_ampl is always set from outside of method 2!
463  if( puCorrMethod == 2 ){
464  std::vector<double> correctedOutput;
465 
466  CaloSamples cs;
467  coder.adc2fC(digi,cs);
468  std::vector<int> capidvec;
469  for(int ip=0; ip<cs.size(); ip++){
470  const int capid = digi[ip].capid();
471  capidvec.push_back(capid);
472  }
473  psFitOOTpuCorr->apply(cs, capidvec, calibs, correctedOutput);
474  if( correctedOutput.back() == 0 && correctedOutput.size() >1 ){
475  time = correctedOutput[1]; ampl = correctedOutput[0];
476  useTriple = correctedOutput[4];
477  }
478  }
479 
480  // S. Brandt - Feb 19th : Adding Section for HLT
481  // Run "Method 3" all the time.
482  {
483  std::vector<double> hltCorrOutput;
484 
485  CaloSamples cs;
486  coder.adc2fC(digi,cs);
487  std::vector<int> capidvec;
488  for(int ip=0; ip<cs.size(); ip++){
489  const int capid = digi[ip].capid();
490  capidvec.push_back(capid);
491  }
492  hltOOTpuCorr->apply(cs, capidvec, calibs, digi, hltCorrOutput);
493  if (hltCorrOutput.size() > 1) {
494  m3_ampl = hltCorrOutput[0];
495  if (puCorrMethod == 3) {
496  time = hltCorrOutput[1]; ampl = hltCorrOutput[0];
497  }
498  }
499  }
500 
501  // Temporary hack to apply energy-dependent corrections to some HB- cells
502  if (runnum > 0) {
503  const HcalDetId& cell = digi.id();
504  if (cell.subdet() == HcalBarrel) {
505  const int ieta = cell.ieta();
506  const int iphi = cell.iphi();
507  uncorr_ampl *= eCorr(ieta, iphi, uncorr_ampl, runnum);
508  ampl *= eCorr(ieta, iphi, ampl, runnum);
509  m3_ampl *= eCorr(ieta, iphi, m3_ampl, runnum);
510  }
511  }
512 
513  // Correction for a leak to pre-sample
514  if(useLeak && !leakCorrApplied) {
515  uncorr_ampl *= leakCorr(uncorr_ampl);
516  if (puCorrMethod < 2)
517  ampl *= leakCorr(ampl);
518  }
519 
520  HBHERecHit rh(digi.id(),ampl,time);
521  if(useTriple)
522  {
524  }
525  setRawEnergy(rh, static_cast<float>(uncorr_ampl));
526  setAuxEnergy(rh, static_cast<float>(m3_ampl));
527  return rh;
528  }
529 }
530 
531 
532 HBHERecHit HcalSimpleRecAlgo::reconstruct(const HBHEDataFrame& digi, int first, int toadd, const HcalCoder& coder, const HcalCalibrations& calibs) const {
533  return HcalSimpleRecAlgoImpl::recoHBHE<HBHEDataFrame,HBHERecHit>(digi,coder,calibs,
535  pulseCorr_->get(digi.id(), toadd, phaseNS_),
538  hbhePileupCorr_.get(),
540 }
541 
542 
543 HORecHit HcalSimpleRecAlgo::reconstruct(const HODataFrame& digi, int first, int toadd, const HcalCoder& coder, const HcalCalibrations& calibs) const {
544  return HcalSimpleRecAlgoImpl::reco<HODataFrame,HORecHit>(digi,coder,calibs,
546  pulseCorr_->get(digi.id(), toadd, phaseNS_),
548  runnum_, false, hoPileupCorr_.get(),
550 }
551 
552 
553 HcalCalibRecHit HcalSimpleRecAlgo::reconstruct(const HcalCalibDataFrame& digi, int first, int toadd, const HcalCoder& coder, const HcalCalibrations& calibs) const {
554  return HcalSimpleRecAlgoImpl::reco<HcalCalibDataFrame,HcalCalibRecHit>(digi,coder,calibs,
556  pulseCorr_->get(digi.id(), toadd, phaseNS_),
558  runnum_, false, 0,
560 }
561 
562 
563 HBHERecHit HcalSimpleRecAlgo::reconstructHBHEUpgrade(const HcalUpgradeDataFrame& digi, int first, int toadd, const HcalCoder& coder, const HcalCalibrations& calibs) const {
564  HBHERecHit result = HcalSimpleRecAlgoImpl::reco<HcalUpgradeDataFrame,HBHERecHit>(digi, coder, calibs,
566  pulseCorr_->get(digi.id(), toadd, phaseNS_),
567  HcalTimeSlew::Medium, 0, false,
568  hbhePileupCorr_.get(),
571  tdcReco.reconstruct(digi, result);
572  return result;
573 }
574 
575 
577  const int first,
578  const int toadd,
579  const HcalCoder& coder,
580  const HcalCalibrations& calibs) const
581 {
582  const HcalPulseContainmentCorrection* corr = pulseCorr_->get(digi.id(), toadd, phaseNS_);
583 
584  double amp_fC, ampl, uncorr_ampl, maxA;
585  int nRead, maxI;
586  bool leakCorrApplied;
587  float t0, t2;
588 
589  HcalSimpleRecAlgoImpl::removePileup(digi, coder, calibs, first, toadd,
590  correctForPulse_, corr, hfPileupCorr_.get(),
592  &maxA, &ampl, &uncorr_ampl, &amp_fC, &nRead,
593  &maxI, &leakCorrApplied, &t0, &t2);
594 
595  float time=-9999.f;
596  if (maxI > 0 && maxI < (nRead - 1))
597  time = HcalSimpleRecAlgoImpl::recoHFTime(digi,maxI,amp_fC,correctForTimeslew_,maxA,t0,t2) -
598  calibs.timecorr();
599 
600  HFRecHit rh(digi.id(),ampl,time);
601  setRawEnergy(rh, static_cast<float>(uncorr_ampl));
602  return rh;
603 }
604 
605 
606 // NB: Upgrade HFRecHit method content is just the same as regular HFRecHit
607 // with one exclusion: double time (second is dummy) in constructor
609  const int first,
610  const int toadd,
611  const HcalCoder& coder,
612  const HcalCalibrations& calibs) const
613 {
614  const HcalPulseContainmentCorrection* corr = pulseCorr_->get(digi.id(), toadd, phaseNS_);
615 
616  double amp_fC, ampl, uncorr_ampl, maxA;
617  int nRead, maxI;
618  bool leakCorrApplied;
619  float t0, t2;
620 
621  HcalSimpleRecAlgoImpl::removePileup(digi, coder, calibs, first, toadd,
622  correctForPulse_, corr, hfPileupCorr_.get(),
624  &maxA, &ampl, &uncorr_ampl, &amp_fC, &nRead,
625  &maxI, &leakCorrApplied, &t0, &t2);
626 
627  float time=-9999.f;
628  if (maxI > 0 && maxI < (nRead - 1))
629  time = HcalSimpleRecAlgoImpl::recoHFTime(digi,maxI,amp_fC,correctForTimeslew_,maxA,t0,t2) -
630  calibs.timecorr();
631 
632  HFRecHit rh(digi.id(),ampl,time); // new RecHit gets second time = 0.
633  setRawEnergy(rh, static_cast<float>(uncorr_ampl));
634  return rh;
635 }
636 
637 
639 float eCorr(int ieta, int iphi, double energy, int runnum) {
640 // return energy correction factor for HBM channels
641 // iphi=6 ieta=(-1,-15) and iphi=32 ieta=(-1,-7)
642 // I.Vodopianov 28 Feb. 2011
643  static const float low32[7] = {0.741,0.721,0.730,0.698,0.708,0.751,0.861};
644  static const float high32[7] = {0.973,0.925,0.900,0.897,0.950,0.935,1};
645  static const float low6[15] = {0.635,0.623,0.670,0.633,0.644,0.648,0.600,
646  0.570,0.595,0.554,0.505,0.513,0.515,0.561,0.579};
647  static const float high6[15] = {0.875,0.937,0.942,0.900,0.922,0.925,0.901,
648  0.850,0.852,0.818,0.731,0.717,0.782,0.853,0.778};
649 
650 
651  double slope, mid, en;
652  double corr = 1.0;
653 
654  if (!(iphi==6 && ieta<0 && ieta>-16) && !(iphi==32 && ieta<0 && ieta>-8))
655  return corr;
656 
657  int jeta = -ieta-1;
658  double xeta = (double) ieta;
659  if (energy > 0.) en=energy;
660  else en = 0.;
661 
662  if (iphi == 32) {
663  slope = 0.2272;
664  mid = 17.14 + 0.7147*xeta;
665  if (en > 100.) corr = high32[jeta];
666  else corr = low32[jeta]+(high32[jeta]-low32[jeta])/(1.0+exp(-(en-mid)*slope));
667  }
668  else if (iphi == 6 && runnum < 216091 ) {
669  slope = 0.1956;
670  mid = 15.96 + 0.3075*xeta;
671  if (en > 100.0) corr = high6[jeta];
672  else corr = low6[jeta]+(high6[jeta]-low6[jeta])/(1.0+exp(-(en-mid)*slope));
673  }
674 
675  // std::cout << "HBHE cell: ieta, iphi = " << ieta << " " << iphi
676  // << " -> energy = " << en << " corr = " << corr << std::endl;
677 
678  return corr;
679 }
680 
681 
682 // Actual leakage (to pre-sample) correction
683 float leakCorr(double energy) {
684  double corr = 1.0;
685  return corr;
686 }
687 
688 
689 // timeshift implementation
690 
691 static const float wpksamp0_hbheho = 0.5;
692 static const int num_bins_hbheho = 61;
693 
694 static const float actual_ns_hbheho[num_bins_hbheho] = {
695 -5.44000, // 0.500, 0.000-0.017
696 -4.84250, // 0.517, 0.017-0.033
697 -4.26500, // 0.533, 0.033-0.050
698 -3.71000, // 0.550, 0.050-0.067
699 -3.18000, // 0.567, 0.067-0.083
700 -2.66250, // 0.583, 0.083-0.100
701 -2.17250, // 0.600, 0.100-0.117
702 -1.69000, // 0.617, 0.117-0.133
703 -1.23000, // 0.633, 0.133-0.150
704 -0.78000, // 0.650, 0.150-0.167
705 -0.34250, // 0.667, 0.167-0.183
706  0.08250, // 0.683, 0.183-0.200
707  0.50250, // 0.700, 0.200-0.217
708  0.90500, // 0.717, 0.217-0.233
709  1.30500, // 0.733, 0.233-0.250
710  1.69500, // 0.750, 0.250-0.267
711  2.07750, // 0.767, 0.267-0.283
712  2.45750, // 0.783, 0.283-0.300
713  2.82500, // 0.800, 0.300-0.317
714  3.19250, // 0.817, 0.317-0.333
715  3.55750, // 0.833, 0.333-0.350
716  3.91750, // 0.850, 0.350-0.367
717  4.27500, // 0.867, 0.367-0.383
718  4.63000, // 0.883, 0.383-0.400
719  4.98500, // 0.900, 0.400-0.417
720  5.33750, // 0.917, 0.417-0.433
721  5.69500, // 0.933, 0.433-0.450
722  6.05000, // 0.950, 0.450-0.467
723  6.40500, // 0.967, 0.467-0.483
724  6.77000, // 0.983, 0.483-0.500
725  7.13500, // 1.000, 0.500-0.517
726  7.50000, // 1.017, 0.517-0.533
727  7.88250, // 1.033, 0.533-0.550
728  8.26500, // 1.050, 0.550-0.567
729  8.66000, // 1.067, 0.567-0.583
730  9.07000, // 1.083, 0.583-0.600
731  9.48250, // 1.100, 0.600-0.617
732  9.92750, // 1.117, 0.617-0.633
733 10.37750, // 1.133, 0.633-0.650
734 10.87500, // 1.150, 0.650-0.667
735 11.38000, // 1.167, 0.667-0.683
736 11.95250, // 1.183, 0.683-0.700
737 12.55000, // 1.200, 0.700-0.717
738 13.22750, // 1.217, 0.717-0.733
739 13.98500, // 1.233, 0.733-0.750
740 14.81500, // 1.250, 0.750-0.767
741 15.71500, // 1.267, 0.767-0.783
742 16.63750, // 1.283, 0.783-0.800
743 17.53750, // 1.300, 0.800-0.817
744 18.38500, // 1.317, 0.817-0.833
745 19.16500, // 1.333, 0.833-0.850
746 19.89750, // 1.350, 0.850-0.867
747 20.59250, // 1.367, 0.867-0.883
748 21.24250, // 1.383, 0.883-0.900
749 21.85250, // 1.400, 0.900-0.917
750 22.44500, // 1.417, 0.917-0.933
751 22.99500, // 1.433, 0.933-0.950
752 23.53250, // 1.450, 0.950-0.967
753 24.03750, // 1.467, 0.967-0.983
754 24.53250, // 1.483, 0.983-1.000
755 25.00000 // 1.500, 1.000-1.017 - keep for interpolation
756 };
757 
758 float timeshift_ns_hbheho(float wpksamp) {
759  float flx = (num_bins_hbheho-1)*(wpksamp - wpksamp0_hbheho);
760  int index = (int)flx;
761  float yval;
762 
763  if (index < 0) return actual_ns_hbheho[0];
764  else if (index >= num_bins_hbheho-1) return actual_ns_hbheho[num_bins_hbheho-1];
765 
766  // else interpolate:
767  float y1 = actual_ns_hbheho[index];
768  float y2 = actual_ns_hbheho[index+1];
769 
770  yval = y1 + (y2-y1)*(flx-(float)index);
771 
772  return yval;
773 }
774 
775 static const int num_bins_hf = 101;
776 static const float wpksamp0_hf = 0.5;
777 
778 static const float actual_ns_hf[num_bins_hf] = {
779  0.00250, // 0.000-0.010
780  0.04500, // 0.010-0.020
781  0.08750, // 0.020-0.030
782  0.13000, // 0.030-0.040
783  0.17250, // 0.040-0.050
784  0.21500, // 0.050-0.060
785  0.26000, // 0.060-0.070
786  0.30250, // 0.070-0.080
787  0.34500, // 0.080-0.090
788  0.38750, // 0.090-0.100
789  0.42750, // 0.100-0.110
790  0.46000, // 0.110-0.120
791  0.49250, // 0.120-0.130
792  0.52500, // 0.130-0.140
793  0.55750, // 0.140-0.150
794  0.59000, // 0.150-0.160
795  0.62250, // 0.160-0.170
796  0.65500, // 0.170-0.180
797  0.68750, // 0.180-0.190
798  0.72000, // 0.190-0.200
799  0.75250, // 0.200-0.210
800  0.78500, // 0.210-0.220
801  0.81750, // 0.220-0.230
802  0.85000, // 0.230-0.240
803  0.88250, // 0.240-0.250
804  0.91500, // 0.250-0.260
805  0.95500, // 0.260-0.270
806  0.99250, // 0.270-0.280
807  1.03250, // 0.280-0.290
808  1.07000, // 0.290-0.300
809  1.10750, // 0.300-0.310
810  1.14750, // 0.310-0.320
811  1.18500, // 0.320-0.330
812  1.22500, // 0.330-0.340
813  1.26250, // 0.340-0.350
814  1.30000, // 0.350-0.360
815  1.34000, // 0.360-0.370
816  1.37750, // 0.370-0.380
817  1.41750, // 0.380-0.390
818  1.48750, // 0.390-0.400
819  1.55750, // 0.400-0.410
820  1.62750, // 0.410-0.420
821  1.69750, // 0.420-0.430
822  1.76750, // 0.430-0.440
823  1.83750, // 0.440-0.450
824  1.90750, // 0.450-0.460
825  2.06750, // 0.460-0.470
826  2.23250, // 0.470-0.480
827  2.40000, // 0.480-0.490
828  2.82250, // 0.490-0.500
829  3.81000, // 0.500-0.510
830  6.90500, // 0.510-0.520
831  8.99250, // 0.520-0.530
832 10.50000, // 0.530-0.540
833 11.68250, // 0.540-0.550
834 12.66250, // 0.550-0.560
835 13.50250, // 0.560-0.570
836 14.23750, // 0.570-0.580
837 14.89750, // 0.580-0.590
838 15.49000, // 0.590-0.600
839 16.03250, // 0.600-0.610
840 16.53250, // 0.610-0.620
841 17.00000, // 0.620-0.630
842 17.44000, // 0.630-0.640
843 17.85250, // 0.640-0.650
844 18.24000, // 0.650-0.660
845 18.61000, // 0.660-0.670
846 18.96750, // 0.670-0.680
847 19.30500, // 0.680-0.690
848 19.63000, // 0.690-0.700
849 19.94500, // 0.700-0.710
850 20.24500, // 0.710-0.720
851 20.54000, // 0.720-0.730
852 20.82250, // 0.730-0.740
853 21.09750, // 0.740-0.750
854 21.37000, // 0.750-0.760
855 21.62750, // 0.760-0.770
856 21.88500, // 0.770-0.780
857 22.13000, // 0.780-0.790
858 22.37250, // 0.790-0.800
859 22.60250, // 0.800-0.810
860 22.83000, // 0.810-0.820
861 23.04250, // 0.820-0.830
862 23.24500, // 0.830-0.840
863 23.44250, // 0.840-0.850
864 23.61000, // 0.850-0.860
865 23.77750, // 0.860-0.870
866 23.93500, // 0.870-0.880
867 24.05500, // 0.880-0.890
868 24.17250, // 0.890-0.900
869 24.29000, // 0.900-0.910
870 24.40750, // 0.910-0.920
871 24.48250, // 0.920-0.930
872 24.55500, // 0.930-0.940
873 24.62500, // 0.940-0.950
874 24.69750, // 0.950-0.960
875 24.77000, // 0.960-0.970
876 24.84000, // 0.970-0.980
877 24.91250, // 0.980-0.990
878 24.95500, // 0.990-1.000
879 24.99750, // 1.000-1.010 - keep for interpolation
880 };
881 
882 float timeshift_ns_hf(float wpksamp) {
883  float flx = (num_bins_hf-1)*(wpksamp-wpksamp0_hf);
884  int index = (int)flx;
885  float yval;
886 
887  if (index < 0) return actual_ns_hf[0];
888  else if (index >= num_bins_hf-1) return actual_ns_hf[num_bins_hf-1];
889 
890  // else interpolate:
891  float y1 = actual_ns_hf[index];
892  float y2 = actual_ns_hf[index+1];
893 
894  // float delta_x = 1/(float)num_bins_hf;
895  // yval = y1 + (y2-y1)*(flx-(float)index)/delta_x;
896 
897  yval = y1 + (y2-y1)*(flx-(float)index);
898  return yval;
899 }
#define LogDebug(id)
void setMeth3Params(int iPedSubMethod, float iPedSubThreshold, int iTimeSlewParsType, std::vector< double > iTimeSlewPars, double irespCorrM3)
virtual bool inputIsEnergy() const =0
const Shape & getShape(int shapeType) const
int i
Definition: DBlmapReader.cc:9
static const float actual_ns_hf[num_bins_hf]
static const TGPicture * info(bool iBackgroundIsBlack)
HBHERecHit reconstruct(const HBHEDataFrame &digi, int first, int toadd, const HcalCoder &coder, const HcalCalibrations &calibs) const
double respcorrgain(int fCapId) const
get response corrected gain for capid=0..3
void apply(const CaloSamples &cs, const std::vector< int > &capidvec, const HcalCalibrations &calibs, const Digi &digi, std::vector< double > &Output) const
auto_ptr< ClusterSequence > cs
double MaximumFractionalError
static float timeshift_ns_hf(float wpksamp)
Same as above, but for the HF PMTs.
HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:49
static const int MAXSAMPLES
Definition: CaloSamples.h:73
double getCorrection(double fc_ampl) const
void beginRun(edm::EventSetup const &es)
static const double slope[3]
assert(m_qm.get())
void setFlagField(uint32_t value, int base, int width=1)
Definition: CaloRecHit.cc:20
void setHOPileupCorrection(boost::shared_ptr< AbsOOTPileupCorrection > corr)
double pedestal(int fCapId) const
get pedestal for capid=0..3
HcalPulseShapes theHcalPulseShapes_
std::tuple< unsigned int, int, int, DigiType, int, int, int, float > Digi
Definition: GenericDigi.h:30
const BunchXParameter * bunchCrossingInfo_
boost::shared_ptr< AbsOOTPileupCorrection > hoPileupCorr_
#define constexpr
static const float wpksamp0_hbheho
void setBXInfo(const BunchXParameter *info, unsigned lenInfo)
boost::shared_ptr< AbsOOTPileupCorrection > hbhePileupCorr_
static const float wpksamp0_hf
int HPDShapev3MCNum
void setRawEnergy(HcalRecHit &h, float e)
Definition: rawEnergy.h:215
static const float actual_ns_hbheho[num_bins_hbheho]
auto const T2 &decltype(t1.eta()) t2
Definition: deltaR.h:16
tuple puCorrMethod
unsigned lenBunchCrossingInfo_
HFRecHit reconstructHFUpgrade(const HcalUpgradeDataFrame &digi, int first, int toadd, const HcalCoder &coder, const HcalCalibrations &calibs) const
tuple result
Definition: query.py:137
int ieta() const
get the cell ieta
Definition: HcalDetId.h:56
virtual void apply(const HcalDetId &id, const double *inputCharge, unsigned lenInputCharge, const BunchXParameter *bcParams, unsigned lenBcParams, unsigned firstTimeSlice, unsigned nTimeSlices, double *correctedCharge, unsigned lenCorrectedCharge, bool *pulseShapeCorrApplied, bool *leakCorrApplied, bool *readjustTiming) const =0
void initPulseCorr(int toadd)
double f[11][100]
HcalSimpleRecAlgo(bool correctForTimeslew, bool correctForContainment, float fixedPhaseNs)
static float leakCorr(double energy)
Leak correction.
T min(T a, T b)
Definition: MathUtil.h:58
std::auto_ptr< HcalPulseContainmentManager > pulseCorr_
std::auto_ptr< PulseShapeFitOOTPileupCorrection > psFitOOTpuCorr_
void removePileup(const Digi &digi, const HcalCoder &coder, const HcalCalibrations &calibs, const int ifirst, const int n, const bool pulseCorrect, const HcalPulseContainmentCorrection *corr, const AbsOOTPileupCorrection *pileupCorrection, const BunchXParameter *bxInfo, const unsigned lenInfo, double *p_maxA, double *p_ampl, double *p_uncorr_ampl, double *p_fc_ampl, int *p_nRead, int *p_maxI, bool *leakCorrApplied, float *p_t0, float *p_t2)
JetCorrectorParameters corr
Definition: classes.h:5
static float timeshift_ns_hbheho(float wpksamp)
int iphi() const
get the cell iphi
Definition: HcalDetId.cc:101
int HPDShapev3DataNum
void setRecoParams(bool correctForTimeslew, bool correctForPulse, bool setLeakCorrection, int pileupCleaningID, float phaseNS)
int size() const
get the size
Definition: CaloSamples.h:24
std::auto_ptr< PedestalSub > pedSubFxn_
double timecorr() const
get time correction factor
void setHBHEPileupCorrection(boost::shared_ptr< AbsOOTPileupCorrection > corr)
void setForData(int runnum)
HBHERecHit reconstructHBHEUpgrade(const HcalUpgradeDataFrame &digi, int first, int toadd, const HcalCoder &coder, const HcalCalibrations &calibs) const
static const int num_bins_hf
boost::shared_ptr< AbsOOTPileupCorrection > hfPileupCorr_
const HcalDetId & id() const
static const int num_bins_hbheho
virtual void adc2fC(const HBHEDataFrame &df, CaloSamples &lf) const =0
void apply(const CaloSamples &cs, const std::vector< int > &capidvec, const HcalCalibrations &calibs, std::vector< double > &correctedOutput) const
RecHit recoHBHE(const Digi &digi, const HcalCoder &coder, const HcalCalibrations &calibs, const int ifirst, const int n, const bool slewCorrect, const bool pulseCorrect, const HcalPulseContainmentCorrection *corr, const HcalTimeSlew::BiasSetting slewFlavor, const int runnum, const bool useLeak, const AbsOOTPileupCorrection *pileupCorrection, const BunchXParameter *bxInfo, const unsigned lenInfo, const int puCorrMethod, const PulseShapeFitOOTPileupCorrection *psFitOOTpuCorr, HcalDeterministicFit *hltOOTpuCorr, PedestalSub *hltPedSub)
float recoHFTime(const Digi &digi, const int maxI, const double amp_fC, const bool slewCorrect, double maxA, float t0, float t2)
tuple runnum
Definition: summaryLumi.py:210
std::auto_ptr< HcalDeterministicFit > hltOOTpuCorr_
volatile std::atomic< bool > shutdown_flag false
const HcalDetId & id() const
Definition: HFDataFrame.h:22
RecHit reco(const Digi &digi, const HcalCoder &coder, const HcalCalibrations &calibs, const int ifirst, const int n, const bool slewCorrect, const bool pulseCorrect, const HcalPulseContainmentCorrection *corr, const HcalTimeSlew::BiasSetting slewFlavor, const int runnum, const bool useLeak, const AbsOOTPileupCorrection *pileupCorrection, const BunchXParameter *bxInfo, const unsigned lenInfo, const int puCorrMethod, const PulseShapeFitOOTPileupCorrection *psFitOOTpuCorr, HcalDeterministicFit *hltOOTpuCorr, PedestalSub *hltPedSub)
static float eCorr(int ieta, int iphi, double ampl, int runnum)
Ugly hack to apply energy corrections to some HB- cells.
void setHFPileupCorrection(boost::shared_ptr< AbsOOTPileupCorrection > corr)
void setAuxEnergy(HcalRecHit &h, float e)
Definition: rawEnergy.h:232
static double delay(double fC, BiasSetting bias=Medium)
Returns the amount (ns) by which a pulse of the given number of fC will be delayed by the timeslew ef...
Definition: HcalTimeSlew.cc:12
void setpuCorrParams(bool iPedestalConstraint, bool iTimeConstraint, bool iAddPulseJitter, bool iUnConstrainedFit, bool iApplyTimeSlew, double iTS4Min, double iTS4Max, double iPulseJitter, double iTimeMean, double iTimeSig, double iPedMean, double iPedSig, double iNoise, double iTMin, double iTMax, double its3Chi2, double its4Chi2, double its345Chi2, double iChargeThreshold, int iFitTimes)