CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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)
 
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 
)
inline

Definition at line 286 of file HcalSimpleRecAlgo.cc.

References HcalTimeSlew::delay(), eCorr(), HcalBarrel, HcalDetId::ieta(), HcalDetId::iphi(), leakCorr(), max(), findQualityFiles::maxI, removePileup(), setRawEnergy(), HcalDetId::subdet(), reco::t2, cond::rpcobgas::time, HcalCalibrations::timecorr(), and timeshift_ns_hbheho().

294  {
295  double fc_ampl, ampl, uncorr_ampl, maxA;
296  int nRead, maxI;
297  bool leakCorrApplied;
298  float t0, t2;
299 
300  removePileup(digi, coder, calibs, ifirst, n,
301  pulseCorrect, corr, pileupCorrection,
302  bxInfo, lenInfo, &maxA, &ampl,
303  &uncorr_ampl, &fc_ampl, &nRead, &maxI,
304  &leakCorrApplied, &t0, &t2);
305 
306  float time = -9999;
307  if (maxI > 0 && maxI < (nRead - 1))
308  {
309  // Handle negative excursions by moving "zero":
310  float minA=t0;
311  if (maxA<minA) minA=maxA;
312  if (t2<minA) minA=t2;
313  if (minA<0) { maxA-=minA; t0-=minA; t2-=minA; } // positivizes all samples
314 
315  float wpksamp = (t0 + maxA + t2);
316  if (wpksamp!=0) wpksamp=(maxA + 2.0*t2) / wpksamp;
317  time = (maxI - digi.presamples())*25.0 + timeshift_ns_hbheho(wpksamp);
318 
319  if (slewCorrect) time-=HcalTimeSlew::delay(std::max(1.0,fc_ampl),slewFlavor);
320 
321  time=time-calibs.timecorr(); // time calibration
322  }
323 
324  // Temporary hack to apply energy-dependent corrections to some HB- cells
325  if (runnum > 0) {
326  const HcalDetId& cell = digi.id();
327  if (cell.subdet() == HcalBarrel) {
328  const int ieta = cell.ieta();
329  const int iphi = cell.iphi();
330  ampl *= eCorr(ieta, iphi, ampl, runnum);
331  uncorr_ampl *= eCorr(ieta, iphi, uncorr_ampl, runnum);
332  }
333  }
334 
335  // Correction for a leak to pre-sample
336  if(useLeak && !leakCorrApplied) {
337  ampl *= leakCorr(ampl);
338  uncorr_ampl *= leakCorr(uncorr_ampl);
339  }
340 
341  RecHit rh(digi.id(),ampl,time);
342  setRawEnergy(rh, static_cast<float>(uncorr_ampl));
343  return rh;
344  }
HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:30
void setRawEnergy(HcalRecHit &h, float e)
Definition: rawEnergy.h:139
const T & max(const T &a, const T &b)
auto const T2 &decltype(t1.eta()) t2
Definition: deltaR.h:18
int ieta() const
get the cell ieta
Definition: HcalDetId.h:36
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)
static float timeshift_ns_hbheho(float wpksamp)
int iphi() const
get the cell iphi
Definition: HcalDetId.h:38
double timecorr() const
get time correction factor
tuple runnum
Definition: summaryLumi.py:210
static float eCorr(int ieta, int iphi, double ampl, int runnum)
Ugly hack to apply energy corrections to some HB- cells.
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:8
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 102 of file HcalSimpleRecAlgo.cc.

References alignCSCRings::e, create_public_lumi_plots::exp, f, bookConverter::min, reco::t2, cond::rpcobgas::time, and timeshift_ns_hf().

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

104  {
105  // Handle negative excursions by moving "zero":
106  float zerocorr=std::min(t0,t2);
107  if (zerocorr<0.f) {
108  t0 -= zerocorr;
109  t2 -= zerocorr;
110  maxA -= zerocorr;
111  }
112 
113  // pair the peak with the larger of the two neighboring time samples
114  float wpksamp=0.f;
115  if (t0>t2) {
116  wpksamp = t0+maxA;
117  if (wpksamp != 0.f) wpksamp = maxA/wpksamp;
118  } else {
119  wpksamp = maxA+t2;
120  if (wpksamp != 0.f) wpksamp = 1.+(t2/wpksamp);
121  }
122 
123  float time = (maxI - digi.presamples())*25.0 + timeshift_ns_hf(wpksamp);
124 
125  if (slewCorrect && amp_fC > 0.0) {
126  // -5.12327 - put in calibs.timecorr()
127  double tslew=exp(0.337681-5.94689e-4*amp_fC)+exp(2.44628-1.34888e-2*amp_fC);
128  time -= (float)tslew;
129  }
130 
131  return time;
132  }
static float timeshift_ns_hf(float wpksamp)
Same as above, but for the HF PMTs.
auto const T2 &decltype(t1.eta()) t2
Definition: deltaR.h:18
double f[11][100]
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 136 of file HcalSimpleRecAlgo.cc.

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

Referenced by reco(), HcalSimpleRecAlgo::reconstruct(), and HcalSimpleRecAlgo::reconstructHFUpgrade().

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