CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10/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: 2012/06/21 13:40:22 $
00011  * $Revision: 1.7 $
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_)
00087     {
00088       this->setup();
00089       this->reset();
00090     }
00091   else if (tevt_ == 0)
00092     {
00093       this->setup();
00094       this->reset();
00095     }
00096 } // beginRun(const edm::Run& run, const edm::EventSetup& c)
00097 
00098 void HcalBaseDQMonitor::endRun(const edm::Run& run, const edm::EventSetup& c)
00099 {
00100   if (debug_>0) std::cout <<"HcalBaseDQMonitor::endRun:  task = "<<subdir_<<std::endl;
00101 } //endRun(...)
00102 
00103 void HcalBaseDQMonitor::reset(void)
00104 {
00105   if (debug_>0) std::cout <<"HcalBaseDQMonitor::reset():  task = "<<subdir_<<std::endl;
00106   if (meIevt_) meIevt_->Fill(-1);
00107   ievt_=0;
00108   if (meLevt_) meLevt_->Fill(-1);
00109   levt_=0;
00110   if (meTevt_) meTevt_->Fill(-1);
00111   tevt_=0;
00112   if (meTevtHist_) meTevtHist_->Reset();
00113   if (ProblemsCurrentLB) ProblemsCurrentLB->Reset();
00114   HBpresent_=false;
00115   HEpresent_=false;
00116   HOpresent_=false;
00117   HFpresent_=false;
00118   currentLS=0;
00119   currenttype_=-1;
00120 } //reset()
00121 
00122 void HcalBaseDQMonitor::cleanup(void)
00123 {
00124 
00125 } //cleanup()
00126 
00127 void HcalBaseDQMonitor::setup(void)
00128 {
00129   if (debug_>3) std::cout <<"<HcalBaseDQMonitor> setup in progress"<<std::endl;
00130   dbe_->setCurrentFolder(subdir_);
00131   meIevt_ = dbe_->bookInt("EventsProcessed");
00132   if (meIevt_) meIevt_->Fill(-1);
00133   meLevt_ = dbe_->bookInt("EventsProcessed_currentLS");
00134   if (meLevt_) meLevt_->Fill(-1);
00135   meTevt_ = dbe_->bookInt("EventsProcessed_Total");
00136   if (meTevt_) meTevt_->Fill(-1);
00137   meTevtHist_=dbe_->book1D("Events_Processed_Task_Histogram","Counter of Events Processed By This Task",1,0.5,1.5);
00138   if (meTevtHist_) meTevtHist_->Reset();
00139   dbe_->setCurrentFolder(subdir_+"LSvalues");
00140   ProblemsCurrentLB=dbe_->book2D("ProblemsThisLS","Problem Channels in current Lumi Section",
00141                                  7,0,7,1,0,1);
00142   if (ProblemsCurrentLB)
00143     {
00144       (ProblemsCurrentLB->getTH2F())->GetXaxis()->SetBinLabel(1,"HB");
00145       (ProblemsCurrentLB->getTH2F())->GetXaxis()->SetBinLabel(2,"HE");
00146       (ProblemsCurrentLB->getTH2F())->GetXaxis()->SetBinLabel(3,"HO");
00147       (ProblemsCurrentLB->getTH2F())->GetXaxis()->SetBinLabel(4,"HF");
00148       (ProblemsCurrentLB->getTH2F())->GetXaxis()->SetBinLabel(5,"HO0");
00149       (ProblemsCurrentLB->getTH2F())->GetXaxis()->SetBinLabel(6,"HO12");
00150       (ProblemsCurrentLB->getTH2F())->GetXaxis()->SetBinLabel(7,"HFlumi");
00151       (ProblemsCurrentLB->getTH2F())->GetYaxis()->SetBinLabel(1,"Status");
00152       ProblemsCurrentLB->Reset();
00153     }
00154 } // setup()
00155 
00156 
00157 void HcalBaseDQMonitor::beginLuminosityBlock(const edm::LuminosityBlock& lumiSeg,
00158                           const edm::EventSetup& c)
00159 {
00160   if (this->LumiInOrder(lumiSeg.luminosityBlock())==false) return;
00161   currentLS=lumiSeg.luminosityBlock();
00162   levt_=0;
00163   if (meLevt_!=0) meLevt_->Fill(-1);
00164   if (ProblemsCurrentLB)
00165     ProblemsCurrentLB->Reset();
00166 }
00167 
00168 void HcalBaseDQMonitor::endLuminosityBlock(const edm::LuminosityBlock& lumiSeg,
00169                         const edm::EventSetup& c)
00170 {
00171   if (this->LumiInOrder(lumiSeg.luminosityBlock())==false) return;
00172   // Inherited classes do end-of-lumi functions here
00173 }
00174 
00175 
00176 bool HcalBaseDQMonitor::LumiInOrder(int lumisec)
00177 {
00178   if (skipOutOfOrderLS_==false) return true; // don't skip out-of-order lumi sections
00179   // check that latest lumi section is >= last processed
00180   if (Online_ && lumisec<currentLS)
00181     return false;
00182   return true;
00183 }
00184 
00185 bool HcalBaseDQMonitor::IsAllowedCalibType()
00186 {
00187   if (debug_>9) std::cout <<"<HcalBaseDQMonitor::IsAllowedCalibType>"<<std::endl;
00188   if (AllowedCalibTypes_.size()==0)
00189     {
00190       if (debug_>9) std::cout <<"\tNo calib types specified by user; All events allowed"<<std::endl;
00191       return true;
00192     }
00193   MonitorElement* me = dbe_->get((prefixME_+"HcalInfo/CURRENT_EVENT_TYPE").c_str());
00194   if (me) currenttype_=me->getIntValue();
00195   else 
00196     {
00197       if (debug_>9) std::cout <<"\tCalib Type cannot be determined from HcalMonitorModule"<<std::endl;
00198       return true; // is current type can't be determined, assume event is allowed
00199     }
00200   if (debug_>9) std::cout <<"\tHcalBaseDQMonitor::IsAllowedCalibType  checking if calibration type = "<<currenttype_<<" is allowed...";
00201   for (std::vector<int>::size_type i=0;i<AllowedCalibTypes_.size();++i)
00202     {
00203       if (AllowedCalibTypes_[i]==currenttype_)
00204         {
00205           if (debug_>9) std::cout <<"\t Type allowed!"<<std::endl;
00206           return true;
00207         }
00208     }
00209   if (debug_>9) std::cout <<"\t Type not allowed!"<<std::endl;
00210   return false;
00211 } // bool HcalBaseDQMonitor::IsAllowedCalibType()
00212 
00213 void HcalBaseDQMonitor::analyze(const edm::Event& e, const edm::EventSetup& c)
00214 {
00215   if (debug_>5) std::cout <<"\t<HcalBaseDQMonitor::analyze>  event = "<<ievt_<<std::endl;
00216   eventAllowed_=true; // assume event is allowed
00217 
00218   // fill with total events seen (this differs from ievent, which is total # of good events)
00219   ++tevt_;
00220   if (meTevt_) meTevt_->Fill(tevt_);
00221   if (meTevtHist_) meTevtHist_->Fill(1);
00222   // skip out of order lumi events
00223   if (this->LumiInOrder(e.luminosityBlock())==false)
00224     {
00225       eventAllowed_=false;
00226       return;
00227     }
00228   // skip events of wrong calibration type
00229   eventAllowed_&=(this->IsAllowedCalibType());
00230   if (!eventAllowed_) return;
00231 
00232   // Event is good; count it
00233   ++ievt_;
00234   ++levt_;
00235   if (meIevt_) meIevt_->Fill(ievt_);
00236   if (meLevt_) meLevt_->Fill(levt_);
00237 
00238 
00239   MonitorElement* me;
00240   if (HBpresent_==false)
00241     {
00242       me = dbe_->get((prefixME_+"HcalInfo/HBpresent"));
00243       if (me==0 || me->getIntValue()>0) HBpresent_=true;
00244     }
00245   if (HEpresent_==false)
00246     {
00247       me = dbe_->get((prefixME_+"HcalInfo/HEpresent"));
00248       if (me==0 || me->getIntValue()>0) HEpresent_=true;
00249     }
00250   if (HOpresent_==false)
00251     {
00252       me = dbe_->get((prefixME_+"HcalInfo/HOpresent"));
00253       if (me==0 || me->getIntValue()>0) HOpresent_=true;
00254     }
00255   if (HFpresent_==false)
00256     {
00257       me = dbe_->get((prefixME_+"HcalInfo/HOpresent"));
00258       if (me ==0 || me->getIntValue()>0) HFpresent_=true;
00259     }
00260 
00261 
00262 } // void HcalBaseDQMonitor::analyze(const edm::Event& e, const edm::EventSetup& c)
00263 
00264 DEFINE_FWK_MODULE(HcalBaseDQMonitor);