Back to home page

Project CMSSW displayed by LXR

 
 

    


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

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