CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch2/src/DQM/L1TMonitorClient/src/L1TTestsSummary.cc

Go to the documentation of this file.
00001 #include "DQM/L1TMonitorClient/interface/L1TTestsSummary.h"
00002 
00003 #include "FWCore/ServiceRegistry/interface/Service.h"
00004 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00005 #include "FWCore/Framework/interface/ESHandle.h"
00006 #include "FWCore/Framework/interface/EventSetup.h"
00007 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00008 #include "DQMServices/Core/interface/QReport.h"
00009 #include "DQMServices/Core/interface/DQMStore.h"
00010 #include "DQMServices/Core/interface/MonitorElement.h"
00011 #include "DQMServices/Core/interface/DQMChannel.h"
00012 #include "DataFormats/Histograms/interface/MEtoEDMFormat.h"
00013 #include <stdio.h>
00014 #include <sstream>
00015 #include <math.h>
00016 #include <vector>
00017 #include <TMath.h>
00018 #include <limits.h>
00019 #include <TFile.h>
00020 #include <TDirectory.h>
00021 #include <TProfile.h>
00022 
00023 using namespace std;
00024 using namespace edm;
00025 
00026 //____________________________________________________________________________
00027 // Function: L1TTestsSummary
00028 // Description: This is the constructor, basic variable initialization
00029 // Inputs: 
00030 // * const edm::ParameterSet& ps = Parameter for this analyzer
00031 //____________________________________________________________________________
00032 L1TTestsSummary::L1TTestsSummary(const edm::ParameterSet& ps){
00033 
00034   if(mVerbose){cout << "[L1TTestsSummary:] Called constructor" << endl;}
00035   
00036   // Get parameters
00037   mParameters          = ps;
00038   mVerbose             = ps.getUntrackedParameter<bool>("verbose"            ,true);
00039   mMonitorL1TRate      = ps.getUntrackedParameter<bool>("MonitorL1TRate"     ,true);
00040   mMonitorL1TSync      = ps.getUntrackedParameter<bool>("MonitorL1TSync"     ,true);
00041   mMonitorL1TOccupancy = ps.getUntrackedParameter<bool>("MonitorL1TOccupancy",true);
00042 
00043   mL1TRatePath      = ps.getUntrackedParameter<string>("L1TRatePath"     ,"L1T/L1TRate/Certification/");
00044   mL1TSyncPath      = ps.getUntrackedParameter<string>("L1TSyncPath"     ,"L1T/L1TSync/Certification/");
00045   mL1TOccupancyPath = ps.getUntrackedParameter<string>("L1TOccupancyPath","L1T/L1TOccupancy/Certification/");
00046   
00047   // Get back-end interface
00048   mDBE  = Service<DQMStore>().operator->();   
00049 
00050 }
00051 
00052 //____________________________________________________________________________
00053 // Function: ~L1TTestsSummary
00054 // Description: This is the destructor, basic variable deletion
00055 //____________________________________________________________________________
00056 L1TTestsSummary::~L1TTestsSummary(){
00057   if(mVerbose){cout << "[L1TTestsSummary:] Called destructor" << endl;}
00058 }
00059 
00060 //____________________________________________________________________________
00061 // Function: beginJob
00062 // Description: This will be run at the beginning of each job
00063 //____________________________________________________________________________
00064 void L1TTestsSummary::beginJob(void){
00065 
00066   if(mVerbose){cout << "[L1TTestsSummary:] Called BeginJob" << endl;}
00067 
00068   // get backend interface  
00069   mDBE = Service<DQMStore>().operator->();
00070 
00071   if (mDBE) {
00072     mDBE->setCurrentFolder("L1T/L1TOccupancy");
00073     mDBE->rmdir("L1T/L1TOccupancy");
00074   }
00075 }
00076 
00077 //____________________________________________________________________________
00078 // Function: endJob
00079 // Description: This will be run at the end of each job
00080 //____________________________________________________________________________
00081 void L1TTestsSummary::endJob(){
00082   if(mVerbose){cout << "[L1TTestsSummary:] Called endJob" << endl;}
00083 }
00084 
00085 //____________________________________________________________________________
00086 // Function: beginRun
00087 // Description: This is will be run at the begining of each run
00088 // Inputs: 
00089 // * const Run&        r       = Run information 
00090 // * const EventSetup& context = Event Setup information
00091 //____________________________________________________________________________
00092 void L1TTestsSummary::beginRun(const Run& r, const EventSetup& context){  
00093   
00094   if(mVerbose){cout << "[L1TTestsSummary:] Called beginRun" << endl;}
00095   
00096   int maxLS = 2500;
00097   
00098   if(mMonitorL1TRate){
00099 
00100     if(mVerbose){cout << "[L1TTestsSummary:] Initializing L1TRate Module Monitoring" << endl;}
00101     
00102     mDBE->setCurrentFolder(mL1TRatePath);
00103     vector<string> histToMonitor = mDBE->getMEs();
00104     int            histLines     = histToMonitor.size()+1;
00105     
00106     mDBE->setCurrentFolder("L1T/L1TTestsSummary/");    
00107     mL1TRateMonitor = mDBE->book2D("RateQualitySummary","L1T Rates Monitor Summary",maxLS,+0.5,double(maxLS)+0.5,histLines,0,histLines); 
00108     mL1TRateMonitor->setAxisTitle("Lumi Section" ,1);
00109 
00110     mL1TRateMonitor->setBinLabel(1,"Summary",2);
00111     for(unsigned int i=0 ; i<histToMonitor.size() ; i++){
00112       string name = mDBE->get(mL1TRatePath+histToMonitor[i])->getTH1()->GetName();
00113       mL1TRateMonitor->setBinLabel(i+2,name,2);  
00114     }
00115   }
00116   if(mMonitorL1TSync){
00117 
00118     if(mVerbose){cout << "[L1TTestsSummary:] Initializing L1TSync Module Monitoring" << endl;}
00119     
00120     mDBE->setCurrentFolder(mL1TSyncPath);
00121     vector<string> histToMonitor = mDBE->getMEs();
00122     int            histLines     = histToMonitor.size()+1;
00123 
00124     mDBE->setCurrentFolder("L1T/L1TTestsSummary/");  
00125     mL1TSyncMonitor = mDBE->book2D("SyncQualitySummary","L1T Synchronization Monitor Summary",maxLS,0.5,double(maxLS)+0.5,histLines,0,histLines); 
00126     mL1TSyncMonitor->setAxisTitle("Lumi Section" ,1);
00127 
00128     mL1TSyncMonitor->setBinLabel(1,"Summary",2);
00129     for(unsigned int i=0 ; i<histToMonitor.size() ; i++){
00130       string name = mDBE->get(mL1TSyncPath+histToMonitor[i])->getTH1()->GetName();
00131       mL1TSyncMonitor->setBinLabel(i+2,name,2);  
00132     }
00133     
00134   }
00135   if(mMonitorL1TOccupancy){    
00136         
00137     if(mVerbose){cout << "[L1TTestsSummary:] Initializing L1TOccupancy Module Monitoring" << endl;}
00138     
00139     mDBE->setCurrentFolder(mL1TOccupancyPath);
00140     vector<string> histToMonitor = mDBE->getMEs();
00141     int            histLines     = histToMonitor.size()+1;
00142 
00143     mDBE->setCurrentFolder("L1T/L1TTestsSummary/");  
00144     mL1TOccupancyMonitor = mDBE->book2D("OccupancySummary","L1T Occupancy Monitor Summary",maxLS,+0.5,double(maxLS)+0.5,histLines,0,histLines);
00145     mL1TOccupancyMonitor->setAxisTitle("Lumi Section" ,1);   
00146 
00147     mL1TOccupancyMonitor->setBinLabel(1,"Summary",2);
00148     for(unsigned int i=0 ; i<histToMonitor.size() ; i++){
00149       string name = mDBE->get(mL1TOccupancyPath+histToMonitor[i])->getTH1()->GetName();
00150       mL1TOccupancyMonitor->setBinLabel(i+2,name,2);  
00151     } 
00152   }
00153   
00154   //-> Making the summary of summaries
00155   int testsToMonitor=1;
00156   if(mMonitorL1TRate)     {testsToMonitor++;}
00157   if(mMonitorL1TSync)     {testsToMonitor++;}
00158   if(mMonitorL1TOccupancy){testsToMonitor++;}
00159 
00160   // Creating
00161   mDBE->setCurrentFolder("L1T/L1TTestsSummary/");  
00162   mL1TSummary = mDBE->book2D("L1TQualitySummary","L1 Tests Summary",maxLS,+0.5,double(maxLS)+0.5,testsToMonitor,0,testsToMonitor);
00163   mL1TSummary->setAxisTitle("Lumi Section" ,1);
00164   mL1TSummary->setBinLabel(1,"L1T Summary",2);
00165   
00166   int it=2;
00167   if(mMonitorL1TRate)     {mL1TSummary->setBinLabel(it,"Rates"          ,2); binYRate    =it; it++;}
00168   if(mMonitorL1TSync)     {mL1TSummary->setBinLabel(it,"Synchronization",2); binYSync    =it; it++;}
00169   if(mMonitorL1TOccupancy){mL1TSummary->setBinLabel(it,"Occupancy"      ,2); binYOccpancy=it;}
00170 
00171 }
00172 
00173 //____________________________________________________________________________
00174 // Function: endRun
00175 // Description: This is will be run at the end of each run
00176 // Inputs: 
00177 // * const Run&        r       = Run information 
00178 // * const EventSetup& context = Event Setup information
00179 //____________________________________________________________________________
00180 void L1TTestsSummary::endRun(const Run& r, const EventSetup& context){
00181   
00182   if(mVerbose){cout << "[L1TTestsSummary:] Called endRun()" << endl;}
00183   
00184   if(mMonitorL1TRate)     {updateL1TRateMonitor();}
00185   if(mMonitorL1TSync)     {updateL1TSyncMonitor();}
00186   if(mMonitorL1TOccupancy){updateL1TOccupancyMonitor();}
00187   updateL1TSummary();
00188   
00189 }
00190 
00191 //____________________________________________________________________________
00192 // Function: beginLuminosityBlock
00193 // Description: This is will be run at the begining of each luminosity block
00194 // Inputs: 
00195 // * const LuminosityBlock& lumiSeg = Luminosity Block information 
00196 // * const EventSetup&      context = Event Setup information
00197 //____________________________________________________________________________
00198 void L1TTestsSummary::beginLuminosityBlock(const LuminosityBlock& lumiSeg, const EventSetup& context) {
00199   if(mVerbose){cout << "[L1TTestsSummary:] Called beginLuminosityBlock()" << endl;}
00200 }
00201 
00202 //____________________________________________________________________________
00203 // Function: endLuminosityBlock
00204 // Description: This is will be run at the end of each luminosity block
00205 // Inputs: 
00206 // * const LuminosityBlock& lumiSeg = Luminosity Block information 
00207 // * const EventSetup&      context = Event Setup information
00208 //____________________________________________________________________________
00209 void L1TTestsSummary::endLuminosityBlock(const edm::LuminosityBlock& lumiSeg, const edm::EventSetup& c){
00210   
00211   int eventLS = lumiSeg.id().luminosityBlock();
00212 
00213   mProcessedLS.push_back(eventLS);
00214   
00215   if(mVerbose) {
00216     cout << "[L1TTestsSummary:] Called endLuminosityBlock()" << endl;
00217     cout << "[L1TTestsSummary:] Lumisection: " << eventLS << endl;
00218   }
00219 
00220   if(mMonitorL1TRate)     {updateL1TRateMonitor();}
00221   if(mMonitorL1TSync)     {updateL1TSyncMonitor();}
00222   if(mMonitorL1TOccupancy){updateL1TOccupancyMonitor();}
00223   updateL1TSummary();
00224   
00225 }
00226 
00227 //____________________________________________________________________________
00228 // Function: analyze
00229 // Description: This is will be run for every event
00230 // Inputs: 
00231 // * const Event&      e       = Event information 
00232 // * const EventSetup& context = Event Setup information
00233 //____________________________________________________________________________
00234 void L1TTestsSummary::analyze(const Event& e, const EventSetup& context){}
00235 
00236 //____________________________________________________________________________
00237 // Function: updateL1TRateMonitor
00238 // Description: 
00239 //____________________________________________________________________________
00240 void L1TTestsSummary::updateL1TRateMonitor(){
00241   
00242   mDBE->setCurrentFolder(mL1TRatePath);
00243   vector<string> histToMonitor = mDBE->getMEs();
00244   
00245   for(unsigned int i=0 ; i<histToMonitor.size() ; i++){
00246 
00247     MonitorElement* me = mDBE->get(mL1TRatePath+histToMonitor[i]);
00248     if(mVerbose) {cout << "[L1TTestsSummary:] Found ME: " << me->getTH1()->GetName() << endl;}
00249 
00250     const QReport * myQReport = me->getQReport("L1TRateTest"); //get QReport associated to your ME  
00251     if(myQReport) {
00252       float  qtresult  = myQReport->getQTresult(); // get QT result value
00253       int    qtstatus  = myQReport->getStatus();   // get QT status value (see table below)
00254       string qtmessage = myQReport->getMessage() ; // get the whole QT result message
00255       vector<DQMChannel> qtBadChannels = myQReport->getBadChannels();
00256 
00257       if(mVerbose) {
00258         cout << "[L1TTestsSummary:] Found QReport for ME: " << me->getTH1()->GetName() << endl;
00259         cout << "[L1TTestsSummary:] Result=" << qtresult << " status=" << qtstatus << " message=" << qtmessage << endl;
00260         cout << "[L1TTestsSummary:] Bad Channels size=" << qtBadChannels.size() << endl;
00261       }
00262 
00263       for(unsigned int i=0 ; i<mProcessedLS.size()-1 ; i++){
00264         int binx = mL1TRateMonitor->getTH2F()->GetXaxis()->FindBin(mProcessedLS[i]);
00265         int biny = mL1TRateMonitor->getTH2F()->GetYaxis()->FindBin(me->getTH1()->GetName());
00266         mL1TRateMonitor->setBinContent(binx,biny,100);
00267       }
00268 
00269       for(unsigned int a=0 ; a<qtBadChannels.size() ; a++){
00270         for(unsigned int b=0 ; b<mProcessedLS.size()-1 ; b++){
00271    
00272           // Converting bin to value
00273           double valueBinBad = me->getTH1()->GetBinCenter(qtBadChannels[a].getBin());
00274      
00275           if(valueBinBad==(mProcessedLS[b])){
00276             int binx = mL1TRateMonitor->getTH2F()->GetXaxis()->FindBin(valueBinBad);
00277             int biny = mL1TRateMonitor->getTH2F()->GetYaxis()->FindBin(me->getTH1()->GetName());
00278             mL1TRateMonitor->setBinContent(binx,biny,300);
00279           }
00280         }
00281       }
00282     }
00283   } 
00284 
00285   //-> Filling the summaries
00286   int nBinX = mL1TRateMonitor->getTH2F()->GetXaxis()->GetNbins();
00287   int nBinY = mL1TRateMonitor->getTH2F()->GetYaxis()->GetNbins();
00288   for(int binx=1; binx<=nBinX ; binx++){
00289     int GlobalStatus=0;
00290     for(int biny=2; biny<=nBinY ; biny++){
00291       double flag = mL1TRateMonitor->getBinContent(binx,biny);
00292       if(GlobalStatus<flag){GlobalStatus=flag;}
00293     }
00294     
00295     // NOTE: Assumes mL1TSummary has same size then mL1TRateMonitor
00296     mL1TRateMonitor->setBinContent(binx,       1,GlobalStatus);
00297     mL1TSummary    ->setBinContent(binx,binYRate,GlobalStatus);
00298   } 
00299 }
00300 
00301 //____________________________________________________________________________
00302 // Function: updateL1TSyncMonitor
00303 // Description: 
00304 // Note: Values certified by L1TRates are always about currentLS-1 since we
00305 //       use rate averages from the previous LS from GT
00306 //____________________________________________________________________________
00307 void L1TTestsSummary::updateL1TSyncMonitor(){
00308   
00309   mDBE->setCurrentFolder(mL1TSyncPath);
00310   vector<string> histToMonitor = mDBE->getMEs();
00311   
00312   for(unsigned int i=0 ; i<histToMonitor.size() ; i++){
00313 
00314     MonitorElement* me = mDBE->get(mL1TSyncPath+histToMonitor[i]);
00315     if(mVerbose) {cout << "[L1TTestsSummary:] Found ME: " << me->getTH1()->GetName() << endl;}
00316 
00317     const QReport * myQReport = me->getQReport("L1TSyncTest"); //get QReport associated to your ME  
00318     if(myQReport) {
00319       float              qtresult      = myQReport->getQTresult();    // get QT result value
00320       int                qtstatus      = myQReport->getStatus();      // get QT status value (see table below)
00321       string             qtmessage     = myQReport->getMessage() ;    // get the whole QT result message
00322       vector<DQMChannel> qtBadChannels = myQReport->getBadChannels();
00323 
00324       if(mVerbose) {
00325         cout << "[L1TTestsSummary:] Found QReport for ME: " << me->getTH1()->GetName() << endl;
00326         cout << "[L1TTestsSummary:] Result=" << qtresult << " status=" << qtstatus << " message=" << qtmessage << endl;
00327         cout << "[L1TTestsSummary:] Bad Channels size=" << qtBadChannels.size() << endl;
00328       }
00329 
00330       for(unsigned int i=0 ; i<mProcessedLS.size() ; i++){
00331         int binx = mL1TSyncMonitor->getTH2F()->GetXaxis()->FindBin(mProcessedLS[i]);
00332         int biny = mL1TSyncMonitor->getTH2F()->GetYaxis()->FindBin(me->getTH1()->GetName());
00333         mL1TSyncMonitor->setBinContent(binx,biny,100);
00334       }
00335 
00336       for(unsigned int a=0 ; a<qtBadChannels.size() ; a++){ 
00337         for(unsigned int b=0 ; b<mProcessedLS.size() ; b++){
00338           
00339           // Converting bin to value
00340           double valueBinBad = me->getTH1()->GetBinCenter(qtBadChannels[a].getBin());
00341           
00342           if(valueBinBad==mProcessedLS[b]){
00343             int binx = mL1TSyncMonitor->getTH2F()->GetXaxis()->FindBin(valueBinBad);
00344             int biny = mL1TSyncMonitor->getTH2F()->GetYaxis()->FindBin(me->getTH1()->GetName());
00345             mL1TSyncMonitor->setBinContent(binx,biny,300);
00346           }
00347         }
00348       }
00349     }
00350   } 
00351 
00352   //-> Filling the summaries
00353   int nBinX = mL1TSyncMonitor->getTH2F()->GetXaxis()->GetNbins();
00354   int nBinY = mL1TSyncMonitor->getTH2F()->GetYaxis()->GetNbins();
00355   for(int binx=1; binx<=nBinX ; binx++){
00356     int GlobalStatus=0;
00357     for(int biny=2; biny<=nBinY ; biny++){
00358       double flag = mL1TSyncMonitor->getBinContent(binx,biny);
00359       if(GlobalStatus<flag){GlobalStatus=flag;}
00360     }
00361     
00362     // NOTE: Assumes mL1TSummary has same size then mL1TSyncMonitor
00363     mL1TSyncMonitor->setBinContent(binx,       1,GlobalStatus);
00364     mL1TSummary    ->setBinContent(binx,binYSync,GlobalStatus);
00365   } 
00366 }
00367 
00368 //____________________________________________________________________________
00369 // Function: updateL1TOccupancyMonitor
00370 // Description: 
00371 //____________________________________________________________________________
00372 void L1TTestsSummary::updateL1TOccupancyMonitor(){
00373   
00374   mDBE->setCurrentFolder(mL1TOccupancyPath);
00375   vector<string> histToMonitor = mDBE->getMEs();
00376   
00377   for(unsigned int i=0 ; i<histToMonitor.size() ; i++){
00378 
00379     MonitorElement* me = mDBE->get(mL1TOccupancyPath+histToMonitor[i]);
00380     if(mVerbose) {cout << "[L1TTestsSummary:] Found ME: " << me->getTH1()->GetName() << endl;}
00381 
00382     const QReport * myQReport = me->getQReport("L1TOccupancyTest"); //get QReport associated to your ME  
00383     if(myQReport) {
00384       float  qtresult  = myQReport->getQTresult();        // get QT result value
00385       int    qtstatus  = myQReport->getStatus();          // get QT status value (see table below)
00386       string qtmessage = myQReport->getMessage() ; // get the whole QT result message
00387       vector<DQMChannel> qtBadChannels = myQReport->getBadChannels();
00388 
00389       if(mVerbose) {
00390         cout << "[L1TTestsSummary:] Found QReport for ME: " << me->getTH1()->GetName() << endl;
00391         cout << "[L1TTestsSummary:] Result=" << qtresult << " status=" << qtstatus << " message=" << qtmessage << endl;
00392         cout << "[L1TTestsSummary:] Bad Channels size=" << qtBadChannels.size() << endl;
00393       }
00394 
00395       for(unsigned int i=0 ; i<mProcessedLS.size() ; i++){
00396         int binx = mL1TOccupancyMonitor->getTH2F()->GetXaxis()->FindBin(mProcessedLS[i]);
00397         int biny = mL1TOccupancyMonitor->getTH2F()->GetYaxis()->FindBin(me->getTH1()->GetName());
00398         mL1TOccupancyMonitor->setBinContent(binx,biny,100);
00399       }
00400 
00401       for(unsigned int a=0 ; a<qtBadChannels.size() ; a++){
00402         for(unsigned int b=0 ; b<mProcessedLS.size() ; b++){
00403           
00404           // Converting bin to value
00405           double valueBinBad = me->getTH1()->GetBinCenter(qtBadChannels[a].getBin());
00406   
00407           if(valueBinBad==mProcessedLS[b]){
00408             int binx = mL1TOccupancyMonitor->getTH2F()->GetXaxis()->FindBin(valueBinBad);
00409             int biny = mL1TOccupancyMonitor->getTH2F()->GetYaxis()->FindBin(me->getTH1()->GetName());
00410             mL1TOccupancyMonitor->setBinContent(binx,biny,300);
00411           }
00412         }
00413       }
00414     }
00415   } 
00416     
00417   //-> Filling the summaries
00418   int nBinX = mL1TOccupancyMonitor->getTH2F()->GetXaxis()->GetNbins();
00419   int nBinY = mL1TOccupancyMonitor->getTH2F()->GetYaxis()->GetNbins();
00420   for(int binx=1; binx<=nBinX ; binx++){
00421     int GlobalStatus=0;
00422     for(int biny=2; biny<=nBinY ; biny++){
00423       double flag = mL1TOccupancyMonitor->getBinContent(binx,biny);
00424       if(GlobalStatus<flag){GlobalStatus=flag;}
00425     }
00426     
00427     // NOTE: Assumes mL1TSummary has same size then mL1TOccupancyMonitor
00428     mL1TOccupancyMonitor->setBinContent(binx,           1,GlobalStatus);
00429     mL1TSummary         ->setBinContent(binx,binYOccpancy,GlobalStatus);
00430   } 
00431 }
00432 
00433 //____________________________________________________________________________
00434 // Function: updateL1TOccupancyMonitor
00435 // Description: 
00436 //____________________________________________________________________________
00437 void L1TTestsSummary::updateL1TSummary(){
00438 
00439   int nBinX = mL1TSummary->getTH2F()->GetXaxis()->GetNbins();  
00440   for(int binx=1; binx<=nBinX ; binx++){
00441     int GlobalStatus=0;
00442     if(mMonitorL1TRate){
00443       if(mL1TSummary->getBinContent(binx,binYRate)>GlobalStatus){
00444         GlobalStatus=mL1TSummary->getBinContent(binx,binYRate);
00445       }
00446     }
00447     if(mMonitorL1TSync)     {
00448       if(mL1TSummary->getBinContent(binx,binYSync)>GlobalStatus){
00449         GlobalStatus=mL1TSummary->getBinContent(binx,binYSync);
00450       }
00451     }
00452     if(mMonitorL1TOccupancy){
00453       if(mL1TSummary->getBinContent(binx,binYOccpancy)>GlobalStatus){
00454         GlobalStatus=mL1TSummary->getBinContent(binx,binYOccpancy);
00455       }
00456     }
00457     mL1TSummary->setBinContent(binx,1,GlobalStatus);
00458   }  
00459 }
00460 
00461 //define this as a plug-in
00462 DEFINE_FWK_MODULE(L1TTestsSummary);