Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 10:55:21

0001 /*
0002  * \file L1TdeStage2CaloLayer1.cc
0003  *
0004  * N. Smith <nick.smith@cern.ch>
0005  */
0006 //Modified by Bhawna Gomber <bhawna.gomber@cern.ch>
0007 
0008 #include "DQM/L1TMonitor/interface/L1TdeStage2CaloLayer1.h"
0009 
0010 #include "FWCore/Framework/interface/ESHandle.h"
0011 #include "FWCore/Framework/interface/LuminosityBlock.h"
0012 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0013 #include "DataFormats/Provenance/interface/EventAuxiliary.h"
0014 
0015 #include "CondFormats/RunInfo/interface/RunInfo.h"
0016 #include "CondFormats/DataRecord/interface/RunSummaryRcd.h"
0017 
0018 using namespace l1t;
0019 
0020 L1TdeStage2CaloLayer1::L1TdeStage2CaloLayer1(const edm::ParameterSet& ps)
0021     : dataLabel_(ps.getParameter<edm::InputTag>("dataSource")),
0022       dataSource_(consumes<CaloTowerBxCollection>(dataLabel_)),
0023       emulLabel_(ps.getParameter<edm::InputTag>("emulSource")),
0024       emulSource_(consumes<CaloTowerBxCollection>(emulLabel_)),
0025       hcalTowers_(consumes<HcalTrigPrimDigiCollection>(edm::InputTag("l1tCaloLayer1Digis"))),
0026       fedRawData_(consumes<FEDRawDataCollection>(ps.getParameter<edm::InputTag>("fedRawDataLabel"))),
0027       histFolder_(ps.getParameter<std::string>("histFolder")),
0028       tpFillThreshold_(ps.getUntrackedParameter<int>("etDistributionsFillThreshold", 0)) {
0029   dataEmulDenominator_ = 0.;
0030   for (size_t i = 0; i < NSummaryColumns; ++i) {
0031     dataEmulNumerator_[i] = 0.;
0032   }
0033 }
0034 
0035 L1TdeStage2CaloLayer1::~L1TdeStage2CaloLayer1() {}
0036 
0037 void L1TdeStage2CaloLayer1::analyze(const edm::Event& event, const edm::EventSetup& es) {
0038   edm::Handle<CaloTowerBxCollection> dataTowers;
0039   event.getByToken(dataSource_, dataTowers);
0040   edm::Handle<CaloTowerBxCollection> emulTowers;
0041   event.getByToken(emulSource_, emulTowers);
0042   edm::Handle<HcalTrigPrimDigiCollection> hcalTowers;
0043   event.getByToken(hcalTowers_, hcalTowers);
0044 
0045   // Best way I know to check if FED in a run
0046   edm::Handle<FEDRawDataCollection> fedRawDataCollection;
0047   event.getByToken(fedRawData_, fedRawDataCollection);
0048   bool caloLayer1OutOfRun{true};
0049   bool caloLayer2OutOfRun{true};
0050   if (fedRawDataCollection.isValid()) {
0051     caloLayer1OutOfRun = false;
0052     caloLayer2OutOfRun = false;
0053     for (int iFed = 1354; iFed < 1360; iFed += 2) {
0054       const FEDRawData& fedRawData = fedRawDataCollection->FEDData(iFed);
0055       if (fedRawData.size() == 0) {
0056         caloLayer1OutOfRun = true;
0057         continue;  // In case one of 3 layer 1 FEDs not in
0058       }
0059     }
0060     const FEDRawData& fedRawData = fedRawDataCollection->FEDData(1360);
0061     if (fedRawData.size() == 0) {
0062       caloLayer2OutOfRun = true;
0063     }
0064   }
0065 
0066   if (caloLayer1OutOfRun or caloLayer2OutOfRun) {
0067     // No point in comparing
0068     return;
0069   }
0070 
0071   // We'll fill sets, compare, and then dissect comparison failures after.
0072   SimpleTowerSet dataTowerSet;
0073   // BXVector::begin(int bx)
0074   for (auto iter = dataTowers->begin(0); iter != dataTowers->end(0); ++iter) {
0075     const auto& tower = *iter;
0076     int eta = tower.hwEta();
0077     if (eta == 29)
0078       eta = 30;
0079     if (eta == -29)
0080       eta = -30;
0081     dataTowerSet.emplace(eta, tower.hwPhi(), tower.hwPt() + (tower.hwEtRatio() << 9) + (tower.hwQual() << 12), true);
0082     if (tower.hwPt() > tpFillThreshold_) {
0083       dataOcc_->Fill(eta, tower.hwPhi());
0084       dataEtDistribution_->Fill(tower.hwPt());
0085     }
0086   }
0087   SimpleTowerSet emulTowerSet;
0088   for (auto iter = emulTowers->begin(0); iter != emulTowers->end(0); ++iter) {
0089     const auto& tower = *iter;
0090     emulTowerSet.emplace(
0091         tower.hwEta(), tower.hwPhi(), tower.hwPt() + (tower.hwEtRatio() << 9) + (tower.hwQual() << 12), false);
0092     if (tower.hwPt() > tpFillThreshold_) {
0093       emulOcc_->Fill(tower.hwEta(), tower.hwPhi());
0094       emulEtDistribution_->Fill(tower.hwPt());
0095     }
0096   }
0097 
0098   bool etMsmThisEvent{false};
0099   bool erMsmThisEvent{false};
0100   bool fbMsmThisEvent{false};
0101   bool towerCountMsmThisEvent{false};
0102 
0103   if (dataTowerSet.size() != emulTowerSet.size()) {
0104     // This will happen if either CaloLayer1 or CaloLayer2 are out of run (in which case we exit early)
0105     // The problematic situation that we are monitoring is when we see both in, but one MP7 card is not sending fat events when it should
0106     towerCountMismatchesPerBx_->Fill(event.bunchCrossing());
0107     towerCountMsmThisEvent = true;
0108   } else {
0109     auto dataIter = dataTowerSet.begin();
0110     auto emulIter = emulTowerSet.begin();
0111     while (dataIter != dataTowerSet.end() && emulIter != emulTowerSet.end()) {
0112       auto dataTower = *(dataIter++);
0113       auto emulTower = *(emulIter++);
0114       assert(dataTower.ieta_ == emulTower.ieta_ && dataTower.iphi_ == emulTower.iphi_);
0115 
0116       etCorrelation_->Fill(dataTower.et(), emulTower.et());
0117 
0118       if (abs(dataTower.ieta_) >= 30) {
0119         fbCorrelationHF_->Fill(dataTower.fb(), emulTower.fb());
0120       } else {
0121         fbCorrelation_->Fill(dataTower.fb(), emulTower.fb());
0122       }
0123 
0124       if (dataTower.data_ == emulTower.data_) {
0125         // Perfect match
0126         if (dataTower.et() > tpFillThreshold_) {
0127           matchOcc_->Fill(dataTower.ieta_, dataTower.iphi_);
0128         }
0129       } else {
0130         // Ok, now dissect the failure
0131         if (dataTower.et() != emulTower.et()) {
0132           if (dataTower.et() == 0)
0133             failureOccEtDataZero_->Fill(dataTower.ieta_, dataTower.iphi_);
0134           else if (emulTower.et() == 0)
0135             failureOccEtEmulZero_->Fill(dataTower.ieta_, dataTower.iphi_);
0136           else
0137             failureOccEtMismatch_->Fill(dataTower.ieta_, dataTower.iphi_);
0138 
0139           etMismatchDiff_->Fill(dataTower.et() - emulTower.et());
0140           etMismatchByLumi_->Fill(event.id().luminosityBlock());
0141           etMismatchesPerBx_->Fill(event.bunchCrossing());
0142           etMsmThisEvent = true;
0143           updateMismatch(event, 0);
0144         }
0145         if (dataTower.er() != emulTower.er()) {
0146           failureOccErMismatch_->Fill(dataTower.ieta_, dataTower.iphi_);
0147           erMismatchByLumi_->Fill(event.id().luminosityBlock());
0148           erMismatchesPerBx_->Fill(event.bunchCrossing());
0149           erMsmThisEvent = true;
0150           updateMismatch(event, 1);
0151         }
0152         if (dataTower.fb() != emulTower.fb()) {
0153           failureOccFbMismatch_->Fill(dataTower.ieta_, dataTower.iphi_);
0154           fbMismatchByLumi_->Fill(event.id().luminosityBlock());
0155           fbMismatchesPerBx_->Fill(event.bunchCrossing());
0156           dataEtDistributionFBMismatch_->Fill(dataTower.et());
0157           fbMsmThisEvent = true;
0158           updateMismatch(event, 2);
0159         }
0160       }
0161     }
0162   }
0163 
0164   dataEmulDenominator_ += 1;
0165   if (etMsmThisEvent)
0166     dataEmulNumerator_[EtMismatch] += 1;
0167   if (erMsmThisEvent)
0168     dataEmulNumerator_[ErMismatch] += 1;
0169   if (fbMsmThisEvent)
0170     dataEmulNumerator_[FbMismatch] += 1;
0171   if (towerCountMsmThisEvent)
0172     dataEmulNumerator_[TowerCountMismatch] += 1;
0173 
0174   for (size_t i = 0; i < NSummaryColumns; ++i) {
0175     dataEmulSummary_->setBinContent(1 + i, dataEmulNumerator_[i] / dataEmulDenominator_);
0176   }
0177   // GetEntries() increments every SetBinContent() call
0178   dataEmulSummary_->getTH1F()->SetEntries(dataEmulDenominator_);
0179 
0180   // To see if problems correlate with TMT cycle (i.e. an MP7-side issue)
0181   if (etMsmThisEvent or erMsmThisEvent or fbMsmThisEvent or towerCountMsmThisEvent) {
0182     mismatchesPerBxMod9_->Fill(event.bunchCrossing() % 9);
0183   }
0184 }
0185 
0186 void L1TdeStage2CaloLayer1::updateMismatch(const edm::Event& e, int mismatchType) {
0187   auto id = e.id();
0188   std::string eventString{std::to_string(id.run()) + ":" + std::to_string(id.luminosityBlock()) + ":" +
0189                           std::to_string(id.event())};
0190   if (last20MismatchArray_.at(lastMismatchIndex_).first == eventString) {
0191     // same event
0192     last20MismatchArray_.at(lastMismatchIndex_).second |= 1 << mismatchType;
0193   } else {
0194     // New event, advance
0195     lastMismatchIndex_ = (lastMismatchIndex_ + 1) % 20;
0196     last20MismatchArray_.at(lastMismatchIndex_) = {eventString, 1 << mismatchType};
0197   }
0198 }
0199 
0200 void L1TdeStage2CaloLayer1::beginLuminosityBlock(const edm::LuminosityBlock&, const edm::EventSetup&) {
0201   // Ugly way to loop backwards through the last 20 mismatches
0202   for (size_t ibin = 1, imatch = lastMismatchIndex_; ibin <= 20; ibin++, imatch = (imatch + 19) % 20) {
0203     last20Mismatches_->setBinLabel(ibin, last20MismatchArray_.at(imatch).first, /* axis */ 2);
0204     for (int itype = 0; itype < 4; ++itype) {
0205       int binContent = (last20MismatchArray_.at(imatch).second >> itype) & 1;
0206       last20Mismatches_->setBinContent(itype + 1, ibin, binContent);
0207     }
0208   }
0209 }
0210 
0211 void L1TdeStage2CaloLayer1::endLuminosityBlock(const edm::LuminosityBlock& lumi, const edm::EventSetup&) {
0212   auto id = static_cast<double>(lumi.id().luminosityBlock());  // uint64_t
0213   // Simple way to embed current lumi to auto-scale axis limits in render plugin
0214   etMismatchByLumi_->setBinContent(0, id);
0215   erMismatchByLumi_->setBinContent(0, id);
0216   fbMismatchByLumi_->setBinContent(0, id);
0217 }
0218 
0219 void L1TdeStage2CaloLayer1::bookHistograms(DQMStore::IBooker& ibooker, const edm::Run& run, const edm::EventSetup& es) {
0220   auto bookEt = [&ibooker](std::string name, std::string title) {
0221     title += ";Raw ET value";
0222     return ibooker.book1D(name, title, 512, -0.5, 511.5);
0223   };
0224   auto bookEtDiff = [&ibooker](std::string name, std::string title) {
0225     title += ";#Delta raw ET value";
0226     return ibooker.book1D(name, title, 1023, -511.5, 511.5);
0227   };
0228   auto bookEtCorrelation = [&ibooker](std::string name, std::string title) {
0229     return ibooker.book2D(name, title, 512, -0.5, 511.5, 512, -0.5, 511.5);
0230   };
0231   auto bookOccupancy = [&ibooker](std::string name, std::string title) {
0232     title += ";i#eta;i#phi";
0233     return ibooker.book2D(name, title, 83, -41.5, 41.5, 72, 0.5, 72.5);
0234   };
0235 
0236   ibooker.setCurrentFolder(histFolder_ + "/");
0237   dataEmulSummary_ = ibooker.book1D("dataEmulSummary",
0238                                     "CaloLayer1 data-emul mismatch frac. (entries=evts processed)",
0239                                     NSummaryColumns,
0240                                     0.,
0241                                     double(NSummaryColumns));
0242   dataEmulSummary_->setAxisTitle("Fraction events with mismatch", /* axis */ 2);
0243   dataEmulSummary_->setBinLabel(1 + EtMismatch, "Et");
0244   dataEmulSummary_->setBinLabel(1 + ErMismatch, "Et ratio");
0245   dataEmulSummary_->setBinLabel(1 + FbMismatch, "Feature bit");
0246   dataEmulSummary_->setBinLabel(1 + TowerCountMismatch, "CaloTower readout");
0247   mismatchesPerBxMod9_ = ibooker.book1D(
0248       "mismatchesPerBxMod9", "CaloLayer1 data-emulator mismatch per bx%9;Bunch crossing mod 9;Counts", 9, -0.5, 8.5);
0249 
0250   ibooker.setCurrentFolder(histFolder_ + "/Occupancies");
0251 
0252   dataOcc_ = bookOccupancy("dataOcc", "Tower Occupancy for Data");
0253   emulOcc_ = bookOccupancy("emulOcc", "Tower Occupancy for Emulator");
0254   matchOcc_ = bookOccupancy("matchOcc", "Tower Occupancy for Data/Emulator Full Matches");
0255   failureOccEtMismatch_ = bookOccupancy("failureOccEtMismatch", "Tower Occupancy for Data/Emulator ET Mismatch");
0256   failureOccEtDataZero_ = bookOccupancy("failureOccEtDataZero", "Tower Occupancy for Data ET Zero, Emul Nonzero");
0257   failureOccEtEmulZero_ = bookOccupancy("failureOccEtEmulZero", "Tower Occupancy for Data ET Nonzero, Emul Zero");
0258   failureOccErMismatch_ = bookOccupancy("failureOccErMismatch", "Tower Occupancy for Data/Emulator ET Ratio Mismatch");
0259   failureOccFbMismatch_ =
0260       bookOccupancy("failureOccFbMismatch", "Tower Occupancy for Data/Emulator Feature Bit Mismatch");
0261 
0262   ibooker.setCurrentFolder(histFolder_ + "/EtDistributions");
0263   dataEtDistribution_ = bookEt("dataEtDistribution", "ET distribution for towers in data");
0264   dataEtDistributionFBMismatch_ =
0265       bookEt("dataEtDistributionFBMismatch", "ET distribution for towers in data when FB Mismatch");
0266   emulEtDistribution_ = bookEt("emulEtDistribution", "ET distribution for towers in emulator");
0267   etCorrelation_ = bookEtCorrelation("EtCorrelation", "Et correlation for all towers;Data tower Et;Emulator tower Et");
0268   matchEtDistribution_ = bookEt("matchEtDistribution", "ET distribution for towers matched between data and emulator");
0269   etMismatchDiff_ = bookEtDiff("etMismatchDiff", "ET difference (data-emulator) for ET mismatches");
0270   fbCorrelation_ =
0271       ibooker.book2D("FbCorrelation", "Feature Bit correlation for BE;Data;Emulator", 16, -0.5, 15.5, 16, -0.5, 15.5);
0272   fbCorrelationHF_ =
0273       ibooker.book2D("FbCorrelationHF", "Feature Bit correlation for HF;Data;Emulator", 16, -0.5, 15.5, 16, -0.5, 15.5);
0274 
0275   ibooker.setCurrentFolder(histFolder_ + "/MismatchDetail");
0276 
0277   const int nLumis = 2000;
0278   etMismatchByLumi_ = ibooker.book1D(
0279       "etMismatchByLumi", "ET Mismatch counts per lumi section;LumiSection;Counts", nLumis, .5, nLumis + 0.5);
0280   erMismatchByLumi_ = ibooker.book1D(
0281       "erMismatchByLumi", "ET Ratio Mismatch counts per lumi section;LumiSection;Counts", nLumis, .5, nLumis + 0.5);
0282   fbMismatchByLumi_ = ibooker.book1D(
0283       "fbMismatchByLumi", "Feature Bit Mismatch counts per lumi section;LumiSection;Counts", nLumis, .5, nLumis + 0.5);
0284 
0285   etMismatchesPerBx_ = ibooker.book1D("etMismatchesPerBx", "ET Mismatch counts per bunch crossing", 3564, -.5, 3563.5);
0286   erMismatchesPerBx_ =
0287       ibooker.book1D("erMismatchesPerBx", "ET Ratio Mismatch counts per bunch crossing", 3564, -.5, 3563.5);
0288   fbMismatchesPerBx_ =
0289       ibooker.book1D("fbMismatchesPerBx", "Feature Bit Mismatch counts per bunch crossing", 3564, -.5, 3563.5);
0290   towerCountMismatchesPerBx_ = ibooker.book1D(
0291       "towerCountMismatchesPerBx", "CaloTower size mismatch counts per bunch crossing", 3564, -.5, 3563.5);
0292 
0293   last20Mismatches_ =
0294       ibooker.book2D("last20Mismatches", "Log of last 20 mismatches (use json tool to copy/paste)", 4, 0, 4, 20, 0, 20);
0295   last20Mismatches_->setBinLabel(1, "Et Mismatch");
0296   last20Mismatches_->setBinLabel(2, "Et ratio Mismatch");
0297   last20Mismatches_->setBinLabel(3, "Feature bit Mismatch");
0298   last20Mismatches_->setBinLabel(4, "-");
0299   for (size_t i = 0; i < last20MismatchArray_.size(); ++i)
0300     last20MismatchArray_[i] = {"-" + std::to_string(i), 0};
0301   for (size_t i = 1; i <= 20; ++i)
0302     last20Mismatches_->setBinLabel(i, "-" + std::to_string(i), /* axis */ 2);
0303 }