Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:10:06

0001 /*
0002  * \file DQMFEDIntegrityClient.cc
0003  * \author M. Marienfeld
0004  * Last Update:
0005  *
0006  * Description: Summing up FED entries from all subdetectors.
0007  *
0008  */
0009 
0010 #include <string>
0011 #include <vector>
0012 
0013 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0014 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0015 #include "DQMServices/Core/interface/DQMStore.h"
0016 #include "FWCore/ServiceRegistry/interface/Service.h"
0017 
0018 //
0019 // class declaration
0020 //
0021 
0022 class DQMFEDIntegrityClient : public edm::one::EDAnalyzer<edm::one::WatchRuns, edm::one::WatchLuminosityBlocks> {
0023 public:
0024   typedef dqm::legacy::DQMStore DQMStore;
0025   typedef dqm::legacy::MonitorElement MonitorElement;
0026   DQMFEDIntegrityClient(const edm::ParameterSet&);
0027   ~DQMFEDIntegrityClient() override = default;
0028 
0029 protected:
0030   void beginJob() override;
0031   void beginRun(const edm::Run& r, const edm::EventSetup& c) override;
0032 
0033   /// Analyze
0034   void analyze(const edm::Event& e, const edm::EventSetup& c) override;
0035 
0036   void beginLuminosityBlock(const edm::LuminosityBlock& l, const edm::EventSetup& c) override;
0037   void endLuminosityBlock(const edm::LuminosityBlock& l, const edm::EventSetup& c) override;
0038 
0039   void endRun(const edm::Run& r, const edm::EventSetup& c) override;
0040   void endJob() override;
0041 
0042 private:
0043   void initialize();
0044   void fillHistograms();
0045 
0046   edm::ParameterSet parameters_;
0047 
0048   DQMStore* dbe_;
0049 
0050   // ---------- member data ----------
0051 
0052   int NBINS;
0053   float XMIN, XMAX;
0054   float SummaryContent[10];
0055 
0056   MonitorElement* FedEntries;
0057   MonitorElement* FedFatal;
0058   MonitorElement* FedNonFatal;
0059 
0060   MonitorElement* reportSummary;
0061   MonitorElement* reportSummaryContent[10];
0062   MonitorElement* reportSummaryMap;
0063 
0064   bool fillInEventloop;
0065   bool fillOnEndRun;
0066   bool fillOnEndJob;
0067   bool fillOnEndLumi;
0068   std::string moduleName;
0069   std::string fedFolderName;
0070 };
0071 
0072 // -----------------------------
0073 //  constructors and destructor
0074 // -----------------------------
0075 
0076 DQMFEDIntegrityClient::DQMFEDIntegrityClient(const edm::ParameterSet& ps) {
0077   parameters_ = ps;
0078   initialize();
0079   fillInEventloop = ps.getUntrackedParameter<bool>("fillInEventloop", false);
0080   fillOnEndRun = ps.getUntrackedParameter<bool>("fillOnEndRun", false);
0081   fillOnEndJob = ps.getUntrackedParameter<bool>("fillOnEndJob", false);
0082   fillOnEndLumi = ps.getUntrackedParameter<bool>("fillOnEndLumi", true);
0083   moduleName = ps.getUntrackedParameter<std::string>("moduleName", "FED");
0084   fedFolderName = ps.getUntrackedParameter<std::string>("fedFolderName", "FEDIntegrity");
0085 }
0086 
0087 void DQMFEDIntegrityClient::initialize() {
0088   // get back-end interface
0089   dbe_ = edm::Service<DQMStore>().operator->();
0090 }
0091 
0092 void DQMFEDIntegrityClient::beginJob() {
0093   NBINS = 850;
0094   XMIN = 0.;
0095   XMAX = 850.;
0096 
0097   dbe_ = edm::Service<DQMStore>().operator->();
0098 
0099   // ----------------------------------------------------------------------------------
0100   std::string currentFolder = moduleName + "/" + fedFolderName;
0101   dbe_->setCurrentFolder(currentFolder);
0102 
0103   FedEntries = dbe_->book1D("FedEntries", "FED Entries", NBINS, XMIN, XMAX);
0104   FedFatal = dbe_->book1D("FedFatal", "FED Fatal Errors", NBINS, XMIN, XMAX);
0105   FedNonFatal = dbe_->book1D("FedNonFatal", "FED Non Fatal Errors", NBINS, XMIN, XMAX);
0106 
0107   FedEntries->setAxisTitle("", 1);
0108   FedFatal->setAxisTitle("", 1);
0109   FedNonFatal->setAxisTitle("", 1);
0110 
0111   FedEntries->setAxisTitle("", 2);
0112   FedFatal->setAxisTitle("", 2);
0113   FedNonFatal->setAxisTitle("", 2);
0114 
0115   FedEntries->setBinLabel(11, "PIXEL", 1);
0116   FedEntries->setBinLabel(221, "SIST", 1);
0117   FedEntries->setBinLabel(606, "EE", 1);
0118   FedEntries->setBinLabel(628, "EB", 1);
0119   FedEntries->setBinLabel(651, "EE", 1);
0120   FedEntries->setBinLabel(550, "ES", 1);
0121   FedEntries->setBinLabel(716, "HCAL", 1);
0122   FedEntries->setBinLabel(754, "CSC", 1);
0123   FedEntries->setBinLabel(772, "DT", 1);
0124   FedEntries->setBinLabel(791, "RPC", 1);
0125   FedEntries->setBinLabel(804, "L1T", 1);
0126 
0127   FedFatal->setBinLabel(11, "PIXEL", 1);
0128   FedFatal->setBinLabel(221, "SIST", 1);
0129   FedFatal->setBinLabel(606, "EE", 1);
0130   FedFatal->setBinLabel(628, "EB", 1);
0131   FedFatal->setBinLabel(651, "EE", 1);
0132   FedFatal->setBinLabel(550, "ES", 1);
0133   FedFatal->setBinLabel(716, "HCAL", 1);
0134   FedFatal->setBinLabel(754, "CSC", 1);
0135   FedFatal->setBinLabel(772, "DT", 1);
0136   FedFatal->setBinLabel(791, "RPC", 1);
0137   FedFatal->setBinLabel(804, "L1T", 1);
0138 
0139   FedNonFatal->setBinLabel(11, "PIXEL", 1);
0140   FedNonFatal->setBinLabel(221, "SIST", 1);
0141   FedNonFatal->setBinLabel(606, "EE", 1);
0142   FedNonFatal->setBinLabel(628, "EB", 1);
0143   FedNonFatal->setBinLabel(651, "EE", 1);
0144   FedNonFatal->setBinLabel(550, "ES", 1);
0145   FedNonFatal->setBinLabel(716, "HCAL", 1);
0146   FedNonFatal->setBinLabel(754, "CSC", 1);
0147   FedNonFatal->setBinLabel(772, "DT", 1);
0148   FedNonFatal->setBinLabel(791, "RPC", 1);
0149   FedNonFatal->setBinLabel(804, "L1T", 1);
0150 
0151   //-----------------------------------------------------------------------------------
0152   currentFolder = moduleName + "/EventInfo";
0153   dbe_->setCurrentFolder(currentFolder);
0154 
0155   reportSummary = dbe_->bookFloat("reportSummary");
0156 
0157   int nSubsystems = 10;
0158 
0159   if (reportSummary)
0160     reportSummary->Fill(1.);
0161 
0162   currentFolder = moduleName + "/EventInfo/reportSummaryContents";
0163   dbe_->setCurrentFolder(currentFolder);
0164 
0165   reportSummaryContent[0] = dbe_->bookFloat("CSC FEDs");
0166   reportSummaryContent[1] = dbe_->bookFloat("DT FEDs");
0167   reportSummaryContent[2] = dbe_->bookFloat("EB FEDs");
0168   reportSummaryContent[3] = dbe_->bookFloat("EE FEDs");
0169   reportSummaryContent[4] = dbe_->bookFloat("ES FEDs");
0170   reportSummaryContent[5] = dbe_->bookFloat("Hcal FEDs");
0171   reportSummaryContent[6] = dbe_->bookFloat("L1T FEDs");
0172   reportSummaryContent[7] = dbe_->bookFloat("Pixel FEDs");
0173   reportSummaryContent[8] = dbe_->bookFloat("RPC FEDs");
0174   reportSummaryContent[9] = dbe_->bookFloat("SiStrip FEDs");
0175 
0176   // initialize reportSummaryContents to 1
0177   for (int i = 0; i < nSubsystems; ++i) {
0178     SummaryContent[i] = 1.;
0179     reportSummaryContent[i]->Fill(1.);
0180   }
0181 
0182   currentFolder = moduleName + "/EventInfo";
0183   dbe_->setCurrentFolder(currentFolder);
0184 
0185   reportSummaryMap = dbe_->book2D("reportSummaryMap", "FED Report Summary Map", 1, 1, 2, 10, 1, 11);
0186 
0187   reportSummaryMap->setAxisTitle("", 1);
0188   reportSummaryMap->setAxisTitle("", 2);
0189 
0190   reportSummaryMap->setBinLabel(1, " ", 1);
0191   reportSummaryMap->setBinLabel(10, "CSC", 2);
0192   reportSummaryMap->setBinLabel(9, "DT", 2);
0193   reportSummaryMap->setBinLabel(8, "EB", 2);
0194   reportSummaryMap->setBinLabel(7, "EE", 2);
0195   reportSummaryMap->setBinLabel(6, "ES", 2);
0196   reportSummaryMap->setBinLabel(5, "Hcal", 2);
0197   reportSummaryMap->setBinLabel(4, "L1T", 2);
0198   reportSummaryMap->setBinLabel(3, "Pixel", 2);
0199   reportSummaryMap->setBinLabel(2, "RPC", 2);
0200   reportSummaryMap->setBinLabel(1, "SiStrip", 2);
0201 }
0202 
0203 void DQMFEDIntegrityClient::beginRun(const edm::Run& r, const edm::EventSetup& context) {}
0204 
0205 void DQMFEDIntegrityClient::analyze(const edm::Event& e, const edm::EventSetup& context) {
0206   if (fillInEventloop)
0207     fillHistograms();
0208 }
0209 
0210 void DQMFEDIntegrityClient::beginLuminosityBlock(const edm::LuminosityBlock& l, const edm::EventSetup& c) {}
0211 
0212 void DQMFEDIntegrityClient::endLuminosityBlock(const edm::LuminosityBlock& lumiBlock, const edm::EventSetup& context) {
0213   if (fillOnEndLumi)
0214     fillHistograms();
0215 }
0216 
0217 void DQMFEDIntegrityClient::fillHistograms() {
0218   // FED Entries
0219 
0220   std::vector<std::string> entries;
0221   entries.push_back("CSC/" + fedFolderName + "/FEDEntries");
0222   entries.push_back("DT/" + fedFolderName + "/FEDEntries");
0223   entries.push_back("EcalBarrel/" + fedFolderName + "/FEDEntries");
0224   entries.push_back("EcalEndcap/" + fedFolderName + "/FEDEntries");
0225   entries.push_back("EcalPreshower/" + fedFolderName + "/FEDEntries");
0226   entries.push_back("Hcal/" + fedFolderName + "/FEDEntries");
0227   entries.push_back("L1T/" + fedFolderName + "/FEDEntries");
0228   entries.push_back("Pixel/" + fedFolderName + "/FEDEntries");
0229   entries.push_back("RPC/" + fedFolderName + "/FEDEntries");
0230   entries.push_back("SiStrip/" + fedFolderName + "/FEDEntries");
0231 
0232   for (auto ent = entries.begin(); ent != entries.end(); ++ent) {
0233     if (!(dbe_->get(*ent))) {
0234       continue;
0235     }
0236 
0237     MonitorElement* me = dbe_->get(*ent);
0238     if (TH1F* rootHisto = me->getTH1F()) {
0239       int Nbins = me->getNbinsX();
0240       float entry = 0.;
0241       int xmin = (int)rootHisto->GetXaxis()->GetXmin();
0242       if (*ent == "L1T/" + fedFolderName + "/FEDEntries")
0243         xmin = xmin + 800;
0244       if (*ent == "DT/" + fedFolderName + "/FEDEntries")
0245         xmin = 770;  //Real DT FEDIDs are 1369-1371
0246 
0247       for (int bin = 1; bin <= Nbins; ++bin) {
0248         int id = xmin + bin;
0249         entry = rootHisto->GetBinContent(bin);
0250         if (entry > 0.)
0251           FedEntries->setBinContent(id, entry);
0252       }
0253     }
0254   }
0255 
0256   // FED Fatal
0257 
0258   int nSubsystems = 10;
0259 
0260   std::vector<std::string> fatal;
0261   fatal.push_back("CSC/" + fedFolderName + "/FEDFatal");
0262   fatal.push_back("DT/" + fedFolderName + "/FEDFatal");
0263   fatal.push_back("EcalBarrel/" + fedFolderName + "/FEDFatal");
0264   fatal.push_back("EcalEndcap/" + fedFolderName + "/FEDFatal");
0265   fatal.push_back("EcalPreshower/" + fedFolderName + "/FEDFatal");
0266   fatal.push_back("Hcal/" + fedFolderName + "/FEDFatal");
0267   fatal.push_back("L1T/" + fedFolderName + "/FEDFatal");
0268   fatal.push_back("Pixel/" + fedFolderName + "/FEDFatal");
0269   fatal.push_back("RPC/" + fedFolderName + "/FEDFatal");
0270   fatal.push_back("SiStrip/" + fedFolderName + "/FEDFatal");
0271 
0272   int k = 0, count = 0;
0273 
0274   float sum = 0.;
0275 
0276   auto ent = entries.begin();
0277   for (auto fat = fatal.begin(); fat != fatal.end(); ++fat) {
0278     if (!(dbe_->get(*fat))) {
0279       reportSummaryContent[k]->Fill(-1);
0280       reportSummaryMap->setBinContent(1, nSubsystems - k, -1);
0281       k++;
0282       ent++;
0283       continue;
0284     }
0285 
0286     MonitorElement* me = dbe_->get(*fat);
0287     MonitorElement* meNorm = dbe_->get(*ent);
0288 
0289     float entry = 0.;
0290     float norm = 0.;
0291 
0292     if (TH1F* rootHisto = me->getTH1F()) {
0293       if (TH1F* rootHistoNorm = meNorm->getTH1F()) {
0294         int Nbins = me->getNbinsX();
0295         int xmin = (int)rootHisto->GetXaxis()->GetXmin();
0296         if (*fat == "L1T/" + fedFolderName + "/FEDFatal")
0297           xmin = xmin + 800;
0298         if (*fat == "DT/" + fedFolderName + "/FEDFatal")
0299           xmin = 770;  //Real DT FED IDs are 1369-1371
0300 
0301         float binentry = 0.;
0302         for (int bin = 1; bin <= Nbins; ++bin) {
0303           int id = xmin + bin;
0304           binentry = rootHisto->GetBinContent(bin);
0305           entry += binentry;
0306           norm += rootHistoNorm->GetBinContent(bin);
0307           FedFatal->setBinContent(id, binentry);
0308         }
0309       }
0310     }
0311 
0312     if (norm > 0)
0313       SummaryContent[k] = 1.0 - entry / norm;
0314     reportSummaryContent[k]->Fill(SummaryContent[k]);
0315     if ((k == 2 || k == 3)  // for EE and EB only show yellow when more than 1% errors.
0316         && SummaryContent[k] >= 0.95 && SummaryContent[k] < 0.99)
0317       SummaryContent[k] = 0.949;
0318     reportSummaryMap->setBinContent(1, nSubsystems - k, SummaryContent[k]);
0319     sum = sum + SummaryContent[k];
0320 
0321     k++;
0322     ent++;
0323     count++;
0324   }
0325 
0326   if (count > 0)
0327     reportSummary->Fill(sum / (float)count);
0328 
0329   // FED Non Fatal
0330 
0331   std::vector<std::string> nonfatal;
0332   nonfatal.push_back("CSC/" + fedFolderName + "/FEDNonFatal");
0333   nonfatal.push_back("DT/" + fedFolderName + "/FEDNonFatal");
0334   nonfatal.push_back("EcalBarrel/" + fedFolderName + "/FEDNonFatal");
0335   nonfatal.push_back("EcalEndcap/" + fedFolderName + "/FEDNonFatal");
0336   nonfatal.push_back("EcalPreshower/" + fedFolderName + "/FEDNonFatal");
0337   nonfatal.push_back("Hcal/" + fedFolderName + "/FEDNonFatal");
0338   nonfatal.push_back("L1T/" + fedFolderName + "/FEDNonFatal");
0339   nonfatal.push_back("Pixel/" + fedFolderName + "/FEDNonFatal");
0340   nonfatal.push_back("RPC/" + fedFolderName + "/FEDNonFatal");
0341   nonfatal.push_back("SiStrip/" + fedFolderName + "/FEDNonFatal");
0342 
0343   for (auto non = nonfatal.begin(); non != nonfatal.end(); ++non) {
0344     if (!(dbe_->get(*non))) {
0345       continue;
0346     }
0347 
0348     MonitorElement* me = dbe_->get(*non);
0349 
0350     if (TH1F* rootHisto = me->getTH1F()) {
0351       int Nbins = me->getNbinsX();
0352       float entry = 0.;
0353       int xmin = (int)rootHisto->GetXaxis()->GetXmin();
0354       if (*non == "L1T/" + fedFolderName + "/FEDNonFatal")
0355         xmin = xmin + 800;
0356 
0357       for (int bin = 1; bin <= Nbins; ++bin) {
0358         int id = xmin + bin;
0359         entry = rootHisto->GetBinContent(bin);
0360         if (entry > 0.)
0361           FedNonFatal->setBinContent(id, entry);
0362       }
0363     }
0364   }
0365 }
0366 
0367 void DQMFEDIntegrityClient::endRun(const edm::Run& r, const edm::EventSetup& context) {
0368   if (fillOnEndRun)
0369     fillHistograms();
0370 }
0371 
0372 void DQMFEDIntegrityClient::endJob() {
0373   if (fillOnEndJob)
0374     fillHistograms();
0375 }
0376 
0377 #include "FWCore/PluginManager/interface/ModuleDef.h"
0378 #include "FWCore/Framework/interface/MakerMacros.h"
0379 DEFINE_FWK_MODULE(DQMFEDIntegrityClient);