Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2022-06-20 00:21:32

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 L1TStage2CaloLayer1::L1TStage2CaloLayer1(const edm::ParameterSet& ps)
0023     : ecalTPSourceRecd_(consumes<EcalTrigPrimDigiCollection>(ps.getParameter<edm::InputTag>("ecalTPSourceRecd"))),
0024       ecalTPSourceRecdLabel_(ps.getParameter<edm::InputTag>("ecalTPSourceRecd").label()),
0025       ecalTPSourceRecdBx1_(consumes<EcalTrigPrimDigiCollection>(ps.getParameter<edm::InputTag>("ecalTPSourceRecdBx1"))),
0026       ecalTPSourceRecdBx1Label_(ps.getParameter<edm::InputTag>("ecalTPSourceRecdBx1").label()),
0027       ecalTPSourceRecdBx2_(consumes<EcalTrigPrimDigiCollection>(ps.getParameter<edm::InputTag>("ecalTPSourceRecdBx2"))),
0028       ecalTPSourceRecdBx2Label_(ps.getParameter<edm::InputTag>("ecalTPSourceRecdBx2").label()),
0029       ecalTPSourceRecdBx3_(consumes<EcalTrigPrimDigiCollection>(ps.getParameter<edm::InputTag>("ecalTPSourceRecdBx3"))),
0030       ecalTPSourceRecdBx3Label_(ps.getParameter<edm::InputTag>("ecalTPSourceRecdBx3").label()),
0031       ecalTPSourceRecdBx4_(consumes<EcalTrigPrimDigiCollection>(ps.getParameter<edm::InputTag>("ecalTPSourceRecdBx4"))),
0032       ecalTPSourceRecdBx4Label_(ps.getParameter<edm::InputTag>("ecalTPSourceRecdBx4").label()),
0033       ecalTPSourceRecdBx5_(consumes<EcalTrigPrimDigiCollection>(ps.getParameter<edm::InputTag>("ecalTPSourceRecdBx5"))),
0034       ecalTPSourceRecdBx5Label_(ps.getParameter<edm::InputTag>("ecalTPSourceRecdBx5").label()),
0035       hcalTPSourceRecd_(consumes<HcalTrigPrimDigiCollection>(ps.getParameter<edm::InputTag>("hcalTPSourceRecd"))),
0036       hcalTPSourceRecdLabel_(ps.getParameter<edm::InputTag>("hcalTPSourceRecd").label()),
0037       ecalTPSourceSent_(consumes<EcalTrigPrimDigiCollection>(ps.getParameter<edm::InputTag>("ecalTPSourceSent"))),
0038       ecalTPSourceSentLabel_(ps.getParameter<edm::InputTag>("ecalTPSourceSent").label()),
0039       hcalTPSourceSent_(consumes<HcalTrigPrimDigiCollection>(ps.getParameter<edm::InputTag>("hcalTPSourceSent"))),
0040       hcalTPSourceSentLabel_(ps.getParameter<edm::InputTag>("hcalTPSourceSent").label()),
0041       fedRawData_(consumes<FEDRawDataCollection>(ps.getParameter<edm::InputTag>("fedRawDataLabel"))),
0042       histFolder_(ps.getParameter<std::string>("histFolder")),
0043       tpFillThreshold_(ps.getUntrackedParameter<int>("etDistributionsFillThreshold", 0)),
0044       tpFillThreshold5Bx_(ps.getUntrackedParameter<int>("etDistributionsFillThreshold5Bx", 4)),
0045       ignoreHFfbs_(ps.getUntrackedParameter<bool>("ignoreHFfbs", false)) {}
0046 
0047 L1TStage2CaloLayer1::~L1TStage2CaloLayer1() {}
0048 
0049 void L1TStage2CaloLayer1::dqmAnalyze(const edm::Event& event,
0050                                      const edm::EventSetup& es,
0051                                      const CaloL1Information::monitoringDataHolder& eventMonitors) const {
0052   //This must be moved here compared to previous version to avoid a const function modifying a class variable
0053   //This unfortunately means that the variable is local and reallocated per event
0054   std::vector<std::pair<EcalTriggerPrimitiveDigi, EcalTriggerPrimitiveDigi>> ecalTPSentRecd_;
0055   std::vector<std::pair<HcalTriggerPrimitiveDigi, HcalTriggerPrimitiveDigi>> hcalTPSentRecd_;
0056   ecalTPSentRecd_.reserve(28 * 2 * 72);
0057   hcalTPSentRecd_.reserve(41 * 2 * 72);
0058 
0059   // Monitorables stored in Layer 1 raw data but
0060   // not accessible from existing persistent data formats
0061   edm::Handle<FEDRawDataCollection> fedRawDataCollection;
0062   event.getByToken(fedRawData_, fedRawDataCollection);
0063   bool caloLayer1OutOfRun{true};
0064   if (fedRawDataCollection.isValid()) {
0065     caloLayer1OutOfRun = false;
0066     for (int iFed = 1354; iFed < 1360; iFed += 2) {
0067       const FEDRawData& fedRawData = fedRawDataCollection->FEDData(iFed);
0068       if (fedRawData.size() == 0) {
0069         caloLayer1OutOfRun = true;
0070         continue;  // In case one of 3 layer 1 FEDs not in
0071       }
0072       const uint64_t* fedRawDataArray = (const uint64_t*)fedRawData.data();
0073       UCTDAQRawData daqData(fedRawDataArray);
0074       for (uint32_t i = 0; i < daqData.nAMCs(); i++) {
0075         UCTAMCRawData amcData(daqData.amcPayload(i));
0076         int lPhi = amcData.layer1Phi();
0077         if (daqData.BXID() != amcData.BXID()) {
0078           eventMonitors.bxidErrors_->Fill(lPhi);
0079         }
0080         if (daqData.L1ID() != amcData.L1ID()) {
0081           eventMonitors.l1idErrors_->Fill(lPhi);
0082         }
0083         // AMC payload header has 16 bit orbit number, AMC13 header is full 32
0084         if ((daqData.orbitNumber() & 0xFFFF) != amcData.orbitNo()) {
0085           eventMonitors.orbitErrors_->Fill(lPhi);
0086         }
0087       }
0088     }
0089   }
0090 
0091   edm::Handle<EcalTrigPrimDigiCollection> ecalTPsSent;
0092   event.getByToken(ecalTPSourceSent_, ecalTPsSent);
0093   edm::Handle<EcalTrigPrimDigiCollection> ecalTPsRecd;
0094   event.getByToken(ecalTPSourceRecd_, ecalTPsRecd);
0095 
0096   ecalTPSentRecd_.clear();
0097 
0098   ComparisonHelper::zip(ecalTPsSent->begin(),
0099                         ecalTPsSent->end(),
0100                         ecalTPsRecd->begin(),
0101                         ecalTPsRecd->end(),
0102                         std::inserter(ecalTPSentRecd_, ecalTPSentRecd_.begin()),
0103                         EcalTrigPrimDigiCollection::key_compare());
0104 
0105   int nEcalLinkErrors{0};
0106   int nEcalMismatch{0};
0107 
0108   for (const auto& tpPair : ecalTPSentRecd_) {
0109     auto sentTp = tpPair.first;
0110     if (sentTp.compressedEt() < 0) {
0111       // ECal zero-suppresses digis, and a default-constructed
0112       // digi has et=-1 apparently, but we know it should be zero
0113       EcalTriggerPrimitiveSample sample(0);
0114       EcalTriggerPrimitiveDigi tpg(sentTp.id());
0115       tpg.setSize(1);
0116       tpg.setSample(0, sample);
0117       swap(sentTp, tpg);
0118     }
0119     const auto& recdTp = tpPair.second;
0120     const int ieta = sentTp.id().ieta();
0121     const int iphi = sentTp.id().iphi();
0122     const bool towerMasked = recdTp.sample(0).raw() & (1 << 13);
0123     const bool linkMasked = recdTp.sample(0).raw() & (1 << 14);
0124     const bool linkError = recdTp.sample(0).raw() & (1 << 15);
0125 
0126     // Link status bits from layer 1
0127     if (towerMasked) {
0128       eventMonitors.ecalOccTowerMasked_->Fill(ieta, iphi);
0129     }
0130     if (linkMasked) {
0131       eventMonitors.ecalOccLinkMasked_->Fill(ieta, iphi);
0132     }
0133 
0134     if (sentTp.compressedEt() > tpFillThreshold_) {
0135       eventMonitors.ecalTPRawEtSent_->Fill(sentTp.compressedEt());
0136       eventMonitors.ecalOccSent_->Fill(ieta, iphi);
0137     }
0138     if (sentTp.fineGrain() == 1) {
0139       eventMonitors.ecalOccSentFgVB_->Fill(ieta, iphi);
0140     }
0141 
0142     if (towerMasked || caloLayer1OutOfRun) {
0143       // Do not compare if we have a mask applied
0144       continue;
0145     }
0146 
0147     if (linkError) {
0148       eventMonitors.ecalLinkError_->Fill(ieta, iphi);
0149       eventMonitors.ecalLinkErrorByLumi_->Fill(event.id().luminosityBlock());
0150       nEcalLinkErrors++;
0151       // Don't compare anymore, we already know its bad
0152       continue;
0153     }
0154 
0155     eventMonitors.ecalTPRawEtCorrelation_->Fill(sentTp.compressedEt(), recdTp.compressedEt());
0156 
0157     if (recdTp.compressedEt() > tpFillThreshold_) {
0158       eventMonitors.ecalTPRawEtRecd_->Fill(recdTp.compressedEt());
0159       eventMonitors.ecalOccupancy_->Fill(ieta, iphi);
0160       eventMonitors.ecalOccRecdEtWgt_->Fill(ieta, iphi, recdTp.compressedEt());
0161     }
0162     if (recdTp.fineGrain() == 1) {
0163       eventMonitors.ecalOccRecdFgVB_->Fill(ieta, iphi);
0164     }
0165 
0166     // Now for comparisons
0167 
0168     const bool EetAgreement = sentTp.compressedEt() == recdTp.compressedEt();
0169     const bool EfbAgreement = sentTp.fineGrain() == recdTp.fineGrain();
0170     if (EetAgreement && EfbAgreement) {
0171       // Full match
0172       if (sentTp.compressedEt() > tpFillThreshold_) {
0173         eventMonitors.ecalOccSentAndRecd_->Fill(ieta, iphi);
0174         eventMonitors.ecalTPRawEtSentAndRecd_->Fill(sentTp.compressedEt());
0175       }
0176     } else {
0177       // There is some issue
0178       eventMonitors.ecalDiscrepancy_->Fill(ieta, iphi);
0179       eventMonitors.ecalMismatchByLumi_->Fill(event.id().luminosityBlock());
0180       eventMonitors.ECALmismatchesPerBx_->Fill(event.bunchCrossing());
0181       nEcalMismatch++;
0182 
0183       if (not EetAgreement) {
0184         eventMonitors.ecalOccEtDiscrepancy_->Fill(ieta, iphi);
0185         eventMonitors.ecalTPRawEtDiffNoMatch_->Fill(recdTp.compressedEt() - sentTp.compressedEt());
0186         updateMismatch(event, 0, streamCache(event.streamID())->streamMismatchList);
0187 
0188         if (sentTp.compressedEt() == 0)
0189           eventMonitors.ecalOccRecdNotSent_->Fill(ieta, iphi);
0190         else if (recdTp.compressedEt() == 0)
0191           eventMonitors.ecalOccSentNotRecd_->Fill(ieta, iphi);
0192         else
0193           eventMonitors.ecalOccNoMatch_->Fill(ieta, iphi);
0194       }
0195       if (not EfbAgreement) {
0196         // occ for fine grain mismatch
0197         eventMonitors.ecalOccFgDiscrepancy_->Fill(ieta, iphi);
0198         updateMismatch(event, 1, streamCache(event.streamID())->streamMismatchList);
0199       }
0200     }
0201   }
0202 
0203   if (nEcalLinkErrors > streamCache(event.streamID())->streamNumMaxEvtLinkErrorsECAL)
0204     streamCache(event.streamID())->streamNumMaxEvtLinkErrorsECAL = nEcalLinkErrors;
0205 
0206   if (nEcalMismatch > streamCache(event.streamID())->streamNumMaxEvtMismatchECAL)
0207     streamCache(event.streamID())->streamNumMaxEvtMismatchECAL = nEcalMismatch;
0208 
0209   edm::Handle<EcalTrigPrimDigiCollection> ecalTPsRecdBx1;
0210   event.getByToken(ecalTPSourceRecdBx1_, ecalTPsRecdBx1);
0211   edm::Handle<EcalTrigPrimDigiCollection> ecalTPsRecdBx2;
0212   event.getByToken(ecalTPSourceRecdBx2_, ecalTPsRecdBx2);
0213   edm::Handle<EcalTrigPrimDigiCollection> ecalTPsRecdBx3;
0214   event.getByToken(ecalTPSourceRecdBx3_, ecalTPsRecdBx3);
0215   edm::Handle<EcalTrigPrimDigiCollection> ecalTPsRecdBx4;
0216   event.getByToken(ecalTPSourceRecdBx4_, ecalTPsRecdBx4);
0217   edm::Handle<EcalTrigPrimDigiCollection> ecalTPsRecdBx5;
0218   event.getByToken(ecalTPSourceRecdBx5_, ecalTPsRecdBx5);
0219 
0220   for (const auto& tp : (*ecalTPsRecdBx1)) {
0221     if (tp.compressedEt() > tpFillThreshold5Bx_) {
0222       const int ieta = tp.id().ieta();
0223       const int iphi = tp.id().iphi();
0224       const bool towerMasked = tp.sample(0).raw() & (1 << 13);
0225       const bool linkError = tp.sample(0).raw() & (1 << 15);
0226       if (towerMasked || caloLayer1OutOfRun || linkError) {
0227         continue;
0228       }
0229       eventMonitors.ecalOccRecdBx1_->Fill(ieta, iphi);
0230       eventMonitors.ecalOccRecd5Bx_->Fill(1);
0231       eventMonitors.ecalOccRecd5BxEtWgt_->Fill(1, tp.compressedEt());
0232     }
0233   }
0234   for (const auto& tp : (*ecalTPsRecdBx2)) {
0235     if (tp.compressedEt() > tpFillThreshold5Bx_) {
0236       const int ieta = tp.id().ieta();
0237       const int iphi = tp.id().iphi();
0238       const bool towerMasked = tp.sample(0).raw() & (1 << 13);
0239       const bool linkError = tp.sample(0).raw() & (1 << 15);
0240       if (towerMasked || caloLayer1OutOfRun || linkError) {
0241         continue;
0242       }
0243       eventMonitors.ecalOccRecdBx2_->Fill(ieta, iphi);
0244       eventMonitors.ecalOccRecd5Bx_->Fill(2);
0245       eventMonitors.ecalOccRecd5BxEtWgt_->Fill(2, tp.compressedEt());
0246     }
0247   }
0248   for (const auto& tp : (*ecalTPsRecdBx3)) {
0249     if (tp.compressedEt() > tpFillThreshold5Bx_) {
0250       const int ieta = tp.id().ieta();
0251       const int iphi = tp.id().iphi();
0252       const bool towerMasked = tp.sample(0).raw() & (1 << 13);
0253       const bool linkError = tp.sample(0).raw() & (1 << 15);
0254       if (towerMasked || caloLayer1OutOfRun || linkError) {
0255         continue;
0256       }
0257       eventMonitors.ecalOccRecdBx3_->Fill(ieta, iphi);
0258       eventMonitors.ecalOccRecd5Bx_->Fill(3);
0259       eventMonitors.ecalOccRecd5BxEtWgt_->Fill(3, tp.compressedEt());
0260     }
0261   }
0262   for (const auto& tp : (*ecalTPsRecdBx4)) {
0263     if (tp.compressedEt() > tpFillThreshold5Bx_) {
0264       const int ieta = tp.id().ieta();
0265       const int iphi = tp.id().iphi();
0266       const bool towerMasked = tp.sample(0).raw() & (1 << 13);
0267       const bool linkError = tp.sample(0).raw() & (1 << 15);
0268       if (towerMasked || caloLayer1OutOfRun || linkError) {
0269         continue;
0270       }
0271       eventMonitors.ecalOccRecdBx4_->Fill(ieta, iphi);
0272       eventMonitors.ecalOccRecd5Bx_->Fill(4);
0273       eventMonitors.ecalOccRecd5BxEtWgt_->Fill(4, tp.compressedEt());
0274     }
0275   }
0276   for (const auto& tp : (*ecalTPsRecdBx5)) {
0277     if (tp.compressedEt() > tpFillThreshold5Bx_) {
0278       const int ieta = tp.id().ieta();
0279       const int iphi = tp.id().iphi();
0280       const bool towerMasked = tp.sample(0).raw() & (1 << 13);
0281       const bool linkError = tp.sample(0).raw() & (1 << 15);
0282       if (towerMasked || caloLayer1OutOfRun || linkError) {
0283         continue;
0284       }
0285       eventMonitors.ecalOccRecdBx5_->Fill(ieta, iphi);
0286       eventMonitors.ecalOccRecd5Bx_->Fill(5);
0287       eventMonitors.ecalOccRecd5BxEtWgt_->Fill(5, tp.compressedEt());
0288     }
0289   }
0290 
0291   edm::Handle<HcalTrigPrimDigiCollection> hcalTPsSent;
0292   event.getByToken(hcalTPSourceSent_, hcalTPsSent);
0293   edm::Handle<HcalTrigPrimDigiCollection> hcalTPsRecd;
0294   event.getByToken(hcalTPSourceRecd_, hcalTPsRecd);
0295 
0296   hcalTPSentRecd_.clear();
0297 
0298   ComparisonHelper::zip(hcalTPsSent->begin(),
0299                         hcalTPsSent->end(),
0300                         hcalTPsRecd->begin(),
0301                         hcalTPsRecd->end(),
0302                         std::inserter(hcalTPSentRecd_, hcalTPSentRecd_.begin()),
0303                         HcalTrigPrimDigiCollection::key_compare());
0304 
0305   int nHcalLinkErrors{0};
0306   int nHcalMismatch{0};
0307 
0308   for (const auto& tpPair : hcalTPSentRecd_) {
0309     const auto& sentTp = tpPair.first;
0310     const auto& recdTp = tpPair.second;
0311     const int ieta = sentTp.id().ieta();
0312     if (abs(ieta) > 28 && sentTp.id().version() != 1)
0313       continue;
0314     const int iphi = sentTp.id().iphi();
0315     const bool towerMasked = recdTp.sample(0).raw() & (1 << 13);
0316     const bool linkMasked = recdTp.sample(0).raw() & (1 << 14);
0317     const bool linkError = recdTp.sample(0).raw() & (1 << 15);
0318 
0319     if (towerMasked) {
0320       eventMonitors.hcalOccTowerMasked_->Fill(ieta, iphi);
0321     }
0322     if (linkMasked) {
0323       eventMonitors.hcalOccLinkMasked_->Fill(ieta, iphi);
0324     }
0325 
0326     if (sentTp.SOI_compressedEt() > tpFillThreshold_) {
0327       eventMonitors.hcalTPRawEtSent_->Fill(sentTp.SOI_compressedEt());
0328       eventMonitors.hcalOccSent_->Fill(ieta, iphi);
0329     }
0330     if (sentTp.SOI_fineGrain() == 1) {
0331       eventMonitors.hcalOccSentFb_->Fill(ieta, iphi);
0332     }
0333     if (sentTp.t0().fineGrain(1) == 1) {
0334       eventMonitors.hcalOccSentFb2_->Fill(ieta, iphi);
0335     }
0336 
0337     if (towerMasked || caloLayer1OutOfRun) {
0338       // Do not compare if we have a mask applied
0339       continue;
0340     }
0341 
0342     if (linkError) {
0343       eventMonitors.hcalLinkError_->Fill(ieta, iphi);
0344       eventMonitors.hcalLinkErrorByLumi_->Fill(event.id().luminosityBlock());
0345       nHcalLinkErrors++;
0346       // Don't compare anymore, we already know its bad
0347       continue;
0348     }
0349 
0350     if (recdTp.SOI_compressedEt() > tpFillThreshold_) {
0351       eventMonitors.hcalTPRawEtRecd_->Fill(recdTp.SOI_compressedEt());
0352       eventMonitors.hcalOccupancy_->Fill(ieta, iphi);
0353       eventMonitors.hcalOccRecdEtWgt_->Fill(ieta, iphi, recdTp.SOI_compressedEt());
0354     }
0355     if (recdTp.SOI_fineGrain()) {
0356       eventMonitors.hcalOccRecdFb_->Fill(ieta, iphi);
0357     }
0358     if (recdTp.t0().fineGrain(1)) {
0359       eventMonitors.hcalOccRecdFb2_->Fill(ieta, iphi);
0360     }
0361 
0362     if (abs(ieta) > 29) {
0363       eventMonitors.hcalTPRawEtCorrelationHF_->Fill(sentTp.SOI_compressedEt(), recdTp.SOI_compressedEt());
0364     } else {
0365       eventMonitors.hcalTPRawEtCorrelationHBHE_->Fill(sentTp.SOI_compressedEt(), recdTp.SOI_compressedEt());
0366     }
0367 
0368     const bool HetAgreement = sentTp.SOI_compressedEt() == recdTp.SOI_compressedEt();
0369     const bool Hfb1Agreement = (abs(ieta) < 29) ? true
0370                                                 : (recdTp.SOI_compressedEt() == 0 ||
0371                                                    (sentTp.SOI_fineGrain() == recdTp.SOI_fineGrain()) || ignoreHFfbs_);
0372     // Ignore minBias (FB2) bit if we receieve 0 ET, which means it is likely zero-suppressed on HCal readout side
0373     const bool Hfb2Agreement =
0374         (abs(ieta) < 29)
0375             ? true
0376             : (recdTp.SOI_compressedEt() == 0 || (sentTp.SOI_fineGrain(1) == recdTp.SOI_fineGrain(1)) || ignoreHFfbs_);
0377     if (HetAgreement && Hfb1Agreement && Hfb2Agreement) {
0378       // Full match
0379       if (sentTp.SOI_compressedEt() > tpFillThreshold_) {
0380         eventMonitors.hcalOccSentAndRecd_->Fill(ieta, iphi);
0381         eventMonitors.hcalTPRawEtSentAndRecd_->Fill(sentTp.SOI_compressedEt());
0382       }
0383     } else {
0384       // There is some issue
0385       eventMonitors.hcalDiscrepancy_->Fill(ieta, iphi);
0386       eventMonitors.hcalMismatchByLumi_->Fill(event.id().luminosityBlock());
0387       nHcalMismatch++;
0388 
0389       if (not HetAgreement) {
0390         if (abs(ieta) > 29) {
0391           eventMonitors.HFmismatchesPerBx_->Fill(event.bunchCrossing());
0392         } else {
0393           eventMonitors.HBHEmismatchesPerBx_->Fill(event.bunchCrossing());
0394         }
0395         eventMonitors.hcalOccEtDiscrepancy_->Fill(ieta, iphi);
0396         eventMonitors.hcalTPRawEtDiffNoMatch_->Fill(recdTp.SOI_compressedEt() - sentTp.SOI_compressedEt());
0397         updateMismatch(event, 2, streamCache(event.streamID())->streamMismatchList);
0398 
0399         // Handle HCal discrepancy debug
0400         if (sentTp.SOI_compressedEt() == 0)
0401           eventMonitors.hcalOccRecdNotSent_->Fill(ieta, iphi);
0402         else if (recdTp.SOI_compressedEt() == 0)
0403           eventMonitors.hcalOccSentNotRecd_->Fill(ieta, iphi);
0404         else
0405           eventMonitors.hcalOccNoMatch_->Fill(ieta, iphi);
0406       }
0407       if (not Hfb1Agreement) {
0408         // Handle fine grain discrepancies
0409         eventMonitors.hcalOccFbDiscrepancy_->Fill(ieta, iphi);
0410         updateMismatch(event, 3, streamCache(event.streamID())->streamMismatchList);
0411       }
0412       if (not Hfb2Agreement) {
0413         // Handle fine grain discrepancies
0414         eventMonitors.hcalOccFb2Discrepancy_->Fill(ieta, iphi);
0415         updateMismatch(event, 3, streamCache(event.streamID())->streamMismatchList);
0416       }
0417     }
0418   }
0419 
0420   if (nHcalLinkErrors > streamCache(event.streamID())->streamNumMaxEvtLinkErrorsHCAL)
0421     streamCache(event.streamID())->streamNumMaxEvtLinkErrorsHCAL = nHcalLinkErrors;
0422   if (nHcalMismatch > streamCache(event.streamID())->streamNumMaxEvtMismatchHCAL)
0423     streamCache(event.streamID())->streamNumMaxEvtMismatchHCAL = nHcalMismatch;
0424 
0425   //fill inclusive link error and mismatch cache values based on whether HCAL or ECAL had more this event
0426   if (nEcalLinkErrors >= nHcalLinkErrors && nEcalLinkErrors > streamCache(event.streamID())->streamNumMaxEvtLinkErrors)
0427     streamCache(event.streamID())->streamNumMaxEvtLinkErrors = nEcalLinkErrors;
0428   else if (nEcalLinkErrors < nHcalLinkErrors &&
0429            nHcalLinkErrors > streamCache(event.streamID())->streamNumMaxEvtLinkErrors)
0430     streamCache(event.streamID())->streamNumMaxEvtLinkErrors = nHcalLinkErrors;
0431 
0432   if (nEcalMismatch >= nHcalMismatch && nEcalMismatch > streamCache(event.streamID())->streamNumMaxEvtMismatch)
0433     streamCache(event.streamID())->streamNumMaxEvtMismatch = nEcalMismatch;
0434   else if (nEcalMismatch < nHcalMismatch && nHcalMismatch > streamCache(event.streamID())->streamNumMaxEvtMismatch)
0435     streamCache(event.streamID())->streamNumMaxEvtMismatch = nHcalMismatch;
0436 }
0437 
0438 //push the mismatch type and identifying information onto the back of the stream-based vector
0439 //will maintain the vector's size at 20
0440 void L1TStage2CaloLayer1::updateMismatch(
0441     const edm::Event& e,
0442     int mismatchType,
0443     std::vector<std::tuple<edm::RunID, edm::LuminosityBlockID, edm::EventID, std::vector<int>>>& streamMismatches)
0444     const {
0445   //check if this combination of Run/Lumi/Event already exists in the stream mismatch list
0446   //if it does, update that entry, otherwise insert a new one with this particular combination.
0447   for (auto mismatchIterator = streamMismatches.begin(); mismatchIterator != streamMismatches.end();
0448        ++mismatchIterator) {
0449     if (e.getRun().id() == std::get<0>(*mismatchIterator) &&
0450         e.getLuminosityBlock().id() == std::get<1>(*mismatchIterator) && e.id() == std::get<2>(*mismatchIterator)) {
0451       //the run, lumi and event exist. Check if this kind of mismatch has been reported before
0452       std::vector<int>& mismatchTypeVector = std::get<3>(*mismatchIterator);
0453       for (auto mismatchTypeIterator = mismatchTypeVector.begin(); mismatchTypeIterator != mismatchTypeVector.end();
0454            ++mismatchTypeIterator) {
0455         if (mismatchType == *mismatchTypeIterator) {
0456           //this has already been reported
0457           return;
0458         }
0459       }
0460       //A mismatch exists, but it is not a type that has been previously reported.
0461       //Insert it into the vector of types of mismatches reported
0462       mismatchTypeVector.push_back(mismatchType);
0463       return;
0464     }
0465   }
0466 
0467   //The run/lumi/event does not exist in the list, construct an entry
0468   std::vector<int> newMismatchTypeVector;
0469   newMismatchTypeVector.push_back(mismatchType);
0470   std::tuple<edm::RunID, edm::LuminosityBlockID, edm::EventID, std::vector<int>> mismatchToInsert = {
0471       e.getRun().id(), e.getLuminosityBlock().id(), e.id(), newMismatchTypeVector};
0472   streamMismatches.push_back(mismatchToInsert);
0473 
0474   //maintain the mismatch vector size at 20.
0475   if (streamMismatches.size() > 20)
0476     streamMismatches.erase(streamMismatches.begin());
0477 }
0478 
0479 void L1TStage2CaloLayer1::dqmBeginRun(const edm::Run&,
0480                                       const edm::EventSetup&,
0481                                       CaloL1Information::monitoringDataHolder& eventMonitors) const {}
0482 
0483 void L1TStage2CaloLayer1::dqmEndRun(const edm::Run& run,
0484                                     const edm::EventSetup& es,
0485                                     const CaloL1Information::monitoringDataHolder& runMonitors,
0486                                     const CaloL1Information::perRunSummaryMonitoringInformation&) const {}
0487 
0488 void L1TStage2CaloLayer1::bookHistograms(DQMStore::IBooker& ibooker,
0489                                          const edm::Run& run,
0490                                          const edm::EventSetup& es,
0491                                          CaloL1Information::monitoringDataHolder& eventMonitors) const {
0492   auto bookEt = [&ibooker](std::string name, std::string title) {
0493     return ibooker.book1D(name, title + ";Raw ET;Counts", 256, -0.5, 255.5);
0494   };
0495   auto bookEtCorrelation = [&ibooker](std::string name, std::string title) {
0496     return ibooker.book2D(name, title, 256, -0.5, 255.5, 256, -0.5, 255.5);
0497   };
0498   auto bookEtDiff = [&ibooker](std::string name, std::string title) {
0499     return ibooker.book1D(name, title + ";#Delta Raw ET;Counts", 511, -255.5, 255.5);
0500   };
0501   auto bookEcalOccupancy = [&ibooker](std::string name, std::string title) {
0502     return ibooker.book2D(name, title + ";iEta;iPhi", 57, -28.5, 28.5, 72, 0.5, 72.5);
0503   };
0504   auto bookHcalOccupancy = [&ibooker](std::string name, std::string title) {
0505     return ibooker.book2D(name, title + ";iEta;iPhi", 83, -41.5, 41.5, 72, 0.5, 72.5);
0506   };
0507 
0508   ibooker.setCurrentFolder(histFolder_);
0509 
0510   eventMonitors.ecalDiscrepancy_ =
0511       bookEcalOccupancy("ecalDiscrepancy", "ECAL Discrepancies between TCC and Layer1 Readout");
0512   eventMonitors.ecalLinkError_ = bookEcalOccupancy("ecalLinkError", "ECAL Link Errors");
0513   eventMonitors.ecalOccupancy_ = bookEcalOccupancy("ecalOccupancy", "ECAL TP Occupancy at Layer1");
0514   eventMonitors.ecalOccRecdEtWgt_ = bookEcalOccupancy("ecalOccRecdEtWgt", "ECal TP ET-weighted Occupancy at Layer1");
0515   eventMonitors.hcalDiscrepancy_ =
0516       bookHcalOccupancy("hcalDiscrepancy", "HCAL Discrepancies between uHTR and Layer1 Readout");
0517   eventMonitors.hcalLinkError_ = bookHcalOccupancy("hcalLinkError", "HCAL Link Errors");
0518   eventMonitors.hcalOccupancy_ = bookHcalOccupancy("hcalOccupancy", "HCAL TP Occupancy at Layer1");
0519   eventMonitors.hcalOccRecdEtWgt_ = bookHcalOccupancy("hcalOccRecdEtWgt", "HCal TP ET-weighted Occupancy at Layer1");
0520 
0521   ibooker.setCurrentFolder(histFolder_ + "/ECalDetail");
0522 
0523   eventMonitors.ecalOccEtDiscrepancy_ = bookEcalOccupancy("ecalOccEtDiscrepancy", "ECal Et Discrepancy Occupancy");
0524   eventMonitors.ecalOccFgDiscrepancy_ =
0525       bookEcalOccupancy("ecalOccFgDiscrepancy", "ECal FG Veto Bit Discrepancy Occupancy");
0526   eventMonitors.ecalOccLinkMasked_ = bookEcalOccupancy("ecalOccLinkMasked", "ECal Masked Links");
0527   eventMonitors.ecalOccRecdFgVB_ = bookEcalOccupancy("ecalOccRecdFgVB", "ECal FineGrain Veto Bit Occupancy at Layer1");
0528   eventMonitors.ecalOccSentAndRecd_ = bookEcalOccupancy("ecalOccSentAndRecd", "ECal TP Occupancy FULL MATCH");
0529   eventMonitors.ecalOccSentFgVB_ = bookEcalOccupancy("ecalOccSentFgVB", "ECal FineGrain Veto Bit Occupancy at TCC");
0530   eventMonitors.ecalOccSent_ = bookEcalOccupancy("ecalOccSent", "ECal TP Occupancy at TCC");
0531   eventMonitors.ecalOccTowerMasked_ = bookEcalOccupancy("ecalOccTowerMasked", "ECal Masked towers");
0532   eventMonitors.ecalTPRawEtCorrelation_ =
0533       bookEtCorrelation("ecalTPRawEtCorrelation", "Raw Et correlation TCC and Layer1;TCC Et;Layer1 Et");
0534   eventMonitors.ecalTPRawEtDiffNoMatch_ = bookEtDiff("ecalTPRawEtDiffNoMatch", "ECal Raw Et Difference Layer1 - TCC");
0535   eventMonitors.ecalTPRawEtRecd_ = bookEt("ecalTPRawEtRecd", "ECal Raw Et Layer1 Readout");
0536   eventMonitors.ecalTPRawEtSentAndRecd_ = bookEt("ecalTPRawEtMatch", "ECal Raw Et FULL MATCH");
0537   eventMonitors.ecalTPRawEtSent_ = bookEt("ecalTPRawEtSent", "ECal Raw Et TCC Readout");
0538   eventMonitors.ecalOccRecd5Bx_ = ibooker.book1D("ecalOccRecd5Bx", "Number of TPs vs BX", 5, 1, 6);
0539   eventMonitors.ecalOccRecd5BxEtWgt_ =
0540       ibooker.book1D("ecalOccRecd5BxEtWgt", "Et-weighted number of TPs vs BX", 5, 1, 6);
0541   eventMonitors.ecalOccRecdBx1_ = bookEcalOccupancy("ecalOccRecdBx1", "ECal TP Occupancy for BX1");
0542   eventMonitors.ecalOccRecdBx2_ = bookEcalOccupancy("ecalOccRecdBx2", "ECal TP Occupancy for BX2");
0543   eventMonitors.ecalOccRecdBx3_ = bookEcalOccupancy("ecalOccRecdBx3", "ECal TP Occupancy for BX3");
0544   eventMonitors.ecalOccRecdBx4_ = bookEcalOccupancy("ecalOccRecdBx4", "ECal TP Occupancy for BX4");
0545   eventMonitors.ecalOccRecdBx5_ = bookEcalOccupancy("ecalOccRecdBx5", "ECal TP Occupancy for BX5");
0546 
0547   ibooker.setCurrentFolder(histFolder_ + "/ECalDetail/TCCDebug");
0548   eventMonitors.ecalOccSentNotRecd_ =
0549       bookHcalOccupancy("ecalOccSentNotRecd", "ECal TP Occupancy sent by TCC, zero at Layer1");
0550   eventMonitors.ecalOccRecdNotSent_ =
0551       bookHcalOccupancy("ecalOccRecdNotSent", "ECal TP Occupancy received by Layer1, zero at TCC");
0552   eventMonitors.ecalOccNoMatch_ =
0553       bookHcalOccupancy("ecalOccNoMatch", "ECal TP Occupancy for TCC and Layer1 nonzero, not matching");
0554 
0555   ibooker.setCurrentFolder(histFolder_ + "/HCalDetail");
0556 
0557   eventMonitors.hcalOccEtDiscrepancy_ = bookHcalOccupancy("hcalOccEtDiscrepancy", "HCal Et Discrepancy Occupancy");
0558   eventMonitors.hcalOccFbDiscrepancy_ =
0559       bookHcalOccupancy("hcalOccFbDiscrepancy", "HCal Feature Bit Discrepancy Occupancy");
0560   eventMonitors.hcalOccFb2Discrepancy_ =
0561       bookHcalOccupancy("hcalOccFb2Discrepancy", "HCal Second Feature Bit Discrepancy Occupancy");
0562   eventMonitors.hcalOccLinkMasked_ = bookHcalOccupancy("hcalOccLinkMasked", "HCal Masked Links");
0563   eventMonitors.hcalOccRecdFb_ = bookHcalOccupancy("hcalOccRecdFb", "HCal Feature Bit Occupancy at Layer1");
0564   eventMonitors.hcalOccRecdFb2_ = bookHcalOccupancy("hcalOccRecdFb2", "HF Second Feature Bit Occupancy at Layer1");
0565   eventMonitors.hcalOccSentAndRecd_ = bookHcalOccupancy("hcalOccSentAndRecd", "HCal TP Occupancy FULL MATCH");
0566   eventMonitors.hcalOccSentFb_ = bookHcalOccupancy("hcalOccSentFb", "HCal Feature Bit Occupancy at uHTR");
0567   eventMonitors.hcalOccSentFb2_ = bookHcalOccupancy("hcalOccSentFb2", "HF Second Feature Bit Occupancy at uHTR");
0568   eventMonitors.hcalOccSent_ = bookHcalOccupancy("hcalOccSent", "HCal TP Occupancy at uHTR");
0569   eventMonitors.hcalOccTowerMasked_ = bookHcalOccupancy("hcalOccTowerMasked", "HCal Masked towers");
0570   eventMonitors.hcalTPRawEtCorrelationHBHE_ =
0571       bookEtCorrelation("hcalTPRawEtCorrelationHBHE", "HBHE Raw Et correlation uHTR and Layer1;uHTR Et;Layer1 Et");
0572   eventMonitors.hcalTPRawEtCorrelationHF_ =
0573       bookEtCorrelation("hcalTPRawEtCorrelationHF", "HF Raw Et correlation uHTR and Layer1;uHTR Et;Layer1 Et");
0574   eventMonitors.hcalTPRawEtDiffNoMatch_ = bookEtDiff("hcalTPRawEtDiffNoMatch", "HCal Raw Et Difference Layer1 - uHTR");
0575   eventMonitors.hcalTPRawEtRecd_ = bookEt("hcalTPRawEtRecd", "HCal Raw Et Layer1 Readout");
0576   eventMonitors.hcalTPRawEtSentAndRecd_ = bookEt("hcalTPRawEtMatch", "HCal Raw Et FULL MATCH");
0577   eventMonitors.hcalTPRawEtSent_ = bookEt("hcalTPRawEtSent", "HCal Raw Et uHTR Readout");
0578 
0579   ibooker.setCurrentFolder(histFolder_ + "/HCalDetail/uHTRDebug");
0580   eventMonitors.hcalOccSentNotRecd_ =
0581       bookHcalOccupancy("hcalOccSentNotRecd", "HCal TP Occupancy sent by uHTR, zero at Layer1");
0582   eventMonitors.hcalOccRecdNotSent_ =
0583       bookHcalOccupancy("hcalOccRecdNotSent", "HCal TP Occupancy received by Layer1, zero at uHTR");
0584   eventMonitors.hcalOccNoMatch_ =
0585       bookHcalOccupancy("hcalOccNoMatch", "HCal TP Occupancy for uHTR and Layer1 nonzero, not matching");
0586 
0587   ibooker.setCurrentFolder(histFolder_ + "/MismatchDetail");
0588 
0589   const int nMismatchTypes = 4;
0590   eventMonitors.last20Mismatches_ = ibooker.book2D("last20Mismatches",
0591                                                    "Log of last 20 mismatches (use json tool to copy/paste)",
0592                                                    nMismatchTypes,
0593                                                    0,
0594                                                    nMismatchTypes,
0595                                                    20,
0596                                                    0,
0597                                                    20);
0598   eventMonitors.last20Mismatches_->setBinLabel(1, "Ecal TP Et Mismatch");
0599   eventMonitors.last20Mismatches_->setBinLabel(2, "Ecal TP Fine Grain Bit Mismatch");
0600   eventMonitors.last20Mismatches_->setBinLabel(3, "Hcal TP Et Mismatch");
0601   eventMonitors.last20Mismatches_->setBinLabel(4, "Hcal TP Feature Bit Mismatch");
0602 
0603   const int nLumis = 2000;
0604   eventMonitors.ecalLinkErrorByLumi_ = ibooker.book1D(
0605       "ecalLinkErrorByLumi", "Link error counts per lumi section for ECAL;LumiSection;Counts", nLumis, .5, nLumis + 0.5);
0606   eventMonitors.ecalMismatchByLumi_ = ibooker.book1D(
0607       "ecalMismatchByLumi", "Mismatch counts per lumi section for ECAL;LumiSection;Counts", nLumis, .5, nLumis + 0.5);
0608   eventMonitors.hcalLinkErrorByLumi_ = ibooker.book1D(
0609       "hcalLinkErrorByLumi", "Link error counts per lumi section for HCAL;LumiSection;Counts", nLumis, .5, nLumis + 0.5);
0610   eventMonitors.hcalMismatchByLumi_ = ibooker.book1D(
0611       "hcalMismatchByLumi", "Mismatch counts per lumi section for HCAL;LumiSection;Counts", nLumis, .5, nLumis + 0.5);
0612 
0613   eventMonitors.ECALmismatchesPerBx_ =
0614       ibooker.book1D("ECALmismatchesPerBx", "Mismatch counts per bunch crossing for ECAL", 3564, -.5, 3563.5);
0615   eventMonitors.HBHEmismatchesPerBx_ =
0616       ibooker.book1D("HBHEmismatchesPerBx", "Mismatch counts per bunch crossing for HBHE", 3564, -.5, 3563.5);
0617   eventMonitors.HFmismatchesPerBx_ =
0618       ibooker.book1D("HFmismatchesPerBx", "Mismatch counts per bunch crossing for HF", 3564, -.5, 3563.5);
0619 
0620   eventMonitors.maxEvtLinkErrorsByLumiECAL_ =
0621       ibooker.book1D("maxEvtLinkErrorsByLumiECAL",
0622                      "Max number of single-event ECAL link errors per lumi section;LumiSection;Counts",
0623                      nLumis,
0624                      .5,
0625                      nLumis + 0.5);
0626   eventMonitors.maxEvtLinkErrorsByLumiHCAL_ =
0627       ibooker.book1D("maxEvtLinkErrorsByLumiHCAL",
0628                      "Max number of single-event HCAL link errors per lumi section;LumiSection;Counts",
0629                      nLumis,
0630                      .5,
0631                      nLumis + 0.5);
0632 
0633   eventMonitors.maxEvtMismatchByLumiECAL_ =
0634       ibooker.book1D("maxEvtMismatchByLumiECAL",
0635                      "Max number of single-event ECAL discrepancies per lumi section;LumiSection;Counts",
0636                      nLumis,
0637                      .5,
0638                      nLumis + 0.5);
0639   eventMonitors.maxEvtMismatchByLumiHCAL_ =
0640       ibooker.book1D("maxEvtMismatchByLumiHCAL",
0641                      "Max number of single-event HCAL discrepancies per lumi section;LumiSection;Counts",
0642                      nLumis,
0643                      .5,
0644                      nLumis + 0.5);
0645 
0646   ibooker.setCurrentFolder(histFolder_);
0647   eventMonitors.maxEvtLinkErrorsByLumi_ =
0648       ibooker.book1D("maxEvtLinkErrorsByLumi",
0649                      "Max number of single-event link errors per lumi section;LumiSection;Counts",
0650                      nLumis,
0651                      .5,
0652                      nLumis + 0.5);
0653   eventMonitors.maxEvtMismatchByLumi_ =
0654       ibooker.book1D("maxEvtMismatchByLumi",
0655                      "Max number of single-event discrepancies per lumi section;LumiSection;Counts",
0656                      nLumis,
0657                      .5,
0658                      nLumis + 0.5);
0659 
0660   ibooker.setCurrentFolder(histFolder_ + "/AMC13ErrorCounters");
0661   eventMonitors.bxidErrors_ =
0662       ibooker.book1D("bxidErrors", "bxid mismatch between AMC13 and CTP Cards;Layer1 Phi;Counts", 18, -.5, 17.5);
0663   eventMonitors.l1idErrors_ =
0664       ibooker.book1D("l1idErrors", "l1id mismatch between AMC13 and CTP Cards;Layer1 Phi;Counts", 18, -.5, 17.5);
0665   eventMonitors.orbitErrors_ =
0666       ibooker.book1D("orbitErrors", "orbit mismatch between AMC13 and CTP Cards;Layer1 Phi;Counts", 18, -.5, 17.5);
0667 }
0668 
0669 void L1TStage2CaloLayer1::streamEndLuminosityBlockSummary(
0670     edm::StreamID theStreamID,
0671     edm::LuminosityBlock const& theLumiBlock,
0672     edm::EventSetup const& theEventSetup,
0673     CaloL1Information::perLumiBlockMonitoringInformation* lumiMonitoringInformation) const {
0674   auto theStreamCache = streamCache(theStreamID);
0675 
0676   lumiMonitoringInformation->lumiNumMaxEvtLinkErrorsECAL =
0677       std::max(lumiMonitoringInformation->lumiNumMaxEvtLinkErrorsECAL, theStreamCache->streamNumMaxEvtLinkErrorsECAL);
0678   lumiMonitoringInformation->lumiNumMaxEvtLinkErrorsHCAL =
0679       std::max(lumiMonitoringInformation->lumiNumMaxEvtLinkErrorsHCAL, theStreamCache->streamNumMaxEvtLinkErrorsHCAL);
0680   lumiMonitoringInformation->lumiNumMaxEvtLinkErrors =
0681       std::max(lumiMonitoringInformation->lumiNumMaxEvtLinkErrors, theStreamCache->streamNumMaxEvtLinkErrors);
0682 
0683   lumiMonitoringInformation->lumiNumMaxEvtMismatchECAL =
0684       std::max(lumiMonitoringInformation->lumiNumMaxEvtMismatchECAL, theStreamCache->streamNumMaxEvtMismatchECAL);
0685   lumiMonitoringInformation->lumiNumMaxEvtMismatchHCAL =
0686       std::max(lumiMonitoringInformation->lumiNumMaxEvtMismatchHCAL, theStreamCache->streamNumMaxEvtMismatchHCAL);
0687   lumiMonitoringInformation->lumiNumMaxEvtMismatch =
0688       std::max(lumiMonitoringInformation->lumiNumMaxEvtMismatch, theStreamCache->streamNumMaxEvtMismatch);
0689 
0690   //reset the stream cache here, since we are done with the luminosity block in this stream
0691   //We don't want the stream to be comparing to previous values in the next lumi-block
0692   theStreamCache->streamNumMaxEvtLinkErrorsECAL = 0;
0693   theStreamCache->streamNumMaxEvtLinkErrorsHCAL = 0;
0694   theStreamCache->streamNumMaxEvtLinkErrors = 0;
0695 
0696   theStreamCache->streamNumMaxEvtMismatchECAL = 0;
0697   theStreamCache->streamNumMaxEvtMismatchHCAL = 0;
0698   theStreamCache->streamNumMaxEvtMismatch = 0;
0699 }
0700 
0701 void L1TStage2CaloLayer1::globalEndLuminosityBlockSummary(
0702     edm::LuminosityBlock const& theLumiBlock,
0703     edm::EventSetup const& theEventSetup,
0704     CaloL1Information::perLumiBlockMonitoringInformation* lumiMonitoringInformation) const {
0705   auto theRunCache = runCache(theLumiBlock.getRun().index());
0706   //read the cache out to the relevant Luminosity histogram bin
0707   theRunCache->maxEvtLinkErrorsByLumiECAL_->Fill(theLumiBlock.luminosityBlock(),
0708                                                  lumiMonitoringInformation->lumiNumMaxEvtLinkErrorsECAL);
0709   theRunCache->maxEvtLinkErrorsByLumiHCAL_->Fill(theLumiBlock.luminosityBlock(),
0710                                                  lumiMonitoringInformation->lumiNumMaxEvtLinkErrorsHCAL);
0711   theRunCache->maxEvtLinkErrorsByLumi_->Fill(theLumiBlock.luminosityBlock(),
0712                                              lumiMonitoringInformation->lumiNumMaxEvtLinkErrors);
0713 
0714   theRunCache->maxEvtMismatchByLumiECAL_->Fill(theLumiBlock.luminosityBlock(),
0715                                                lumiMonitoringInformation->lumiNumMaxEvtMismatchECAL);
0716   theRunCache->maxEvtMismatchByLumiHCAL_->Fill(theLumiBlock.luminosityBlock(),
0717                                                lumiMonitoringInformation->lumiNumMaxEvtMismatchHCAL);
0718   theRunCache->maxEvtMismatchByLumi_->Fill(theLumiBlock.luminosityBlock(),
0719                                            lumiMonitoringInformation->lumiNumMaxEvtMismatch);
0720 
0721   // Simple way to embed current lumi to auto-scale axis limits in render plugin
0722   theRunCache->ecalLinkErrorByLumi_->setBinContent(0, theLumiBlock.luminosityBlock());
0723   theRunCache->ecalMismatchByLumi_->setBinContent(0, theLumiBlock.luminosityBlock());
0724   theRunCache->hcalLinkErrorByLumi_->setBinContent(0, theLumiBlock.luminosityBlock());
0725   theRunCache->hcalMismatchByLumi_->setBinContent(0, theLumiBlock.luminosityBlock());
0726   theRunCache->maxEvtLinkErrorsByLumiECAL_->setBinContent(0, theLumiBlock.luminosityBlock());
0727   theRunCache->maxEvtLinkErrorsByLumiHCAL_->setBinContent(0, theLumiBlock.luminosityBlock());
0728   theRunCache->maxEvtLinkErrorsByLumi_->setBinContent(0, theLumiBlock.luminosityBlock());
0729   theRunCache->maxEvtMismatchByLumiECAL_->setBinContent(0, theLumiBlock.luminosityBlock());
0730   theRunCache->maxEvtMismatchByLumiHCAL_->setBinContent(0, theLumiBlock.luminosityBlock());
0731   theRunCache->maxEvtMismatchByLumi_->setBinContent(0, theLumiBlock.luminosityBlock());
0732 }
0733 
0734 //returns true if the new candidate mismatch is from a later mismatch than the comparison mismatch
0735 // based on Run, Lumisection, and event number
0736 //false otherwise
0737 bool L1TStage2CaloLayer1::isLaterMismatch(
0738     std::tuple<edm::RunID, edm::LuminosityBlockID, edm::EventID, std::vector<int>>& candidateMismatch,
0739     std::tuple<edm::RunID, edm::LuminosityBlockID, edm::EventID, std::vector<int>>& comparisonMismatch) const {
0740   //check the run. If the run ID of the candidate mismatch is less than the run ID of the comparison mismatch, it is earlier,
0741   if (std::get<0>(candidateMismatch) < std::get<0>(comparisonMismatch))
0742     return false;
0743   //if it is greater, then it is a later mismatch
0744   else if (std::get<0>(candidateMismatch) > std::get<0>(comparisonMismatch))
0745     return true;
0746   //if it is even, then we need to repeat this comparison on the luminosity block
0747   else {
0748     if (std::get<1>(candidateMismatch) < std::get<1>(comparisonMismatch))
0749       return false;  //if the lumi block is less, it is an earlier mismatch
0750     else if (std::get<1>(candidateMismatch) > std::get<1>(comparisonMismatch))
0751       return true;  // if the lumi block is greater, than it is a later mismatch
0752     // if these are equivalent, then we repeat the comparison on the event
0753     else {
0754       if (std::get<2>(candidateMismatch) < std::get<2>(comparisonMismatch))
0755         return false;
0756       else
0757         return true;  //in the case of even events here, we consider the event later.
0758     }
0759   }
0760   //default return, should never be called.
0761   return false;
0762 }
0763 
0764 //binary search like algorithm for trying to find the appropriate place to put our new mismatch
0765 //will find an interger we add to the iterator to get the proper location to find the insertion location
0766 //will return -1 if the mismatch should not be inserted into the list
0767 int L1TStage2CaloLayer1::findIndex(
0768     std::tuple<edm::RunID, edm::LuminosityBlockID, edm::EventID, std::vector<int>> candidateMismatch,
0769     std::vector<std::tuple<edm::RunID, edm::LuminosityBlockID, edm::EventID, std::vector<int>>> comparisonList,
0770     int lowerIndexToSearch,
0771     int upperIndexToSearch) const {
0772   //Start by getting the spot in the the vector to start searching
0773   int searchLocation = ((upperIndexToSearch + lowerIndexToSearch) / 2);
0774   auto searchIterator = comparisonList.begin() + searchLocation;
0775   //Multiple possible cases:
0776   //case one. Greater than the search element
0777   if (this->isLaterMismatch(candidateMismatch, *searchIterator)) {
0778     //subcase one, there exists a mismatch to its right,
0779     if (searchIterator + 1 != comparisonList.end()) {
0780       //subsubcase one, the candidate is earlier than the one to the right, insert this element between these two
0781       if (not this->isLaterMismatch(candidateMismatch, *(searchIterator + 1))) {
0782         return searchLocation + 1;  //vector insert inserts *before* the position, so we return +1
0783       }
0784       //subsubcase two, the candidate is later than it, in this case we refine the search area.
0785       else {
0786         return this->findIndex(candidateMismatch, comparisonList, searchLocation + 1, upperIndexToSearch);
0787       }
0788     }
0789     //subcase two, there exists no mismatch to it's right (end of the vector), in which case this is the latest mismatch
0790     else {
0791       return searchLocation +
0792              1;  //vector insert inserts *before* the position, so we return +1 (should be synonymous with .end())
0793     }
0794   }
0795   //case two. we are earlier than the current mismatch
0796   else {
0797     //subcase one, these exists a mismatch to its left,
0798     if (searchIterator != comparisonList.begin()) {
0799       //subsubcase one, the candidate mismatch is earlier than the one to the left. in this case we refine the search area
0800       if (not this->isLaterMismatch(candidateMismatch, *(searchIterator - 1))) {
0801         return this->findIndex(candidateMismatch, comparisonList, lowerIndexToSearch, searchLocation - 1);
0802       }
0803       //subsubcase two the candidate is later than than the one to it's left, in which case, we insert the element between these two.
0804       else {
0805         return searchLocation;  //vector insert inserts *before* the position, so we return without addition here.
0806       }
0807     }
0808     //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.
0809     else {
0810       //subsubcase one. It is possible we still insert this if the capacity of the vector is below 20.
0811       //this should probably be rewritten to have the capacity reserved before-hand so that 20 is not potentially a hard coded magic number
0812       //defined by what we know to be true about a non-obvious data type elsewhere.
0813       if (comparisonList.size() < 20) {
0814         return searchLocation;
0815       }
0816       //subsubcase two. The size is 20 or above, and we are earlier than all of them. This should not be inserted into the list.
0817       else {
0818         return -1;
0819       }
0820     }
0821   }
0822   //default return. Should never be called
0823   return -1;
0824 }
0825 
0826 //will shuffle the candidate mismatch list into the comparison mismatch list.
0827 void L1TStage2CaloLayer1::mergeMismatchVectors(
0828     std::vector<std::tuple<edm::RunID, edm::LuminosityBlockID, edm::EventID, std::vector<int>>>& candidateMismatchList,
0829     std::vector<std::tuple<edm::RunID, edm::LuminosityBlockID, edm::EventID, std::vector<int>>>& comparisonMismatchList)
0830     const {
0831   //okay now we loop over our candidate mismatches
0832   for (auto candidateIterator = candidateMismatchList.begin(); candidateIterator != candidateMismatchList.end();
0833        ++candidateIterator) {
0834     int insertionIndex = this->findIndex(*candidateIterator, comparisonMismatchList, 0, comparisonMismatchList.size());
0835     if (insertionIndex < 0) {
0836       continue;  //if we didn't find anywhere to put this mismatch, we move on
0837     }
0838     auto insertionIterator = comparisonMismatchList.begin() + insertionIndex;
0839     comparisonMismatchList.insert(insertionIterator, *candidateIterator);
0840     //now if we have more than 20 mismatches in the list, we erase the earliest one (the beginning)
0841     //this should probably be rewritten to have the capacity reserved before-hand so that 20 is not potentially a hard coded magic number
0842     //defined by what we know to be true about a non-obvious data type elsewhere.
0843     if (comparisonMismatchList.size() > 20) {
0844       comparisonMismatchList.erase(comparisonMismatchList.begin());
0845     }
0846   }
0847 }
0848 
0849 void L1TStage2CaloLayer1::streamEndRunSummary(
0850     edm::StreamID theStreamID,
0851     edm::Run const& theRun,
0852     edm::EventSetup const& theEventSetup,
0853     CaloL1Information::perRunSummaryMonitoringInformation* theRunSummaryMonitoringInformation) const {
0854   auto theStreamCache = streamCache(theStreamID);
0855 
0856   if (theRunSummaryMonitoringInformation->runMismatchList.empty())
0857     theRunSummaryMonitoringInformation->runMismatchList = theStreamCache->streamMismatchList;
0858   else
0859     this->mergeMismatchVectors(theStreamCache->streamMismatchList, theRunSummaryMonitoringInformation->runMismatchList);
0860 
0861   //clear the stream cache so that the next run does not try to compare to the current one.
0862   theStreamCache->streamMismatchList.clear();
0863 }
0864 
0865 void L1TStage2CaloLayer1::globalEndRunSummary(
0866     edm::Run const& theRun,
0867     edm::EventSetup const& theEventSetup,
0868     CaloL1Information::perRunSummaryMonitoringInformation* theRunSummaryInformation) const {
0869   //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.
0870   auto theRunCache = runCache(theRun.index());
0871   //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.
0872   int ibin = 1;
0873   for (auto mismatchIterator = theRunSummaryInformation->runMismatchList.begin();
0874        mismatchIterator != theRunSummaryInformation->runMismatchList.end();
0875        ++mismatchIterator) {
0876     std::string binLabel = std::to_string(std::get<0>(*mismatchIterator).run()) + ":" +
0877                            std::to_string(std::get<1>(*mismatchIterator).luminosityBlock()) + ":" +
0878                            std::to_string(std::get<2>(*mismatchIterator).event());
0879     theRunCache->last20Mismatches_->setBinLabel(ibin, binLabel, 2);
0880     //Get the vector of mismatches for this particular event and iterate through it.
0881     //Set the bin content to 1 for each type of mismatch seen
0882     std::vector<int> mismatchTypeVector = std::get<3>(*mismatchIterator);
0883     for (auto mismatchTypeIterator = mismatchTypeVector.begin(); mismatchTypeIterator != mismatchTypeVector.end();
0884          ++mismatchTypeIterator) {
0885       theRunCache->last20Mismatches_->setBinContent(*mismatchTypeIterator + 1, ibin, 1);
0886     }
0887     ++ibin;
0888   }
0889   //remove the remaining empty string labels to prevent overlap
0890   for (int emptyBinIndex = ibin; emptyBinIndex <= 20; ++emptyBinIndex) {
0891     theRunCache->last20Mismatches_->setBinLabel(emptyBinIndex, std::to_string(emptyBinIndex), 2);
0892   }
0893 }