CMS 3D CMS Logo

HBHEPhase1Reconstructor.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: RecoLocalCalo/HcalRecProducers
4 // Class: HBHEPhase1Reconstructor
5 //
13 //
14 // Original Author: Igor Volobouev
15 // Created: Tue, 21 Jun 2016 00:56:40 GMT
16 //
17 //
18 
19 
20 // system include files
21 #include <cmath>
22 #include <utility>
23 #include <algorithm>
24 
25 // user include files
32 
35 
40 
44 
46 
48 
53 
56 
57 // Parser for Phase 1 HB/HE reco algorithms
59 
60 // Fetcher for reco algorithm data
62 
63 // Some helper functions
64 namespace {
65  // Class for making SiPM/QIE11 look like HPD/QIE8. HPD/QIE8
66  // needs only pedestal and gain to convert charge into energy.
67  // Due to nonlinearities, response of SiPM/QIE11 is substantially
68  // more complicated. It is possible to calculate all necessary
69  // quantities from the charge and the info stored in the DB every
70  // time the raw charge is needed. However, it does not make sense
71  // to retrieve DB contents stored by channel for every time slice.
72  // Because of this, we look things up only once, in the constructor.
73  template<class DFrame>
74  class RawChargeFromSample
75  {
76  public:
77  inline RawChargeFromSample(const int sipmQTSShift,
78  const int sipmQNTStoSum,
79  const bool saveEffectivePedestal,
80  const HcalDbService& cond,
81  const HcalDetId id,
82  const CaloSamples& cs,
83  const int soi,
84  const DFrame& frame,
85  const int maxTS) {}
86 
87  inline double getRawCharge(const double decodedCharge,
88  const double pedestal) const
89  {return decodedCharge;}
90  };
91 
92  template<>
93  class RawChargeFromSample<QIE11DataFrame>
94  {
95  public:
96  inline RawChargeFromSample(const int sipmQTSShift,
97  const int sipmQNTStoSum,
98  const bool saveEffectivePedestal,
99  const HcalDbService& cond,
100  const HcalDetId id,
101  const CaloSamples& cs,
102  const int soi,
103  const QIE11DataFrame& frame,
104  const int maxTS)
105  : siPMParameter_(*cond.getHcalSiPMParameter(id)),
106  fcByPE_(siPMParameter_.getFCByPE()),
107  corr_(cond.getHcalSiPMCharacteristics()->getNonLinearities(siPMParameter_.getType()))
108  {
109  if (fcByPE_ <= 0.0)
110  throw cms::Exception("HBHEPhase1BadDB")
111  << "Invalid fC/PE conversion factor for SiPM " << id
112  << std::endl;
113 
114  const HcalCalibrations& calib = cond.getHcalCalibrations(id);
115  const double darkCurrent = siPMParameter_.getDarkCurrent();
116  const double lambda = cond.getHcalSiPMCharacteristics()->getCrossTalk(siPMParameter_.getType());
117  const int firstTS = std::max(soi + sipmQTSShift, 0);
118  const int lastTS = std::min(firstTS + sipmQNTStoSum, maxTS);
119  double sipmQ = 0.0;
120 
121  for (int ts = firstTS; ts < lastTS; ++ts)
122  {
123  const double pedestal = calib.pedestal(frame[ts].capid()) +
124  (saveEffectivePedestal ? darkCurrent * 25. / (1. - lambda) : 0.);
125  sipmQ += (cs[ts] - pedestal);
126  }
127 
128  const double effectivePixelsFired = sipmQ/fcByPE_;
129  factor_ = corr_.getRecoCorrectionFactor(effectivePixelsFired);
130  }
131 
132  inline double getRawCharge(const double decodedCharge,
133  const double pedestal) const
134  {
135  return (decodedCharge - pedestal)*factor_ + pedestal;
136 
137  // Old version of TS-by-TS corrections looked as follows:
138  // const double sipmQ = decodedCharge - pedestal;
139  // const double nPixelsFired = sipmQ/fcByPE_;
140  // return sipmQ*corr_.getRecoCorrectionFactor(nPixelsFired) + pedestal;
141  }
142 
143  private:
144  const HcalSiPMParameter& siPMParameter_;
145  double fcByPE_;
146  HcalSiPMnonlinearity corr_;
147  double factor_;
148  };
149 
150  float getTDCTimeFromSample(const QIE11DataFrame::Sample& s)
151  {
152  return HcalSpecialTimes::getTDCTime(s.tdc());
153  }
154 
155  float getTDCTimeFromSample(const HcalQIESample&)
156  {
158  }
159 
160  float getDifferentialChargeGain(const HcalQIECoder& coder,
161  const HcalQIEShape& shape,
162  const unsigned adc,
163  const unsigned capid,
164  const bool isQIE11)
165  {
166  // We have 5-bit ADC mantissa in QIE8 and 6-bit in QIE11
167  static const unsigned mantissaMaskQIE8 = 0x1f;
168  static const unsigned mantissaMaskQIE11 = 0x3f;
169 
170  const float q = coder.charge(shape, adc, capid);
171  const unsigned mantissaMask = isQIE11 ? mantissaMaskQIE11 : mantissaMaskQIE8;
172  const unsigned mantissa = adc & mantissaMask;
173 
174  // First, check if we are in the two lowest or two highest ADC
175  // values for this range. Assume that they have the lowest and
176  // the highest gain in the range, respectively.
177  if (mantissa == 0U || mantissa == mantissaMask - 1U)
178  return coder.charge(shape, adc+1U, capid) - q;
179  else if (mantissa == 1U || mantissa == mantissaMask)
180  return q - coder.charge(shape, adc-1U, capid);
181  else
182  {
183  const float qup = coder.charge(shape, adc+1U, capid);
184  const float qdown = coder.charge(shape, adc-1U, capid);
185  const float upGain = qup - q;
186  const float downGain = q - qdown;
187  const float averageGain = (qup - qdown)/2.f;
188  if (std::abs(upGain - downGain) < 0.01f*averageGain)
189  return averageGain;
190  else
191  {
192  // We are in the gain transition region.
193  // Need to determine if we are in the lower
194  // gain ADC count or in the higher one.
195  // This can be done by figuring out if the
196  // "up" gain is more consistent then the
197  // "down" gain.
198  const float q2up = coder.charge(shape, adc+2U, capid);
199  const float q2down = coder.charge(shape, adc-2U, capid);
200  const float upGain2 = q2up - qup;
201  const float downGain2 = qdown - q2down;
202  if (std::abs(upGain2 - upGain) < std::abs(downGain2 - downGain))
203  return upGain;
204  else
205  return downGain;
206  }
207  }
208  }
209 
210  // The first element of the pair indicates presence of optical
211  // link errors. The second indicated presence of capid errors.
212  std::pair<bool,bool> findHWErrors(const HBHEDataFrame& df,
213  const unsigned len)
214  {
215  bool linkErr = false;
216  bool capidErr = false;
217  if (len)
218  {
219  int expectedCapid = df[0].capid();
220  for (unsigned i=0; i<len; ++i)
221  {
222  if (df[i].er())
223  linkErr = true;
224  if (df[i].capid() != expectedCapid)
225  capidErr = true;
226  expectedCapid = (expectedCapid + 1) % 4;
227  }
228  }
229  return std::pair<bool,bool>(linkErr, capidErr);
230  }
231 
232  std::pair<bool,bool> findHWErrors(const QIE11DataFrame& df,
233  const unsigned /* len */)
234  {
235  return std::pair<bool,bool>(df.linkError(), df.capidError());
236  }
237 
238  std::unique_ptr<HBHEStatusBitSetter> parse_HBHEStatusBitSetter(
239  const edm::ParameterSet& psdigi)
240  {
241  return std::make_unique<HBHEStatusBitSetter>(
242  psdigi.getParameter<double>("nominalPedestal"),
243  psdigi.getParameter<double>("hitEnergyMinimum"),
244  psdigi.getParameter<int>("hitMultiplicityThreshold"),
245  psdigi.getParameter<std::vector<edm::ParameterSet> >("pulseShapeParameterSets"));
246  }
247 
248  std::unique_ptr<HBHEPulseShapeFlagSetter> parse_HBHEPulseShapeFlagSetter(
249  const edm::ParameterSet& psPulseShape, const bool setLegacyFlags)
250  {
251  return std::make_unique<HBHEPulseShapeFlagSetter>(
252  psPulseShape.getParameter<double>("MinimumChargeThreshold"),
253  psPulseShape.getParameter<double>("TS4TS5ChargeThreshold"),
254  psPulseShape.getParameter<double>("TS3TS4ChargeThreshold"),
255  psPulseShape.getParameter<double>("TS3TS4UpperChargeThreshold"),
256  psPulseShape.getParameter<double>("TS5TS6ChargeThreshold"),
257  psPulseShape.getParameter<double>("TS5TS6UpperChargeThreshold"),
258  psPulseShape.getParameter<double>("R45PlusOneRange"),
259  psPulseShape.getParameter<double>("R45MinusOneRange"),
260  psPulseShape.getParameter<unsigned int>("TrianglePeakTS"),
261  psPulseShape.getParameter<std::vector<double> >("LinearThreshold"),
262  psPulseShape.getParameter<std::vector<double> >("LinearCut"),
263  psPulseShape.getParameter<std::vector<double> >("RMS8MaxThreshold"),
264  psPulseShape.getParameter<std::vector<double> >("RMS8MaxCut"),
265  psPulseShape.getParameter<std::vector<double> >("LeftSlopeThreshold"),
266  psPulseShape.getParameter<std::vector<double> >("LeftSlopeCut"),
267  psPulseShape.getParameter<std::vector<double> >("RightSlopeThreshold"),
268  psPulseShape.getParameter<std::vector<double> >("RightSlopeCut"),
269  psPulseShape.getParameter<std::vector<double> >("RightSlopeSmallThreshold"),
270  psPulseShape.getParameter<std::vector<double> >("RightSlopeSmallCut"),
271  psPulseShape.getParameter<std::vector<double> >("TS4TS5LowerThreshold"),
272  psPulseShape.getParameter<std::vector<double> >("TS4TS5LowerCut"),
273  psPulseShape.getParameter<std::vector<double> >("TS4TS5UpperThreshold"),
274  psPulseShape.getParameter<std::vector<double> >("TS4TS5UpperCut"),
275  psPulseShape.getParameter<bool>("UseDualFit"),
276  psPulseShape.getParameter<bool>("TriangleIgnoreSlow"),
277  setLegacyFlags);
278  }
279 }
280 
281 
282 //
283 // class declaration
284 //
286 {
287 public:
289  ~HBHEPhase1Reconstructor() override;
290 
291  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
292 
293 private:
294  void beginRun(edm::Run const&, edm::EventSetup const&) override;
295  void endRun(edm::Run const&, edm::EventSetup const&) override;
296  void produce(edm::Event&, const edm::EventSetup&) override;
297 
298  // Configuration parameters
306  bool tsFromDB_;
311 
312  // Parameters for turning status bit setters on/off
319 
320  // Other members
323  std::unique_ptr<AbsHBHEPhase1Algo> reco_;
324  std::unique_ptr<AbsHcalAlgoData> recoConfig_;
325  std::unique_ptr<HcalRecoParams> paramTS_;
326 
327  // Status bit setters
328  const HBHENegativeEFilter* negEFilter_; // We don't manage this pointer
329  std::unique_ptr<HBHEStatusBitSetter> hbheFlagSetterQIE8_;
330  std::unique_ptr<HBHEStatusBitSetter> hbheFlagSetterQIE11_;
331  std::unique_ptr<HBHEPulseShapeFlagSetter> hbhePulseShapeFlagSetterQIE8_;
332  std::unique_ptr<HBHEPulseShapeFlagSetter> hbhePulseShapeFlagSetterQIE11_;
333 
334  // For the function below, arguments "infoColl" and/or "rechits"
335  // are allowed to be null.
336  template<class DataFrame, class Collection>
337  void processData(const Collection& coll,
338  const HcalDbService& cond,
339  const HcalChannelQuality& qual,
341  const bool isRealData,
343  HBHEChannelInfoCollection* infoColl,
345 
346  // Methods for setting rechit status bits
347  void setAsicSpecificBits(const HBHEDataFrame& frame, const HcalCoder& coder,
348  const HBHEChannelInfo& info, const HcalCalibrations& calib,
349  HBHERecHit* rh);
350  void setAsicSpecificBits(const QIE11DataFrame& frame, const HcalCoder& coder,
351  const HBHEChannelInfo& info, const HcalCalibrations& calib,
352  HBHERecHit* rh);
353  void setCommonStatusBits(const HBHEChannelInfo& info, const HcalCalibrations& calib,
354  HBHERecHit* rh);
355 
356  void runHBHENegativeEFilter(const HBHEChannelInfo& info, HBHERecHit* rh);
357 };
358 
359 //
360 // constructors and destructor
361 //
363  : algoConfigClass_(conf.getParameter<std::string>("algoConfigClass")),
364  processQIE8_(conf.getParameter<bool>("processQIE8")),
365  processQIE11_(conf.getParameter<bool>("processQIE11")),
366  saveInfos_(conf.getParameter<bool>("saveInfos")),
367  saveDroppedInfos_(conf.getParameter<bool>("saveDroppedInfos")),
368  makeRecHits_(conf.getParameter<bool>("makeRecHits")),
369  dropZSmarkedPassed_(conf.getParameter<bool>("dropZSmarkedPassed")),
370  tsFromDB_(conf.getParameter<bool>("tsFromDB")),
371  recoParamsFromDB_(conf.getParameter<bool>("recoParamsFromDB")),
372  saveEffectivePedestal_(conf.getParameter<bool>("saveEffectivePedestal")),
373  sipmQTSShift_(conf.getParameter<int>("sipmQTSShift")),
374  sipmQNTStoSum_(conf.getParameter<int>("sipmQNTStoSum")),
375  setNegativeFlagsQIE8_(conf.getParameter<bool>("setNegativeFlagsQIE8")),
376  setNegativeFlagsQIE11_(conf.getParameter<bool>("setNegativeFlagsQIE11")),
377  setNoiseFlagsQIE8_(conf.getParameter<bool>("setNoiseFlagsQIE8")),
378  setNoiseFlagsQIE11_(conf.getParameter<bool>("setNoiseFlagsQIE11")),
379  setPulseShapeFlagsQIE8_(conf.getParameter<bool>("setPulseShapeFlagsQIE8")),
380  setPulseShapeFlagsQIE11_(conf.getParameter<bool>("setPulseShapeFlagsQIE11")),
381  reco_(parseHBHEPhase1AlgoDescription(conf.getParameter<edm::ParameterSet>("algorithm"))),
382  negEFilter_(nullptr)
383 {
384  // Check that the reco algorithm has been successfully configured
385  if (!reco_.get())
386  throw cms::Exception("HBHEPhase1BadConfig")
387  << "Invalid HBHEPhase1Algo algorithm configuration"
388  << std::endl;
389 
390  // Configure the status bit setters that have been turned on
391  if (setNoiseFlagsQIE8_)
392  hbheFlagSetterQIE8_ = parse_HBHEStatusBitSetter(
393  conf.getParameter<edm::ParameterSet>("flagParametersQIE8"));
394 
396  hbheFlagSetterQIE11_ = parse_HBHEStatusBitSetter(
397  conf.getParameter<edm::ParameterSet>("flagParametersQIE11"));
398 
400  hbhePulseShapeFlagSetterQIE8_ = parse_HBHEPulseShapeFlagSetter(
401  conf.getParameter<edm::ParameterSet>("pulseShapeParametersQIE8"),
402  conf.getParameter<bool>("setLegacyFlagsQIE8"));
403 
405  hbhePulseShapeFlagSetterQIE11_ = parse_HBHEPulseShapeFlagSetter(
406  conf.getParameter<edm::ParameterSet>("pulseShapeParametersQIE11"),
407  conf.getParameter<bool>("setLegacyFlagsQIE11"));
408 
409  // Consumes and produces statements
410  if (processQIE8_)
411  tok_qie8_ = consumes<HBHEDigiCollection>(
412  conf.getParameter<edm::InputTag>("digiLabelQIE8"));
413 
414  if (processQIE11_)
415  tok_qie11_ = consumes<QIE11DigiCollection>(
416  conf.getParameter<edm::InputTag>("digiLabelQIE11"));
417 
418  if (saveInfos_)
419  produces<HBHEChannelInfoCollection>();
420 
421  if (makeRecHits_)
422  produces<HBHERecHitCollection>();
423 }
424 
425 
427 {
428  // do anything here that needs to be done at destruction time
429  // (e.g. close files, deallocate resources etc.)
430 }
431 
432 
433 //
434 // member functions
435 //
436 template<class DFrame, class Collection>
438  const HcalDbService& cond,
439  const HcalChannelQuality& qual,
441  const bool isRealData,
442  HBHEChannelInfo* channelInfo,
445 {
446  // If "saveDroppedInfos_" flag is set, fill the info with something
447  // meaningful even if the database tells us to drop this channel.
448  // Note that this flag affects only "infos", the rechits are still
449  // not going to be constructed from such channels.
450  const bool skipDroppedChannels = !(infos && saveDroppedInfos_);
451 
452  // Iterate over the input collection
453  for (typename Collection::const_iterator it = coll.begin();
454  it != coll.end(); ++it)
455  {
456  const DFrame& frame(*it);
457  const HcalDetId cell(frame.id());
458 
459  // Protection against calibration channels which are not
460  // in the database but can still come in the QIE11DataFrame
461  // in the laser calibs, etc.
462  const HcalSubdetector subdet = cell.subdet();
463  if (!(subdet == HcalSubdetector::HcalBarrel ||
464  subdet == HcalSubdetector::HcalEndcap ||
465  subdet == HcalSubdetector::HcalOuter))
466  continue;
467 
468  // Check if the database tells us to drop this channel
469  const HcalChannelStatus* mydigistatus = qual.getValues(cell.rawId());
470  const bool taggedBadByDb = severity.dropChannel(mydigistatus->getValue());
471  if (taggedBadByDb && skipDroppedChannels)
472  continue;
473 
474  // Check if the channel is zero suppressed
475  bool dropByZS = false;
477  if (frame.zsMarkAndPass())
478  dropByZS = true;
479  if (dropByZS && skipDroppedChannels)
480  continue;
481 
482  // Basic ADC decoding tools
483  const HcalRecoParam* param_ts = paramTS_->getValues(cell.rawId());
484  const HcalCalibrations& calib = cond.getHcalCalibrations(cell);
485  const HcalCalibrationWidths& calibWidth = cond.getHcalCalibrationWidths(cell);
486  const HcalQIECoder* channelCoder = cond.getHcalCoder(cell);
487  const HcalQIEShape* shape = cond.getHcalShape(channelCoder);
488  const HcalCoderDb coder(*channelCoder, *shape);
489 
490  // needed for the dark current in the M2
491  const HcalSiPMParameter& siPMParameter(*cond.getHcalSiPMParameter(cell));
492  const double darkCurrent = siPMParameter.getDarkCurrent();
493  const double fcByPE = siPMParameter.getFCByPE();
494  const double lambda = cond.getHcalSiPMCharacteristics()->getCrossTalk(siPMParameter.getType());
495 
496  // ADC to fC conversion
497  CaloSamples cs;
498  coder.adc2fC(frame, cs);
499 
500  // Prepare to iterate over time slices
501  const int nRead = cs.size();
502  const int maxTS = std::min(nRead, static_cast<int>(HBHEChannelInfo::MAXSAMPLES));
503  const int soi = tsFromDB_ ? param_ts->firstSample() : frame.presamples();
504  const RawChargeFromSample<DFrame> rcfs(sipmQTSShift_, sipmQNTStoSum_, saveEffectivePedestal_,
505  cond, cell, cs, soi, frame, maxTS);
506  int soiCapid = 4;
507 
508  // Go over time slices and fill the samples
509  for (int ts = 0; ts < maxTS; ++ts)
510  {
511  auto s(frame[ts]);
512  const uint8_t adc = s.adc();
513  const int capid = s.capid();
514  //optionally store "effective" pedestal = QIE contribution (default, from calib.pedestal()) + SiPM contribution (dark current + crosstalk)
515  //only done for pedestal mean, to be used for pedestal subtraction downstream
516  const double pedestal = calib.pedestal(capid) + (saveEffectivePedestal_ ? darkCurrent * 25. / (1. - lambda) : 0.);
517  const double pedestalWidth = calibWidth.pedestal(capid);
518  const double gain = calib.respcorrgain(capid);
519  const double gainWidth = calibWidth.gain(capid);
520  const double rawCharge = rcfs.getRawCharge(cs[ts], pedestal);
521  const float t = getTDCTimeFromSample(s);
522  const float dfc = getDifferentialChargeGain(*channelCoder, *shape, adc,
523  capid, channelInfo->hasTimeInfo());
524  channelInfo->setSample(ts, adc, dfc, rawCharge,
525  pedestal, pedestalWidth,
526  gain, gainWidth, t);
527  if (ts == soi)
528  soiCapid = capid;
529  }
530 
531  // Fill the overall channel info items
532  const int pulseShapeID = param_ts->pulseShapeID();
533  const std::pair<bool,bool> hwerr = findHWErrors(frame, maxTS);
534  channelInfo->setChannelInfo(cell, pulseShapeID, maxTS, soi, soiCapid,
535  darkCurrent, fcByPE, lambda,
536  hwerr.first, hwerr.second,
537  taggedBadByDb || dropByZS);
538 
539  // If needed, add the channel info to the output collection
540  const bool makeThisRechit = !channelInfo->isDropped();
541  if (infos && (saveDroppedInfos_ || makeThisRechit))
542  infos->push_back(*channelInfo);
543 
544  // Reconstruct the rechit
545  if (rechits && makeThisRechit)
546  {
547  const HcalRecoParam* pptr = nullptr;
548  if (recoParamsFromDB_)
549  pptr = param_ts;
550  HBHERecHit rh = reco_->reconstruct(*channelInfo, pptr, calib, isRealData);
551  if (rh.id().rawId())
552  {
553  setAsicSpecificBits(frame, coder, *channelInfo, calib, &rh);
554  setCommonStatusBits(*channelInfo, calib, &rh);
555  rechits->push_back(rh);
556  }
557  }
558  }
559 }
560 
562  const HBHEChannelInfo& /* info */, const HcalCalibrations& /* calib */,
563  HBHERecHit* /* rh */)
564 {
565 }
566 
568  const HBHEDataFrame& frame, const HcalCoder& coder,
569  const HBHEChannelInfo& info, const HcalCalibrations& calib,
570  HBHERecHit* rh)
571 {
572  if (setNoiseFlagsQIE8_)
573  hbheFlagSetterQIE8_->rememberHit(*rh);
574 
576  hbhePulseShapeFlagSetterQIE8_->SetPulseShapeFlags(*rh, frame, coder, calib);
577 
579  runHBHENegativeEFilter(info, rh);
580 }
581 
583  const QIE11DataFrame& frame, const HcalCoder& coder,
584  const HBHEChannelInfo& info, const HcalCalibrations& calib,
585  HBHERecHit* rh)
586 {
588  hbheFlagSetterQIE11_->rememberHit(*rh);
589 
591  hbhePulseShapeFlagSetterQIE11_->SetPulseShapeFlags(*rh, frame, coder, calib);
592 
594  runHBHENegativeEFilter(info, rh);
595 }
596 
598  HBHERecHit* rh)
599 {
600  double ts[HBHEChannelInfo::MAXSAMPLES];
601  const unsigned nRead = info.nSamples();
602  for (unsigned i=0; i<nRead; ++i)
603  ts[i] = info.tsCharge(i);
604  const bool passes = negEFilter_->checkPassFilter(info.id(), &ts[0], nRead);
605  if (!passes)
607 }
608 
609 // ------------ method called to produce the data ------------
610 void
612 {
613  using namespace edm;
614 
615  // Get the Hcal topology
617  eventSetup.get<HcalRecNumberingRecord>().get(htopo);
618  paramTS_->setTopo(htopo.product());
619 
620  // Fetch the calibrations
621  ESHandle<HcalDbService> conditions;
622  eventSetup.get<HcalDbRecord>().get(conditions);
623 
625  eventSetup.get<HcalChannelQualityRcd>().get("withTopo", p);
626 
628  eventSetup.get<HcalSeverityLevelComputerRcd>().get(mycomputer);
629 
630  // Configure the negative energy filter
633  {
634  eventSetup.get<HBHENegativeEFilterRcd>().get(negEHandle);
635  negEFilter_ = negEHandle.product();
636  }
637 
638  // Find the input data
639  unsigned maxOutputSize = 0;
641  if (processQIE8_)
642  {
643  e.getByToken(tok_qie8_, hbDigis);
644  maxOutputSize += hbDigis->size();
645  }
646 
648  if (processQIE11_)
649  {
650  e.getByToken(tok_qie11_, heDigis);
651  maxOutputSize += heDigis->size();
652  }
653 
654  // Create new output collections
655  std::unique_ptr<HBHEChannelInfoCollection> infos;
656  if (saveInfos_)
657  {
658  infos = std::make_unique<HBHEChannelInfoCollection>();
659  infos->reserve(maxOutputSize);
660  }
661 
662  std::unique_ptr<HBHERecHitCollection> out;
663  if (makeRecHits_)
664  {
665  out = std::make_unique<HBHERecHitCollection>();
666  out->reserve(maxOutputSize);
667  }
668 
669  // Process the input collections, filling the output ones
670  const bool isData = e.isRealData();
671  if (processQIE8_)
672  {
673  if (setNoiseFlagsQIE8_)
674  hbheFlagSetterQIE8_->Clear();
675 
676  HBHEChannelInfo channelInfo(false);
677  processData<HBHEDataFrame>(*hbDigis, *conditions, *p, *mycomputer,
678  isData, &channelInfo, infos.get(), out.get());
679  if (setNoiseFlagsQIE8_)
680  hbheFlagSetterQIE8_->SetFlagsFromRecHits(*out);
681  }
682 
683  if (processQIE11_)
684  {
686  hbheFlagSetterQIE11_->Clear();
687 
688  HBHEChannelInfo channelInfo(true);
689  processData<QIE11DataFrame>(*heDigis, *conditions, *p, *mycomputer,
690  isData, &channelInfo, infos.get(), out.get());
692  hbheFlagSetterQIE11_->SetFlagsFromRecHits(*out);
693  }
694 
695  // Add the output collections to the event record
696  if (saveInfos_)
697  e.put(std::move(infos));
698  if (makeRecHits_)
699  e.put(std::move(out));
700 }
701 
702 // ------------ method called when starting to processes a run ------------
703 void
705 {
707  es.get<HcalRecoParamsRcd>().get(p);
708  paramTS_ = std::make_unique<HcalRecoParams>(*p.product());
709 
710  if (reco_->isConfigurable())
711  {
713  if (!recoConfig_.get())
714  throw cms::Exception("HBHEPhase1BadConfig")
715  << "Invalid HBHEPhase1Reconstructor \"algoConfigClass\" parameter value \""
716  << algoConfigClass_ << '"' << std::endl;
717  if (!reco_->configure(recoConfig_.get()))
718  throw cms::Exception("HBHEPhase1BadConfig")
719  << "Failed to configure HBHEPhase1Algo algorithm from EventSetup"
720  << std::endl;
721  }
722 
724  {
726  es.get<HcalFrontEndMapRcd>().get(hfemap);
727  if (hfemap.isValid())
728  {
729  if (setNoiseFlagsQIE8_)
730  hbheFlagSetterQIE8_->SetFrontEndMap(hfemap.product());
732  hbheFlagSetterQIE11_->SetFrontEndMap(hfemap.product());
733  }
734  else
735  edm::LogWarning("EventSetup") <<
736  "HBHEPhase1Reconstructor failed to get HcalFrontEndMap!" << std::endl;
737  }
738 
739  reco_->beginRun(r, es);
740 }
741 
742 void
744 {
745  reco_->endRun();
746 }
747 
748 #define add_param_set(name) \
749  edm::ParameterSetDescription name; \
750  name.setAllowAnything(); \
751  desc.add<edm::ParameterSetDescription>(#name, name)
752 
753 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
754 void
756 {
758 
759  desc.add<edm::InputTag>("digiLabelQIE8");
760  desc.add<edm::InputTag>("digiLabelQIE11");
761  desc.add<std::string>("algoConfigClass");
762  desc.add<bool>("processQIE8");
763  desc.add<bool>("processQIE11");
764  desc.add<bool>("saveInfos");
765  desc.add<bool>("saveDroppedInfos");
766  desc.add<bool>("makeRecHits");
767  desc.add<bool>("dropZSmarkedPassed");
768  desc.add<bool>("tsFromDB");
769  desc.add<bool>("recoParamsFromDB");
770  desc.add<bool>("saveEffectivePedestal", false);
771  desc.add<int>("sipmQTSShift", 0);
772  desc.add<int>("sipmQNTStoSum", 3);
773  desc.add<bool>("setNegativeFlagsQIE8");
774  desc.add<bool>("setNegativeFlagsQIE11");
775  desc.add<bool>("setNoiseFlagsQIE8");
776  desc.add<bool>("setNoiseFlagsQIE11");
777  desc.add<bool>("setPulseShapeFlagsQIE8");
778  desc.add<bool>("setPulseShapeFlagsQIE11");
779  desc.add<bool>("setLegacyFlagsQIE8");
780  desc.add<bool>("setLegacyFlagsQIE11");
781 
787 
788  descriptions.addDefault(desc);
789 }
790 
791 //define this as a plug-in
int adc(sample_type sample)
get the ADC sample (12 bits)
HBHEPhase1Reconstructor(const edm::ParameterSet &)
std::unique_ptr< HBHEPulseShapeFlagSetter > hbhePulseShapeFlagSetterQIE8_
unsigned int firstSample() const
Definition: HcalRecoParam.h:32
T getParameter(std::string const &) const
const HBHENegativeEFilter * negEFilter_
float getTDCTime(const int tdc)
double gain(int fCapId) const
get gain width for capid=0..3
static const TGPicture * info(bool iBackgroundIsBlack)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:127
double respcorrgain(int fCapId) const
get response corrected gain for capid=0..3
auto_ptr< ClusterSequence > cs
int maxTS(DIGI const &digi, double ped=0)
Definition: Utilities.h:94
HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:49
bool hasTimeInfo() const
void setCommonStatusBits(const HBHEChannelInfo &info, const HcalCalibrations &calib, HBHERecHit *rh)
std::unique_ptr< HBHEStatusBitSetter > hbheFlagSetterQIE11_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:508
std::unique_ptr< HBHEStatusBitSetter > hbheFlagSetterQIE8_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void beginRun(edm::Run const &, edm::EventSetup const &) override
HcalDetId id() const
get the id
Definition: HBHERecHit.h:25
void setFlagField(uint32_t value, int base, int width=1)
Definition: CaloRecHit.cc:20
void push_back(T const &t)
double pedestal(int fCapId) const
get pedestal for capid=0..3
const Item * getValues(DetId fId, bool throwOnFail=true) const
#define nullptr
std::unique_ptr< AbsHcalAlgoData > fetchHcalAlgoData(const std::string &className, const edm::EventSetup &es)
void setSample(const unsigned ts, const uint8_t rawADC, const float differentialChargeGain, const double q, const double ped, const double pedWidth, const double g, const double gainWidth, const float t)
bool isRealData() const
Definition: EventBase.h:64
HcalDetId id() const
float getCrossTalk(int type) const
get cross talk
bool checkPassFilter(const HcalDetId &id, const double *ts, unsigned lenTS) const
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
void adc2fC(const HBHEDataFrame &df, CaloSamples &lf) const override
Definition: HcalCoderDb.cc:68
std::unique_ptr< HcalRecoParams > paramTS_
void addDefault(ParameterSetDescription const &psetDescription)
static const unsigned MAXSAMPLES
double pedestal(int fCapId) const
get pedestal width for capid=0..3
void endRun(edm::Run const &, edm::EventSetup const &) override
bool linkError() const
void processData(const Collection &coll, const HcalDbService &cond, const HcalChannelQuality &qual, const HcalSeverityLevelComputer &severity, const bool isRealData, HBHEChannelInfo *info, HBHEChannelInfoCollection *infoColl, HBHERecHitCollection *rechits)
edm::EDGetTokenT< QIE11DigiCollection > tok_qie11_
edm::EDGetTokenT< HBHEDigiCollection > tok_qie8_
std::unique_ptr< AbsHBHEPhase1Algo > parseHBHEPhase1AlgoDescription(const edm::ParameterSet &ps)
HcalSubdetector
Definition: HcalAssistant.h:31
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool dropChannel(const uint32_t &mystatus) const
double f[11][100]
double tsCharge(const unsigned ts) const
std::unique_ptr< AbsHBHEPhase1Algo > reco_
T min(T a, T b)
Definition: MathUtil.h:58
float getDarkCurrent() const
get dark current
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void runHBHENegativeEFilter(const HBHEChannelInfo &info, HBHERecHit *rh)
const HcalCalibrationWidths & getHcalCalibrationWidths(const HcalGenericDetId &fId) const
bool capidError() const
unsigned int pulseShapeID() const
Definition: HcalRecoParam.h:34
std::unique_ptr< HBHEPulseShapeFlagSetter > hbhePulseShapeFlagSetterQIE11_
JetCorrectorParametersCollection coll
Definition: classes.h:10
int size() const
get the size
Definition: CaloSamples.h:24
const T & get() const
Definition: EventSetup.h:55
const HcalSiPMCharacteristics * getHcalSiPMCharacteristics() const
const HcalQIECoder * getHcalCoder(const HcalGenericDetId &fId) const
void setAsicSpecificBits(const HBHEDataFrame &frame, const HcalCoder &coder, const HBHEChannelInfo &info, const HcalCalibrations &calib, HBHERecHit *rh)
const HcalQIEShape * getHcalShape(const HcalGenericDetId &fId) const
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
constexpr float UNKNOWN_T_NOTDC
Definition: plugin.cc:24
HLT enums.
size_type size() const
void setChannelInfo(const HcalDetId &detId, const int recoShape, const unsigned nSamp, const unsigned iSoi, const int iCapid, const double darkCurrent, const double fcByPE, const double lambda, const bool linkError, const bool capidError, const bool dropThisChannel)
bool isDropped() const
std::unique_ptr< AbsHcalAlgoData > recoConfig_
unsigned nSamples() const
uint32_t getValue() const
bool isValid() const
Definition: ESHandle.h:47
const HcalCalibrations & getHcalCalibrations(const HcalGenericDetId &fId) const
const HcalSiPMParameter * getHcalSiPMParameter(const HcalGenericDetId &fId) const
T const * product() const
Definition: ESHandle.h:86
def move(src, dest)
Definition: eostools.py:510
void produce(edm::Event &, const edm::EventSetup &) override
Definition: Run.h:43
#define add_param_set(name)
float charge(const HcalQIEShape &fShape, unsigned fAdc, unsigned fCapId) const
ADC [0..127] + capid [0..3] -> fC conversion.
Definition: HcalQIECoder.cc:22