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