CMS 3D CMS Logo

CMSSW_4_4_3_patch1/src/DQM/HcalMonitorClient/src/HcalDAQInfo.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    DQMO/HcalMonitorClient/HcalDAQInfo
00004 // Class:      HcalDAQInfo
00005 // 
00013 //
00014 // Original Author:  "Igor Vodopiyanov"
00015 //         Created:  Feb-21 2009
00016 //
00017 //
00018 
00019 // system include files
00020 #include <memory>
00021 #include <stdio.h>
00022 #include <math.h>
00023 #include <sstream>
00024 #include <fstream>
00025 #include <exception>
00026 
00027 #include "FWCore/Framework/interface/Frameworkfwd.h"
00028 #include "FWCore/Framework/interface/EDAnalyzer.h"
00029 #include "FWCore/Framework/interface/Event.h"
00030 #include "FWCore/Framework/interface/ESHandle.h"
00031 #include "FWCore/Framework/interface/EventSetup.h"
00032 #include "FWCore/Framework/interface/MakerMacros.h"
00033 
00034 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00035 
00036 #include "CondFormats/DataRecord/interface/RunSummaryRcd.h"
00037 #include "CondFormats/RunInfo/interface/RunSummary.h"
00038 #include "CondFormats/RunInfo/interface/RunInfo.h"
00039 
00040 #include "DataFormats/FEDRawData/interface/FEDNumbering.h"
00041 #include "DataFormats/FEDRawData/interface/FEDTrailer.h"
00042 #include "DQMServices/Core/interface/DQMStore.h"
00043 #include "DQMServices/Core/interface/MonitorElement.h"
00044 #include "FWCore/ServiceRegistry/interface/Service.h"
00045 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00046 
00047 //
00048 // class declaration
00049 //
00050 
00051 class HcalDAQInfo : public edm::EDAnalyzer {
00052    public:
00053       explicit HcalDAQInfo(const edm::ParameterSet&);
00054       ~HcalDAQInfo();
00055 
00056    private:
00057       virtual void beginJob();
00058       virtual void analyze(const edm::Event&, const edm::EventSetup&);
00059       virtual void endJob() ;
00060       virtual void beginLuminosityBlock(const edm::LuminosityBlock&, const edm::EventSetup&) ;
00061       virtual void endLuminosityBlock(const edm::LuminosityBlock&, const edm::EventSetup&) ;
00062 
00063    // ----------member data ---------------------------
00064 
00065    edm::ParameterSet conf_;
00066    DQMStore * dbe_;
00067    MonitorElement* HcalDaqFraction;
00068    MonitorElement* DAQSummaryMap;
00069    MonitorElement* HBDaqFraction;
00070    MonitorElement* HEDaqFraction;
00071    MonitorElement* HODaqFraction;
00072    MonitorElement* HFDaqFraction;
00073    MonitorElement* HO0DaqFraction;
00074    MonitorElement* HO12DaqFraction;
00075    MonitorElement* HFlumiDaqFraction;
00076    int debug_;
00077    std::string rootFolder_;
00078 };
00079 
00080 //
00081 // constants, enums and typedefs
00082 //
00083 
00084 //
00085 // static data member definitions
00086 //
00087 
00088 //
00089 // constructors and destructor
00090 //
00091 
00092 HcalDAQInfo::HcalDAQInfo(const edm::ParameterSet& iConfig)
00093 {
00094   // now do what ever initialization is needed
00095   debug_=iConfig.getUntrackedParameter<int>("debug",0);
00096   rootFolder_ = iConfig.getUntrackedParameter<std::string>("subSystemFolder","Hcal");
00097   dbe_ = edm::Service<DQMStore>().operator->();  
00098 }
00099 
00100 HcalDAQInfo::~HcalDAQInfo()
00101 { 
00102    // do anything here that needs to be done at destruction time
00103    // (e.g. close files, deallocate resources etc.)
00104 }
00105 
00106 //
00107 // member functions
00108 //
00109 
00110 // ------------ method called to for each event  ------------
00111 void
00112 HcalDAQInfo::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {}
00113 
00114 // ------------ method called once each job just before starting event loop  ------------
00115 void 
00116 HcalDAQInfo::beginJob()
00117 {
00118   if (debug_>0) std::cout<<"<HcalDAQInfo::beginJob>"<< std::endl;
00119 
00120   dbe_->setCurrentFolder(rootFolder_);
00121   std::string currDir = dbe_->pwd();
00122   if (debug_>0) std::cout << "--- Current Directory " << currDir << std::endl;
00123   std::vector<MonitorElement*> mes = dbe_->getAllContents("");
00124   if (debug_>0) std::cout << "found " << mes.size() << " monitoring elements:" << std::endl;
00125 
00126   dbe_->setCurrentFolder(rootFolder_+"/EventInfo/");
00127 
00128   HcalDaqFraction = dbe_->bookFloat("DAQSummary");
00129 
00130   DAQSummaryMap = dbe_->book2D("DAQSummaryMap","HcalDAQSummaryMap",7,0.,7.,1,0.,1.);
00131   DAQSummaryMap->setAxisRange(-1,1,3);
00132   DAQSummaryMap->setBinLabel(1,"HB");
00133   DAQSummaryMap->setBinLabel(2,"HE");
00134   DAQSummaryMap->setBinLabel(3,"HO");
00135   DAQSummaryMap->setBinLabel(4,"HF");
00136   DAQSummaryMap->setBinLabel(5,"H00");
00137   DAQSummaryMap->setBinLabel(6,"H012");
00138   DAQSummaryMap->setBinLabel(7,"HFlumi");
00139   DAQSummaryMap->setBinLabel(1,"Status",2);
00140 
00141   dbe_->setCurrentFolder(rootFolder_+"/EventInfo/DAQContents/");
00142   HBDaqFraction  = dbe_->bookFloat("Hcal_HB");
00143   HEDaqFraction  = dbe_->bookFloat("Hcal_HE");
00144   HODaqFraction  = dbe_->bookFloat("Hcal_HO");
00145   HFDaqFraction  = dbe_->bookFloat("Hcal_HF");
00146   HO0DaqFraction = dbe_->bookFloat("Hcal_HO0");
00147   HO12DaqFraction   = dbe_->bookFloat("Hcal_HO12");
00148   HFlumiDaqFraction = dbe_->bookFloat("Hcal_HFlumi");
00149 
00150 }
00151 
00152 // ------------ method called once each job just after ending the event loop  ------------
00153 void 
00154 HcalDAQInfo::endJob() 
00155 {
00156   if (debug_>0) std::cout << "<HcalDAQInfo::endJob> " << std::endl;
00157 }
00158 
00159 // ------------ method called just before starting a new run  ------------
00160 void 
00161 HcalDAQInfo::beginLuminosityBlock(const edm::LuminosityBlock& run, const edm::EventSetup& c)
00162 {
00163   if (debug_>0) std::cout<<"<HcalDAQInfo::beginLuminosityBlock>"<<std::endl;
00164 }
00165 
00166 // ------------ method called right after a run ends ------------
00167 void 
00168 HcalDAQInfo::endLuminosityBlock(const edm::LuminosityBlock& run, const edm::EventSetup& iSetup)
00169 {
00170   if (debug_>0) {
00171     std::cout <<"<HcalDAQInfo::endLuminosityBlock> "<<std::endl;
00172     dbe_->setCurrentFolder(rootFolder_);
00173     std::string currDir = dbe_->pwd();
00174     std::cout << "--- Current Directory " << currDir << std::endl;
00175     std::vector<MonitorElement*> mes = dbe_->getAllContents("");
00176     std::cout << "found " << mes.size() << " monitoring elements:" << std::endl;
00177   }
00178 
00179   HcalDaqFraction->Fill(-1);
00180 
00181   for (int ii=0; ii<7; ii++) DAQSummaryMap->setBinContent(ii+1,1,-1);
00182 
00183   HBDaqFraction->Fill(-1);
00184   HEDaqFraction->Fill(-1);
00185   HODaqFraction->Fill(-1);
00186   HFDaqFraction->Fill(-1);
00187   HO0DaqFraction->Fill(-1);
00188   HO12DaqFraction->Fill(-1);
00189   HFlumiDaqFraction->Fill(-1);
00190 
00191   edm::eventsetup::EventSetupRecordKey recordKey(edm::eventsetup::EventSetupRecordKey::TypeTag::findType("RunInfoRcd"));
00192 
00193   if( iSetup.find( recordKey ) ) {
00194 
00195     edm::ESHandle<RunInfo> sumFED;
00196     iSetup.get<RunInfoRcd>().get(sumFED);    
00197    
00198     std::vector<int> FedsInIds= sumFED->m_fed_in;   
00199 
00200     float HcalFedCount   = 0.;
00201     float HBFedCount     = 0.;
00202     float HEFedCount     = 0.;
00203     float HFFedCount     = 0.;
00204     float HOFedCount     = 0.;
00205     float HO0FedCount    = 0.;
00206     float HO12FedCount   = 0.;
00207     float HFlumiFedCount = 0.;
00208 
00209     // By FED taking into account Nchannels per FED
00210 
00211     for( unsigned int fedItr=0; fedItr<FedsInIds.size(); ++fedItr ) {
00212 
00213       int fedID=FedsInIds[fedItr];
00214 
00215       if (fedID >= 700 && fedID <= 731) {
00216         HcalFedCount++;
00217         if (fedID >= 700 && fedID <= 717)  {
00218           HBFedCount++;
00219           HEFedCount++;
00220         }
00221         else if (fedID >= 718 && fedID <= 723)  {
00222           HFFedCount++;
00223           HFlumiFedCount++;
00224         }
00225         else if (fedID >= 724 && fedID <= 731)  {
00226           if (fedID%2 == 0) {
00227             HOFedCount += 276;
00228             HO0FedCount += 84;
00229             HO12FedCount += 192;
00230           }
00231           else {
00232             HOFedCount += 264;
00233             HO0FedCount += 60;
00234             HO12FedCount += 204;
00235           }
00236         }
00237       }
00238 
00239       //else if ( fedID == 735 ) std::cout<<fedID<<" -- LumiScaler"<<std::endl;   
00240     }
00241 
00242     HcalFedCount = (HBFedCount*144+HEFedCount*144+HFFedCount*288+HOFedCount)/9072;
00243     HBFedCount /= 18;
00244     HEFedCount /= 18;
00245     HFFedCount /= 6;
00246     HFlumiFedCount /= 6;
00247     HOFedCount /= 2160;
00248     HO0FedCount /= 576;
00249     HO12FedCount /= 1584;
00250 
00251     DAQSummaryMap->setBinContent(1,1,HBFedCount);
00252     DAQSummaryMap->setBinContent(2,1,HEFedCount);
00253     DAQSummaryMap->setBinContent(3,1,HOFedCount);
00254     DAQSummaryMap->setBinContent(4,1,HFFedCount);
00255     DAQSummaryMap->setBinContent(5,1,HO0FedCount);
00256     DAQSummaryMap->setBinContent(6,1,HO12FedCount);
00257     DAQSummaryMap->setBinContent(7,1,HFlumiFedCount);
00258 
00259     HcalDaqFraction->Fill(HcalFedCount);
00260     HBDaqFraction->Fill(HBFedCount);
00261     HEDaqFraction->Fill(HEFedCount);
00262     HFDaqFraction->Fill(HFFedCount);
00263     HODaqFraction->Fill(HOFedCount);
00264     HO0DaqFraction->Fill(HO0FedCount);
00265     HO12DaqFraction->Fill(HO12FedCount);
00266     HFlumiDaqFraction->Fill(HFlumiFedCount);
00267 
00268     if (debug_>0) {
00269       std::cout<<" HcalFedCount= "<<HcalFedCount<<std::endl;
00270       std::cout<<" HBFedCount= "<<HBFedCount<<std::endl;
00271       std::cout<<" HEFedCount= "<<HEFedCount<<std::endl;
00272       std::cout<<" HFFedCount= "<<HFFedCount<<std::endl;
00273       std::cout<<" HOFedCount= "<<HOFedCount<<std::endl;
00274       std::cout<<" HO0FedCount= "<<HO0FedCount<<std::endl;
00275       std::cout<<" HO12FedCount= "<<HO12FedCount<<std::endl;
00276       std::cout<<" HFlumiFedCount= "<<HFlumiFedCount<<std::endl;
00277     }
00278   }
00279   else edm::LogInfo(rootFolder_+"/EventInfo/")<<"No RunInfoRcd"<<std::endl;
00280 
00281 // ---------------------- end of DAQ Info
00282   if (debug_>0) std::cout << "HcalDAQInfo::MEfilled " << std::endl;
00283 
00284 }
00285 
00286 //define this as a plug-in
00287 DEFINE_FWK_MODULE(HcalDAQInfo);