Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:08:01

0001 #include "DQM/L1TMonitorClient/interface/L1TTestsSummary.h"
0002 
0003 #include "FWCore/ServiceRegistry/interface/Service.h"
0004 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0005 #include "FWCore/Framework/interface/ESHandle.h"
0006 #include "FWCore/Framework/interface/EventSetup.h"
0007 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0008 #include "DQMServices/Core/interface/DQMStore.h"
0009 #include <cstdio>
0010 #include <sstream>
0011 #include <cmath>
0012 #include <vector>
0013 #include <TMath.h>
0014 #include <climits>
0015 #include <TFile.h>
0016 #include <TDirectory.h>
0017 #include <TProfile.h>
0018 
0019 using namespace std;
0020 using namespace edm;
0021 
0022 //____________________________________________________________________________
0023 // Function: L1TTestsSummary
0024 // Description: This is the constructor, basic variable initialization
0025 // Inputs:
0026 // * const edm::ParameterSet& ps = Parameter for this analyzer
0027 //____________________________________________________________________________
0028 L1TTestsSummary::L1TTestsSummary(const edm::ParameterSet &ps) {
0029   if (mVerbose) {
0030     cout << "[L1TTestsSummary:] Called constructor" << endl;
0031   }
0032 
0033   // Get parameters
0034   mParameters = ps;
0035   mVerbose = ps.getUntrackedParameter<bool>("verbose", true);
0036   mMonitorL1TRate = ps.getUntrackedParameter<bool>("MonitorL1TRate", true);
0037   mMonitorL1TSync = ps.getUntrackedParameter<bool>("MonitorL1TSync", true);
0038   mMonitorL1TOccupancy = ps.getUntrackedParameter<bool>("MonitorL1TOccupancy", true);
0039 
0040   mL1TRatePath = ps.getUntrackedParameter<string>("L1TRatePath", "L1T/L1TRate/Certification/");
0041   mL1TSyncPath = ps.getUntrackedParameter<string>("L1TSyncPath", "L1T/L1TSync/Certification/");
0042   mL1TOccupancyPath = ps.getUntrackedParameter<string>("L1TOccupancyPath", "L1T/L1TOccupancy/Certification/");
0043 }
0044 
0045 //____________________________________________________________________________
0046 // Function: ~L1TTestsSummary
0047 // Description: This is the destructor, basic variable deletion
0048 //____________________________________________________________________________
0049 L1TTestsSummary::~L1TTestsSummary() {
0050   if (mVerbose) {
0051     cout << "[L1TTestsSummary:] Called destructor" << endl;
0052   }
0053 }
0054 
0055 //____________________________________________________________________________
0056 // Function: beginRun
0057 // Description: This is will be run at the begining of each run
0058 // Inputs:
0059 // * const Run&        r       = Run information
0060 // * const EventSetup& context = Event Setup information
0061 //____________________________________________________________________________
0062 void L1TTestsSummary::book(DQMStore::IBooker &ibooker, DQMStore::IGetter &igetter) {
0063   if (mVerbose) {
0064     cout << "[L1TTestsSummary:] Called beginRun" << endl;
0065   }
0066 
0067   int maxLS = 2500;
0068 
0069   if (mMonitorL1TRate) {
0070     if (mVerbose) {
0071       cout << "[L1TTestsSummary:] Initializing L1TRate Module Monitoring" << endl;
0072     }
0073 
0074     igetter.setCurrentFolder(mL1TRatePath);
0075     vector<string> histToMonitor = igetter.getMEs();
0076     int histLines = histToMonitor.size() + 1;
0077 
0078     ibooker.setCurrentFolder("L1T/L1TTestsSummary/");
0079     mL1TRateMonitor = ibooker.book2D(
0080         "RateQualitySummary", "L1T Rates Monitor Summary", maxLS, +0.5, double(maxLS) + 0.5, histLines, 0, histLines);
0081     mL1TRateMonitor->setAxisTitle("Lumi Section", 1);
0082 
0083     mL1TRateMonitor->setBinLabel(1, "Summary", 2);
0084     for (unsigned int i = 0; i < histToMonitor.size(); i++) {
0085       string name = igetter.get(mL1TRatePath + histToMonitor[i])->getTH1()->GetName();
0086       mL1TRateMonitor->setBinLabel(i + 2, name, 2);
0087     }
0088   }
0089   if (mMonitorL1TSync) {
0090     if (mVerbose) {
0091       cout << "[L1TTestsSummary:] Initializing L1TSync Module Monitoring" << endl;
0092     }
0093 
0094     igetter.setCurrentFolder(mL1TSyncPath);
0095     vector<string> histToMonitor = igetter.getMEs();
0096     int histLines = histToMonitor.size() + 1;
0097 
0098     ibooker.setCurrentFolder("L1T/L1TTestsSummary/");
0099     mL1TSyncMonitor = ibooker.book2D("SyncQualitySummary",
0100                                      "L1T Synchronization Monitor Summary",
0101                                      maxLS,
0102                                      0.5,
0103                                      double(maxLS) + 0.5,
0104                                      histLines,
0105                                      0,
0106                                      histLines);
0107     mL1TSyncMonitor->setAxisTitle("Lumi Section", 1);
0108 
0109     mL1TSyncMonitor->setBinLabel(1, "Summary", 2);
0110     for (unsigned int i = 0; i < histToMonitor.size(); i++) {
0111       string name = igetter.get(mL1TSyncPath + histToMonitor[i])->getTH1()->GetName();
0112       mL1TSyncMonitor->setBinLabel(i + 2, name, 2);
0113     }
0114   }
0115   if (mMonitorL1TOccupancy) {
0116     if (mVerbose) {
0117       cout << "[L1TTestsSummary:] Initializing L1TOccupancy Module Monitoring" << endl;
0118     }
0119 
0120     igetter.setCurrentFolder(mL1TOccupancyPath);
0121     vector<string> histToMonitor = igetter.getMEs();
0122     int histLines = histToMonitor.size() + 1;
0123 
0124     ibooker.setCurrentFolder("L1T/L1TTestsSummary/");
0125     mL1TOccupancyMonitor = ibooker.book2D(
0126         "OccupancySummary", "L1T Occupancy Monitor Summary", maxLS, +0.5, double(maxLS) + 0.5, histLines, 0, histLines);
0127     mL1TOccupancyMonitor->setAxisTitle("Lumi Section", 1);
0128 
0129     mL1TOccupancyMonitor->setBinLabel(1, "Summary", 2);
0130     for (unsigned int i = 0; i < histToMonitor.size(); i++) {
0131       string name = igetter.get(mL1TOccupancyPath + histToMonitor[i])->getTH1()->GetName();
0132       mL1TOccupancyMonitor->setBinLabel(i + 2, name, 2);
0133     }
0134   }
0135 
0136   //-> Making the summary of summaries
0137   int testsToMonitor = 1;
0138   if (mMonitorL1TRate) {
0139     testsToMonitor++;
0140   }
0141   if (mMonitorL1TSync) {
0142     testsToMonitor++;
0143   }
0144   if (mMonitorL1TOccupancy) {
0145     testsToMonitor++;
0146   }
0147 
0148   // Creating
0149   ibooker.setCurrentFolder("L1T/L1TTestsSummary/");
0150   mL1TSummary = ibooker.book2D(
0151       "L1TQualitySummary", "L1 Tests Summary", maxLS, +0.5, double(maxLS) + 0.5, testsToMonitor, 0, testsToMonitor);
0152   mL1TSummary->setAxisTitle("Lumi Section", 1);
0153   mL1TSummary->setBinLabel(1, "L1T Summary", 2);
0154 
0155   int it = 2;
0156   if (mMonitorL1TRate) {
0157     mL1TSummary->setBinLabel(it, "Rates", 2);
0158     binYRate = it;
0159     it++;
0160   }
0161   if (mMonitorL1TSync) {
0162     mL1TSummary->setBinLabel(it, "Synchronization", 2);
0163     binYSync = it;
0164     it++;
0165   }
0166   if (mMonitorL1TOccupancy) {
0167     mL1TSummary->setBinLabel(it, "Occupancy", 2);
0168     binYOccpancy = it;
0169   }
0170 }
0171 
0172 //____________________________________________________________________________
0173 // Function: endRun
0174 // Description: This is will be run at the end of each run
0175 // Inputs:
0176 // * const Run&        r       = Run information
0177 // * const EventSetup& context = Event Setup information
0178 //____________________________________________________________________________
0179 void L1TTestsSummary::dqmEndJob(DQMStore::IBooker &ibooker, DQMStore::IGetter &igetter) {
0180   book(ibooker, igetter);
0181 
0182   if (mVerbose) {
0183     cout << "[L1TTestsSummary:] Called endRun()" << endl;
0184   }
0185 
0186   if (mMonitorL1TRate) {
0187     updateL1TRateMonitor(ibooker, igetter);
0188   }
0189   if (mMonitorL1TSync) {
0190     updateL1TSyncMonitor(ibooker, igetter);
0191   }
0192   if (mMonitorL1TOccupancy) {
0193     updateL1TOccupancyMonitor(ibooker, igetter);
0194   }
0195   updateL1TSummary(ibooker, igetter);
0196 }
0197 
0198 //____________________________________________________________________________
0199 // Function: endLuminosityBlock
0200 // Description: This is will be run at the end of each luminosity block
0201 // Inputs:
0202 // * const LuminosityBlock& lumiSeg = Luminosity Block information
0203 // * const EventSetup&      context = Event Setup information
0204 //____________________________________________________________________________
0205 void L1TTestsSummary::dqmEndLuminosityBlock(DQMStore::IBooker &ibooker,
0206                                             DQMStore::IGetter &igetter,
0207                                             const edm::LuminosityBlock &lumiSeg,
0208                                             const edm::EventSetup &c) {
0209   int eventLS = lumiSeg.id().luminosityBlock();
0210 
0211   book(ibooker, igetter);
0212 
0213   mProcessedLS.push_back(eventLS);
0214 
0215   if (mVerbose) {
0216     cout << "[L1TTestsSummary:] Called endLuminosityBlock()" << endl;
0217     cout << "[L1TTestsSummary:] Lumisection: " << eventLS << endl;
0218   }
0219 
0220   if (mMonitorL1TRate) {
0221     updateL1TRateMonitor(ibooker, igetter);
0222   }
0223   if (mMonitorL1TSync) {
0224     updateL1TSyncMonitor(ibooker, igetter);
0225   }
0226   if (mMonitorL1TOccupancy) {
0227     updateL1TOccupancyMonitor(ibooker, igetter);
0228   }
0229   updateL1TSummary(ibooker, igetter);
0230 }
0231 
0232 //____________________________________________________________________________
0233 // Function: updateL1TRateMonitor
0234 // Description:
0235 //____________________________________________________________________________
0236 void L1TTestsSummary::updateL1TRateMonitor(DQMStore::IBooker &ibooker, DQMStore::IGetter &igetter) {
0237   igetter.setCurrentFolder(mL1TRatePath);
0238   vector<string> histToMonitor = igetter.getMEs();
0239 
0240   for (unsigned int i = 0; i < histToMonitor.size(); i++) {
0241     MonitorElement *me = igetter.get(mL1TRatePath + histToMonitor[i]);
0242     if (mVerbose) {
0243       cout << "[L1TTestsSummary:] Found ME: " << me->getTH1()->GetName() << endl;
0244     }
0245 
0246     const QReport *myQReport = me->getQReport("L1TRateTest");  //get QReport associated to your ME
0247     if (myQReport) {
0248       float qtresult = myQReport->getQTresult();          // get QT result value
0249       int qtstatus = myQReport->getStatus();              // get QT status value (see table below)
0250       const string &qtmessage = myQReport->getMessage();  // get the whole QT result message
0251       vector<DQMChannel> qtBadChannels = myQReport->getBadChannels();
0252 
0253       if (mVerbose) {
0254         cout << "[L1TTestsSummary:] Found QReport for ME: " << me->getTH1()->GetName() << endl;
0255         cout << "[L1TTestsSummary:] Result=" << qtresult << " status=" << qtstatus << " message=" << qtmessage << endl;
0256         cout << "[L1TTestsSummary:] Bad Channels size=" << qtBadChannels.size() << endl;
0257       }
0258 
0259       for (unsigned int i = 0; i < mProcessedLS.size() - 1; i++) {
0260         int binx = mL1TRateMonitor->getTH2F()->GetXaxis()->FindBin(mProcessedLS[i]);
0261         int biny = mL1TRateMonitor->getTH2F()->GetYaxis()->FindBin(me->getTH1()->GetName());
0262         mL1TRateMonitor->setBinContent(binx, biny, 100);
0263       }
0264 
0265       for (unsigned int a = 0; a < qtBadChannels.size(); a++) {
0266         for (unsigned int b = 0; b < mProcessedLS.size() - 1; b++) {
0267           // Converting bin to value
0268           double valueBinBad = me->getTH1()->GetBinCenter(qtBadChannels[a].getBin());
0269 
0270           if (valueBinBad == (mProcessedLS[b])) {
0271             int binx = mL1TRateMonitor->getTH2F()->GetXaxis()->FindBin(valueBinBad);
0272             int biny = mL1TRateMonitor->getTH2F()->GetYaxis()->FindBin(me->getTH1()->GetName());
0273             mL1TRateMonitor->setBinContent(binx, biny, 300);
0274           }
0275         }
0276       }
0277     }
0278   }
0279 
0280   //-> Filling the summaries
0281   int nBinX = mL1TRateMonitor->getTH2F()->GetXaxis()->GetNbins();
0282   int nBinY = mL1TRateMonitor->getTH2F()->GetYaxis()->GetNbins();
0283   for (int binx = 1; binx <= nBinX; binx++) {
0284     int GlobalStatus = 0;
0285     for (int biny = 2; biny <= nBinY; biny++) {
0286       double flag = mL1TRateMonitor->getBinContent(binx, biny);
0287       if (GlobalStatus < flag) {
0288         GlobalStatus = flag;
0289       }
0290     }
0291 
0292     // NOTE: Assumes mL1TSummary has same size then mL1TRateMonitor
0293     mL1TRateMonitor->setBinContent(binx, 1, GlobalStatus);
0294     mL1TSummary->setBinContent(binx, binYRate, GlobalStatus);
0295   }
0296 }
0297 
0298 //____________________________________________________________________________
0299 // Function: updateL1TSyncMonitor
0300 // Description:
0301 // Note: Values certified by L1TRates are always about currentLS-1 since we
0302 //       use rate averages from the previous LS from GT
0303 //____________________________________________________________________________
0304 void L1TTestsSummary::updateL1TSyncMonitor(DQMStore::IBooker &ibooker, DQMStore::IGetter &igetter) {
0305   igetter.setCurrentFolder(mL1TSyncPath);
0306   vector<string> histToMonitor = igetter.getMEs();
0307 
0308   for (unsigned int i = 0; i < histToMonitor.size(); i++) {
0309     MonitorElement *me = igetter.get(mL1TSyncPath + histToMonitor[i]);
0310     if (mVerbose) {
0311       cout << "[L1TTestsSummary:] Found ME: " << me->getTH1()->GetName() << endl;
0312     }
0313 
0314     const QReport *myQReport = me->getQReport("L1TSyncTest");  //get QReport associated to your ME
0315     if (myQReport) {
0316       float qtresult = myQReport->getQTresult();          // get QT result value
0317       int qtstatus = myQReport->getStatus();              // get QT status value (see table below)
0318       const string &qtmessage = myQReport->getMessage();  // get the whole QT result message
0319       vector<DQMChannel> qtBadChannels = myQReport->getBadChannels();
0320 
0321       if (mVerbose) {
0322         cout << "[L1TTestsSummary:] Found QReport for ME: " << me->getTH1()->GetName() << endl;
0323         cout << "[L1TTestsSummary:] Result=" << qtresult << " status=" << qtstatus << " message=" << qtmessage << endl;
0324         cout << "[L1TTestsSummary:] Bad Channels size=" << qtBadChannels.size() << endl;
0325       }
0326 
0327       for (unsigned int i = 0; i < mProcessedLS.size(); i++) {
0328         int binx = mL1TSyncMonitor->getTH2F()->GetXaxis()->FindBin(mProcessedLS[i]);
0329         int biny = mL1TSyncMonitor->getTH2F()->GetYaxis()->FindBin(me->getTH1()->GetName());
0330         mL1TSyncMonitor->setBinContent(binx, biny, 100);
0331       }
0332 
0333       for (unsigned int a = 0; a < qtBadChannels.size(); a++) {
0334         for (unsigned int b = 0; b < mProcessedLS.size(); b++) {
0335           // Converting bin to value
0336           double valueBinBad = me->getTH1()->GetBinCenter(qtBadChannels[a].getBin());
0337 
0338           if (valueBinBad == mProcessedLS[b]) {
0339             int binx = mL1TSyncMonitor->getTH2F()->GetXaxis()->FindBin(valueBinBad);
0340             int biny = mL1TSyncMonitor->getTH2F()->GetYaxis()->FindBin(me->getTH1()->GetName());
0341             mL1TSyncMonitor->setBinContent(binx, biny, 300);
0342           }
0343         }
0344       }
0345     }
0346   }
0347 
0348   //-> Filling the summaries
0349   int nBinX = mL1TSyncMonitor->getTH2F()->GetXaxis()->GetNbins();
0350   int nBinY = mL1TSyncMonitor->getTH2F()->GetYaxis()->GetNbins();
0351   for (int binx = 1; binx <= nBinX; binx++) {
0352     int GlobalStatus = 0;
0353     for (int biny = 2; biny <= nBinY; biny++) {
0354       double flag = mL1TSyncMonitor->getBinContent(binx, biny);
0355       if (GlobalStatus < flag) {
0356         GlobalStatus = flag;
0357       }
0358     }
0359 
0360     // NOTE: Assumes mL1TSummary has same size then mL1TSyncMonitor
0361     mL1TSyncMonitor->setBinContent(binx, 1, GlobalStatus);
0362     mL1TSummary->setBinContent(binx, binYSync, GlobalStatus);
0363   }
0364 }
0365 
0366 //____________________________________________________________________________
0367 // Function: updateL1TOccupancyMonitor
0368 // Description:
0369 //____________________________________________________________________________
0370 void L1TTestsSummary::updateL1TOccupancyMonitor(DQMStore::IBooker &ibooker, DQMStore::IGetter &igetter) {
0371   igetter.setCurrentFolder(mL1TOccupancyPath);
0372   vector<string> histToMonitor = igetter.getMEs();
0373 
0374   for (unsigned int i = 0; i < histToMonitor.size(); i++) {
0375     MonitorElement *me = igetter.get(mL1TOccupancyPath + histToMonitor[i]);
0376     if (mVerbose) {
0377       cout << "[L1TTestsSummary:] Found ME: " << me->getTH1()->GetName() << endl;
0378     }
0379 
0380     const QReport *myQReport = me->getQReport("L1TOccupancyTest");  //get QReport associated to your ME
0381     if (myQReport) {
0382       float qtresult = myQReport->getQTresult();          // get QT result value
0383       int qtstatus = myQReport->getStatus();              // get QT status value (see table below)
0384       const string &qtmessage = myQReport->getMessage();  // get the whole QT result message
0385       vector<DQMChannel> qtBadChannels = myQReport->getBadChannels();
0386 
0387       if (mVerbose) {
0388         cout << "[L1TTestsSummary:] Found QReport for ME: " << me->getTH1()->GetName() << endl;
0389         cout << "[L1TTestsSummary:] Result=" << qtresult << " status=" << qtstatus << " message=" << qtmessage << endl;
0390         cout << "[L1TTestsSummary:] Bad Channels size=" << qtBadChannels.size() << endl;
0391       }
0392 
0393       for (unsigned int i = 0; i < mProcessedLS.size(); i++) {
0394         int binx = mL1TOccupancyMonitor->getTH2F()->GetXaxis()->FindBin(mProcessedLS[i]);
0395         int biny = mL1TOccupancyMonitor->getTH2F()->GetYaxis()->FindBin(me->getTH1()->GetName());
0396         mL1TOccupancyMonitor->setBinContent(binx, biny, 100);
0397       }
0398 
0399       for (unsigned int a = 0; a < qtBadChannels.size(); a++) {
0400         for (unsigned int b = 0; b < mProcessedLS.size(); b++) {
0401           // Converting bin to value
0402           double valueBinBad = me->getTH1()->GetBinCenter(qtBadChannels[a].getBin());
0403 
0404           if (valueBinBad == mProcessedLS[b]) {
0405             int binx = mL1TOccupancyMonitor->getTH2F()->GetXaxis()->FindBin(valueBinBad);
0406             int biny = mL1TOccupancyMonitor->getTH2F()->GetYaxis()->FindBin(me->getTH1()->GetName());
0407             mL1TOccupancyMonitor->setBinContent(binx, biny, 300);
0408           }
0409         }
0410       }
0411     }
0412   }
0413 
0414   //-> Filling the summaries
0415   int nBinX = mL1TOccupancyMonitor->getTH2F()->GetXaxis()->GetNbins();
0416   int nBinY = mL1TOccupancyMonitor->getTH2F()->GetYaxis()->GetNbins();
0417   for (int binx = 1; binx <= nBinX; binx++) {
0418     int GlobalStatus = 0;
0419     for (int biny = 2; biny <= nBinY; biny++) {
0420       double flag = mL1TOccupancyMonitor->getBinContent(binx, biny);
0421       if (GlobalStatus < flag) {
0422         GlobalStatus = flag;
0423       }
0424     }
0425 
0426     // NOTE: Assumes mL1TSummary has same size then mL1TOccupancyMonitor
0427     mL1TOccupancyMonitor->setBinContent(binx, 1, GlobalStatus);
0428     mL1TSummary->setBinContent(binx, binYOccpancy, GlobalStatus);
0429   }
0430 }
0431 
0432 //____________________________________________________________________________
0433 // Function: updateL1TOccupancyMonitor
0434 // Description:
0435 //____________________________________________________________________________
0436 void L1TTestsSummary::updateL1TSummary(DQMStore::IBooker &ibooker, DQMStore::IGetter &igetter) {
0437   int nBinX = mL1TSummary->getTH2F()->GetXaxis()->GetNbins();
0438   for (int binx = 1; binx <= nBinX; binx++) {
0439     int GlobalStatus = 0;
0440     if (mMonitorL1TRate) {
0441       if (mL1TSummary->getBinContent(binx, binYRate) > GlobalStatus) {
0442         GlobalStatus = mL1TSummary->getBinContent(binx, binYRate);
0443       }
0444     }
0445     if (mMonitorL1TSync) {
0446       if (mL1TSummary->getBinContent(binx, binYSync) > GlobalStatus) {
0447         GlobalStatus = mL1TSummary->getBinContent(binx, binYSync);
0448       }
0449     }
0450     if (mMonitorL1TOccupancy) {
0451       if (mL1TSummary->getBinContent(binx, binYOccpancy) > GlobalStatus) {
0452         GlobalStatus = mL1TSummary->getBinContent(binx, binYOccpancy);
0453       }
0454     }
0455     mL1TSummary->setBinContent(binx, 1, GlobalStatus);
0456   }
0457 }
0458 
0459 //define this as a plug-in
0460 DEFINE_FWK_MODULE(L1TTestsSummary);