CMS 3D CMS Logo

Functions
HcalSimpleRecAlgoImpl Namespace Reference

Functions

template<class Digi , class RecHit >
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)
 
template<class Digi , class RecHit >
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)
 
template<class Digi >
float recoHFTime (const Digi &digi, const int maxI, const double amp_fC, const bool slewCorrect, double maxA, float t0, float t2)
 
template<class Digi >
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)
 

Function Documentation

template<class Digi , class RecHit >
RecHit HcalSimpleRecAlgoImpl::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 
)
inline

Definition at line 312 of file HcalSimpleRecAlgo.cc.

References HcalCoder::adc2fC(), HcalDeterministicFit::apply(), PulseShapeFitOOTPileupCorrection::apply(), HiEvtPlane_cfi::chi2, fwrapper::cs, HcalTimeSlew::delay(), hbminus_special_ecorr(), HcalBarrel, HcalDetId::ieta(), HcalDetId::iphi(), leakCorr(), hpstanc_transforms::max, findQualityFiles::maxI, removePileup(), setAuxEnergy(), setRawEnergy(), CaloSamples::size(), HcalDetId::subdet(), cscNeutronWriter_cfi::t0, reco::t2, ntuplemaker::time, HcalCalibrations::timecorr(), and timeshift_ns_hbheho().

320  {
321  double fc_ampl =0, ampl =0, uncorr_ampl =0, m3_ampl =0, maxA = -1.e300;
322  int nRead = 0, maxI = -1;
323  bool leakCorrApplied = false;
324  float t0 =0, t2 =0;
325  float time = -9999;
326  float m3_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 
358  bool useTriple=false;
359  float chi2=-1;
360 
361  CaloSamples cs;
362  coder.adc2fC(digi,cs);
363  std::vector<int> capidvec;
364  for(int ip=0; ip<cs.size(); ip++){
365  const int capid = digi[ip].capid();
366  capidvec.push_back(capid);
367  }
368  psFitOOTpuCorr->apply(cs, capidvec, calibs, ampl, time, useTriple,chi2);
369  }
370 
371  // S. Brandt - Feb 19th : Adding Section for HLT
372  // Run "Method 3" all the time.
373  {
374 
375  CaloSamples cs;
376  coder.adc2fC(digi,cs);
377  std::vector<int> capidvec;
378  for(int ip=0; ip<cs.size(); ip++){
379  const int capid = digi[ip].capid();
380  capidvec.push_back(capid);
381  }
382  hltOOTpuCorr->apply(cs, capidvec, calibs, digi, m3_ampl,m3_time);
383  if (puCorrMethod == 3) {ampl = m3_ampl; time=m3_time;}
384  }
385 
386  // Temporary hack to apply energy-dependent corrections to some HB- cells
387  if (runnum > 0) {
388  const HcalDetId& cell = digi.id();
389  if (cell.subdet() == HcalBarrel) {
390  const int ieta = cell.ieta();
391  const int iphi = cell.iphi();
392  uncorr_ampl *= hbminus_special_ecorr(ieta, iphi, uncorr_ampl, runnum);
393  ampl *= hbminus_special_ecorr(ieta, iphi, ampl, runnum);
394  m3_ampl *= hbminus_special_ecorr(ieta, iphi, m3_ampl, runnum);
395  }
396  }
397 
398  // Correction for a leak to pre-sample
399  if(useLeak && !leakCorrApplied) {
400  uncorr_ampl *= leakCorr(uncorr_ampl);
401  if (puCorrMethod < 2)
402  ampl *= leakCorr(ampl);
403  }
404 
405  RecHit rh(digi.id(),ampl,time);
406  setRawEnergy(rh, static_cast<float>(uncorr_ampl));
407  setAuxEnergy(rh, static_cast<float>(m3_ampl));
408  return rh;
409  }
void apply(const CaloSamples &cs, const std::vector< int > &capidvec, const HcalCalibrations &calibs, double &reconstructedEnergy, float &reconstructedTime, bool &useTriple, float &chi2) const
auto_ptr< ClusterSequence > cs
HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:49
void setRawEnergy(HcalRecHit &h, float e)
Definition: rawEnergy.h:215
auto const T2 &decltype(t1.eta()) t2
Definition: deltaR.h:16
int ieta() const
get the cell ieta
Definition: HcalDetId.h:56
static float leakCorr(double energy)
Leak correction.
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)
virtual void adc2fC(const HBHEDataFrame &df, CaloSamples &lf) const =0
float timeshift_ns_hbheho(float wpksamp)
int iphi() const
get the cell iphi
Definition: HcalDetId.cc:103
void apply(const CaloSamples &cs, const std::vector< int > &capidvec, const HcalCalibrations &calibs, const Digi &digi, double &ampl, float &time) const
int size() const
get the size
Definition: CaloSamples.h:24
double timecorr() const
get time correction factor
float hbminus_special_ecorr(int ieta, int iphi, double energy, int runnum)
Ugly hack to apply energy corrections to some HB- cells.
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:16
template<class Digi , class RecHit >
RecHit HcalSimpleRecAlgoImpl::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 
)
inline

Definition at line 412 of file HcalSimpleRecAlgo.cc.

References HcalCoder::adc2fC(), HcalDeterministicFit::apply(), PulseShapeFitOOTPileupCorrection::apply(), HiEvtPlane_cfi::chi2, fwrapper::cs, HcalTimeSlew::delay(), HcalCaloFlagLabels::HBHEPulseFitBit, hbminus_special_ecorr(), HcalBarrel, HcalDetId::ieta(), HcalDetId::iphi(), leakCorr(), hpstanc_transforms::max, findQualityFiles::maxI, removePileup(), setAuxEnergy(), CaloRecHit::setFlagField(), setRawEnergy(), CaloSamples::size(), HcalDetId::subdet(), cscNeutronWriter_cfi::t0, reco::t2, ntuplemaker::time, HcalCalibrations::timecorr(), and timeshift_ns_hbheho().

420  {
421  double fc_ampl =0, ampl =0, uncorr_ampl =0, m3_ampl =0, maxA = -1.e300;
422  int nRead = 0, maxI = -1;
423  bool leakCorrApplied = false;
424  float t0 =0, t2 =0;
425  float time = -9999;
426  bool useTriple = false;
427  float chi2 = -1;
428 
429  // Disable method 1 inside the removePileup function this way!
430  // Some code in removePileup does NOT do pileup correction & to make sure maximum share of code
431  const AbsOOTPileupCorrection * inputAbsOOTpuCorr = ( puCorrMethod == 1 ? pileupCorrection: 0 );
432 
433  removePileup(digi, coder, calibs, ifirst, n,
434  pulseCorrect, corr, inputAbsOOTpuCorr,
435  bxInfo, lenInfo, &maxA, &ampl,
436  &uncorr_ampl, &fc_ampl, &nRead, &maxI,
437  &leakCorrApplied, &t0, &t2);
438 
439  if (maxI > 0 && maxI < (nRead - 1))
440  {
441  // Handle negative excursions by moving "zero":
442  float minA=t0;
443  if (maxA<minA) minA=maxA;
444  if (t2<minA) minA=t2;
445  if (minA<0) { maxA-=minA; t0-=minA; t2-=minA; } // positivizes all samples
446 
447  float wpksamp = (t0 + maxA + t2);
448  if (wpksamp!=0) wpksamp=(maxA + 2.0*t2) / wpksamp;
449  time = (maxI - digi.presamples())*25.0 + timeshift_ns_hbheho(wpksamp);
450 
451  if (slewCorrect) time-=HcalTimeSlew::delay(std::max(1.0,fc_ampl),slewFlavor);
452 
453  time=time-calibs.timecorr(); // time calibration
454  }
455 
456  // Note that uncorr_ampl is always set from outside of method 2!
457  if( puCorrMethod == 2 ){
458 
459  CaloSamples cs;
460  coder.adc2fC(digi,cs);
461  std::vector<int> capidvec;
462  for(int ip=0; ip<cs.size(); ip++){
463  const int capid = digi[ip].capid();
464  capidvec.push_back(capid);
465  }
466  psFitOOTpuCorr->apply(cs, capidvec, calibs, ampl, time, useTriple,chi2);
467  }
468 
469  // S. Brandt - Feb 19th : Adding Section for HLT
470  // Run "Method 3" all the time.
471  {
472 
473  CaloSamples cs;
474  coder.adc2fC(digi,cs);
475  std::vector<int> capidvec;
476  for(int ip=0; ip<cs.size(); ip++){
477  const int capid = digi[ip].capid();
478  capidvec.push_back(capid);
479  }
480 
481  float m3_time=0;
482 
483  hltOOTpuCorr->apply(cs, capidvec, calibs, digi, m3_ampl, m3_time);
484  if (puCorrMethod == 3) { ampl = m3_ampl; time = m3_time; }
485 
486  }
487 
488  // Temporary hack to apply energy-dependent corrections to some HB- cells
489  if (runnum > 0) {
490  const HcalDetId& cell = digi.id();
491  if (cell.subdet() == HcalBarrel) {
492  const int ieta = cell.ieta();
493  const int iphi = cell.iphi();
494  uncorr_ampl *= hbminus_special_ecorr(ieta, iphi, uncorr_ampl, runnum);
495  ampl *= hbminus_special_ecorr(ieta, iphi, ampl, runnum);
496  m3_ampl *= hbminus_special_ecorr(ieta, iphi, m3_ampl, runnum);
497  }
498  }
499 
500  // Correction for a leak to pre-sample
501  if(useLeak && !leakCorrApplied) {
502  uncorr_ampl *= leakCorr(uncorr_ampl);
503  if (puCorrMethod < 2)
504  ampl *= leakCorr(ampl);
505  }
506 
507  HBHERecHit rh(digi.id(),ampl,time);
508  if(useTriple)
509  {
511  }
512  rh.setChiSquared(chi2);
513  setRawEnergy(rh, static_cast<float>(uncorr_ampl));
514  setAuxEnergy(rh, static_cast<float>(m3_ampl));
515  return rh;
516  }
void apply(const CaloSamples &cs, const std::vector< int > &capidvec, const HcalCalibrations &calibs, double &reconstructedEnergy, float &reconstructedTime, bool &useTriple, float &chi2) const
auto_ptr< ClusterSequence > cs
HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:49
void setFlagField(uint32_t value, int base, int width=1)
Definition: CaloRecHit.cc:20
void setRawEnergy(HcalRecHit &h, float e)
Definition: rawEnergy.h:215
auto const T2 &decltype(t1.eta()) t2
Definition: deltaR.h:16
int ieta() const
get the cell ieta
Definition: HcalDetId.h:56
static float leakCorr(double energy)
Leak correction.
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)
virtual void adc2fC(const HBHEDataFrame &df, CaloSamples &lf) const =0
float timeshift_ns_hbheho(float wpksamp)
int iphi() const
get the cell iphi
Definition: HcalDetId.cc:103
void apply(const CaloSamples &cs, const std::vector< int > &capidvec, const HcalCalibrations &calibs, const Digi &digi, double &ampl, float &time) const
int size() const
get the size
Definition: CaloSamples.h:24
double timecorr() const
get time correction factor
float hbminus_special_ecorr(int ieta, int iphi, double energy, int runnum)
Ugly hack to apply energy corrections to some HB- cells.
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:16
template<class Digi >
float HcalSimpleRecAlgoImpl::recoHFTime ( const Digi &  digi,
const int  maxI,
const double  amp_fC,
const bool  slewCorrect,
double  maxA,
float  t0,
float  t2 
)
inline

Definition at line 128 of file HcalSimpleRecAlgo.cc.

References MillePedeFileConverter_cfg::e, JetChargeProducer_cfi::exp, f, objects.autophobj::float, min(), reco::t2, ntuplemaker::time, and timeshift_ns_hf().

Referenced by HcalSimpleRecAlgo::reconstruct(), and HcalSimpleRecAlgo::reconstructQIE10().

130  {
131  // Handle negative excursions by moving "zero":
132  float zerocorr=std::min(t0,t2);
133  if (zerocorr<0.f) {
134  t0 -= zerocorr;
135  t2 -= zerocorr;
136  maxA -= zerocorr;
137  }
138 
139  // pair the peak with the larger of the two neighboring time samples
140  float wpksamp=0.f;
141  if (t0>t2) {
142  wpksamp = t0+maxA;
143  if (wpksamp != 0.f) wpksamp = maxA/wpksamp;
144  } else {
145  wpksamp = maxA+t2;
146  if (wpksamp != 0.f) wpksamp = 1.+(t2/wpksamp);
147  }
148 
149  float time = (maxI - digi.presamples())*25.0 + timeshift_ns_hf(wpksamp);
150 
151  if (slewCorrect && amp_fC > 0.0) {
152  // -5.12327 - put in calibs.timecorr()
153  double tslew=exp(0.337681-5.94689e-4*amp_fC)+exp(2.44628-1.34888e-2*amp_fC);
154  time -= (float)tslew;
155  }
156 
157  return time;
158  }
static float timeshift_ns_hf(float wpksamp)
Timeshift correction for the HF PMTs.
auto const T2 &decltype(t1.eta()) t2
Definition: deltaR.h:16
double f[11][100]
T min(T a, T b)
Definition: MathUtil.h:58
template<class Digi >
void HcalSimpleRecAlgoImpl::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 
)
inline

Definition at line 162 of file HcalSimpleRecAlgo.cc.

References HcalCoder::adc2fC(), AbsOOTPileupCorrection::apply(), fwrapper::cs, HcalPulseContainmentCorrection::getCorrection(), mps_fire::i, AbsOOTPileupCorrection::inputIsEnergy(), LogDebug, hpstanc_transforms::max, findQualityFiles::maxI, CaloSamples::MAXSAMPLES, min(), gen::n, HcalCalibrations::pedestal(), HcalCalibrations::respcorrgain(), and CaloSamples::size().

Referenced by reco(), recoHBHE(), HcalSimpleRecAlgo::reconstruct(), and HcalSimpleRecAlgo::reconstructQIE10().

172  {
173  CaloSamples cs;
174  coder.adc2fC(digi,cs);
175  const int nRead = cs.size();
176  const int iStop = std::min(nRead, n + ifirst);
177 
178  // Signal energy will be calculated both with
179  // and without OOT pileup corrections. Try to
180  // arrange the calculations so that we do not
181  // repeat them.
182  double uncorrectedEnergy[CaloSamples::MAXSAMPLES] {}, buf[CaloSamples::MAXSAMPLES] {};
183  double* correctedEnergy = 0;
184  double fc_ampl = 0.0, corr_fc_ampl = 0.0;
185  bool pulseShapeCorrApplied = false, readjustTiming = false;
186  *leakCorrApplied = false;
187 
188  if (pileupCorrection)
189  {
190  correctedEnergy = &buf[0];
191 
192  double correctionInput[CaloSamples::MAXSAMPLES] {};
193  double gains[CaloSamples::MAXSAMPLES] {};
194 
195  for (int i=0; i<nRead; ++i)
196  {
197  const int capid = digi[i].capid();
198  correctionInput[i] = cs[i] - calibs.pedestal(capid);
199  gains[i] = calibs.respcorrgain(capid);
200  }
201 
202  for (int i=ifirst; i<iStop; ++i)
203  fc_ampl += correctionInput[i];
204 
205  const bool useGain = pileupCorrection->inputIsEnergy();
206  for (int i=0; i<nRead; ++i)
207  {
208  uncorrectedEnergy[i] = correctionInput[i]*gains[i];
209  if (useGain)
210  correctionInput[i] = uncorrectedEnergy[i];
211  }
212 
213  pileupCorrection->apply(digi.id(), correctionInput, nRead,
214  bxInfo, lenInfo, ifirst, n,
215  correctedEnergy, CaloSamples::MAXSAMPLES,
216  &pulseShapeCorrApplied, leakCorrApplied,
217  &readjustTiming);
218  if (useGain)
219  {
220  // Gain factors have been already applied.
221  // Divide by them for accumulating corr_fc_ampl.
222  for (int i=ifirst; i<iStop; ++i)
223  if (gains[i])
224  corr_fc_ampl += correctedEnergy[i]/gains[i];
225  }
226  else
227  {
228  for (int i=ifirst; i<iStop; ++i)
229  corr_fc_ampl += correctedEnergy[i];
230  for (int i=0; i<nRead; ++i)
231  correctedEnergy[i] *= gains[i];
232  }
233  }
234  else
235  {
236  correctedEnergy = &uncorrectedEnergy[0];
237 
238  // In this situation, we do not need to process all time slices
239  const int istart = std::max(ifirst - 1, 0);
240  const int iend = std::min(n + ifirst + 1, nRead);
241  for (int i=istart; i<iend; ++i)
242  {
243  const int capid = digi[i].capid();
244  float ta = cs[i] - calibs.pedestal(capid);
245  if (i >= ifirst && i < iStop)
246  fc_ampl += ta;
247  ta *= calibs.respcorrgain(capid);
248  uncorrectedEnergy[i] = ta;
249  }
250  corr_fc_ampl = fc_ampl;
251  }
252 
253  // Uncorrected and corrected energies
254  double ampl = 0.0, corr_ampl = 0.0;
255  for (int i=ifirst; i<iStop; ++i)
256  {
257  ampl += uncorrectedEnergy[i];
258  corr_ampl += correctedEnergy[i];
259  }
260 
261  // Apply phase-based amplitude correction:
262  if (corr && pulseCorrect)
263  {
264  ampl *= corr->getCorrection(fc_ampl);
265  if (pileupCorrection)
266  {
267  if (!pulseShapeCorrApplied)
268  corr_ampl *= corr->getCorrection(corr_fc_ampl);
269  }
270  else
271  corr_ampl = ampl;
272  }
273 
274  // Which energies we want to use for timing?
275  const double *etime = readjustTiming ? &correctedEnergy[0] : &uncorrectedEnergy[0];
276  int maxI = -1; double maxA = -1.e300;
277  for (int i=ifirst; i<iStop; ++i)
278  if (etime[i] > maxA)
279  {
280  maxA = etime[i];
281  maxI = i;
282  }
283 
284  // Fill out the output
285  *p_maxA = maxA;
286  *p_ampl = corr_ampl;
287  *p_uncorr_ampl = ampl;
288  *p_fc_ampl = readjustTiming ? corr_fc_ampl : fc_ampl;
289  *p_nRead = nRead;
290  *p_maxI = maxI;
291 
292  if (maxI <= 0 || maxI >= (nRead-1))
293  {
294  LogDebug("HCAL Pulse") << "HcalSimpleRecAlgoImpl::removePileup :"
295  << " Invalid max amplitude position, "
296  << " max Amplitude: " << maxI
297  << " first: " << ifirst
298  << " last: " << ifirst + n
299  << std::endl;
300  *p_t0 = 0.f;
301  *p_t2 = 0.f;
302  }
303  else
304  {
305  *p_t0 = etime[maxI - 1];
306  *p_t2 = etime[maxI + 1];
307  }
308  }
#define LogDebug(id)
double respcorrgain(int fCapId) const
get response corrected gain for capid=0..3
auto_ptr< ClusterSequence > cs
static const int MAXSAMPLES
Definition: CaloSamples.h:74
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
double getCorrection(double fc_ampl) const
double pedestal(int fCapId) const
get pedestal for capid=0..3
T min(T a, T b)
Definition: MathUtil.h:58
virtual void adc2fC(const HBHEDataFrame &df, CaloSamples &lf) const =0
int size() const
get the size
Definition: CaloSamples.h:24
virtual bool inputIsEnergy() const =0