Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:07:54

0001 /*
0002  * \file L1TStage2CaloLayer1.cc
0003  *
0004  * N. Smith <nick.smith@cern.ch>
0005  */
0006 //Modified by Bhawna Gomber <bhawna.gomber@cern.ch>
0007 //Modified by Andrew Loeliger <andrew.loeliger@cern.ch>
0008 //Modified by Ho-Fung Tsoi <ho.fung.tsoi@cern.ch>
0009 
0010 #include "DQM/L1TMonitor/interface/L1TStage2CaloLayer1.h"
0011 
0012 #include "FWCore/Framework/interface/ESHandle.h"
0013 #include "FWCore/Framework/interface/LuminosityBlock.h"
0014 #include "DataFormats/Provenance/interface/EventAuxiliary.h"
0015 
0016 #include "CondFormats/RunInfo/interface/RunInfo.h"
0017 #include "CondFormats/DataRecord/interface/RunSummaryRcd.h"
0018 
0019 #include "EventFilter/L1TXRawToDigi/interface/UCTDAQRawData.h"
0020 #include "EventFilter/L1TXRawToDigi/interface/UCTAMCRawData.h"
0021 
0022 using namespace l1t;
0023 
0024 L1TStage2CaloLayer1::L1TStage2CaloLayer1(const edm::ParameterSet& ps)
0025     : ecalTPSourceRecd_(consumes<EcalTrigPrimDigiCollection>(ps.getParameter<edm::InputTag>("ecalTPSourceRecd"))),
0026       ecalTPSourceRecdLabel_(ps.getParameter<edm::InputTag>("ecalTPSourceRecd").label()),
0027       ecalTPSourceRecdBx1_(consumes<EcalTrigPrimDigiCollection>(ps.getParameter<edm::InputTag>("ecalTPSourceRecdBx1"))),
0028       ecalTPSourceRecdBx1Label_(ps.getParameter<edm::InputTag>("ecalTPSourceRecdBx1").label()),
0029       ecalTPSourceRecdBx2_(consumes<EcalTrigPrimDigiCollection>(ps.getParameter<edm::InputTag>("ecalTPSourceRecdBx2"))),
0030       ecalTPSourceRecdBx2Label_(ps.getParameter<edm::InputTag>("ecalTPSourceRecdBx2").label()),
0031       ecalTPSourceRecdBx3_(consumes<EcalTrigPrimDigiCollection>(ps.getParameter<edm::InputTag>("ecalTPSourceRecdBx3"))),
0032       ecalTPSourceRecdBx3Label_(ps.getParameter<edm::InputTag>("ecalTPSourceRecdBx3").label()),
0033       ecalTPSourceRecdBx4_(consumes<EcalTrigPrimDigiCollection>(ps.getParameter<edm::InputTag>("ecalTPSourceRecdBx4"))),
0034       ecalTPSourceRecdBx4Label_(ps.getParameter<edm::InputTag>("ecalTPSourceRecdBx4").label()),
0035       ecalTPSourceRecdBx5_(consumes<EcalTrigPrimDigiCollection>(ps.getParameter<edm::InputTag>("ecalTPSourceRecdBx5"))),
0036       ecalTPSourceRecdBx5Label_(ps.getParameter<edm::InputTag>("ecalTPSourceRecdBx5").label()),
0037       hcalTPSourceRecd_(consumes<HcalTrigPrimDigiCollection>(ps.getParameter<edm::InputTag>("hcalTPSourceRecd"))),
0038       hcalTPSourceRecdLabel_(ps.getParameter<edm::InputTag>("hcalTPSourceRecd").label()),
0039       ecalTPSourceSent_(consumes<EcalTrigPrimDigiCollection>(ps.getParameter<edm::InputTag>("ecalTPSourceSent"))),
0040       ecalTPSourceSentLabel_(ps.getParameter<edm::InputTag>("ecalTPSourceSent").label()),
0041       hcalTPSourceSent_(consumes<HcalTrigPrimDigiCollection>(ps.getParameter<edm::InputTag>("hcalTPSourceSent"))),
0042       hcalTPSourceSentLabel_(ps.getParameter<edm::InputTag>("hcalTPSourceSent").label()),
0043       CaloTowerCollectionData_(
0044           consumes<l1t::CaloTowerBxCollection>(ps.getParameter<edm::InputTag>("CaloTowerCollectionData"))),
0045       CaloTowerCollectionDataLabel_(ps.getParameter<edm::InputTag>("CaloTowerCollectionData").label()),
0046       fedRawData_(consumes<FEDRawDataCollection>(ps.getParameter<edm::InputTag>("fedRawDataLabel"))),
0047       histFolder_(ps.getParameter<std::string>("histFolder")),
0048       tpFillThreshold_(ps.getUntrackedParameter<int>("etDistributionsFillThreshold", 0)),
0049       tpFillThreshold5Bx_(ps.getUntrackedParameter<int>("etDistributionsFillThreshold5Bx", 4)),
0050       ignoreHFfbs_(ps.getUntrackedParameter<bool>("ignoreHFfbs", false)) {}
0051 
0052 L1TStage2CaloLayer1::~L1TStage2CaloLayer1() {}
0053 
0054 void L1TStage2CaloLayer1::dqmAnalyze(const edm::Event& event,
0055                                      const edm::EventSetup& es,
0056                                      const CaloL1Information::monitoringDataHolder& eventMonitors) const {
0057   //This must be moved here compared to previous version to avoid a const function modifying a class variable
0058   //This unfortunately means that the variable is local and reallocated per event
0059   std::vector<std::pair<EcalTriggerPrimitiveDigi, EcalTriggerPrimitiveDigi>> ecalTPSentRecd_;
0060   std::vector<std::pair<HcalTriggerPrimitiveDigi, HcalTriggerPrimitiveDigi>> hcalTPSentRecd_;
0061   ecalTPSentRecd_.reserve(28 * 2 * 72);
0062   hcalTPSentRecd_.reserve(41 * 2 * 72);
0063 
0064   // Monitorables stored in Layer 1 raw data but
0065   // not accessible from existing persistent data formats
0066   edm::Handle<FEDRawDataCollection> fedRawDataCollection;
0067   event.getByToken(fedRawData_, fedRawDataCollection);
0068   bool caloLayer1OutOfRun{true};
0069   bool FATevent{false};
0070   bool additionalFB{false};
0071   bool card7flag1354{false};
0072   bool card7flag1356{false};
0073   bool card7flag1358{false};
0074   std::vector<uint32_t> FED1354_slot7_bits;
0075   std::vector<uint32_t> FED1356_slot7_bits;
0076   std::vector<uint32_t> FED1358_slot7_bits;
0077   for (int i = 0; i < 6; i++) {
0078     FED1354_slot7_bits.push_back(0);
0079     FED1356_slot7_bits.push_back(0);
0080     FED1358_slot7_bits.push_back(0);
0081   }
0082   if (fedRawDataCollection.isValid()) {
0083     caloLayer1OutOfRun = false;
0084     for (int iFed = 1354; iFed < 1360; iFed += 2) {
0085       const FEDRawData& fedRawData = fedRawDataCollection->FEDData(iFed);
0086       if (fedRawData.size() == 0) {
0087         caloLayer1OutOfRun = true;
0088         continue;  // In case one of 3 layer 1 FEDs not in
0089       }
0090       const uint64_t* fedRawDataArray = (const uint64_t*)fedRawData.data();
0091       UCTDAQRawData daqData(fedRawDataArray);
0092       for (uint32_t i = 0; i < daqData.nAMCs(); i++) {
0093         UCTAMCRawData amcData(daqData.amcPayload(i));
0094         const uint32_t* amcPtr = amcData.dataPtr();
0095         FATevent = ((amcPtr[5] >> 16) & 0xf) == 5;
0096         additionalFB = (amcPtr[5] >> 15) & 0x1;
0097         if (((amcPtr[5] >> 14) & 0x1) == 0) {
0098           int lPhi = amcData.layer1Phi();
0099           if (daqData.BXID() != amcData.BXID()) {
0100             eventMonitors.bxidErrors_->Fill(lPhi);
0101           }
0102           if (daqData.L1ID() != amcData.L1ID()) {
0103             eventMonitors.l1idErrors_->Fill(lPhi);
0104           }
0105           // AMC payload header has 16 bit orbit number, AMC13 header is full 32
0106           if ((daqData.orbitNumber() & 0xFFFF) != amcData.orbitNo()) {
0107             eventMonitors.orbitErrors_->Fill(lPhi);
0108           }
0109         }
0110         if (iFed == 1354 && daqData.nAMCs() == 7 && i == 3 && amcData.amcNo() == 7 && ((amcPtr[5] >> 14) & 0x1) == 1) {
0111           card7flag1354 = true;
0112           for (int j = 0; j < 6; j++) {
0113             FED1354_slot7_bits[j] = amcPtr[j + 6] & 0xFFFFFFFF;
0114           }
0115         } else if (iFed == 1356 && daqData.nAMCs() == 7 && i == 3 && amcData.amcNo() == 7 &&
0116                    ((amcPtr[5] >> 14) & 0x1) == 1) {
0117           card7flag1356 = true;
0118           for (int j = 0; j < 6; j++) {
0119             FED1356_slot7_bits[j] = amcPtr[j + 6] & 0xFFFFFFFF;
0120           }
0121         } else if (iFed == 1358 && daqData.nAMCs() == 7 && i == 3 && amcData.amcNo() == 7 &&
0122                    ((amcPtr[5] >> 14) & 0x1) == 1) {
0123           card7flag1358 = true;
0124           for (int j = 0; j < 6; j++) {
0125             FED1358_slot7_bits[j] = amcPtr[j + 6] & 0xFFFFFFFF;
0126           }
0127         }
0128       }
0129     }
0130   }
0131 
0132   if (not caloLayer1OutOfRun) {
0133     for (int iword = 0; iword < 6; iword++) {
0134       for (int ibit = 0; ibit < 32; ibit++) {
0135         if (card7flag1354 && ((FED1354_slot7_bits[iword] >> ibit) & 0x1) == 1) {
0136           eventMonitors.slot7bit_->Fill(ibit, iword);
0137         }
0138         if (card7flag1356 && ((FED1356_slot7_bits[iword] >> ibit) & 0x1) == 1) {
0139           eventMonitors.slot7bit_->Fill(ibit, iword + 7);
0140         }
0141         if (card7flag1358 && ((FED1358_slot7_bits[iword] >> ibit) & 0x1) == 1) {
0142           eventMonitors.slot7bit_->Fill(ibit, iword + 14);
0143         }
0144       }
0145     }
0146   }
0147 
0148   edm::Handle<EcalTrigPrimDigiCollection> ecalTPsSent;
0149   event.getByToken(ecalTPSourceSent_, ecalTPsSent);
0150   edm::Handle<EcalTrigPrimDigiCollection> ecalTPsRecd;
0151   event.getByToken(ecalTPSourceRecd_, ecalTPsRecd);
0152 
0153   ecalTPSentRecd_.clear();
0154 
0155   ComparisonHelper::zip(ecalTPsSent->begin(),
0156                         ecalTPsSent->end(),
0157                         ecalTPsRecd->begin(),
0158                         ecalTPsRecd->end(),
0159                         std::inserter(ecalTPSentRecd_, ecalTPSentRecd_.begin()),
0160                         EcalTrigPrimDigiCollection::key_compare());
0161 
0162   int nEcalLinkErrors{0};
0163   int nEcalMismatch{0};
0164 
0165   for (const auto& tpPair : ecalTPSentRecd_) {
0166     auto sentTp = tpPair.first;
0167     if (sentTp.compressedEt() < 0) {
0168       // ECal zero-suppresses digis, and a default-constructed
0169       // digi has et=-1 apparently, but we know it should be zero
0170       EcalTriggerPrimitiveSample sample(0);
0171       EcalTriggerPrimitiveDigi tpg(sentTp.id());
0172       tpg.setSize(1);
0173       tpg.setSample(0, sample);
0174       swap(sentTp, tpg);
0175     }
0176     const auto& recdTp = tpPair.second;
0177     const int ieta = sentTp.id().ieta();
0178     const int iphi = sentTp.id().iphi();
0179     const bool towerMasked = recdTp.sample(0).raw() & (1 << 13);
0180     const bool linkMasked = recdTp.sample(0).raw() & (1 << 14);
0181     const bool linkError = recdTp.sample(0).raw() & (1 << 15);
0182 
0183     // Link status bits from layer 1
0184     if (towerMasked) {
0185       eventMonitors.ecalOccTowerMasked_->Fill(ieta, iphi);
0186     }
0187     if (linkMasked) {
0188       eventMonitors.ecalOccLinkMasked_->Fill(ieta, iphi);
0189     }
0190 
0191     if (sentTp.compressedEt() > tpFillThreshold_) {
0192       eventMonitors.ecalTPRawEtSent_->Fill(sentTp.compressedEt());
0193       eventMonitors.ecalOccSent_->Fill(ieta, iphi);
0194     }
0195     if (sentTp.fineGrain() == 1) {
0196       eventMonitors.ecalOccSentFgVB_->Fill(ieta, iphi);
0197     }
0198 
0199     if (towerMasked || caloLayer1OutOfRun) {
0200       // Do not compare if we have a mask applied
0201       continue;
0202     }
0203 
0204     if (linkError) {
0205       eventMonitors.ecalLinkError_->Fill(ieta, iphi);
0206       eventMonitors.ecalLinkErrorByLumi_->Fill(event.id().luminosityBlock());
0207       nEcalLinkErrors++;
0208       // Don't compare anymore, we already know its bad
0209       continue;
0210     }
0211 
0212     eventMonitors.ecalTPRawEtCorrelation_->Fill(sentTp.compressedEt(), recdTp.compressedEt());
0213 
0214     if (recdTp.compressedEt() > tpFillThreshold_) {
0215       eventMonitors.ecalTPRawEtRecd_->Fill(recdTp.compressedEt());
0216       eventMonitors.ecalOccupancy_->Fill(ieta, iphi);
0217       eventMonitors.ecalOccRecdEtWgt_->Fill(ieta, iphi, recdTp.compressedEt());
0218     }
0219     if (recdTp.fineGrain() == 1) {
0220       eventMonitors.ecalOccRecdFgVB_->Fill(ieta, iphi);
0221     }
0222 
0223     // Now for comparisons
0224 
0225     const bool EetAgreement = sentTp.compressedEt() == recdTp.compressedEt();
0226     const bool EfbAgreement = sentTp.fineGrain() == recdTp.fineGrain();
0227     if (EetAgreement && EfbAgreement) {
0228       // Full match
0229       if (sentTp.compressedEt() > tpFillThreshold_) {
0230         eventMonitors.ecalOccSentAndRecd_->Fill(ieta, iphi);
0231         eventMonitors.ecalTPRawEtSentAndRecd_->Fill(sentTp.compressedEt());
0232       }
0233     } else {
0234       // There is some issue
0235       eventMonitors.ecalDiscrepancy_->Fill(ieta, iphi);
0236       eventMonitors.ecalMismatchByLumi_->Fill(event.id().luminosityBlock());
0237       eventMonitors.ECALmismatchesPerBx_->Fill(event.bunchCrossing());
0238       nEcalMismatch++;
0239 
0240       if (not EetAgreement) {
0241         eventMonitors.ecalOccEtDiscrepancy_->Fill(ieta, iphi);
0242         eventMonitors.ecalTPRawEtDiffNoMatch_->Fill(recdTp.compressedEt() - sentTp.compressedEt());
0243         updateMismatch(event, 0, streamCache(event.streamID())->streamMismatchList);
0244 
0245         if (sentTp.compressedEt() == 0)
0246           eventMonitors.ecalOccRecdNotSent_->Fill(ieta, iphi);
0247         else if (recdTp.compressedEt() == 0)
0248           eventMonitors.ecalOccSentNotRecd_->Fill(ieta, iphi);
0249         else
0250           eventMonitors.ecalOccNoMatch_->Fill(ieta, iphi);
0251       }
0252       if (not EfbAgreement) {
0253         // occ for fine grain mismatch
0254         eventMonitors.ecalOccFgDiscrepancy_->Fill(ieta, iphi);
0255         updateMismatch(event, 1, streamCache(event.streamID())->streamMismatchList);
0256       }
0257     }
0258   }
0259 
0260   if (nEcalLinkErrors > streamCache(event.streamID())->streamNumMaxEvtLinkErrorsECAL)
0261     streamCache(event.streamID())->streamNumMaxEvtLinkErrorsECAL = nEcalLinkErrors;
0262 
0263   if (nEcalMismatch > streamCache(event.streamID())->streamNumMaxEvtMismatchECAL)
0264     streamCache(event.streamID())->streamNumMaxEvtMismatchECAL = nEcalMismatch;
0265 
0266   edm::Handle<EcalTrigPrimDigiCollection> ecalTPsRecdBx1;
0267   event.getByToken(ecalTPSourceRecdBx1_, ecalTPsRecdBx1);
0268   edm::Handle<EcalTrigPrimDigiCollection> ecalTPsRecdBx2;
0269   event.getByToken(ecalTPSourceRecdBx2_, ecalTPsRecdBx2);
0270   edm::Handle<EcalTrigPrimDigiCollection> ecalTPsRecdBx3;
0271   event.getByToken(ecalTPSourceRecdBx3_, ecalTPsRecdBx3);
0272   edm::Handle<EcalTrigPrimDigiCollection> ecalTPsRecdBx4;
0273   event.getByToken(ecalTPSourceRecdBx4_, ecalTPsRecdBx4);
0274   edm::Handle<EcalTrigPrimDigiCollection> ecalTPsRecdBx5;
0275   event.getByToken(ecalTPSourceRecdBx5_, ecalTPsRecdBx5);
0276 
0277   for (const auto& tp : (*ecalTPsRecdBx1)) {
0278     if (tp.compressedEt() > tpFillThreshold5Bx_) {
0279       const int ieta = tp.id().ieta();
0280       const int iphi = tp.id().iphi();
0281       const bool towerMasked = tp.sample(0).raw() & (1 << 13);
0282       const bool linkError = tp.sample(0).raw() & (1 << 15);
0283       if (towerMasked || caloLayer1OutOfRun || linkError) {
0284         continue;
0285       }
0286       eventMonitors.ecalOccRecdBx1_->Fill(ieta, iphi);
0287       eventMonitors.ecalOccRecd5Bx_->Fill(1);
0288       eventMonitors.ecalOccRecd5BxEtWgt_->Fill(1, tp.compressedEt());
0289     }
0290   }
0291   for (const auto& tp : (*ecalTPsRecdBx2)) {
0292     if (tp.compressedEt() > tpFillThreshold5Bx_) {
0293       const int ieta = tp.id().ieta();
0294       const int iphi = tp.id().iphi();
0295       const bool towerMasked = tp.sample(0).raw() & (1 << 13);
0296       const bool linkError = tp.sample(0).raw() & (1 << 15);
0297       if (towerMasked || caloLayer1OutOfRun || linkError) {
0298         continue;
0299       }
0300       eventMonitors.ecalOccRecdBx2_->Fill(ieta, iphi);
0301       eventMonitors.ecalOccRecd5Bx_->Fill(2);
0302       eventMonitors.ecalOccRecd5BxEtWgt_->Fill(2, tp.compressedEt());
0303     }
0304   }
0305   for (const auto& tp : (*ecalTPsRecdBx3)) {
0306     if (tp.compressedEt() > tpFillThreshold5Bx_) {
0307       const int ieta = tp.id().ieta();
0308       const int iphi = tp.id().iphi();
0309       const bool towerMasked = tp.sample(0).raw() & (1 << 13);
0310       const bool linkError = tp.sample(0).raw() & (1 << 15);
0311       if (towerMasked || caloLayer1OutOfRun || linkError) {
0312         continue;
0313       }
0314       eventMonitors.ecalOccRecdBx3_->Fill(ieta, iphi);
0315       eventMonitors.ecalOccRecd5Bx_->Fill(3);
0316       eventMonitors.ecalOccRecd5BxEtWgt_->Fill(3, tp.compressedEt());
0317     }
0318   }
0319   for (const auto& tp : (*ecalTPsRecdBx4)) {
0320     if (tp.compressedEt() > tpFillThreshold5Bx_) {
0321       const int ieta = tp.id().ieta();
0322       const int iphi = tp.id().iphi();
0323       const bool towerMasked = tp.sample(0).raw() & (1 << 13);
0324       const bool linkError = tp.sample(0).raw() & (1 << 15);
0325       if (towerMasked || caloLayer1OutOfRun || linkError) {
0326         continue;
0327       }
0328       eventMonitors.ecalOccRecdBx4_->Fill(ieta, iphi);
0329       eventMonitors.ecalOccRecd5Bx_->Fill(4);
0330       eventMonitors.ecalOccRecd5BxEtWgt_->Fill(4, tp.compressedEt());
0331     }
0332   }
0333   for (const auto& tp : (*ecalTPsRecdBx5)) {
0334     if (tp.compressedEt() > tpFillThreshold5Bx_) {
0335       const int ieta = tp.id().ieta();
0336       const int iphi = tp.id().iphi();
0337       const bool towerMasked = tp.sample(0).raw() & (1 << 13);
0338       const bool linkError = tp.sample(0).raw() & (1 << 15);
0339       if (towerMasked || caloLayer1OutOfRun || linkError) {
0340         continue;
0341       }
0342       eventMonitors.ecalOccRecdBx5_->Fill(ieta, iphi);
0343       eventMonitors.ecalOccRecd5Bx_->Fill(5);
0344       eventMonitors.ecalOccRecd5BxEtWgt_->Fill(5, tp.compressedEt());
0345     }
0346   }
0347 
0348   edm::Handle<HcalTrigPrimDigiCollection> hcalTPsSent;
0349   event.getByToken(hcalTPSourceSent_, hcalTPsSent);
0350   edm::Handle<HcalTrigPrimDigiCollection> hcalTPsRecd;
0351   event.getByToken(hcalTPSourceRecd_, hcalTPsRecd);
0352   edm::Handle<l1t::CaloTowerBxCollection> caloTowerDataCol;
0353   event.getByToken(CaloTowerCollectionData_, caloTowerDataCol);
0354 
0355   hcalTPSentRecd_.clear();
0356 
0357   ComparisonHelper::zip(hcalTPsSent->begin(),
0358                         hcalTPsSent->end(),
0359                         hcalTPsRecd->begin(),
0360                         hcalTPsRecd->end(),
0361                         std::inserter(hcalTPSentRecd_, hcalTPSentRecd_.begin()),
0362                         HcalTrigPrimDigiCollection::key_compare());
0363 
0364   int nHcalLinkErrors{0};
0365   int nHcalMismatch{0};
0366 
0367   for (const auto& tpPair : hcalTPSentRecd_) {
0368     const auto& sentTp = tpPair.first;
0369     const auto& recdTp = tpPair.second;
0370     const int ieta = sentTp.id().ieta();
0371     if (abs(ieta) > 28 && sentTp.id().version() != 1)
0372       continue;
0373     const int iphi = sentTp.id().iphi();
0374     const bool towerMasked = recdTp.sample(0).raw() & (1 << 13);
0375     const bool linkMasked = recdTp.sample(0).raw() & (1 << 14);
0376     const bool linkError = recdTp.sample(0).raw() & (1 << 15);
0377 
0378     if (towerMasked) {
0379       eventMonitors.hcalOccTowerMasked_->Fill(ieta, iphi);
0380     }
0381     if (linkMasked) {
0382       eventMonitors.hcalOccLinkMasked_->Fill(ieta, iphi);
0383     }
0384 
0385     if (sentTp.SOI_compressedEt() > tpFillThreshold_) {
0386       eventMonitors.hcalTPRawEtSent_->Fill(sentTp.SOI_compressedEt());
0387       eventMonitors.hcalOccSent_->Fill(ieta, iphi);
0388     }
0389 
0390     // 6 HCAL fine grain bits from uHTR readout
0391     bool uHTRfg0 = sentTp.SOI_fineGrain(0);
0392     bool uHTRfg1 = sentTp.SOI_fineGrain(1);
0393     bool uHTRfg2 = sentTp.SOI_fineGrain(2);
0394     bool uHTRfg3 = sentTp.SOI_fineGrain(3);
0395     bool uHTRfg4 = sentTp.SOI_fineGrain(4);
0396     bool uHTRfg5 = sentTp.SOI_fineGrain(5);
0397 
0398     if (uHTRfg0) {
0399       eventMonitors.hcalOccSentFg0_->Fill(ieta, iphi);
0400     }
0401     if (uHTRfg1) {
0402       eventMonitors.hcalOccSentFg1_->Fill(ieta, iphi);
0403     }
0404     if (uHTRfg2) {
0405       eventMonitors.hcalOccSentFg2_->Fill(ieta, iphi);
0406     }
0407     if (uHTRfg3) {
0408       eventMonitors.hcalOccSentFg3_->Fill(ieta, iphi);
0409     }
0410     if (uHTRfg4) {
0411       eventMonitors.hcalOccSentFg4_->Fill(ieta, iphi);
0412     }
0413     if (uHTRfg5) {
0414       eventMonitors.hcalOccSentFg5_->Fill(ieta, iphi);
0415     }
0416 
0417     if (towerMasked || caloLayer1OutOfRun) {
0418       // Do not compare if we have a mask applied
0419       continue;
0420     }
0421 
0422     if (linkError) {
0423       eventMonitors.hcalLinkError_->Fill(ieta, iphi);
0424       eventMonitors.hcalLinkErrorByLumi_->Fill(event.id().luminosityBlock());
0425       nHcalLinkErrors++;
0426       // Don't compare anymore, we already know its bad
0427       continue;
0428     }
0429 
0430     // 6 HCAL fine grain bits from Layer1 readout
0431     // Fg 0 is always there in ctp7 payload, set as bit 8 after Et bits in unpacked towerDatum
0432     // When additionalFB flag is set:
0433     // Fg 1-5 are unpacked and set as bits 0-5 in another unpacked towerDatum2 and packed in HCALTP sample(1)
0434     // since standard sample size is 16bit and there is no room in sample(0) which contains already Et and link status
0435     // Otherwise:
0436     // Fg 1-5 are all zero
0437     bool layer1fg0 = recdTp.SOI_fineGrain(0);
0438     bool layer1fg1 = false;
0439     bool layer1fg2 = false;
0440     bool layer1fg3 = false;
0441     bool layer1fg4 = false;
0442     bool layer1fg5 = false;
0443     if (additionalFB && (abs(ieta) < 29)) {
0444       for (const auto& tp : (*hcalTPsRecd)) {
0445         if (not(tp.id().ieta() == ieta && tp.id().iphi() == iphi)) {
0446           continue;
0447         }
0448         layer1fg1 = tp.sample(1).raw() & (1 << 0);
0449         layer1fg2 = tp.sample(1).raw() & (1 << 1);
0450         layer1fg3 = tp.sample(1).raw() & (1 << 2);
0451         layer1fg4 = tp.sample(1).raw() & (1 << 3);
0452         layer1fg5 = tp.sample(1).raw() & (1 << 4);
0453       }
0454     }
0455 
0456     // Check mismatches only for HBHE
0457     const bool Hfg0Agreement = (abs(ieta) < 29) ? (layer1fg0 == uHTRfg0) : true;
0458     const bool Hfg1Agreement = (abs(ieta) < 29) ? (layer1fg1 == uHTRfg1) : true;
0459     const bool Hfg2Agreement = (abs(ieta) < 29) ? (layer1fg2 == uHTRfg2) : true;
0460     const bool Hfg3Agreement = (abs(ieta) < 29) ? (layer1fg3 == uHTRfg3) : true;
0461     //const bool Hfg4Agreement = (abs(ieta) < 29) ? (layer1fg4 == uHTRfg4) : true;
0462     //const bool Hfg5Agreement = (abs(ieta) < 29) ? (layer1fg5 == uHTRfg5) : true;
0463     // Mute fg4 and fg5 for now (reserved bits not used anyway)
0464     const bool HfgAgreement = (Hfg0Agreement && Hfg1Agreement && Hfg2Agreement && Hfg3Agreement);
0465 
0466     // Construct an 6-bit integer from the layer1 fine grain readout (input to 6:1 logic emulation)
0467     uint64_t fg_bits = 0;
0468     if (layer1fg0) {
0469       fg_bits |= 0x1;
0470     }
0471     if (layer1fg1) {
0472       fg_bits |= 0x1 << 1;
0473     }
0474     if (layer1fg2) {
0475       fg_bits |= 0x1 << 2;
0476     }
0477     if (layer1fg3) {
0478       fg_bits |= 0x1 << 3;
0479     }
0480     if (layer1fg4) {
0481       fg_bits |= 0x1 << 4;
0482     }
0483     if (layer1fg5) {
0484       fg_bits |= 0x1 << 5;
0485     }
0486 
0487     // Current 6:1 LUT in fw
0488     const uint64_t HCalFbLUT = 0xBBBABBBABBBABBBA;
0489     // Expected LLP bit output (mute emulation for normal events, since layer2 only reads out FAT events)
0490     const bool LLPfb_Expd = (FATevent == 1) ? ((HCalFbLUT >> fg_bits) & 1) : false;
0491     // Actual LLP bit output in layer2 data collection
0492     uint32_t tower_hwqual = 0;
0493     for (auto tower = caloTowerDataCol->begin(0); tower != caloTowerDataCol->end(0); ++tower) {
0494       if (not(tower->hwEta() == ieta && tower->hwPhi() == iphi)) {
0495         continue;
0496       }
0497       tower_hwqual = tower->hwQual();
0498     }
0499     // CaloTower hwQual is 4-bit long, LLP output bit is set at the 2nd bit (counting from 0)
0500     const bool LLPfb_Data = ((tower_hwqual & 0b0100) >> 2) & 1;
0501 
0502     const bool LLPfbAgreement = (abs(ieta) < 29) ? (LLPfb_Expd == LLPfb_Data) : true;
0503 
0504     // Fill feature bits Occ for HBHE only
0505     if (abs(ieta) < 29) {
0506       if (layer1fg0) {
0507         eventMonitors.hcalOccRecdFg0_->Fill(ieta, iphi);
0508       }
0509       if (layer1fg1) {
0510         eventMonitors.hcalOccRecdFg1_->Fill(ieta, iphi);
0511       }
0512       if (layer1fg2) {
0513         eventMonitors.hcalOccRecdFg2_->Fill(ieta, iphi);
0514       }
0515       if (layer1fg3) {
0516         eventMonitors.hcalOccRecdFg3_->Fill(ieta, iphi);
0517       }
0518       if (layer1fg4) {
0519         eventMonitors.hcalOccRecdFg4_->Fill(ieta, iphi);
0520       }
0521       if (layer1fg5) {
0522         eventMonitors.hcalOccRecdFg5_->Fill(ieta, iphi);
0523       }
0524       // fg4-5 are reserved bits and not used
0525       // so compare here and not stream to mismatch list for now
0526       //if (not Hfg4Agreement) {
0527       //  eventMonitors.hcalOccFg4Discrepancy_->Fill(ieta, iphi);
0528       //}
0529       //if (not Hfg5Agreement) {
0530       //  eventMonitors.hcalOccFg5Discrepancy_->Fill(ieta, iphi);
0531       //}
0532       // Fill Fb Occ and compare between layer1 emulated and layer2 data readout
0533       // FAT events only!!
0534       if (LLPfb_Expd) {
0535         eventMonitors.hcalOccLLPFbExpd_->Fill(ieta, iphi);
0536       }
0537       if (LLPfb_Data) {
0538         eventMonitors.hcalOccLLPFbData_->Fill(ieta, iphi);
0539       }
0540     }
0541 
0542     if (recdTp.SOI_compressedEt() > tpFillThreshold_) {
0543       eventMonitors.hcalTPRawEtRecd_->Fill(recdTp.SOI_compressedEt());
0544       eventMonitors.hcalOccupancy_->Fill(ieta, iphi);
0545       eventMonitors.hcalOccRecdEtWgt_->Fill(ieta, iphi, recdTp.SOI_compressedEt());
0546     }
0547 
0548     if (abs(ieta) > 29) {
0549       eventMonitors.hcalTPRawEtCorrelationHF_->Fill(sentTp.SOI_compressedEt(), recdTp.SOI_compressedEt());
0550     } else {
0551       eventMonitors.hcalTPRawEtCorrelationHBHE_->Fill(sentTp.SOI_compressedEt(), recdTp.SOI_compressedEt());
0552     }
0553 
0554     const bool HetAgreement = sentTp.SOI_compressedEt() == recdTp.SOI_compressedEt();
0555     if (HetAgreement && HfgAgreement && LLPfbAgreement) {
0556       // Full match
0557       if (sentTp.SOI_compressedEt() > tpFillThreshold_) {
0558         eventMonitors.hcalOccSentAndRecd_->Fill(ieta, iphi);
0559         eventMonitors.hcalTPRawEtSentAndRecd_->Fill(sentTp.SOI_compressedEt());
0560       }
0561     } else {
0562       // There is some issue
0563       eventMonitors.hcalDiscrepancy_->Fill(ieta, iphi);
0564       eventMonitors.hcalMismatchByLumi_->Fill(event.id().luminosityBlock());
0565       nHcalMismatch++;
0566 
0567       if (not HetAgreement) {
0568         if (abs(ieta) > 29) {
0569           eventMonitors.HFmismatchesPerBx_->Fill(event.bunchCrossing());
0570         } else {
0571           eventMonitors.HBHEmismatchesPerBx_->Fill(event.bunchCrossing());
0572         }
0573         eventMonitors.hcalOccEtDiscrepancy_->Fill(ieta, iphi);
0574         eventMonitors.hcalTPRawEtDiffNoMatch_->Fill(recdTp.SOI_compressedEt() - sentTp.SOI_compressedEt());
0575         updateMismatch(event, 2, streamCache(event.streamID())->streamMismatchList);
0576 
0577         // Handle HCal discrepancy debug
0578         if (sentTp.SOI_compressedEt() == 0)
0579           eventMonitors.hcalOccRecdNotSent_->Fill(ieta, iphi);
0580         else if (recdTp.SOI_compressedEt() == 0)
0581           eventMonitors.hcalOccSentNotRecd_->Fill(ieta, iphi);
0582         else
0583           eventMonitors.hcalOccNoMatch_->Fill(ieta, iphi);
0584       }
0585       if (not(HfgAgreement && LLPfbAgreement)) {
0586         if (not Hfg0Agreement) {
0587           eventMonitors.hcalOccFg0Discrepancy_->Fill(ieta, iphi);
0588         }
0589         if (not Hfg1Agreement) {
0590           eventMonitors.hcalOccFg1Discrepancy_->Fill(ieta, iphi);
0591         }
0592         if (not Hfg2Agreement) {
0593           eventMonitors.hcalOccFg2Discrepancy_->Fill(ieta, iphi);
0594         }
0595         if (not Hfg3Agreement) {
0596           eventMonitors.hcalOccFg3Discrepancy_->Fill(ieta, iphi);
0597         }
0598         if (not LLPfbAgreement) {
0599           eventMonitors.hcalOccLLPFbDiscrepancy_->Fill(ieta, iphi);
0600         }
0601         updateMismatch(event, 3, streamCache(event.streamID())->streamMismatchList);
0602       }
0603     }
0604   }
0605 
0606   if (nHcalLinkErrors > streamCache(event.streamID())->streamNumMaxEvtLinkErrorsHCAL)
0607     streamCache(event.streamID())->streamNumMaxEvtLinkErrorsHCAL = nHcalLinkErrors;
0608   if (nHcalMismatch > streamCache(event.streamID())->streamNumMaxEvtMismatchHCAL)
0609     streamCache(event.streamID())->streamNumMaxEvtMismatchHCAL = nHcalMismatch;
0610 
0611   //fill inclusive link error and mismatch cache values based on whether HCAL or ECAL had more this event
0612   if (nEcalLinkErrors >= nHcalLinkErrors && nEcalLinkErrors > streamCache(event.streamID())->streamNumMaxEvtLinkErrors)
0613     streamCache(event.streamID())->streamNumMaxEvtLinkErrors = nEcalLinkErrors;
0614   else if (nEcalLinkErrors < nHcalLinkErrors &&
0615            nHcalLinkErrors > streamCache(event.streamID())->streamNumMaxEvtLinkErrors)
0616     streamCache(event.streamID())->streamNumMaxEvtLinkErrors = nHcalLinkErrors;
0617 
0618   if (nEcalMismatch >= nHcalMismatch && nEcalMismatch > streamCache(event.streamID())->streamNumMaxEvtMismatch)
0619     streamCache(event.streamID())->streamNumMaxEvtMismatch = nEcalMismatch;
0620   else if (nEcalMismatch < nHcalMismatch && nHcalMismatch > streamCache(event.streamID())->streamNumMaxEvtMismatch)
0621     streamCache(event.streamID())->streamNumMaxEvtMismatch = nHcalMismatch;
0622 }
0623 
0624 //push the mismatch type and identifying information onto the back of the stream-based vector
0625 //will maintain the vector's size at 20
0626 void L1TStage2CaloLayer1::updateMismatch(
0627     const edm::Event& e,
0628     int mismatchType,
0629     std::vector<std::tuple<edm::RunID, edm::LuminosityBlockID, edm::EventID, std::vector<int>>>& streamMismatches)
0630     const {
0631   //check if this combination of Run/Lumi/Event already exists in the stream mismatch list
0632   //if it does, update that entry, otherwise insert a new one with this particular combination.
0633   for (auto mismatchIterator = streamMismatches.begin(); mismatchIterator != streamMismatches.end();
0634        ++mismatchIterator) {
0635     if (e.getRun().id() == std::get<0>(*mismatchIterator) &&
0636         e.getLuminosityBlock().id() == std::get<1>(*mismatchIterator) && e.id() == std::get<2>(*mismatchIterator)) {
0637       //the run, lumi and event exist. Check if this kind of mismatch has been reported before
0638       std::vector<int>& mismatchTypeVector = std::get<3>(*mismatchIterator);
0639       for (auto mismatchTypeIterator = mismatchTypeVector.begin(); mismatchTypeIterator != mismatchTypeVector.end();
0640            ++mismatchTypeIterator) {
0641         if (mismatchType == *mismatchTypeIterator) {
0642           //this has already been reported
0643           return;
0644         }
0645       }
0646       //A mismatch exists, but it is not a type that has been previously reported.
0647       //Insert it into the vector of types of mismatches reported
0648       mismatchTypeVector.push_back(mismatchType);
0649       return;
0650     }
0651   }
0652 
0653   //The run/lumi/event does not exist in the list, construct an entry
0654   std::vector<int> newMismatchTypeVector;
0655   newMismatchTypeVector.push_back(mismatchType);
0656   std::tuple<edm::RunID, edm::LuminosityBlockID, edm::EventID, std::vector<int>> mismatchToInsert = {
0657       e.getRun().id(), e.getLuminosityBlock().id(), e.id(), newMismatchTypeVector};
0658   streamMismatches.push_back(mismatchToInsert);
0659 
0660   //maintain the mismatch vector size at 20.
0661   if (streamMismatches.size() > 20)
0662     streamMismatches.erase(streamMismatches.begin());
0663 }
0664 
0665 void L1TStage2CaloLayer1::dqmBeginRun(const edm::Run&,
0666                                       const edm::EventSetup&,
0667                                       CaloL1Information::monitoringDataHolder& eventMonitors) const {}
0668 
0669 void L1TStage2CaloLayer1::dqmEndRun(const edm::Run& run,
0670                                     const edm::EventSetup& es,
0671                                     const CaloL1Information::monitoringDataHolder& runMonitors,
0672                                     const CaloL1Information::perRunSummaryMonitoringInformation&) const {}
0673 
0674 void L1TStage2CaloLayer1::bookHistograms(DQMStore::IBooker& ibooker,
0675                                          const edm::Run& run,
0676                                          const edm::EventSetup& es,
0677                                          CaloL1Information::monitoringDataHolder& eventMonitors) const {
0678   auto bookEt = [&ibooker](std::string name, std::string title) {
0679     return ibooker.book1D(name, title + ";Raw ET;Counts", 256, -0.5, 255.5);
0680   };
0681   auto bookEtCorrelation = [&ibooker](std::string name, std::string title) {
0682     return ibooker.book2D(name, title, 256, -0.5, 255.5, 256, -0.5, 255.5);
0683   };
0684   auto bookEtDiff = [&ibooker](std::string name, std::string title) {
0685     return ibooker.book1D(name, title + ";#Delta Raw ET;Counts", 511, -255.5, 255.5);
0686   };
0687   auto bookEcalOccupancy = [&ibooker](std::string name, std::string title) {
0688     return ibooker.book2D(name, title + ";iEta;iPhi", 57, -28.5, 28.5, 72, 0.5, 72.5);
0689   };
0690   auto bookHcalOccupancy = [&ibooker](std::string name, std::string title) {
0691     return ibooker.book2D(name, title + ";iEta;iPhi", 83, -41.5, 41.5, 72, 0.5, 72.5);
0692   };
0693   auto bookSlot7Occupancy = [&ibooker](std::string name, std::string title) {
0694     return ibooker.book2D(
0695         name, title + ";bit;slot-7 word in FED13(54[0-5],56[7-12],58[14-19])", 32, -0.5, 31.5, 20, -0.5, 19.5);
0696   };
0697 
0698   ibooker.setCurrentFolder(histFolder_);
0699 
0700   eventMonitors.ecalDiscrepancy_ =
0701       bookEcalOccupancy("ecalDiscrepancy", "ECAL Discrepancies between TCC and Layer1 Readout");
0702   eventMonitors.ecalLinkError_ = bookEcalOccupancy("ecalLinkError", "ECAL Link Errors");
0703   eventMonitors.ecalOccupancy_ = bookEcalOccupancy("ecalOccupancy", "ECAL TP Occupancy at Layer1");
0704   eventMonitors.ecalOccRecdEtWgt_ = bookEcalOccupancy("ecalOccRecdEtWgt", "ECal TP ET-weighted Occupancy at Layer1");
0705   eventMonitors.hcalDiscrepancy_ =
0706       bookHcalOccupancy("hcalDiscrepancy", "HCAL Discrepancies between uHTR and Layer1 Readout");
0707   eventMonitors.hcalLinkError_ = bookHcalOccupancy("hcalLinkError", "HCAL Link Errors");
0708   eventMonitors.hcalOccupancy_ = bookHcalOccupancy("hcalOccupancy", "HCAL TP Occupancy at Layer1");
0709   eventMonitors.hcalOccRecdEtWgt_ = bookHcalOccupancy("hcalOccRecdEtWgt", "HCal TP ET-weighted Occupancy at Layer1");
0710   eventMonitors.slot7bit_ = bookSlot7Occupancy("slot7bitOcc", "Bit Occupancy of the Three Slot-7 Cards");
0711 
0712   ibooker.setCurrentFolder(histFolder_ + "/ECalDetail");
0713 
0714   eventMonitors.ecalOccEtDiscrepancy_ = bookEcalOccupancy("ecalOccEtDiscrepancy", "ECal Et Discrepancy Occupancy");
0715   eventMonitors.ecalOccFgDiscrepancy_ =
0716       bookEcalOccupancy("ecalOccFgDiscrepancy", "ECal FG Veto Bit Discrepancy Occupancy");
0717   eventMonitors.ecalOccLinkMasked_ = bookEcalOccupancy("ecalOccLinkMasked", "ECal Masked Links");
0718   eventMonitors.ecalOccRecdFgVB_ = bookEcalOccupancy("ecalOccRecdFgVB", "ECal FineGrain Veto Bit Occupancy at Layer1");
0719   eventMonitors.ecalOccSentAndRecd_ = bookEcalOccupancy("ecalOccSentAndRecd", "ECal TP Occupancy FULL MATCH");
0720   eventMonitors.ecalOccSentFgVB_ = bookEcalOccupancy("ecalOccSentFgVB", "ECal FineGrain Veto Bit Occupancy at TCC");
0721   eventMonitors.ecalOccSent_ = bookEcalOccupancy("ecalOccSent", "ECal TP Occupancy at TCC");
0722   eventMonitors.ecalOccTowerMasked_ = bookEcalOccupancy("ecalOccTowerMasked", "ECal Masked towers");
0723   eventMonitors.ecalTPRawEtCorrelation_ =
0724       bookEtCorrelation("ecalTPRawEtCorrelation", "Raw Et correlation TCC and Layer1;TCC Et;Layer1 Et");
0725   eventMonitors.ecalTPRawEtDiffNoMatch_ = bookEtDiff("ecalTPRawEtDiffNoMatch", "ECal Raw Et Difference Layer1 - TCC");
0726   eventMonitors.ecalTPRawEtRecd_ = bookEt("ecalTPRawEtRecd", "ECal Raw Et Layer1 Readout");
0727   eventMonitors.ecalTPRawEtSentAndRecd_ = bookEt("ecalTPRawEtMatch", "ECal Raw Et FULL MATCH");
0728   eventMonitors.ecalTPRawEtSent_ = bookEt("ecalTPRawEtSent", "ECal Raw Et TCC Readout");
0729   eventMonitors.ecalOccRecd5Bx_ = ibooker.book1D("ecalOccRecd5Bx", "Number of TPs vs BX", 5, 1, 6);
0730   eventMonitors.ecalOccRecd5BxEtWgt_ =
0731       ibooker.book1D("ecalOccRecd5BxEtWgt", "Et-weighted number of TPs vs BX", 5, 1, 6);
0732   eventMonitors.ecalOccRecdBx1_ = bookEcalOccupancy("ecalOccRecdBx1", "ECal TP Occupancy for BX1");
0733   eventMonitors.ecalOccRecdBx2_ = bookEcalOccupancy("ecalOccRecdBx2", "ECal TP Occupancy for BX2");
0734   eventMonitors.ecalOccRecdBx3_ = bookEcalOccupancy("ecalOccRecdBx3", "ECal TP Occupancy for BX3");
0735   eventMonitors.ecalOccRecdBx4_ = bookEcalOccupancy("ecalOccRecdBx4", "ECal TP Occupancy for BX4");
0736   eventMonitors.ecalOccRecdBx5_ = bookEcalOccupancy("ecalOccRecdBx5", "ECal TP Occupancy for BX5");
0737 
0738   ibooker.setCurrentFolder(histFolder_ + "/ECalDetail/TCCDebug");
0739   eventMonitors.ecalOccSentNotRecd_ =
0740       bookHcalOccupancy("ecalOccSentNotRecd", "ECal TP Occupancy sent by TCC, zero at Layer1");
0741   eventMonitors.ecalOccRecdNotSent_ =
0742       bookHcalOccupancy("ecalOccRecdNotSent", "ECal TP Occupancy received by Layer1, zero at TCC");
0743   eventMonitors.ecalOccNoMatch_ =
0744       bookHcalOccupancy("ecalOccNoMatch", "ECal TP Occupancy for TCC and Layer1 nonzero, not matching");
0745 
0746   ibooker.setCurrentFolder(histFolder_ + "/HCalDetail");
0747 
0748   eventMonitors.hcalOccEtDiscrepancy_ = bookHcalOccupancy("hcalOccEtDiscrepancy", "HCal Et Discrepancy Occupancy");
0749   eventMonitors.hcalOccLLPFbDiscrepancy_ =
0750       bookHcalOccupancy("hcalOccLLPFbDiscrepancy", "HCal LLP Feature Bit Discrepancy between Expected and Data");
0751   eventMonitors.hcalOccLLPFbExpd_ = bookHcalOccupancy("hcalOccLLPFbExpd", "HCal LLP Feature Bit Occupancy Expected");
0752   eventMonitors.hcalOccLLPFbData_ = bookHcalOccupancy("hcalOccLLPFbData", "HCal LLP Feature Bit Occupancy in Data");
0753   eventMonitors.hcalOccFg0Discrepancy_ =
0754       bookHcalOccupancy("hcalOccFg0Discrepancy", "HCal Fine Grain 0 Discrepancy between uHTR and Layer1");
0755   eventMonitors.hcalOccFg1Discrepancy_ =
0756       bookHcalOccupancy("hcalOccFg1Discrepancy", "HCal Fine Grain 1 Discrepancy between uHTR and Layer1");
0757   eventMonitors.hcalOccFg2Discrepancy_ =
0758       bookHcalOccupancy("hcalOccFg2Discrepancy", "HCal Fine Grain 2 Discrepancy between uHTR and Layer1");
0759   eventMonitors.hcalOccFg3Discrepancy_ =
0760       bookHcalOccupancy("hcalOccFg3Discrepancy", "HCal Fine Grain 3 Discrepancy between uHTR and Layer1");
0761   //eventMonitors.hcalOccFg4Discrepancy_ =
0762   //    bookHcalOccupancy("hcalOccFg4Discrepancy", "HCal Fine Grain 4 Discrepancy between uHTR and Layer1");
0763   //eventMonitors.hcalOccFg5Discrepancy_ =
0764   //    bookHcalOccupancy("hcalOccFg5Discrepancy", "HCal Fine Grain 5 Discrepancy between uHTR and Layer1");
0765   eventMonitors.hcalOccSentFg0_ = bookHcalOccupancy("hcalOccSentFg0", "HCal Fine Grain 0 Occupancy at uHTR");
0766   eventMonitors.hcalOccSentFg1_ = bookHcalOccupancy("hcalOccSentFg1", "HCal Fine Grain 1 Occupancy at uHTR");
0767   eventMonitors.hcalOccSentFg2_ = bookHcalOccupancy("hcalOccSentFg2", "HCal Fine Grain 2 Occupancy at uHTR");
0768   eventMonitors.hcalOccSentFg3_ = bookHcalOccupancy("hcalOccSentFg3", "HCal Fine Grain 3 Occupancy at uHTR");
0769   eventMonitors.hcalOccSentFg4_ = bookHcalOccupancy("hcalOccSentFg4", "HCal Fine Grain 4 Occupancy at uHTR");
0770   eventMonitors.hcalOccSentFg5_ = bookHcalOccupancy("hcalOccSentFg5", "HCal Fine Grain 5 Occupancy at uHTR");
0771   eventMonitors.hcalOccRecdFg0_ = bookHcalOccupancy("hcalOccRecdFg0", "HCal Fine Grain 0 Occupancy at Layer1");
0772   eventMonitors.hcalOccRecdFg1_ = bookHcalOccupancy("hcalOccRecdFg1", "HCal Fine Grain 1 Occupancy at Layer1");
0773   eventMonitors.hcalOccRecdFg2_ = bookHcalOccupancy("hcalOccRecdFg2", "HCal Fine Grain 2 Occupancy at Layer1");
0774   eventMonitors.hcalOccRecdFg3_ = bookHcalOccupancy("hcalOccRecdFg3", "HCal Fine Grain 3 Occupancy at Layer1");
0775   eventMonitors.hcalOccRecdFg4_ = bookHcalOccupancy("hcalOccRecdFg4", "HCal Fine Grain 4 Occupancy at Layer1");
0776   eventMonitors.hcalOccRecdFg5_ = bookHcalOccupancy("hcalOccRecdFg5", "HCal Fine Grain 5 Occupancy at Layer1");
0777   eventMonitors.hcalOccLinkMasked_ = bookHcalOccupancy("hcalOccLinkMasked", "HCal Masked Links");
0778   eventMonitors.hcalOccSentAndRecd_ = bookHcalOccupancy("hcalOccSentAndRecd", "HCal TP Occupancy FULL MATCH");
0779   eventMonitors.hcalOccSent_ = bookHcalOccupancy("hcalOccSent", "HCal TP Occupancy at uHTR");
0780   eventMonitors.hcalOccTowerMasked_ = bookHcalOccupancy("hcalOccTowerMasked", "HCal Masked towers");
0781   eventMonitors.hcalTPRawEtCorrelationHBHE_ =
0782       bookEtCorrelation("hcalTPRawEtCorrelationHBHE", "HBHE Raw Et correlation uHTR and Layer1;uHTR Et;Layer1 Et");
0783   eventMonitors.hcalTPRawEtCorrelationHF_ =
0784       bookEtCorrelation("hcalTPRawEtCorrelationHF", "HF Raw Et correlation uHTR and Layer1;uHTR Et;Layer1 Et");
0785   eventMonitors.hcalTPRawEtDiffNoMatch_ = bookEtDiff("hcalTPRawEtDiffNoMatch", "HCal Raw Et Difference Layer1 - uHTR");
0786   eventMonitors.hcalTPRawEtRecd_ = bookEt("hcalTPRawEtRecd", "HCal Raw Et Layer1 Readout");
0787   eventMonitors.hcalTPRawEtSentAndRecd_ = bookEt("hcalTPRawEtMatch", "HCal Raw Et FULL MATCH");
0788   eventMonitors.hcalTPRawEtSent_ = bookEt("hcalTPRawEtSent", "HCal Raw Et uHTR Readout");
0789 
0790   ibooker.setCurrentFolder(histFolder_ + "/HCalDetail/uHTRDebug");
0791   eventMonitors.hcalOccSentNotRecd_ =
0792       bookHcalOccupancy("hcalOccSentNotRecd", "HCal TP Occupancy sent by uHTR, zero at Layer1");
0793   eventMonitors.hcalOccRecdNotSent_ =
0794       bookHcalOccupancy("hcalOccRecdNotSent", "HCal TP Occupancy received by Layer1, zero at uHTR");
0795   eventMonitors.hcalOccNoMatch_ =
0796       bookHcalOccupancy("hcalOccNoMatch", "HCal TP Occupancy for uHTR and Layer1 nonzero, not matching");
0797 
0798   ibooker.setCurrentFolder(histFolder_ + "/MismatchDetail");
0799 
0800   const int nMismatchTypes = 4;
0801   eventMonitors.last20Mismatches_ = ibooker.book2D("last20Mismatches",
0802                                                    "Log of last 20 mismatches (use json tool to copy/paste)",
0803                                                    nMismatchTypes,
0804                                                    0,
0805                                                    nMismatchTypes,
0806                                                    20,
0807                                                    0,
0808                                                    20);
0809   eventMonitors.last20Mismatches_->setBinLabel(1, "Ecal TP Et Mismatch");
0810   eventMonitors.last20Mismatches_->setBinLabel(2, "Ecal TP Fine Grain Bit Mismatch");
0811   eventMonitors.last20Mismatches_->setBinLabel(3, "Hcal TP Et Mismatch");
0812   eventMonitors.last20Mismatches_->setBinLabel(4, "Hcal TP Feature Bit Mismatch");
0813 
0814   const int nLumis = 2000;
0815   eventMonitors.ecalLinkErrorByLumi_ = ibooker.book1D(
0816       "ecalLinkErrorByLumi", "Link error counts per lumi section for ECAL;LumiSection;Counts", nLumis, .5, nLumis + 0.5);
0817   eventMonitors.ecalMismatchByLumi_ = ibooker.book1D(
0818       "ecalMismatchByLumi", "Mismatch counts per lumi section for ECAL;LumiSection;Counts", nLumis, .5, nLumis + 0.5);
0819   eventMonitors.hcalLinkErrorByLumi_ = ibooker.book1D(
0820       "hcalLinkErrorByLumi", "Link error counts per lumi section for HCAL;LumiSection;Counts", nLumis, .5, nLumis + 0.5);
0821   eventMonitors.hcalMismatchByLumi_ = ibooker.book1D(
0822       "hcalMismatchByLumi", "Mismatch counts per lumi section for HCAL;LumiSection;Counts", nLumis, .5, nLumis + 0.5);
0823 
0824   eventMonitors.ECALmismatchesPerBx_ =
0825       ibooker.book1D("ECALmismatchesPerBx", "Mismatch counts per bunch crossing for ECAL", 3564, -.5, 3563.5);
0826   eventMonitors.HBHEmismatchesPerBx_ =
0827       ibooker.book1D("HBHEmismatchesPerBx", "Mismatch counts per bunch crossing for HBHE", 3564, -.5, 3563.5);
0828   eventMonitors.HFmismatchesPerBx_ =
0829       ibooker.book1D("HFmismatchesPerBx", "Mismatch counts per bunch crossing for HF", 3564, -.5, 3563.5);
0830 
0831   eventMonitors.maxEvtLinkErrorsByLumiECAL_ =
0832       ibooker.book1D("maxEvtLinkErrorsByLumiECAL",
0833                      "Max number of single-event ECAL link errors per lumi section;LumiSection;Counts",
0834                      nLumis,
0835                      .5,
0836                      nLumis + 0.5);
0837   eventMonitors.maxEvtLinkErrorsByLumiHCAL_ =
0838       ibooker.book1D("maxEvtLinkErrorsByLumiHCAL",
0839                      "Max number of single-event HCAL link errors per lumi section;LumiSection;Counts",
0840                      nLumis,
0841                      .5,
0842                      nLumis + 0.5);
0843 
0844   eventMonitors.maxEvtMismatchByLumiECAL_ =
0845       ibooker.book1D("maxEvtMismatchByLumiECAL",
0846                      "Max number of single-event ECAL discrepancies per lumi section;LumiSection;Counts",
0847                      nLumis,
0848                      .5,
0849                      nLumis + 0.5);
0850   eventMonitors.maxEvtMismatchByLumiHCAL_ =
0851       ibooker.book1D("maxEvtMismatchByLumiHCAL",
0852                      "Max number of single-event HCAL discrepancies per lumi section;LumiSection;Counts",
0853                      nLumis,
0854                      .5,
0855                      nLumis + 0.5);
0856 
0857   ibooker.setCurrentFolder(histFolder_);
0858   eventMonitors.maxEvtLinkErrorsByLumi_ =
0859       ibooker.book1D("maxEvtLinkErrorsByLumi",
0860                      "Max number of single-event link errors per lumi section;LumiSection;Counts",
0861                      nLumis,
0862                      .5,
0863                      nLumis + 0.5);
0864   eventMonitors.maxEvtMismatchByLumi_ =
0865       ibooker.book1D("maxEvtMismatchByLumi",
0866                      "Max number of single-event discrepancies per lumi section;LumiSection;Counts",
0867                      nLumis,
0868                      .5,
0869                      nLumis + 0.5);
0870 
0871   ibooker.setCurrentFolder(histFolder_ + "/AMC13ErrorCounters");
0872   eventMonitors.bxidErrors_ =
0873       ibooker.book1D("bxidErrors", "bxid mismatch between AMC13 and CTP Cards;Layer1 Phi;Counts", 18, -.5, 17.5);
0874   eventMonitors.l1idErrors_ =
0875       ibooker.book1D("l1idErrors", "l1id mismatch between AMC13 and CTP Cards;Layer1 Phi;Counts", 18, -.5, 17.5);
0876   eventMonitors.orbitErrors_ =
0877       ibooker.book1D("orbitErrors", "orbit mismatch between AMC13 and CTP Cards;Layer1 Phi;Counts", 18, -.5, 17.5);
0878 }
0879 
0880 void L1TStage2CaloLayer1::streamEndLuminosityBlockSummary(
0881     edm::StreamID theStreamID,
0882     edm::LuminosityBlock const& theLumiBlock,
0883     edm::EventSetup const& theEventSetup,
0884     CaloL1Information::perLumiBlockMonitoringInformation* lumiMonitoringInformation) const {
0885   auto theStreamCache = streamCache(theStreamID);
0886 
0887   lumiMonitoringInformation->lumiNumMaxEvtLinkErrorsECAL =
0888       std::max(lumiMonitoringInformation->lumiNumMaxEvtLinkErrorsECAL, theStreamCache->streamNumMaxEvtLinkErrorsECAL);
0889   lumiMonitoringInformation->lumiNumMaxEvtLinkErrorsHCAL =
0890       std::max(lumiMonitoringInformation->lumiNumMaxEvtLinkErrorsHCAL, theStreamCache->streamNumMaxEvtLinkErrorsHCAL);
0891   lumiMonitoringInformation->lumiNumMaxEvtLinkErrors =
0892       std::max(lumiMonitoringInformation->lumiNumMaxEvtLinkErrors, theStreamCache->streamNumMaxEvtLinkErrors);
0893 
0894   lumiMonitoringInformation->lumiNumMaxEvtMismatchECAL =
0895       std::max(lumiMonitoringInformation->lumiNumMaxEvtMismatchECAL, theStreamCache->streamNumMaxEvtMismatchECAL);
0896   lumiMonitoringInformation->lumiNumMaxEvtMismatchHCAL =
0897       std::max(lumiMonitoringInformation->lumiNumMaxEvtMismatchHCAL, theStreamCache->streamNumMaxEvtMismatchHCAL);
0898   lumiMonitoringInformation->lumiNumMaxEvtMismatch =
0899       std::max(lumiMonitoringInformation->lumiNumMaxEvtMismatch, theStreamCache->streamNumMaxEvtMismatch);
0900 
0901   //reset the stream cache here, since we are done with the luminosity block in this stream
0902   //We don't want the stream to be comparing to previous values in the next lumi-block
0903   theStreamCache->streamNumMaxEvtLinkErrorsECAL = 0;
0904   theStreamCache->streamNumMaxEvtLinkErrorsHCAL = 0;
0905   theStreamCache->streamNumMaxEvtLinkErrors = 0;
0906 
0907   theStreamCache->streamNumMaxEvtMismatchECAL = 0;
0908   theStreamCache->streamNumMaxEvtMismatchHCAL = 0;
0909   theStreamCache->streamNumMaxEvtMismatch = 0;
0910 }
0911 
0912 void L1TStage2CaloLayer1::globalEndLuminosityBlockSummary(
0913     edm::LuminosityBlock const& theLumiBlock,
0914     edm::EventSetup const& theEventSetup,
0915     CaloL1Information::perLumiBlockMonitoringInformation* lumiMonitoringInformation) const {
0916   auto theRunCache = runCache(theLumiBlock.getRun().index());
0917   //read the cache out to the relevant Luminosity histogram bin
0918   theRunCache->maxEvtLinkErrorsByLumiECAL_->Fill(theLumiBlock.luminosityBlock(),
0919                                                  lumiMonitoringInformation->lumiNumMaxEvtLinkErrorsECAL);
0920   theRunCache->maxEvtLinkErrorsByLumiHCAL_->Fill(theLumiBlock.luminosityBlock(),
0921                                                  lumiMonitoringInformation->lumiNumMaxEvtLinkErrorsHCAL);
0922   theRunCache->maxEvtLinkErrorsByLumi_->Fill(theLumiBlock.luminosityBlock(),
0923                                              lumiMonitoringInformation->lumiNumMaxEvtLinkErrors);
0924 
0925   theRunCache->maxEvtMismatchByLumiECAL_->Fill(theLumiBlock.luminosityBlock(),
0926                                                lumiMonitoringInformation->lumiNumMaxEvtMismatchECAL);
0927   theRunCache->maxEvtMismatchByLumiHCAL_->Fill(theLumiBlock.luminosityBlock(),
0928                                                lumiMonitoringInformation->lumiNumMaxEvtMismatchHCAL);
0929   theRunCache->maxEvtMismatchByLumi_->Fill(theLumiBlock.luminosityBlock(),
0930                                            lumiMonitoringInformation->lumiNumMaxEvtMismatch);
0931 
0932   // Simple way to embed current lumi to auto-scale axis limits in render plugin
0933   theRunCache->ecalLinkErrorByLumi_->setBinContent(0, theLumiBlock.luminosityBlock());
0934   theRunCache->ecalMismatchByLumi_->setBinContent(0, theLumiBlock.luminosityBlock());
0935   theRunCache->hcalLinkErrorByLumi_->setBinContent(0, theLumiBlock.luminosityBlock());
0936   theRunCache->hcalMismatchByLumi_->setBinContent(0, theLumiBlock.luminosityBlock());
0937   theRunCache->maxEvtLinkErrorsByLumiECAL_->setBinContent(0, theLumiBlock.luminosityBlock());
0938   theRunCache->maxEvtLinkErrorsByLumiHCAL_->setBinContent(0, theLumiBlock.luminosityBlock());
0939   theRunCache->maxEvtLinkErrorsByLumi_->setBinContent(0, theLumiBlock.luminosityBlock());
0940   theRunCache->maxEvtMismatchByLumiECAL_->setBinContent(0, theLumiBlock.luminosityBlock());
0941   theRunCache->maxEvtMismatchByLumiHCAL_->setBinContent(0, theLumiBlock.luminosityBlock());
0942   theRunCache->maxEvtMismatchByLumi_->setBinContent(0, theLumiBlock.luminosityBlock());
0943 }
0944 
0945 //returns true if the new candidate mismatch is from a later mismatch than the comparison mismatch
0946 // based on Run, Lumisection, and event number
0947 //false otherwise
0948 bool L1TStage2CaloLayer1::isLaterMismatch(
0949     std::tuple<edm::RunID, edm::LuminosityBlockID, edm::EventID, std::vector<int>>& candidateMismatch,
0950     std::tuple<edm::RunID, edm::LuminosityBlockID, edm::EventID, std::vector<int>>& comparisonMismatch) const {
0951   //check the run. If the run ID of the candidate mismatch is less than the run ID of the comparison mismatch, it is earlier,
0952   if (std::get<0>(candidateMismatch) < std::get<0>(comparisonMismatch))
0953     return false;
0954   //if it is greater, then it is a later mismatch
0955   else if (std::get<0>(candidateMismatch) > std::get<0>(comparisonMismatch))
0956     return true;
0957   //if it is even, then we need to repeat this comparison on the luminosity block
0958   else {
0959     if (std::get<1>(candidateMismatch) < std::get<1>(comparisonMismatch))
0960       return false;  //if the lumi block is less, it is an earlier mismatch
0961     else if (std::get<1>(candidateMismatch) > std::get<1>(comparisonMismatch))
0962       return true;  // if the lumi block is greater, than it is a later mismatch
0963     // if these are equivalent, then we repeat the comparison on the event
0964     else {
0965       if (std::get<2>(candidateMismatch) < std::get<2>(comparisonMismatch))
0966         return false;
0967       else
0968         return true;  //in the case of even events here, we consider the event later.
0969     }
0970   }
0971   //default return, should never be called.
0972   return false;
0973 }
0974 
0975 //binary search like algorithm for trying to find the appropriate place to put our new mismatch
0976 //will find an interger we add to the iterator to get the proper location to find the insertion location
0977 //will return -1 if the mismatch should not be inserted into the list
0978 int L1TStage2CaloLayer1::findIndex(
0979     std::tuple<edm::RunID, edm::LuminosityBlockID, edm::EventID, std::vector<int>> candidateMismatch,
0980     std::vector<std::tuple<edm::RunID, edm::LuminosityBlockID, edm::EventID, std::vector<int>>> comparisonList,
0981     int lowerIndexToSearch,
0982     int upperIndexToSearch) const {
0983   //Start by getting the spot in the the vector to start searching
0984   int searchLocation = ((upperIndexToSearch + lowerIndexToSearch) / 2);
0985   auto searchIterator = comparisonList.begin() + searchLocation;
0986   //Multiple possible cases:
0987   //case one. Greater than the search element
0988   if (this->isLaterMismatch(candidateMismatch, *searchIterator)) {
0989     //subcase one, there exists a mismatch to its right,
0990     if (searchIterator + 1 != comparisonList.end()) {
0991       //subsubcase one, the candidate is earlier than the one to the right, insert this element between these two
0992       if (not this->isLaterMismatch(candidateMismatch, *(searchIterator + 1))) {
0993         return searchLocation + 1;  //vector insert inserts *before* the position, so we return +1
0994       }
0995       //subsubcase two, the candidate is later than it, in this case we refine the search area.
0996       else {
0997         return this->findIndex(candidateMismatch, comparisonList, searchLocation + 1, upperIndexToSearch);
0998       }
0999     }
1000     //subcase two, there exists no mismatch to it's right (end of the vector), in which case this is the latest mismatch
1001     else {
1002       return searchLocation +
1003              1;  //vector insert inserts *before* the position, so we return +1 (should be synonymous with .end())
1004     }
1005   }
1006   //case two. we are earlier than the current mismatch
1007   else {
1008     //subcase one, these exists a mismatch to its left,
1009     if (searchIterator != comparisonList.begin()) {
1010       //subsubcase one, the candidate mismatch is earlier than the one to the left. in this case we refine the search area
1011       if (not this->isLaterMismatch(candidateMismatch, *(searchIterator - 1))) {
1012         return this->findIndex(candidateMismatch, comparisonList, lowerIndexToSearch, searchLocation - 1);
1013       }
1014       //subsubcase two the candidate is later than than the one to it's left, in which case, we insert the element between these two.
1015       else {
1016         return searchLocation;  //vector insert inserts *before* the position, so we return without addition here.
1017       }
1018     }
1019     //subcase two, there exists no mismatch to its left (beginning of vector), in which case this is the is earlier than all mismatches on the list.
1020     else {
1021       //subsubcase one. It is possible we still insert this if the capacity of the vector is below 20.
1022       //this should probably be rewritten to have the capacity reserved before-hand so that 20 is not potentially a hard coded magic number
1023       //defined by what we know to be true about a non-obvious data type elsewhere.
1024       if (comparisonList.size() < 20) {
1025         return searchLocation;
1026       }
1027       //subsubcase two. The size is 20 or above, and we are earlier than all of them. This should not be inserted into the list.
1028       else {
1029         return -1;
1030       }
1031     }
1032   }
1033   //default return. Should never be called
1034   return -1;
1035 }
1036 
1037 //will shuffle the candidate mismatch list into the comparison mismatch list.
1038 void L1TStage2CaloLayer1::mergeMismatchVectors(
1039     std::vector<std::tuple<edm::RunID, edm::LuminosityBlockID, edm::EventID, std::vector<int>>>& candidateMismatchList,
1040     std::vector<std::tuple<edm::RunID, edm::LuminosityBlockID, edm::EventID, std::vector<int>>>& comparisonMismatchList)
1041     const {
1042   //okay now we loop over our candidate mismatches
1043   for (auto candidateIterator = candidateMismatchList.begin(); candidateIterator != candidateMismatchList.end();
1044        ++candidateIterator) {
1045     int insertionIndex = this->findIndex(*candidateIterator, comparisonMismatchList, 0, comparisonMismatchList.size());
1046     if (insertionIndex < 0) {
1047       continue;  //if we didn't find anywhere to put this mismatch, we move on
1048     }
1049     auto insertionIterator = comparisonMismatchList.begin() + insertionIndex;
1050     comparisonMismatchList.insert(insertionIterator, *candidateIterator);
1051     //now if we have more than 20 mismatches in the list, we erase the earliest one (the beginning)
1052     //this should probably be rewritten to have the capacity reserved before-hand so that 20 is not potentially a hard coded magic number
1053     //defined by what we know to be true about a non-obvious data type elsewhere.
1054     if (comparisonMismatchList.size() > 20) {
1055       comparisonMismatchList.erase(comparisonMismatchList.begin());
1056     }
1057   }
1058 }
1059 
1060 void L1TStage2CaloLayer1::streamEndRunSummary(
1061     edm::StreamID theStreamID,
1062     edm::Run const& theRun,
1063     edm::EventSetup const& theEventSetup,
1064     CaloL1Information::perRunSummaryMonitoringInformation* theRunSummaryMonitoringInformation) const {
1065   auto theStreamCache = streamCache(theStreamID);
1066 
1067   if (theRunSummaryMonitoringInformation->runMismatchList.empty())
1068     theRunSummaryMonitoringInformation->runMismatchList = theStreamCache->streamMismatchList;
1069   else
1070     this->mergeMismatchVectors(theStreamCache->streamMismatchList, theRunSummaryMonitoringInformation->runMismatchList);
1071 
1072   //clear the stream cache so that the next run does not try to compare to the current one.
1073   theStreamCache->streamMismatchList.clear();
1074 }
1075 
1076 void L1TStage2CaloLayer1::globalEndRunSummary(
1077     edm::Run const& theRun,
1078     edm::EventSetup const& theEventSetup,
1079     CaloL1Information::perRunSummaryMonitoringInformation* theRunSummaryInformation) const {
1080   //The mismatch vector should be properly ordered and sized by other functions, so all we need to do is read it into a monitoring element.
1081   auto theRunCache = runCache(theRun.index());
1082   //loop over the accepted mismatch list in the run, and based on the run/lumi/event create the bin label, the mismatch type will define which bin gets content.
1083   int ibin = 1;
1084   for (auto mismatchIterator = theRunSummaryInformation->runMismatchList.begin();
1085        mismatchIterator != theRunSummaryInformation->runMismatchList.end();
1086        ++mismatchIterator) {
1087     std::string binLabel = std::to_string(std::get<0>(*mismatchIterator).run()) + ":" +
1088                            std::to_string(std::get<1>(*mismatchIterator).luminosityBlock()) + ":" +
1089                            std::to_string(std::get<2>(*mismatchIterator).event());
1090     theRunCache->last20Mismatches_->setBinLabel(ibin, binLabel, 2);
1091     //Get the vector of mismatches for this particular event and iterate through it.
1092     //Set the bin content to 1 for each type of mismatch seen
1093     std::vector<int> mismatchTypeVector = std::get<3>(*mismatchIterator);
1094     for (auto mismatchTypeIterator = mismatchTypeVector.begin(); mismatchTypeIterator != mismatchTypeVector.end();
1095          ++mismatchTypeIterator) {
1096       theRunCache->last20Mismatches_->setBinContent(*mismatchTypeIterator + 1, ibin, 1);
1097     }
1098     ++ibin;
1099   }
1100   //remove the remaining empty string labels to prevent overlap
1101   for (int emptyBinIndex = ibin; emptyBinIndex <= 20; ++emptyBinIndex) {
1102     theRunCache->last20Mismatches_->setBinLabel(emptyBinIndex, std::to_string(emptyBinIndex), 2);
1103   }
1104 }