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