CMS 3D CMS Logo

HcaluLUTTPGCoder.cc
Go to the documentation of this file.
1 #include <iostream>
2 #include <fstream>
3 #include <cmath>
4 #include <string>
31 
32 const float HcaluLUTTPGCoder::lsb_ = 1. / 16;
33 
37 
38 constexpr double MaximumFractionalError = 0.002; // 0.2% error allowed from this source
39 
41  : topo_{},
42  delay_{},
45  bitToMask_{},
46  firstHBEta_{},
47  lastHBEta_{},
48  nHBEta_{},
49  maxDepthHB_{},
50  sizeHB_{},
51  firstHEEta_{},
52  lastHEEta_{},
53  nHEEta_{},
54  maxDepthHE_{},
55  sizeHE_{},
56  firstHFEta_{},
57  lastHFEta_{},
58  nHFEta_{},
59  maxDepthHF_{},
60  sizeHF_{},
64  allLinear_{},
68 
69 HcaluLUTTPGCoder::HcaluLUTTPGCoder(const HcalTopology* top, const HcalTimeSlew* delay) { init(top, delay); }
70 
71 void HcaluLUTTPGCoder::init(const HcalTopology* top, const HcalTimeSlew* delay) {
72  topo_ = top;
73  delay_ = delay;
74  LUTGenerationMode_ = true;
75  FG_HF_thresholds_ = {0, 0};
76  bitToMask_ = 0;
77  allLinear_ = false;
78  linearLSB_QIE8_ = 1.;
79  linearLSB_QIE11_ = 1.;
81  pulseCorr_ = std::make_unique<HcalPulseContainmentManager>(MaximumFractionalError);
84  nHBEta_ = (lastHBEta_ - firstHBEta_ + 1);
86  sizeHB_ = 2 * nHBEta_ * nFi_ * maxDepthHB_;
89  nHEEta_ = (lastHEEta_ - firstHEEta_ + 1);
91  sizeHE_ = 2 * nHEEta_ * nFi_ * maxDepthHE_;
94  nHFEta_ = (lastHFEta_ - firstHFEta_ + 1);
96  sizeHF_ = 2 * nHFEta_ * nFi_ * maxDepthHF_;
97  size_t nluts = (size_t)(sizeHB_ + sizeHE_ + sizeHF_ + 1);
98  inputLUT_ = std::vector<HcaluLUTTPGCoder::Lut>(nluts);
99  gain_ = std::vector<float>(nluts, 0.);
100  ped_ = std::vector<float>(nluts, 0.);
102 }
103 
105  const std::vector<bool>& featureBits,
106  HcalTriggerPrimitiveDigi& tp) const {
107  throw cms::Exception("PROBLEM: This method should never be invoked!");
108 }
109 
111 
113  int retval(0);
114  if (id == HcalBarrel) {
115  retval = (depth - 1) + maxDepthHB_ * (iphi - 1);
116  if (ieta > 0)
117  retval += maxDepthHB_ * nFi_ * (ieta - firstHBEta_);
118  else
119  retval += maxDepthHB_ * nFi_ * (ieta + lastHBEta_ + nHBEta_);
120  } else if (id == HcalEndcap) {
121  retval = sizeHB_;
122  retval += (depth - 1) + maxDepthHE_ * (iphi - 1);
123  if (ieta > 0)
124  retval += maxDepthHE_ * nFi_ * (ieta - firstHEEta_);
125  else
126  retval += maxDepthHE_ * nFi_ * (ieta + lastHEEta_ + nHEEta_);
127  } else if (id == HcalForward) {
128  retval = sizeHB_ + sizeHE_;
129  retval += (depth - 1) + maxDepthHF_ * (iphi - 1);
130  if (ieta > 0)
131  retval += maxDepthHF_ * nFi_ * (ieta - firstHFEta_);
132  else
133  retval += maxDepthHF_ * nFi_ * (ieta + lastHFEta_ + nHFEta_);
134  }
135  return retval;
136 }
137 
138 int HcaluLUTTPGCoder::getLUTId(uint32_t rawid) const {
139  HcalDetId detid(rawid);
140  return getLUTId(detid.subdet(), detid.ieta(), detid.iphi(), detid.depth());
141 }
142 
143 int HcaluLUTTPGCoder::getLUTId(const HcalDetId& detid) const {
144  return getLUTId(detid.subdet(), detid.ieta(), detid.iphi(), detid.depth());
145 }
146 
147 void HcaluLUTTPGCoder::update(const char* filename, bool appendMSB) {
148  std::ifstream file(filename, std::ios::in);
149  assert(file.is_open());
150 
151  std::vector<HcalSubdetector> subdet;
153 
154  // Drop first (comment) line
155  std::getline(file, buffer);
156  std::getline(file, buffer);
157 
158  unsigned int index = buffer.find("H", 0);
159  while (index < buffer.length()) {
160  std::string subdetStr = buffer.substr(index, 2);
161  if (subdetStr == "HB")
162  subdet.push_back(HcalBarrel);
163  else if (subdetStr == "HE")
164  subdet.push_back(HcalEndcap);
165  else if (subdetStr == "HF")
166  subdet.push_back(HcalForward);
167  //TODO Check subdet
168  //else exception
169  index += 2;
170  index = buffer.find("H", index);
171  }
172 
173  // Get upper/lower ranges for ieta/iphi/depth
174  size_t nCol = subdet.size();
175  assert(nCol > 0);
176 
177  std::vector<int> ietaU;
178  std::vector<int> ietaL;
179  std::vector<int> iphiU;
180  std::vector<int> iphiL;
181  std::vector<int> depU;
182  std::vector<int> depL;
183  std::vector<Lut> lutFromFile(nCol);
184  LutElement lutValue;
185 
186  for (size_t i = 0; i < nCol; ++i) {
187  int ieta;
188  file >> ieta;
189  ietaL.push_back(ieta);
190  }
191 
192  for (size_t i = 0; i < nCol; ++i) {
193  int ieta;
194  file >> ieta;
195  ietaU.push_back(ieta);
196  }
197 
198  for (size_t i = 0; i < nCol; ++i) {
199  int iphi;
200  file >> iphi;
201  iphiL.push_back(iphi);
202  }
203 
204  for (size_t i = 0; i < nCol; ++i) {
205  int iphi;
206  file >> iphi;
207  iphiU.push_back(iphi);
208  }
209 
210  for (size_t i = 0; i < nCol; ++i) {
211  int dep;
212  file >> dep;
213  depL.push_back(dep);
214  }
215 
216  for (size_t i = 0; i < nCol; ++i) {
217  int dep;
218  file >> dep;
219  depU.push_back(dep);
220  }
221 
222  // Read Lut Entry
223  for (size_t i = 0; file >> lutValue; i = (i + 1) % nCol) {
224  lutFromFile[i].push_back(lutValue);
225  }
226 
227  // Check lut size
228  for (size_t i = 0; i < nCol; ++i)
229  assert(lutFromFile[i].size() == INPUT_LUT_SIZE);
230 
231  for (size_t i = 0; i < nCol; ++i) {
232  for (int ieta = ietaL[i]; ieta <= ietaU[i]; ++ieta) {
233  for (int iphi = iphiL[i]; iphi <= iphiU[i]; ++iphi) {
234  for (int depth = depL[i]; depth <= depU[i]; ++depth) {
235  HcalDetId id(subdet[i], ieta, iphi, depth);
236  if (!topo_->valid(id))
237  continue;
238 
239  int lutId = getLUTId(id);
240  for (size_t adc = 0; adc < INPUT_LUT_SIZE; ++adc) {
241  if (appendMSB) {
242  // Append FG bit LUT to MSB
243  // MSB = Most Significant Bit = bit 10
244  // Overwrite bit 10
245  LutElement msb = (lutFromFile[i][adc] != 0 ? QIE8_LUT_MSB : 0);
246  inputLUT_[lutId][adc] = (msb | (inputLUT_[lutId][adc] & QIE8_LUT_BITMASK));
247  } else
248  inputLUT_[lutId][adc] = lutFromFile[i][adc];
249  } // for adc
250  } // for depth
251  } // for iphi
252  } // for ieta
253  } // for nCol
254 }
255 
257  LutXml* _xml = new LutXml(filename);
258  _xml->create_lut_map();
261  for (unsigned int iphi = 0; iphi <= HcalDetId::kHcalPhiMask2; ++iphi) {
262  for (unsigned int depth = 1; depth < HcalDetId::kHcalDepthMask2; ++depth) {
263  for (int isub = 0; isub < 3; ++isub) {
264  HcalDetId detid(subdet[isub], ieta, iphi, depth);
265  if (!topo_->valid(detid))
266  continue;
267  int id = getLUTId(subdet[isub], ieta, iphi, depth);
268  std::vector<unsigned int>* lut = _xml->getLutFast(detid);
269  if (lut == nullptr)
270  throw cms::Exception("PROBLEM: No inputLUT_ in xml file for ") << detid << std::endl;
271  if (lut->size() != INPUT_LUT_SIZE)
272  throw cms::Exception("PROBLEM: Wrong inputLUT_ size in xml file for ") << detid << std::endl;
273  for (unsigned int i = 0; i < INPUT_LUT_SIZE; ++i)
274  inputLUT_[id][i] = (LutElement)lut->at(i);
275  }
276  }
277  }
278  }
279  delete _xml;
281 }
282 
284  // ieta = 28 and 29 are both associated with trigger tower 28
285  // so special handling is required. HF ieta=29 channels included in TT30
286  // are already handled correctly in cosh_ieta_
287  if (abs(ieta) >= 28 && subdet == HcalEndcap && allLinear_) {
288  if (abs(ieta) == 29)
289  return cosh_ieta_29_HE_;
290  if (abs(ieta) == 28) {
291  if (depth <= 3)
293  else
295  }
296  }
297 
298  return cosh_ieta_[ieta];
299 }
300 
302  cosh_ieta_ = std::vector<double>(lastHFEta_ + 1, -1.0);
303 
304  HcalTrigTowerGeometry triggeo(topo_);
305 
306  for (int i = 1; i <= firstHFEta_; ++i) {
307  double eta_low = 0., eta_high = 0.;
308  triggeo.towerEtaBounds(i, 0, eta_low, eta_high);
309  cosh_ieta_[i] = cosh((eta_low + eta_high) / 2.);
310  }
311  for (int i = firstHFEta_; i <= lastHFEta_; ++i) {
312  std::pair<double, double> etas = topo_->etaRange(HcalForward, i);
313  double eta1 = etas.first;
314  double eta2 = etas.second;
315  cosh_ieta_[i] = cosh((eta1 + eta2) / 2.);
316  }
317 
318  // trigger tower 28 in HE has a more complicated geometry
319  std::pair<double, double> eta28 = topo_->etaRange(HcalEndcap, 28);
320  std::pair<double, double> eta29 = topo_->etaRange(HcalEndcap, 29);
321  cosh_ieta_29_HE_ = cosh((eta29.first + eta29.second) / 2.);
322  cosh_ieta_28_HE_low_depths_ = cosh((eta28.first + eta28.second) / 2.);
323  // for higher depths in ieta = 28, the trigger tower extends past
324  // the ieta = 29 channels
325  cosh_ieta_28_HE_high_depths_ = cosh((eta28.first + eta29.second) / 2.);
326 }
327 
328 void HcaluLUTTPGCoder::update(const HcalDbService& conditions) {
329  const HcalLutMetadata* metadata = conditions.getHcalLutMetadata();
330  assert(metadata != nullptr);
331  float nominalgain_ = metadata->getNominalGain();
332 
333  pulseCorr_->beginRun(&conditions, delay_);
334 
336 
337  for (const auto& id : metadata->getAllChannels()) {
338  if (not(id.det() == DetId::Hcal and topo_->valid(id)))
339  continue;
340 
341  HcalDetId cell(id);
342  HcalSubdetector subdet = cell.subdet();
343 
344  if (subdet != HcalBarrel and subdet != HcalEndcap and subdet != HcalForward)
345  continue;
346 
347  const HcalQIECoder* channelCoder = conditions.getHcalCoder(cell);
348  const HcalQIEShape* shape = conditions.getHcalShape(cell);
349  HcalCoderDb coder(*channelCoder, *shape);
350  const HcalLutMetadatum* meta = metadata->getValues(cell);
351 
352  unsigned int mipMax = 0;
353  unsigned int mipMin = 0;
354 
355  bool is2018OrLater = topo_->triggerMode() >= HcalTopologyMode::TriggerMode_2018 or
357  if (is2018OrLater or topo_->dddConstants()->isPlan1(cell)) {
358  const HcalTPChannelParameter* channelParameters = conditions.getHcalTPChannelParameter(cell);
359  mipMax = channelParameters->getFGBitInfo() >> 16;
360  mipMin = channelParameters->getFGBitInfo() & 0xFFFF;
361  }
362 
363  int lutId = getLUTId(cell);
364  Lut& lut = inputLUT_[lutId];
365  float ped = 0;
366  float gain = 0;
367  uint32_t status = 0;
368 
369  if (LUTGenerationMode_) {
370  const HcalCalibrations& calibrations = conditions.getHcalCalibrations(cell);
371  for (auto capId : {0, 1, 2, 3}) {
372  ped += calibrations.effpedestal(capId);
373  gain += calibrations.LUTrespcorrgain(capId);
374  }
375  ped /= 4.0;
376  gain /= 4.0;
377 
378  //Get Channel Quality
379  const HcalChannelStatus* channelStatus = conditions.getHcalChannelStatus(cell);
380  status = channelStatus->getValue();
381 
382  } else {
383  const HcalL1TriggerObject* myL1TObj = conditions.getHcalL1TriggerObject(cell);
384  ped = myL1TObj->getPedestal();
385  gain = myL1TObj->getRespGain();
386  status = myL1TObj->getFlag();
387  } // LUTGenerationMode_
388 
389  ped_[lutId] = ped;
390  gain_[lutId] = gain;
391  bool isMasked = ((status & bitToMask_) > 0);
392  float rcalib = meta->getRCalib();
393 
394  auto adc2fC = [channelCoder, shape](unsigned int adc) {
395  float fC = 0;
396  for (auto capId : {0, 1, 2, 3})
397  fC += channelCoder->charge(*shape, adc, capId);
398  return fC / 4;
399  };
400 
401  int qieType = conditions.getHcalQIEType(cell)->getValue();
402 
403  const size_t SIZE = qieType == QIE8 ? INPUT_LUT_SIZE : UPGRADE_LUT_SIZE;
404  const int MASK = qieType == QIE8 ? QIE8_LUT_BITMASK : qieType == QIE10 ? QIE10_LUT_BITMASK : QIE11_LUT_BITMASK;
405  double linearLSB = linearLSB_QIE8_;
406  if (qieType == QIE11 and cell.ietaAbs() == topo_->lastHBRing())
407  linearLSB = linearLSB_QIE11Overlap_;
408  else if (qieType == QIE11)
409  linearLSB = linearLSB_QIE11_;
410 
411  lut.resize(SIZE, 0);
412 
413  // Input LUT for HB/HE/HF
414  if (subdet == HcalBarrel || subdet == HcalEndcap) {
415  int granularity = meta->getLutGranularity();
416 
417  double correctionPhaseNS = conditions.getHcalRecoParam(cell)->correctionPhaseNS();
418  for (unsigned int adc = 0; adc < SIZE; ++adc) {
419  if (isMasked)
420  lut[adc] = 0;
421  else {
422  double nonlinearityCorrection = 1.0;
423  double containmentCorrection2TSCorrected = 1.0;
424  // SiPM nonlinearity was not corrected in 2017
425  // and containment corrections were not
426  // ET-dependent prior to 2018
427  if (is2018OrLater) {
428  double containmentCorrection1TS = pulseCorr_->correction(cell, 1, correctionPhaseNS, adc2fC(adc));
429  // Use the 1-TS containment correction to estimate the charge of the pulse
430  // from the individual samples
431  double correctedCharge = containmentCorrection1TS * adc2fC(adc);
432  containmentCorrection2TSCorrected = pulseCorr_->correction(cell, 2, correctionPhaseNS, correctedCharge);
433  if (qieType == QIE11) {
434  const HcalSiPMParameter& siPMParameter(*conditions.getHcalSiPMParameter(cell));
436  conditions.getHcalSiPMCharacteristics()->getNonLinearities(siPMParameter.getType()));
437  const double fcByPE = siPMParameter.getFCByPE();
438  const double effectivePixelsFired = correctedCharge / fcByPE;
439  nonlinearityCorrection = corr.getRecoCorrectionFactor(effectivePixelsFired);
440  }
441  }
442  if (allLinear_)
443  lut[adc] = (LutElement)std::min(std::max(0,
444  int((adc2fC(adc) - ped) * gain * rcalib * nonlinearityCorrection *
445  containmentCorrection2TSCorrected / linearLSB /
446  cosh_ieta(cell.ietaAbs(), cell.depth(), HcalEndcap))),
447  MASK);
448  else
449  lut[adc] =
451  int((adc2fC(adc) - ped) * gain * rcalib * nonlinearityCorrection *
452  containmentCorrection2TSCorrected / nominalgain_ / granularity)),
453  MASK);
454 
455  if (qieType == QIE11) {
456  if (adc >= mipMin and adc < mipMax)
457  lut[adc] |= QIE11_LUT_MSB0;
458  else if (adc >= mipMax)
459  lut[adc] |= QIE11_LUT_MSB1;
460  }
461  }
462  }
463  } else if (subdet == HcalForward) {
464  for (unsigned int adc = 0; adc < SIZE; ++adc) {
465  if (isMasked)
466  lut[adc] = 0;
467  else {
468  lut[adc] =
469  std::min(std::max(0, int((adc2fC(adc) - ped) * gain * rcalib / lsb_ / cosh_ieta_[cell.ietaAbs()])), MASK);
470  if (adc > FG_HF_thresholds_[0])
471  lut[adc] |= QIE10_LUT_MSB0;
472  if (adc > FG_HF_thresholds_[1])
473  lut[adc] |= QIE10_LUT_MSB1;
474  }
475  }
476  }
477  }
478 }
479 
481  int lutId = getLUTId(df.id());
482  const Lut& lut = inputLUT_.at(lutId);
483  for (int i = 0; i < df.size(); i++) {
484  ics[i] = (lut.at(df[i].adc()) & QIE8_LUT_BITMASK);
485  }
486 }
487 
489  int lutId = getLUTId(df.id());
490  const Lut& lut = inputLUT_.at(lutId);
491  for (int i = 0; i < df.size(); i++) {
492  ics[i] = (lut.at(df[i].adc()) & QIE8_LUT_BITMASK);
493  }
494 }
495 
497  int lutId = getLUTId(HcalDetId(df.id()));
498  const Lut& lut = inputLUT_.at(lutId);
499  for (int i = 0; i < df.samples(); i++) {
500  ics[i] = (lut.at(df[i].adc()) & QIE10_LUT_BITMASK);
501  }
502 }
503 
505  int lutId = getLUTId(HcalDetId(df.id()));
506  const Lut& lut = inputLUT_.at(lutId);
507  for (int i = 0; i < df.samples(); i++) {
508  ics[i] = (lut.at(df[i].adc()) & QIE11_LUT_BITMASK);
509  }
510 }
511 
513  int lutId = getLUTId(id);
514  return ((inputLUT_.at(lutId)).at(sample.adc()) & QIE8_LUT_BITMASK);
515 }
516 
518  int lutId = getLUTId(id);
519  return ped_.at(lutId);
520 }
521 
523  int lutId = getLUTId(id);
524  return gain_.at(lutId);
525 }
526 
527 std::vector<unsigned short> HcaluLUTTPGCoder::getLinearizationLUT(HcalDetId id) const {
528  int lutId = getLUTId(id);
529  return inputLUT_.at(lutId);
530 }
531 
532 void HcaluLUTTPGCoder::lookupMSB(const HBHEDataFrame& df, std::vector<bool>& msb) const {
533  msb.resize(df.size());
534  for (int i = 0; i < df.size(); ++i)
535  msb[i] = getMSB(df.id(), df.sample(i).adc());
536 }
537 
538 bool HcaluLUTTPGCoder::getMSB(const HcalDetId& id, int adc) const {
539  int lutId = getLUTId(id);
540  const Lut& lut = inputLUT_.at(lutId);
541  return (lut.at(adc) & QIE8_LUT_MSB);
542 }
543 
544 void HcaluLUTTPGCoder::lookupMSB(const QIE10DataFrame& df, std::vector<std::bitset<2>>& msb) const {
545  msb.resize(df.samples());
546  int lutId = getLUTId(HcalDetId(df.id()));
547  const Lut& lut = inputLUT_.at(lutId);
548  for (int i = 0; i < df.samples(); ++i) {
549  msb[i][0] = lut.at(df[i].adc()) & QIE10_LUT_MSB0;
550  msb[i][1] = lut.at(df[i].adc()) & QIE10_LUT_MSB1;
551  }
552 }
553 
554 void HcaluLUTTPGCoder::lookupMSB(const QIE11DataFrame& df, std::vector<std::bitset<2>>& msb) const {
555  int lutId = getLUTId(HcalDetId(df.id()));
556  const Lut& lut = inputLUT_.at(lutId);
557  for (int i = 0; i < df.samples(); ++i) {
558  msb[i][0] = lut.at(df[i].adc()) & QIE11_LUT_MSB0;
559  msb[i][1] = lut.at(df[i].adc()) & QIE11_LUT_MSB1;
560  }
561 }
std::vector< uint32_t > FG_HF_thresholds_
int samples() const
total number of samples in the digi
size
Write out results.
uint32_t getFlag() const
int maxDepth(void) const
int firstHFRing() const
Definition: HcalTopology.h:96
float getPedestal() const
static const int nFi_
void adc2Linear(const HBHEDataFrame &df, IntegerCaloSamples &ics) const override
const HcalDDDRecConstants * dddConstants() const
Definition: HcalTopology.h:164
double cosh_ieta_28_HE_high_depths_
int getValue() const
Definition: HcalQIEType.h:19
Definition: LutXml.h:27
HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:138
const HcalTPChannelParameter * getHcalTPChannelParameter(const HcalGenericDetId &fId) const
bool valid(const DetId &id) const override
double linearLSB_QIE11Overlap_
static uint32_t kHcalEtaMask2
Definition: HcalDetId.h:19
const HcalRecoParam * getHcalRecoParam(const HcalGenericDetId &fId) const
static const int QIE11_LUT_MSB1
static const int QIE11_LUT_BITMASK
const HcalChannelStatus * getHcalChannelStatus(const HcalGenericDetId &fId) const
Definition: HcalQIENum.h:4
float getLUTGain(HcalDetId id) const override
int size() const
total number of samples in the digi
Definition: HBHEDataFrame.h:27
int create_lut_map(void)
Definition: LutXml.cc:343
const HcalTopology * topo_
int firstHBRing() const
Definition: HcalTopology.h:91
edm::DataFrame::id_type id() const
int lastHBRing() const
Definition: HcalTopology.h:92
static const int QIE11_LUT_MSB0
static const float lsb_
void towerEtaBounds(int ieta, int version, double &eta1, double &eta2) const
where this tower begins and ends in eta
const Item * getValues(DetId fId, bool throwOnFail=true) const
double cosh_ieta_28_HE_low_depths_
void init(const HcalTopology *top, const HcalTimeSlew *delay)
static const int QIE8_LUT_BITMASK
uint8_t getLutGranularity() const
HcalTopologyMode::TriggerMode triggerMode() const
Definition: HcalTopology.h:35
unsigned short LutElement
void update(const HcalDbService &conditions)
void updateXML(const char *filename)
static const int QIE10_LUT_MSB0
static const size_t UPGRADE_LUT_SIZE
int depth() const
get the tower depth
Definition: HcalDetId.h:164
std::vector< Lut > inputLUT_
void lookupMSB(const HBHEDataFrame &df, std::vector< bool > &msb) const
int getLUTId(HcalSubdetector id, int ieta, int iphi, int depth) const
float getRespGain() const
double MaximumFractionalError
HcalDetId const & id() const
Definition: HFDataFrame.h:23
int terminate(void)
std::vector< DetId > getAllChannels() const
constexpr float correctionPhaseNS() const
Definition: HcalRecoParam.h:30
float getNominalGain() const
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
static uint32_t kHcalDepthMask2
Definition: HcalDetId.h:25
const HcalTimeSlew * delay_
float getRCalib() const
const HcalL1TriggerObject * getHcalL1TriggerObject(const HcalGenericDetId &fId) const
constexpr double effpedestal(int fCapId) const
get effective pedestal for capid=0..3
int ieta() const
get the cell ieta
Definition: HcalDetId.h:155
int lastHFRing() const
Definition: HcalTopology.h:97
edm::DataFrame::id_type id() const
HcalSubdetector
Definition: HcalAssistant.h:31
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const HcalLutMetadata * getHcalLutMetadata() const
T min(T a, T b)
Definition: MathUtil.h:58
constexpr int adc(sample_type sample)
get the ADC sample (12 bits)
static const int QIE8_LUT_MSB
constexpr int adc() const
get the ADC sample
Definition: HcalQIESample.h:43
JetCorrectorParameters corr
Definition: classes.h:5
static const int QIE10_LUT_BITMASK
int ietaAbs() const
get the absolute value of the cell ieta
Definition: HcalDetId.h:148
const HcalQIEType * getHcalQIEType(const HcalGenericDetId &fId) const
std::vector< double > cosh_ieta_
double const adc2fC[256]
Definition: Constants.h:240
int iphi() const
get the cell iphi
Definition: HcalDetId.h:157
std::vector< float > ped_
double cosh_ieta(int ieta, int depth, HcalSubdetector subdet)
std::vector< float > gain_
int size() const
total number of samples in the digi
Definition: HFDataFrame.h:27
HcalQIESample const & sample(int i) const
access a sample
Definition: HBHEDataFrame.h:40
static uint32_t kHcalPhiMask2
Definition: HcalDetId.h:15
~HcaluLUTTPGCoder() override
int firstHERing() const
Definition: HcalTopology.h:93
constexpr double LUTrespcorrgain(int fCapId) const
get LUT corrected and response corrected gain for capid=0..3
const HcalSiPMCharacteristics * getHcalSiPMCharacteristics() const
const HcalQIECoder * getHcalCoder(const HcalGenericDetId &fId) const
std::vector< float > getNonLinearities(int type) const
get nonlinearity constants
const HcalQIEShape * getHcalShape(const HcalGenericDetId &fId) const
std::vector< unsigned int > * getLutFast(uint32_t det_id)
Definition: LutXml.cc:70
static const size_t INPUT_LUT_SIZE
std::pair< double, double > etaRange(HcalSubdetector subdet, int ieta) const
bool getMSB(const HcalDetId &id, int adc) const
float getLUTPedestal(HcalDetId id) const override
uint32_t getFGBitInfo() const
get FG bit information
const HcalDetId & id() const
Definition: HBHEDataFrame.h:23
bool isPlan1(const HcalDetId &id) const
void make_cosh_ieta_map(void)
uint32_t getValue() const
const HcalCalibrations & getHcalCalibrations(const HcalGenericDetId &fId) const
int samples() const
total number of samples in the digi
std::unique_ptr< HcalPulseContainmentManager > pulseCorr_
const HcalSiPMParameter * getHcalSiPMParameter(const HcalGenericDetId &fId) const
#define constexpr
static const int QIE10_LUT_MSB1
static XMLProcessor * getInstance()
Definition: XMLProcessor.h:134
int lastHERing() const
Definition: HcalTopology.h:94
Definition: Lut.h:31
void compress(const IntegerCaloSamples &ics, const std::vector< bool > &featureBits, HcalTriggerPrimitiveDigi &tp) const override
std::vector< unsigned short > getLinearizationLUT(HcalDetId id) const override
Get the full linearization LUT (128 elements). Default implementation just uses adc2Linear to get all...
float charge(const HcalQIEShape &fShape, unsigned fAdc, unsigned fCapId) const
ADC [0..127] + capid [0..3] -> fC conversion.
Definition: HcalQIECoder.cc:20