CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/DQM/HcalMonitorTasks/src/HcalBaseDQMonitor.cc

Go to the documentation of this file.
00001 #include <DQM/HcalMonitorTasks/interface/HcalBaseDQMonitor.h>
00002 #include "FWCore/Framework/interface/LuminosityBlock.h"
00003 
00004 #include <iostream>
00005 #include <vector>
00006 
00007 /*
00008  * \file HcalBaseDQMonitor.cc
00009  *
00010  * $Date: 2010/11/24 18:55:12 $
00011  * $Revision: 1.6 $
00012  * \author J Temple
00013  *
00014  * Base class for all Hcal DQM analyzers
00015  *
00016 
00017 */
00018 
00019 // constructor
00020 
00021 HcalBaseDQMonitor::HcalBaseDQMonitor(const edm::ParameterSet& ps)
00022 {
00023   Online_                = ps.getUntrackedParameter<bool>("online",false);
00024   mergeRuns_             = ps.getUntrackedParameter<bool>("mergeRuns",false);
00025   enableCleanup_         = ps.getUntrackedParameter<bool>("enableCleanup",false);
00026   debug_                 = ps.getUntrackedParameter<int>("debug",0);
00027   prefixME_              = ps.getUntrackedParameter<std::string>("subSystemFolder","Hcal/"); 
00028   if (prefixME_.substr(prefixME_.size()-1,prefixME_.size())!="/")
00029     prefixME_.append("/");
00030   subdir_                = ps.getUntrackedParameter<std::string>("TaskFolder","Test/"); 
00031   if (subdir_.size()>0 && subdir_.substr(subdir_.size()-1,subdir_.size())!="/")
00032     subdir_.append("/");
00033   subdir_=prefixME_+subdir_;
00034   AllowedCalibTypes_     = ps.getUntrackedParameter<std::vector<int> > ("AllowedCalibTypes");
00035   skipOutOfOrderLS_      = ps.getUntrackedParameter<bool>("skipOutOfOrderLS",false);
00036   NLumiBlocks_           = ps.getUntrackedParameter<int>("NLumiBlocks",4000);
00037   makeDiagnostics_       = ps.getUntrackedParameter<bool>("makeDiagnostics",false);
00038   
00039   meIevt_=0;
00040   meLevt_=0;
00041   meTevtHist_=0;
00042   ProblemsVsLB=0;
00043   ProblemsVsLB_HB=0;
00044   ProblemsVsLB_HE=0;
00045   ProblemsVsLB_HF=0;
00046   ProblemsVsLB_HBHEHF=0;
00047   ProblemsVsLB_HO=0;
00048   ProblemsCurrentLB=0;
00049 } //HcalBaseDQMonitor::HcalBaseDQMonitor(const ParameterSet& ps)
00050 
00051 
00052 // destructor
00053 
00054 HcalBaseDQMonitor::~HcalBaseDQMonitor()
00055 {
00056 
00057 }
00058 
00059 void HcalBaseDQMonitor::beginJob(void)
00060 {
00061 
00062   if (debug_>0) std::cout <<"HcalBaseDQMonitor::beginJob():  task =  '"<<subdir_<<"'"<<std::endl;
00063   dbe_ = edm::Service<DQMStore>().operator->();
00064 
00065   ievt_=0;
00066   levt_=0;
00067   tevt_=0;
00068   currenttype_=-1;
00069   HBpresent_=false;
00070   HEpresent_=false;
00071   HOpresent_=false;
00072   HFpresent_=false;
00073 
00074 
00075 } // beginJob()
00076 
00077 void HcalBaseDQMonitor::endJob(void)
00078 {
00079   if (enableCleanup_)
00080     cleanup();
00081 } // endJob()
00082 
00083 void HcalBaseDQMonitor::beginRun(const edm::Run& run, const edm::EventSetup& c)
00084 {
00085   if (debug_>0) std::cout <<"HcalBaseDQMonitor::beginRun():  task =  '"<<subdir_<<"'"<<std::endl;
00086   if (mergeRuns_ && tevt_>0) return;
00087   this->setup();
00088   this->reset();
00089 } // beginRun(const edm::Run& run, const edm::EventSetup& c)
00090 
00091 void HcalBaseDQMonitor::endRun(const edm::Run& run, const edm::EventSetup& c)
00092 {
00093   if (debug_>0) std::cout <<"HcalBaseDQMonitor::endRun:  task = "<<subdir_<<std::endl;
00094 } //endRun(...)
00095 
00096 void HcalBaseDQMonitor::reset(void)
00097 {
00098   if (debug_>0) std::cout <<"HcalBaseDQMonitor::reset():  task = "<<subdir_<<std::endl;
00099   if (meIevt_) meIevt_->Fill(-1);
00100   ievt_=0;
00101   if (meLevt_) meLevt_->Fill(-1);
00102   levt_=0;
00103   if (meTevt_) meTevt_->Fill(-1);
00104   tevt_=0;
00105   if (meTevtHist_) meTevtHist_->Reset();
00106   if (ProblemsCurrentLB) ProblemsCurrentLB->Reset();
00107   HBpresent_=false;
00108   HEpresent_=false;
00109   HOpresent_=false;
00110   HFpresent_=false;
00111   currentLS=0;
00112   currenttype_=-1;
00113 } //reset()
00114 
00115 void HcalBaseDQMonitor::cleanup(void)
00116 {
00117 
00118 } //cleanup()
00119 
00120 void HcalBaseDQMonitor::setup(void)
00121 {
00122   if (debug_>3) std::cout <<"<HcalBaseDQMonitor> setup in progress"<<std::endl;
00123   dbe_->setCurrentFolder(subdir_);
00124   meIevt_ = dbe_->bookInt("EventsProcessed");
00125   if (meIevt_) meIevt_->Fill(-1);
00126   meLevt_ = dbe_->bookInt("EventsProcessed_currentLS");
00127   if (meLevt_) meLevt_->Fill(-1);
00128   meTevt_ = dbe_->bookInt("EventsProcessed_Total");
00129   if (meTevt_) meTevt_->Fill(-1);
00130   meTevtHist_=dbe_->book1D("Events_Processed_Task_Histogram","Counter of Events Processed By This Task",1,0.5,1.5);
00131   if (meTevtHist_) meTevtHist_->Reset();
00132   dbe_->setCurrentFolder(subdir_+"LSvalues");
00133   ProblemsCurrentLB=dbe_->book2D("ProblemsThisLS","Problem Channels in current Lumi Section",
00134                                  7,0,7,1,0,1);
00135   if (ProblemsCurrentLB)
00136     {
00137       (ProblemsCurrentLB->getTH2F())->GetXaxis()->SetBinLabel(1,"HB");
00138       (ProblemsCurrentLB->getTH2F())->GetXaxis()->SetBinLabel(2,"HE");
00139       (ProblemsCurrentLB->getTH2F())->GetXaxis()->SetBinLabel(3,"HO");
00140       (ProblemsCurrentLB->getTH2F())->GetXaxis()->SetBinLabel(4,"HF");
00141       (ProblemsCurrentLB->getTH2F())->GetXaxis()->SetBinLabel(5,"HO0");
00142       (ProblemsCurrentLB->getTH2F())->GetXaxis()->SetBinLabel(6,"HO12");
00143       (ProblemsCurrentLB->getTH2F())->GetXaxis()->SetBinLabel(7,"HFlumi");
00144       (ProblemsCurrentLB->getTH2F())->GetYaxis()->SetBinLabel(1,"Status");
00145       ProblemsCurrentLB->Reset();
00146     }
00147 } // setup()
00148 
00149 
00150 void HcalBaseDQMonitor::beginLuminosityBlock(const edm::LuminosityBlock& lumiSeg,
00151                           const edm::EventSetup& c)
00152 {
00153   if (this->LumiInOrder(lumiSeg.luminosityBlock())==false) return;
00154   currentLS=lumiSeg.luminosityBlock();
00155   levt_=0;
00156   if (meLevt_!=0) meLevt_->Fill(-1);
00157   if (ProblemsCurrentLB)
00158     ProblemsCurrentLB->Reset();
00159 }
00160 
00161 void HcalBaseDQMonitor::endLuminosityBlock(const edm::LuminosityBlock& lumiSeg,
00162                         const edm::EventSetup& c)
00163 {
00164   if (this->LumiInOrder(lumiSeg.luminosityBlock())==false) return;
00165   // Inherited classes do end-of-lumi functions here
00166 }
00167 
00168 
00169 bool HcalBaseDQMonitor::LumiInOrder(int lumisec)
00170 {
00171   if (skipOutOfOrderLS_==false) return true; // don't skip out-of-order lumi sections
00172   // check that latest lumi section is >= last processed
00173   if (Online_ && lumisec<currentLS)
00174     return false;
00175   return true;
00176 }
00177 
00178 bool HcalBaseDQMonitor::IsAllowedCalibType()
00179 {
00180   if (debug_>9) std::cout <<"<HcalBaseDQMonitor::IsAllowedCalibType>"<<std::endl;
00181   if (AllowedCalibTypes_.size()==0)
00182     {
00183       if (debug_>9) std::cout <<"\tNo calib types specified by user; All events allowed"<<std::endl;
00184       return true;
00185     }
00186   MonitorElement* me = dbe_->get((prefixME_+"HcalInfo/CURRENT_EVENT_TYPE").c_str());
00187   if (me) currenttype_=me->getIntValue();
00188   else 
00189     {
00190       if (debug_>9) std::cout <<"\tCalib Type cannot be determined from HcalMonitorModule"<<std::endl;
00191       return true; // is current type can't be determined, assume event is allowed
00192     }
00193   if (debug_>9) std::cout <<"\tHcalBaseDQMonitor::IsAllowedCalibType  checking if calibration type = "<<currenttype_<<" is allowed...";
00194   for (std::vector<int>::size_type i=0;i<AllowedCalibTypes_.size();++i)
00195     {
00196       if (AllowedCalibTypes_[i]==currenttype_)
00197         {
00198           if (debug_>9) std::cout <<"\t Type allowed!"<<std::endl;
00199           return true;
00200         }
00201     }
00202   if (debug_>9) std::cout <<"\t Type not allowed!"<<std::endl;
00203   return false;
00204 } // bool HcalBaseDQMonitor::IsAllowedCalibType()
00205 
00206 void HcalBaseDQMonitor::analyze(const edm::Event& e, const edm::EventSetup& c)
00207 {
00208   if (debug_>5) std::cout <<"\t<HcalBaseDQMonitor::analyze>  event = "<<ievt_<<std::endl;
00209   eventAllowed_=true; // assume event is allowed
00210 
00211   // fill with total events seen (this differs from ievent, which is total # of good events)
00212   ++tevt_;
00213   if (meTevt_) meTevt_->Fill(tevt_);
00214   if (meTevtHist_) meTevtHist_->Fill(1);
00215   // skip out of order lumi events
00216   if (this->LumiInOrder(e.luminosityBlock())==false)
00217     {
00218       eventAllowed_=false;
00219       return;
00220     }
00221   // skip events of wrong calibration type
00222   eventAllowed_&=(this->IsAllowedCalibType());
00223   if (!eventAllowed_) return;
00224 
00225   // Event is good; count it
00226   ++ievt_;
00227   ++levt_;
00228   if (meIevt_) meIevt_->Fill(ievt_);
00229   if (meLevt_) meLevt_->Fill(levt_);
00230 
00231 
00232   MonitorElement* me;
00233   if (HBpresent_==false)
00234     {
00235       me = dbe_->get((prefixME_+"HcalInfo/HBpresent"));
00236       if (me==0 || me->getIntValue()>0) HBpresent_=true;
00237     }
00238   if (HEpresent_==false)
00239     {
00240       me = dbe_->get((prefixME_+"HcalInfo/HEpresent"));
00241       if (me==0 || me->getIntValue()>0) HEpresent_=true;
00242     }
00243   if (HOpresent_==false)
00244     {
00245       me = dbe_->get((prefixME_+"HcalInfo/HOpresent"));
00246       if (me==0 || me->getIntValue()>0) HOpresent_=true;
00247     }
00248   if (HFpresent_==false)
00249     {
00250       me = dbe_->get((prefixME_+"HcalInfo/HOpresent"));
00251       if (me ==0 || me->getIntValue()>0) HFpresent_=true;
00252     }
00253 
00254 
00255 } // void HcalBaseDQMonitor::analyze(const edm::Event& e, const edm::EventSetup& c)
00256 
00257 DEFINE_FWK_MODULE(HcalBaseDQMonitor);