File indexing completed on 2022-04-26 22:29:45
0001 #include "SimCalorimetry/HcalTrigPrimAlgos/interface/HcalTriggerPrimitiveAlgo.h"
0002
0003 #include "CalibFormats/CaloObjects/interface/IntegerCaloSamples.h"
0004 #include "CondFormats/HcalObjects/interface/HcalTPParameters.h"
0005 #include "CondFormats/HcalObjects/interface/HcalTPChannelParameters.h"
0006
0007 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
0008 #include "DataFormats/FEDRawData/interface/FEDNumbering.h"
0009 #include "DataFormats/HcalDetId/interface/HcalElectronicsId.h"
0010
0011 #include "EventFilter/HcalRawToDigi/interface/HcalDCCHeader.h"
0012 #include "EventFilter/HcalRawToDigi/interface/HcalHTRData.h"
0013
0014 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0015 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0016
0017 #include "Geometry/HcalTowerAlgo/interface/HcalTrigTowerGeometry.h"
0018
0019 #include <iostream>
0020
0021 using namespace std;
0022
0023 HcalTriggerPrimitiveAlgo::HcalTriggerPrimitiveAlgo(bool pf,
0024 const std::vector<double>& w,
0025 int latency,
0026 uint32_t FG_threshold,
0027 const std::vector<uint32_t>& FG_HF_thresholds,
0028 uint32_t ZS_threshold,
0029 int numberOfSamples,
0030 int numberOfPresamples,
0031 int numberOfFilterPresamplesHBQIE11,
0032 int numberOfFilterPresamplesHEQIE11,
0033 int numberOfSamplesHF,
0034 int numberOfPresamplesHF,
0035 bool useTDCInMinBiasBits,
0036 uint32_t minSignalThreshold,
0037 uint32_t PMT_NoiseThreshold)
0038 : incoder_(nullptr),
0039 outcoder_(nullptr),
0040 theThreshold(0),
0041 peakfind_(pf),
0042 weights_(w),
0043 latency_(latency),
0044 FG_threshold_(FG_threshold),
0045 FG_HF_thresholds_(FG_HF_thresholds),
0046 ZS_threshold_(ZS_threshold),
0047 numberOfSamples_(numberOfSamples),
0048 numberOfPresamples_(numberOfPresamples),
0049 numberOfFilterPresamplesHBQIE11_(numberOfFilterPresamplesHBQIE11),
0050 numberOfFilterPresamplesHEQIE11_(numberOfFilterPresamplesHEQIE11),
0051 numberOfSamplesHF_(numberOfSamplesHF),
0052 numberOfPresamplesHF_(numberOfPresamplesHF),
0053 useTDCInMinBiasBits_(useTDCInMinBiasBits),
0054 minSignalThreshold_(minSignalThreshold),
0055 PMT_NoiseThreshold_(PMT_NoiseThreshold),
0056 NCTScaleShift(0),
0057 RCTScaleShift(0),
0058 peak_finder_algorithm_(2),
0059 override_parameters_() {
0060
0061 if (!peakfind_) {
0062 numberOfSamples_ = 1;
0063 numberOfPresamples_ = 0;
0064 numberOfSamplesHF_ = 1;
0065 numberOfPresamplesHF_ = 0;
0066 }
0067
0068 ZS_threshold_I_ = ZS_threshold_;
0069 }
0070
0071 HcalTriggerPrimitiveAlgo::~HcalTriggerPrimitiveAlgo() {}
0072
0073 void HcalTriggerPrimitiveAlgo::setUpgradeFlags(bool hb, bool he, bool hf) {
0074 upgrade_hb_ = hb;
0075 upgrade_he_ = he;
0076 upgrade_hf_ = hf;
0077 }
0078
0079 void HcalTriggerPrimitiveAlgo::setFixSaturationFlag(bool fix_saturation) { fix_saturation_ = fix_saturation; }
0080
0081 void HcalTriggerPrimitiveAlgo::overrideParameters(const edm::ParameterSet& ps) {
0082 override_parameters_ = ps;
0083
0084 if (override_parameters_.exists("ADCThresholdHF")) {
0085 override_adc_hf_ = true;
0086 override_adc_hf_value_ = override_parameters_.getParameter<uint32_t>("ADCThresholdHF");
0087 }
0088 if (override_parameters_.exists("TDCMaskHF")) {
0089 override_tdc_hf_ = true;
0090 override_tdc_hf_value_ = override_parameters_.getParameter<unsigned long long>("TDCMaskHF");
0091 }
0092 }
0093
0094 void HcalTriggerPrimitiveAlgo::addSignal(const HBHEDataFrame& frame) {
0095
0096
0097 if (frame.id().depth() == 5)
0098 return;
0099
0100 std::vector<HcalTrigTowerDetId> ids = theTrigTowerGeometry->towerIds(frame.id());
0101 assert(ids.size() == 1 || ids.size() == 2);
0102 IntegerCaloSamples samples1(ids[0], int(frame.size()));
0103
0104 samples1.setPresamples(frame.presamples());
0105 incoder_->adc2Linear(frame, samples1);
0106
0107 std::vector<bool> msb;
0108 incoder_->lookupMSB(frame, msb);
0109
0110 if (ids.size() == 2) {
0111
0112 IntegerCaloSamples samples2(ids[1], samples1.size());
0113 for (int i = 0; i < samples1.size(); ++i) {
0114 samples1[i] = uint32_t(samples1[i] * 0.5);
0115 samples2[i] = samples1[i];
0116 }
0117 samples2.setPresamples(frame.presamples());
0118 addSignal(samples2);
0119 addFG(ids[1], msb);
0120 }
0121 addSignal(samples1);
0122 addFG(ids[0], msb);
0123 }
0124
0125 void HcalTriggerPrimitiveAlgo::addSignal(const HFDataFrame& frame) {
0126 if (frame.id().depth() == 1 || frame.id().depth() == 2) {
0127 std::vector<HcalTrigTowerDetId> ids = theTrigTowerGeometry->towerIds(frame.id());
0128 std::vector<HcalTrigTowerDetId>::const_iterator it;
0129 for (it = ids.begin(); it != ids.end(); ++it) {
0130 HcalTrigTowerDetId trig_tower_id = *it;
0131 IntegerCaloSamples samples(trig_tower_id, frame.size());
0132 samples.setPresamples(frame.presamples());
0133 incoder_->adc2Linear(frame, samples);
0134
0135
0136
0137 IntegerCaloSamples zero_samples(trig_tower_id, frame.size());
0138 zero_samples.setPresamples(frame.presamples());
0139 addSignal(zero_samples);
0140
0141
0142 if (trig_tower_id.version() == 0) {
0143
0144 uint32_t fgid = (frame.id().maskDepth());
0145
0146 if (theTowerMapFGSum.find(trig_tower_id) == theTowerMapFGSum.end()) {
0147 SumFGContainer sumFG;
0148 theTowerMapFGSum.insert(std::pair<HcalTrigTowerDetId, SumFGContainer>(trig_tower_id, sumFG));
0149 }
0150
0151 SumFGContainer& sumFG = theTowerMapFGSum[trig_tower_id];
0152 SumFGContainer::iterator sumFGItr;
0153 for (sumFGItr = sumFG.begin(); sumFGItr != sumFG.end(); ++sumFGItr) {
0154 if (sumFGItr->id() == fgid) {
0155 break;
0156 }
0157 }
0158
0159 if (sumFGItr != sumFG.end()) {
0160 for (int i = 0; i < samples.size(); ++i) {
0161 (*sumFGItr)[i] += samples[i];
0162 }
0163 } else {
0164
0165 IntegerCaloSamples sumFGSamples(DetId(fgid), samples.size());
0166 sumFGSamples.setPresamples(samples.presamples());
0167 for (int i = 0; i < samples.size(); ++i) {
0168 sumFGSamples[i] = samples[i];
0169 }
0170 sumFG.push_back(sumFGSamples);
0171 }
0172
0173
0174 if (HF_Veto.find(fgid) == HF_Veto.end()) {
0175 vector<bool> vetoBits(samples.size(), false);
0176 HF_Veto[fgid] = vetoBits;
0177 }
0178 for (int i = 0; i < samples.size(); ++i) {
0179 if (samples[i] < minSignalThreshold_) {
0180 HF_Veto[fgid][i] = true;
0181 }
0182 }
0183 }
0184
0185 else if (trig_tower_id.version() == 1) {
0186 uint32_t fgid = (frame.id().maskDepth());
0187 HFDetails& details = theHFDetailMap[trig_tower_id][fgid];
0188
0189 if (frame.id().depth() == 1) {
0190 details.long_fiber = samples;
0191 details.LongDigi = frame;
0192 } else if (frame.id().depth() == 2) {
0193 details.short_fiber = samples;
0194 details.ShortDigi = frame;
0195 } else {
0196
0197 edm::LogWarning("HcalTPAlgo") << "Unable to figure out what to do with data frame for " << frame.id();
0198 return;
0199 }
0200 }
0201
0202 else {
0203 return;
0204 }
0205 }
0206 }
0207 }
0208
0209 void HcalTriggerPrimitiveAlgo::addSignal(const QIE10DataFrame& frame) {
0210 HcalDetId detId = frame.detid();
0211
0212 if (detId.subdet() != HcalForward)
0213 return;
0214
0215 auto ids = theTrigTowerGeometry->towerIds(frame.id());
0216 for (const auto& id : ids) {
0217 if (id.version() == 0) {
0218 edm::LogError("HcalTPAlgo") << "Encountered QIE10 data frame mapped to TP version 0:" << id;
0219 continue;
0220 }
0221
0222 int nsamples = frame.samples();
0223
0224 IntegerCaloSamples samples(id, nsamples);
0225 samples.setPresamples(frame.presamples());
0226 incoder_->adc2Linear(frame, samples);
0227
0228
0229
0230 IntegerCaloSamples zero_samples(id, nsamples);
0231 zero_samples.setPresamples(frame.presamples());
0232 addSignal(zero_samples);
0233
0234 auto fid = HcalDetId(frame.id());
0235 auto& details = theHFUpgradeDetailMap[id][fid.maskDepth()];
0236 auto& detail = details[fid.depth() - 1];
0237 detail.samples = samples;
0238 detail.digi = frame;
0239 detail.validity.resize(nsamples);
0240 detail.passTDC.resize(nsamples);
0241 incoder_->lookupMSB(frame, detail.fgbits);
0242 for (int idx = 0; idx < nsamples; ++idx) {
0243 detail.validity[idx] = validChannel(frame, idx);
0244 detail.passTDC[idx] = passTDC(frame, idx);
0245 }
0246 }
0247 }
0248
0249 void HcalTriggerPrimitiveAlgo::addSignal(const QIE11DataFrame& frame) {
0250 HcalDetId detId(frame.id());
0251
0252 if (detId.subdet() != HcalEndcap && detId.subdet() != HcalBarrel)
0253 return;
0254
0255 std::vector<HcalTrigTowerDetId> ids = theTrigTowerGeometry->towerIds(detId);
0256 assert(ids.size() == 1 || ids.size() == 2);
0257 IntegerCaloSamples samples1(ids[0], int(frame.samples()));
0258
0259 samples1.setPresamples(frame.presamples());
0260 incoder_->adc2Linear(frame, samples1);
0261
0262 std::vector<std::bitset<2>> msb(frame.samples(), 0);
0263 incoder_->lookupMSB(frame, msb);
0264
0265 if (ids.size() == 2) {
0266
0267 IntegerCaloSamples samples2(ids[1], samples1.size());
0268 for (int i = 0; i < samples1.size(); ++i) {
0269 samples1[i] = uint32_t(samples1[i]);
0270 samples2[i] = samples1[i];
0271 }
0272 samples2.setPresamples(frame.presamples());
0273 addSignal(samples2);
0274 addUpgradeFG(ids[1], detId.depth(), msb);
0275 addUpgradeTDCFG(ids[1], frame);
0276 }
0277 addSignal(samples1);
0278 addUpgradeFG(ids[0], detId.depth(), msb);
0279 addUpgradeTDCFG(ids[0], frame);
0280 }
0281
0282 void HcalTriggerPrimitiveAlgo::addSignal(const IntegerCaloSamples& samples) {
0283 HcalTrigTowerDetId id(samples.id());
0284 SumMap::iterator itr = theSumMap.find(id);
0285
0286 if (itr == theSumMap.end()) {
0287 theSumMap.insert(std::make_pair(id, samples));
0288 } else {
0289
0290 for (int i = 0; i < samples.size(); ++i) {
0291 (itr->second)[i] += samples[i];
0292 }
0293 }
0294
0295
0296 if (fix_saturation_) {
0297 SatMap::iterator itr_sat = theSatMap.find(id);
0298
0299 assert((itr == theSumMap.end()) == (itr_sat == theSatMap.end()));
0300
0301 if (itr_sat == theSatMap.end()) {
0302 vector<bool> check_sat;
0303 for (int i = 0; i < samples.size(); ++i) {
0304 if (!(samples[i] < QIE11_LINEARIZATION_ET)) {
0305 check_sat.push_back(true);
0306 } else
0307 check_sat.push_back(false);
0308 }
0309 theSatMap.insert(std::make_pair(id, check_sat));
0310 } else {
0311 for (int i = 0; i < samples.size(); ++i) {
0312 if (!(samples[i] < QIE11_LINEARIZATION_ET))
0313 (itr_sat->second)[i] = true;
0314 }
0315 }
0316 }
0317 }
0318
0319 void HcalTriggerPrimitiveAlgo::analyze(IntegerCaloSamples& samples, HcalTriggerPrimitiveDigi& result) {
0320 int shrink = weights_.size() - 1;
0321 std::vector<bool>& msb = fgMap_[samples.id()];
0322 IntegerCaloSamples sum(samples.id(), samples.size());
0323
0324
0325 for (int ibin = 0; ibin < int(samples.size()) - shrink; ++ibin) {
0326 int algosumvalue = 0;
0327 for (unsigned int i = 0; i < weights_.size(); i++) {
0328
0329 algosumvalue += int(samples[ibin + i] * weights_[i]);
0330 }
0331 if (algosumvalue < 0)
0332 sum[ibin] = 0;
0333
0334
0335 else
0336 sum[ibin] = algosumvalue;
0337 }
0338
0339
0340 int dgPresamples = samples.presamples();
0341 int tpPresamples = numberOfPresamples_;
0342 int shift = dgPresamples - tpPresamples;
0343 int dgSamples = samples.size();
0344 int tpSamples = numberOfSamples_;
0345 if (peakfind_) {
0346 if ((shift < shrink) || (shift + tpSamples + shrink > dgSamples - (peak_finder_algorithm_ - 1))) {
0347 edm::LogInfo("HcalTriggerPrimitiveAlgo::analyze")
0348 << "TP presample or size from the configuration file is out of the accessible range. Using digi values from "
0349 "data instead...";
0350 shift = shrink;
0351 tpPresamples = dgPresamples - shrink;
0352 tpSamples = dgSamples - (peak_finder_algorithm_ - 1) - shrink - shift;
0353 }
0354 }
0355
0356 std::vector<int> finegrain(tpSamples, false);
0357
0358 IntegerCaloSamples output(samples.id(), tpSamples);
0359 output.setPresamples(tpPresamples);
0360
0361 for (int ibin = 0; ibin < tpSamples; ++ibin) {
0362
0363
0364 int idx = ibin + shift;
0365
0366
0367 if (peakfind_) {
0368 bool isPeak = false;
0369 switch (peak_finder_algorithm_) {
0370 case 1:
0371 isPeak = (samples[idx] > samples[idx - 1] && samples[idx] >= samples[idx + 1] && samples[idx] > theThreshold);
0372 break;
0373 case 2:
0374 isPeak = (sum[idx] > sum[idx - 1] && sum[idx] >= sum[idx + 1] && sum[idx] > theThreshold);
0375 break;
0376 default:
0377 break;
0378 }
0379
0380 if (isPeak) {
0381 output[ibin] = std::min<unsigned int>(sum[idx], QIE8_LINEARIZATION_ET);
0382 finegrain[ibin] = msb[idx];
0383 }
0384
0385 else
0386 output[ibin] = 0;
0387 } else {
0388 output[ibin] = std::min<unsigned int>(sum[idx], QIE8_LINEARIZATION_ET);
0389 finegrain[ibin] = msb[idx];
0390 }
0391
0392
0393 if (peak_finder_algorithm_ == 1) {
0394 if (samples[idx] >= QIE8_LINEARIZATION_ET)
0395 output[ibin] = QIE8_LINEARIZATION_ET;
0396 }
0397 }
0398 outcoder_->compress(output, finegrain, result);
0399 }
0400
0401 void HcalTriggerPrimitiveAlgo::analyzeQIE11(IntegerCaloSamples& samples,
0402 vector<bool> sample_saturation,
0403 HcalTriggerPrimitiveDigi& result,
0404 const HcalFinegrainBit& fg_algo) {
0405 HcalDetId detId(samples.id());
0406
0407
0408 int theIeta = detId.ietaAbs();
0409
0410 unsigned int dgSamples = samples.size();
0411 unsigned int dgPresamples = samples.presamples();
0412
0413 unsigned int tpSamples = numberOfSamples_;
0414 unsigned int tpPresamples = numberOfPresamples_;
0415
0416 unsigned int filterSamples = weightsQIE11_[theIeta].size();
0417 unsigned int filterPresamples = theIeta > theTrigTowerGeometry->topology().lastHBRing()
0418 ? numberOfFilterPresamplesHEQIE11_
0419 : numberOfFilterPresamplesHBQIE11_;
0420
0421 unsigned int shift = dgPresamples - tpPresamples;
0422
0423
0424 unsigned int shrink = filterSamples - 1;
0425
0426 auto& msb = fgUpgradeMap_[samples.id()];
0427 auto& timingTDC = fgUpgradeTDCMap_[samples.id()];
0428 IntegerCaloSamples sum(samples.id(), samples.size());
0429
0430 std::vector<HcalTrigTowerDetId> ids = theTrigTowerGeometry->towerIds(detId);
0431
0432
0433 bool force_saturation[samples.size()];
0434 for (int i = 0; i < samples.size(); i++) {
0435 force_saturation[i] = false;
0436 }
0437
0438
0439 for (unsigned int ibin = 0; ibin < dgSamples - shrink; ++ibin) {
0440 int algosumvalue = 0;
0441 bool check_sat = false;
0442 for (unsigned int i = 0; i < filterSamples; i++) {
0443
0444
0445
0446 unsigned int sample = samples[ibin + i];
0447
0448 if (fix_saturation_ && (sample_saturation.size() > ibin + i))
0449 check_sat = (check_sat | sample_saturation[ibin + i] | (sample > QIE11_MAX_LINEARIZATION_ET));
0450
0451 if (sample > QIE11_MAX_LINEARIZATION_ET)
0452 sample = QIE11_MAX_LINEARIZATION_ET;
0453
0454
0455 double segmentationFactor = 1.0;
0456 if (ids.size() == 2) {
0457 segmentationFactor = 0.5;
0458 }
0459
0460
0461 double theWeight = weightsQIE11_[theIeta][i];
0462
0463 algosumvalue += int(sample * segmentationFactor * theWeight);
0464 }
0465 if (algosumvalue < 0)
0466 sum[ibin] = 0;
0467
0468
0469 else
0470 sum[ibin] = algosumvalue;
0471
0472 if (check_sat)
0473 force_saturation[ibin] = true;
0474 }
0475
0476 std::vector<int> finegrain(tpSamples, false);
0477
0478 IntegerCaloSamples output(samples.id(), tpSamples);
0479 output.setPresamples(tpPresamples);
0480
0481 for (unsigned int ibin = 0; ibin < tpSamples; ++ibin) {
0482
0483
0484 int idx = ibin + shift - filterPresamples;
0485
0486
0487 if (idx <= 0) {
0488 output[ibin] = 0;
0489 continue;
0490 }
0491
0492 bool isPeak = (sum[idx] > sum[idx - 1] && sum[idx] >= sum[idx + 1] && sum[idx] > theThreshold);
0493
0494 if (isPeak) {
0495 output[ibin] = std::min<unsigned int>(sum[idx], QIE11_MAX_LINEARIZATION_ET);
0496
0497 if (fix_saturation_ && force_saturation[idx] && ids.size() == 2)
0498 output[ibin] = QIE11_MAX_LINEARIZATION_ET * 0.5;
0499 else if (fix_saturation_ && force_saturation[idx])
0500 output[ibin] = QIE11_MAX_LINEARIZATION_ET;
0501
0502 } else {
0503
0504 output[ibin] = 0;
0505 }
0506
0507
0508 finegrain[ibin] = fg_algo.compute(timingTDC[idx + filterPresamples], ids[0]).to_ulong() |
0509 fg_algo.compute(msb[idx + filterPresamples]).to_ulong() << 4;
0510 if (ibin == tpPresamples && (idx + filterPresamples) != dgPresamples)
0511 edm::LogError("HcalTriggerPritimveAlgo")
0512 << "TP SOI (tpPresamples = " << tpPresamples
0513 << ") is not aligned with digi SOI (dgPresamples = " << dgPresamples << ")";
0514 }
0515 outcoder_->compress(output, finegrain, result);
0516 }
0517
0518 void HcalTriggerPrimitiveAlgo::analyzeHF(IntegerCaloSamples& samples,
0519 HcalTriggerPrimitiveDigi& result,
0520 const int hf_lumi_shift) {
0521 HcalTrigTowerDetId detId(samples.id());
0522
0523
0524 int dgPresamples = samples.presamples();
0525 int tpPresamples = numberOfPresamplesHF_;
0526 int shift = dgPresamples - tpPresamples;
0527 int dgSamples = samples.size();
0528 int tpSamples = numberOfSamplesHF_;
0529 if (shift < 0 || shift + tpSamples > dgSamples) {
0530 edm::LogInfo("HcalTriggerPrimitiveAlgo::analyzeHF")
0531 << "TP presample or size from the configuration file is out of the accessible range. Using digi values from "
0532 "data instead...";
0533 tpPresamples = dgPresamples;
0534 shift = 0;
0535 tpSamples = dgSamples;
0536 }
0537
0538 std::vector<int> finegrain(tpSamples, false);
0539
0540 TowerMapFGSum::const_iterator tower2fg = theTowerMapFGSum.find(detId);
0541 assert(tower2fg != theTowerMapFGSum.end());
0542
0543 const SumFGContainer& sumFG = tower2fg->second;
0544
0545
0546 for (SumFGContainer::const_iterator sumFGItr = sumFG.begin(); sumFGItr != sumFG.end(); ++sumFGItr) {
0547 const std::vector<bool>& veto = HF_Veto[sumFGItr->id().rawId()];
0548 for (int ibin = 0; ibin < tpSamples; ++ibin) {
0549 int idx = ibin + shift;
0550
0551 bool vetoed = idx < int(veto.size()) && veto[idx];
0552 if (!(vetoed && (*sumFGItr)[idx] > PMT_NoiseThreshold_)) {
0553 samples[idx] += (*sumFGItr)[idx];
0554 finegrain[ibin] = (finegrain[ibin] || (*sumFGItr)[idx] >= FG_threshold_);
0555 }
0556 }
0557 }
0558
0559 IntegerCaloSamples output(samples.id(), tpSamples);
0560 output.setPresamples(tpPresamples);
0561
0562 for (int ibin = 0; ibin < tpSamples; ++ibin) {
0563 int idx = ibin + shift;
0564 output[ibin] = samples[idx] >> hf_lumi_shift;
0565 static const int MAX_OUTPUT = QIE8_LINEARIZATION_ET;
0566 if (output[ibin] > MAX_OUTPUT)
0567 output[ibin] = MAX_OUTPUT;
0568 }
0569 outcoder_->compress(output, finegrain, result);
0570 }
0571
0572 void HcalTriggerPrimitiveAlgo::analyzeHF2016(const IntegerCaloSamples& samples,
0573 HcalTriggerPrimitiveDigi& result,
0574 const int hf_lumi_shift,
0575 const HcalFeatureBit* embit) {
0576
0577 const int SHIFT = samples.presamples() - numberOfPresamplesHF_;
0578 assert(SHIFT >= 0);
0579 assert((SHIFT + numberOfSamplesHF_) <= samples.size());
0580
0581
0582 const HcalTrigTowerDetId detId(samples.id());
0583 HFDetailMap::const_iterator it = theHFDetailMap.find(detId);
0584
0585 if (it == theHFDetailMap.end()) {
0586 return;
0587 }
0588
0589 std::vector<std::bitset<2>> finegrain(numberOfSamplesHF_, false);
0590
0591
0592 IntegerCaloSamples output(samples.id(), numberOfSamplesHF_);
0593 output.setPresamples(numberOfPresamplesHF_);
0594
0595 for (const auto& item : it->second) {
0596 auto& details = item.second;
0597 for (int ibin = 0; ibin < numberOfSamplesHF_; ++ibin) {
0598 const int IDX = ibin + SHIFT;
0599 int long_fiber_val = 0;
0600 if (IDX < details.long_fiber.size()) {
0601 long_fiber_val = details.long_fiber[IDX];
0602 }
0603 int short_fiber_val = 0;
0604 if (IDX < details.short_fiber.size()) {
0605 short_fiber_val = details.short_fiber[IDX];
0606 }
0607 output[ibin] += (long_fiber_val + short_fiber_val);
0608
0609 uint32_t ADCLong = details.LongDigi[ibin].adc();
0610 uint32_t ADCShort = details.ShortDigi[ibin].adc();
0611
0612 if (details.LongDigi.id().ietaAbs() >= FIRST_FINEGRAIN_TOWER) {
0613 finegrain[ibin][1] = (ADCLong > FG_HF_thresholds_[0] || ADCShort > FG_HF_thresholds_[0]);
0614
0615 if (embit != nullptr)
0616 finegrain[ibin][0] = embit->fineGrainbit(details.ShortDigi, details.LongDigi, ibin);
0617 }
0618 }
0619 }
0620
0621 for (int bin = 0; bin < numberOfSamplesHF_; ++bin) {
0622 static const unsigned int MAX_OUTPUT = QIE8_LINEARIZATION_ET;
0623 output[bin] = min({MAX_OUTPUT, output[bin] >> hf_lumi_shift});
0624 }
0625
0626 std::vector<int> finegrain_converted;
0627 finegrain_converted.reserve(finegrain.size());
0628 for (const auto& fg : finegrain)
0629 finegrain_converted.push_back(fg.to_ulong());
0630 outcoder_->compress(output, finegrain_converted, result);
0631 }
0632
0633 bool HcalTriggerPrimitiveAlgo::passTDC(const QIE10DataFrame& digi, int ts) const {
0634 auto parameters = conditions_->getHcalTPParameters();
0635 auto adc_threshold = parameters->getADCThresholdHF();
0636 auto tdc_mask = parameters->getTDCMaskHF();
0637
0638 if (override_adc_hf_)
0639 adc_threshold = override_adc_hf_value_;
0640 if (override_tdc_hf_)
0641 tdc_mask = override_tdc_hf_value_;
0642
0643 if (digi[ts].adc() < adc_threshold)
0644 return true;
0645
0646 return (1ul << digi[ts].le_tdc()) & tdc_mask;
0647 }
0648
0649 bool HcalTriggerPrimitiveAlgo::validChannel(const QIE10DataFrame& digi, int ts) const {
0650
0651 if (digi.linkError() || ts >= digi.samples() || !digi[ts].ok())
0652 return false;
0653
0654 auto mask = conditions_->getHcalTPChannelParameter(HcalDetId(digi.id()))->getMask();
0655 if (mask)
0656 return false;
0657
0658 return true;
0659 }
0660
0661 void HcalTriggerPrimitiveAlgo::analyzeHFQIE10(const IntegerCaloSamples& samples,
0662 HcalTriggerPrimitiveDigi& result,
0663 const int hf_lumi_shift,
0664 const HcalFeatureBit* embit) {
0665
0666 const int shift = samples.presamples() - numberOfPresamplesHF_;
0667 assert(shift >= 0);
0668 assert((shift + numberOfSamplesHF_) <= samples.size());
0669 assert(hf_lumi_shift >= 2);
0670
0671
0672 const HcalTrigTowerDetId detId(samples.id());
0673 auto it = theHFUpgradeDetailMap.find(detId);
0674
0675 if (it == theHFUpgradeDetailMap.end()) {
0676 return;
0677 }
0678
0679 std::vector<std::bitset<2>> finegrain(numberOfSamplesHF_, false);
0680
0681
0682 IntegerCaloSamples output(samples.id(), numberOfSamplesHF_);
0683 output.setPresamples(numberOfPresamplesHF_);
0684
0685 for (const auto& item : it->second) {
0686 auto& details = item.second;
0687 for (int ibin = 0; ibin < numberOfSamplesHF_; ++ibin) {
0688 const int idx = ibin + shift;
0689
0690 int long_fiber_val = 0;
0691 int long_fiber_count = 0;
0692 int short_fiber_val = 0;
0693 int short_fiber_count = 0;
0694
0695 bool saturated = false;
0696
0697 for (auto i : {0, 2}) {
0698 if (idx < details[i].samples.size() and details[i].validity[idx] and details[i].passTDC[idx]) {
0699 long_fiber_val += details[i].samples[idx];
0700 saturated = saturated || (details[i].samples[idx] == QIE10_LINEARIZATION_ET);
0701 ++long_fiber_count;
0702 }
0703 }
0704 for (auto i : {1, 3}) {
0705 if (idx < details[i].samples.size() and details[i].validity[idx] and details[i].passTDC[idx]) {
0706 short_fiber_val += details[i].samples[idx];
0707 saturated = saturated || (details[i].samples[idx] == QIE10_LINEARIZATION_ET);
0708 ++short_fiber_count;
0709 }
0710 }
0711
0712 if (saturated) {
0713 output[ibin] = QIE10_MAX_LINEARIZATION_ET;
0714 } else {
0715
0716
0717
0718
0719 if (long_fiber_count == 2)
0720 long_fiber_val >>= 1;
0721 if (short_fiber_count == 2)
0722 short_fiber_val >>= 1;
0723
0724 auto sum = long_fiber_val + short_fiber_val;
0725
0726
0727
0728 if (long_fiber_count > 0 and short_fiber_count > 0)
0729 sum >>= 1;
0730
0731 output[ibin] += sum;
0732 }
0733
0734 for (const auto& detail : details) {
0735 if (idx < int(detail.digi.size()) and detail.validity[idx] and
0736 HcalDetId(detail.digi.id()).ietaAbs() >= FIRST_FINEGRAIN_TOWER) {
0737 if (useTDCInMinBiasBits_ && !detail.passTDC[idx])
0738 continue;
0739 finegrain[ibin][1] = finegrain[ibin][1] or detail.fgbits[idx][0];
0740
0741
0742
0743 finegrain[ibin][0] = finegrain[ibin][0] or detail.fgbits[idx][1];
0744 }
0745 }
0746
0747 if (embit != nullptr and FG_HF_thresholds_.at(1) != 255) {
0748 finegrain[ibin][0] = embit->fineGrainbit(details[1].digi,
0749 details[3].digi,
0750 details[0].digi,
0751 details[2].digi,
0752 details[1].validity[idx],
0753 details[3].validity[idx],
0754 details[0].validity[idx],
0755 details[2].validity[idx],
0756 idx);
0757 }
0758 }
0759 }
0760
0761 for (int bin = 0; bin < numberOfSamplesHF_; ++bin) {
0762 output[bin] = min({(unsigned int)QIE10_MAX_LINEARIZATION_ET, output[bin] >> (hf_lumi_shift - 2)});
0763 }
0764 std::vector<int> finegrain_converted;
0765 finegrain_converted.reserve(finegrain.size());
0766 for (const auto& fg : finegrain)
0767 finegrain_converted.push_back(fg.to_ulong());
0768 outcoder_->compress(output, finegrain_converted, result);
0769 }
0770
0771 void HcalTriggerPrimitiveAlgo::runZS(HcalTrigPrimDigiCollection& result) {
0772 for (HcalTrigPrimDigiCollection::iterator tp = result.begin(); tp != result.end(); ++tp) {
0773 bool ZS = true;
0774 for (int i = 0; i < tp->size(); ++i) {
0775 if (tp->sample(i).compressedEt() > ZS_threshold_I_) {
0776 ZS = false;
0777 break;
0778 }
0779 }
0780 if (ZS)
0781 tp->setZSInfo(false, true);
0782 else
0783 tp->setZSInfo(true, false);
0784 }
0785 }
0786
0787 void HcalTriggerPrimitiveAlgo::runFEFormatError(const FEDRawDataCollection* rawraw,
0788 const HcalElectronicsMap* emap,
0789 HcalTrigPrimDigiCollection& result) {
0790 std::set<uint32_t> FrontEndErrors;
0791
0792 for (int i = FEDNumbering::MINHCALFEDID; i <= FEDNumbering::MAXHCALFEDID; ++i) {
0793 const FEDRawData& raw = rawraw->FEDData(i);
0794 if (raw.size() < 12)
0795 continue;
0796 const HcalDCCHeader* dccHeader = (const HcalDCCHeader*)(raw.data());
0797 if (!dccHeader)
0798 continue;
0799 HcalHTRData htr;
0800 for (int spigot = 0; spigot < HcalDCCHeader::SPIGOT_COUNT; spigot++) {
0801 if (!dccHeader->getSpigotPresent(spigot))
0802 continue;
0803 dccHeader->getSpigotData(spigot, htr, raw.size());
0804 int dccid = dccHeader->getSourceId();
0805 int errWord = htr.getErrorsWord() & 0x1FFFF;
0806 bool HTRError = (!htr.check() || htr.isHistogramEvent() || (errWord & 0x800) != 0);
0807
0808 if (HTRError) {
0809 bool valid = false;
0810 for (int fchan = 0; fchan < 3 && !valid; fchan++) {
0811 for (int fib = 0; fib < 9 && !valid; fib++) {
0812 HcalElectronicsId eid(fchan, fib, spigot, dccid - FEDNumbering::MINHCALFEDID);
0813 eid.setHTR(htr.readoutVMECrateId(), htr.htrSlot(), htr.htrTopBottom());
0814 DetId detId = emap->lookup(eid);
0815 if (detId.null())
0816 continue;
0817 HcalSubdetector subdet = (HcalSubdetector(detId.subdetId()));
0818 if (detId.det() != 4 || (subdet != HcalBarrel && subdet != HcalEndcap && subdet != HcalForward))
0819 continue;
0820 std::vector<HcalTrigTowerDetId> ids = theTrigTowerGeometry->towerIds(detId);
0821 for (std::vector<HcalTrigTowerDetId>::const_iterator triggerId = ids.begin(); triggerId != ids.end();
0822 ++triggerId) {
0823 FrontEndErrors.insert(triggerId->rawId());
0824 }
0825
0826 }
0827 }
0828 }
0829 }
0830 }
0831
0832
0833
0834 HcalTriggerPrimitiveSample zeroSample(0);
0835 for (HcalTrigPrimDigiCollection::iterator tp = result.begin(); tp != result.end(); ++tp) {
0836 if (FrontEndErrors.find(tp->id().rawId()) != FrontEndErrors.end()) {
0837 for (int i = 0; i < tp->size(); ++i)
0838 tp->setSample(i, zeroSample);
0839 }
0840 }
0841 }
0842
0843 void HcalTriggerPrimitiveAlgo::addFG(const HcalTrigTowerDetId& id, std::vector<bool>& msb) {
0844 FGbitMap::iterator itr = fgMap_.find(id);
0845 if (itr != fgMap_.end()) {
0846 std::vector<bool>& _msb = itr->second;
0847 for (size_t i = 0; i < msb.size(); ++i)
0848 _msb[i] = _msb[i] || msb[i];
0849 } else
0850 fgMap_[id] = msb;
0851 }
0852
0853 bool HcalTriggerPrimitiveAlgo::validUpgradeFG(const HcalTrigTowerDetId& id, int depth) const {
0854 if (depth > LAST_FINEGRAIN_DEPTH)
0855 return false;
0856 if (id.ietaAbs() > LAST_FINEGRAIN_TOWER)
0857 return false;
0858 if (id.ietaAbs() == HBHE_OVERLAP_TOWER and not upgrade_hb_)
0859 return false;
0860 return true;
0861 }
0862
0863 bool HcalTriggerPrimitiveAlgo::needLegacyFG(const HcalTrigTowerDetId& id) const {
0864
0865
0866
0867 if (id.ietaAbs() == HBHE_OVERLAP_TOWER and not upgrade_hb_)
0868 return true;
0869 return false;
0870 }
0871
0872 bool HcalTriggerPrimitiveAlgo::needUpgradeID(const HcalTrigTowerDetId& id, int depth) const {
0873
0874
0875
0876 unsigned int aieta = id.ietaAbs();
0877 if (aieta >= FIRST_DEPTH7_TOWER and aieta <= LAST_FINEGRAIN_TOWER and depth > LAST_FINEGRAIN_DEPTH)
0878 return true;
0879 return false;
0880 }
0881
0882 void HcalTriggerPrimitiveAlgo::addUpgradeFG(const HcalTrigTowerDetId& id,
0883 int depth,
0884 const std::vector<std::bitset<2>>& bits) {
0885 if (not validUpgradeFG(id, depth)) {
0886 if (needLegacyFG(id)) {
0887 std::vector<bool> pseudo(bits.size(), false);
0888 addFG(id, pseudo);
0889 } else if (needUpgradeID(id, depth)) {
0890
0891
0892
0893 auto it = fgUpgradeMap_.find(id);
0894 if (it == fgUpgradeMap_.end()) {
0895 FGUpgradeContainer element;
0896 element.resize(bits.size());
0897 fgUpgradeMap_.insert(std::make_pair(id, element));
0898 }
0899 }
0900
0901 return;
0902 }
0903
0904 auto it = fgUpgradeMap_.find(id);
0905 if (it == fgUpgradeMap_.end()) {
0906 FGUpgradeContainer element;
0907 element.resize(bits.size());
0908 it = fgUpgradeMap_.insert(std::make_pair(id, element)).first;
0909 }
0910 for (unsigned int i = 0; i < bits.size(); ++i) {
0911 it->second[i][0][depth - 1] = bits[i][0];
0912 it->second[i][1][depth - 1] = bits[i][1];
0913 }
0914 }
0915
0916 void HcalTriggerPrimitiveAlgo::addUpgradeTDCFG(const HcalTrigTowerDetId& id, const QIE11DataFrame& frame) {
0917 HcalDetId detId(frame.id());
0918 if (detId.subdet() != HcalEndcap && detId.subdet() != HcalBarrel)
0919 return;
0920
0921 std::vector<HcalTrigTowerDetId> ids = theTrigTowerGeometry->towerIds(detId);
0922 assert(ids.size() == 1 || ids.size() == 2);
0923 IntegerCaloSamples samples1(ids[0], int(frame.samples()));
0924 samples1.setPresamples(frame.presamples());
0925 incoder_->adc2Linear(frame, samples1);
0926 std::vector<unsigned short> bits12_15 = incoder_->group0FGbits(frame);
0927
0928 bool is_compressed = false;
0929 if (detId.subdet() == HcalBarrel) {
0930 is_compressed = (frame.flavor() == 3);
0931
0932 }
0933
0934 auto it = fgUpgradeTDCMap_.find(id);
0935 if (it == fgUpgradeTDCMap_.end()) {
0936 FGUpgradeTDCContainer element;
0937 element.resize(frame.samples());
0938 it = fgUpgradeTDCMap_.insert(std::make_pair(id, element)).first;
0939 }
0940 for (int i = 0; i < frame.samples(); i++) {
0941 it->second[i][detId.depth() - 1] =
0942 std::make_pair(std::make_pair(bits12_15[i], is_compressed), std::make_pair(frame[i].tdc(), samples1[i]));
0943 }
0944 }
0945
0946 void HcalTriggerPrimitiveAlgo::setWeightsQIE11(const edm::ParameterSet& weightsQIE11) {
0947
0948 std::vector<std::string> ietaStrs = weightsQIE11.getParameterNames();
0949 for (auto& ietaStr : ietaStrs) {
0950
0951 auto const& v = weightsQIE11.getParameter<std::vector<double>>(ietaStr);
0952 weightsQIE11_[std::stoi(ietaStr.substr(4))] = {{v[0], v[1]}};
0953 }
0954 }
0955
0956 void HcalTriggerPrimitiveAlgo::setWeightQIE11(int aieta, double weight) {
0957
0958
0959 weightsQIE11_[aieta] = {{weight, 1.0}};
0960 }
0961
0962 void HcalTriggerPrimitiveAlgo::setPeakFinderAlgorithm(int algo) {
0963 if (algo <= 0 || algo > 2)
0964 throw cms::Exception("ERROR: Only algo 1 & 2 are supported.") << std::endl;
0965 peak_finder_algorithm_ = algo;
0966 }
0967
0968 void HcalTriggerPrimitiveAlgo::setNCTScaleShift(int shift) { NCTScaleShift = shift; }
0969
0970 void HcalTriggerPrimitiveAlgo::setRCTScaleShift(int shift) { RCTScaleShift = shift; }