CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_8_patch3/src/DQM/L1TMonitorClient/src/L1TGMTClient.cc

Go to the documentation of this file.
00001 #include "DQM/L1TMonitorClient/interface/L1TGMTClient.h"
00002 
00003 #include "FWCore/ServiceRegistry/interface/Service.h"
00004 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00006 #include "DQMServices/Core/interface/DQMStore.h"
00007 #include "DQMServices/Core/interface/MonitorElement.h"
00008 
00009 #include <sstream>
00010 
00011 L1TGMTClient::L1TGMTClient(const edm::ParameterSet& ps) {
00012     parameters_ = ps;
00013     initialize();
00014 }
00015 
00016 L1TGMTClient::~L1TGMTClient() {
00017     LogDebug("TriggerDQM") << "[TriggerDQM]: ending... ";
00018 }
00019 
00020 //--------------------------------------------------------
00021 void L1TGMTClient::initialize() {
00022 
00023     // get back-end interface
00024     dbe_ = edm::Service<DQMStore>().operator->();
00025 
00026     // base folder for the contents of this job
00027     monitorName_ = parameters_.getUntrackedParameter<std::string> (
00028             "monitorName", "");
00029     LogDebug("TriggerDQM") << "Monitor name = " << monitorName_ << std::endl;
00030 
00031     output_dir_ = parameters_.getUntrackedParameter<std::string> ("output_dir",
00032             "");
00033     LogDebug("TriggerDQM") << "DQM output dir = " << output_dir_ << std::endl;
00034 
00035     input_dir_ = parameters_.getUntrackedParameter<std::string> ("input_dir",
00036             "");
00037     LogDebug("TriggerDQM") << "DQM input dir = " << input_dir_ << std::endl;
00038 
00039     m_runInEventLoop = parameters_.getUntrackedParameter<bool> (
00040             "runInEventLoop", false);
00041     m_runInEndLumi = parameters_.getUntrackedParameter<bool> ("runInEndLumi",
00042             false);
00043     m_runInEndRun = parameters_.getUntrackedParameter<bool> ("runInEndRun",
00044             false);
00045     m_runInEndJob = parameters_.getUntrackedParameter<bool> ("runInEndJob",
00046             false);
00047 
00048 }
00049 
00050 //--------------------------------------------------------
00051 void L1TGMTClient::beginJob() {
00052 
00053     LogDebug("TriggerDQM") << "[TriggerDQM]: Begin Job";
00054 
00055 }
00056 
00057 //--------------------------------------------------------
00058 void L1TGMTClient::beginRun(const edm::Run& r, const edm::EventSetup& evSetup) {
00059 
00060     // booking histograms in the output_dir_
00061 
00062     dbe_->setCurrentFolder(output_dir_);
00063 
00064     eff_eta_dtcsc = bookClone1DVB("eff_eta_dtcsc", "efficiency DTCSC vs eta",
00065             "eta_DTCSC_and_RPC");
00066 
00067     if (eff_eta_dtcsc != 0) {
00068         eff_eta_dtcsc->setAxisTitle("eta", 1);
00069         eff_eta_dtcsc->getTH1F()->Sumw2();
00070 
00071     }
00072 
00073     eff_eta_rpc = bookClone1DVB("eff_eta_rpc", "efficiency RPC vs eta",
00074             "eta_DTCSC_and_RPC");
00075 
00076     if (eff_eta_rpc != 0) {
00077         eff_eta_rpc->setAxisTitle("eta", 1);
00078         eff_eta_rpc->getTH1F()->Sumw2();
00079 
00080     }
00081 
00082     eff_phi_dtcsc = bookClone1D("eff_phi_dtcsc", "efficiency DTCSC vs phi",
00083             "phi_DTCSC_and_RPC");
00084 
00085     if (eff_phi_dtcsc != 0) {
00086         eff_phi_dtcsc->setAxisTitle("phi (deg)", 1);
00087         eff_phi_dtcsc->getTH1F()->Sumw2();
00088 
00089     }
00090 
00091     eff_phi_rpc = bookClone1D("eff_phi_rpc", "efficiency RPC vs phi",
00092             "phi_DTCSC_and_RPC");
00093 
00094     if (eff_phi_rpc != 0) {
00095         eff_phi_rpc->setAxisTitle("phi (deg)", 1);
00096         eff_phi_rpc->getTH1F()->Sumw2();
00097 
00098     }
00099 
00100     eff_etaphi_dtcsc = bookClone2D("eff_etaphi_dtcsc",
00101             "efficiency DTCSC vs eta and phi", "etaphi_DTCSC_and_RPC");
00102 
00103     if (eff_etaphi_dtcsc != 0) {
00104         eff_etaphi_dtcsc->setAxisTitle("eta", 1);
00105         eff_etaphi_dtcsc->setAxisTitle("phi (deg)", 2);
00106         eff_etaphi_dtcsc->getTH2F()->Sumw2();
00107 
00108     }
00109 
00110     eff_etaphi_rpc = bookClone2D("eff_etaphi_rpc",
00111             "efficiency RPC vs eta and phi", "etaphi_DTCSC_and_RPC");
00112 
00113     if (eff_etaphi_rpc != 0) {
00114         eff_etaphi_rpc->setAxisTitle("eta", 1);
00115         eff_etaphi_rpc->setAxisTitle("phi (deg)", 2);
00116         eff_etaphi_rpc->getTH2F()->Sumw2();
00117 
00118     }
00119 }
00120 
00121 //--------------------------------------------------------
00122 void L1TGMTClient::beginLuminosityBlock(const edm::LuminosityBlock& lumiSeg,
00123         const edm::EventSetup& evSetup) {
00124 
00125     // empty
00126 
00127 }
00128 //--------------------------------------------------------
00129 
00130 void L1TGMTClient::endLuminosityBlock(const edm::LuminosityBlock& lumiSeg,
00131         const edm::EventSetup& evSetup) {
00132 
00133     if (m_runInEndLumi) {
00134 
00135         processHistograms();
00136     }
00137 
00138 }
00139 
00140 //--------------------------------------------------------
00141 void L1TGMTClient::analyze(const edm::Event& iEvent,
00142         const edm::EventSetup& evSetup) {
00143 
00144     // there is no loop on events in the offline harvesting step
00145     // code here will not be executed offline
00146 
00147     if (m_runInEventLoop) {
00148 
00149         processHistograms();
00150     }
00151 
00152 }
00153 
00154 //--------------------------------------------------------
00155 void L1TGMTClient::processHistograms() {
00156 
00157     LogDebug("TriggerDQM") << "L1TGMTClient: processing..." << std::endl;
00158 
00159     makeEfficiency1D(eff_eta_dtcsc, "eta_DTCSC_and_RPC", "eta_RPC_only");
00160     makeEfficiency1D(eff_eta_rpc, "eta_DTCSC_and_RPC", "eta_DTCSC_only");
00161 
00162     makeEfficiency1D(eff_phi_dtcsc, "phi_DTCSC_and_RPC", "phi_RPC_only");
00163     makeEfficiency1D(eff_phi_rpc, "phi_DTCSC_and_RPC", "phi_DTCSC_only");
00164 
00165     makeEfficiency2D(eff_etaphi_dtcsc, "etaphi_DTCSC_and_RPC",
00166             "etaphi_RPC_only");
00167     makeEfficiency2D(eff_etaphi_rpc, "etaphi_DTCSC_and_RPC",
00168             "etaphi_DTCSC_only");
00169 
00170 }
00171 
00172 //--------------------------------------------------------
00173 void L1TGMTClient::endRun(const edm::Run& r, const edm::EventSetup& evSetup) {
00174 
00175     if (m_runInEndRun) {
00176 
00177         processHistograms();
00178     }
00179 
00180 }
00181 
00182 //--------------------------------------------------------
00183 void L1TGMTClient::endJob() {
00184 
00185     if (m_runInEndJob) {
00186 
00187         processHistograms();
00188     }
00189 
00190 }
00191 
00193 void L1TGMTClient::makeRatio1D(MonitorElement* mer, std::string h1Name,
00194         std::string h2Name) {
00195 
00196     dbe_->setCurrentFolder(output_dir_);
00197 
00198     TH1F* h1 = get1DHisto(input_dir_ + "/" + h1Name, dbe_);
00199     TH1F* h2 = get1DHisto(input_dir_ + "/" + h2Name, dbe_);
00200 
00201     if (mer == 0) {
00202         LogDebug("TriggerDQM")
00203                 << "\nL1TGMTClient::makeRatio1D: monitoring element zero, not able to retrieve histogram"
00204                 << std::endl;
00205         return;
00206     }
00207 
00208     TH1F* hr = mer->getTH1F();
00209 
00210     if (hr && h1 && h2) {
00211         hr->Divide(h1, h2, 1., 1., " ");
00212     }
00213 }
00214 
00216 void L1TGMTClient::makeEfficiency1D(MonitorElement* meeff, std::string heName,
00217         std::string hiName) {
00218 
00219     dbe_->setCurrentFolder(output_dir_);
00220 
00221     TH1F* he = get1DHisto(input_dir_ + "/" + heName, dbe_);
00222     TH1F* hi = get1DHisto(input_dir_ + "/" + hiName, dbe_);
00223 
00224     if (meeff == 0) {
00225         LogDebug("TriggerDQM")
00226                 << "L1TGMTClient::makeEfficiency1D: monitoring element zero, not able to retrieve histogram"
00227                 << std::endl;
00228         return;
00229     }
00230 
00231     TH1F* heff = meeff->getTH1F();
00232 
00233     if (heff && he && hi) {
00234         TH1F* hall = (TH1F*) he->Clone("hall");
00235         hall->Add(hi);
00236         heff->Divide(he, hall, 1., 1., "B");
00237         delete hall;
00238     }
00239 }
00240 
00242 void L1TGMTClient::makeEfficiency2D(MonitorElement* meeff, std::string heName,
00243         std::string hiName) {
00244 
00245     dbe_->setCurrentFolder(output_dir_);
00246 
00247     TH2F* he = get2DHisto(input_dir_ + "/" + heName, dbe_);
00248     TH2F* hi = get2DHisto(input_dir_ + "/" + hiName, dbe_);
00249 
00250     if (meeff == 0) {
00251         LogDebug("TriggerDQM")
00252                 << "\nL1TGMTClient::makeEfficiency2D: monitoring element zero, not able to retrieve histogram"
00253                 << std::endl;
00254         return;
00255     }
00256 
00257     TH2F* heff = meeff->getTH2F();
00258 
00259     if (heff && he && hi) {
00260         TH2F* hall = (TH2F*) he->Clone("hall");
00261         hall->Add(hi);
00262         heff->Divide(he, hall, 1., 1., "B");
00263         delete hall;
00264     }
00265 }
00266 
00268 TH1F* L1TGMTClient::get1DHisto(std::string meName, DQMStore* dbi) {
00269 
00270     MonitorElement * me_ = dbi->get(meName);
00271 
00272     if (!me_) {
00273         LogDebug("TriggerDQM") << "\nL1TGMTClient: " << meName << " NOT FOUND.";
00274         return 0;
00275     }
00276 
00277     return me_->getTH1F();
00278 }
00279 
00281 TH2F* L1TGMTClient::get2DHisto(std::string meName, DQMStore* dbi) {
00282     MonitorElement * me_ = dbi->get(meName);
00283 
00284     if (!me_) {
00285         LogDebug("TriggerDQM") << "\nL1TGMTClient: " << meName << " NOT FOUND.";
00286         return 0;
00287     }
00288     return me_->getTH2F();
00289 }
00290 
00292 MonitorElement* L1TGMTClient::bookClone1D(const std::string& name,
00293         const std::string& title, const std::string& hrefName) {
00294 
00295     MonitorElement* me;
00296 
00297     TH1F* href = get1DHisto(input_dir_ + "/" + hrefName, dbe_);
00298 
00299     if (href) {
00300         LogDebug("TriggerDQM")
00301                 << "\nL1TGMTClient::bookClone1D: booking histogram "
00302                 << hrefName << std::endl;
00303         const unsigned nbx = href->GetNbinsX();
00304         const double xmin = href->GetXaxis()->GetXmin();
00305         const double xmax = href->GetXaxis()->GetXmax();
00306         me = dbe_->book1D(name, title, nbx, xmin, xmax);
00307     } else {
00308         LogDebug("TriggerDQM")
00309                 << "\nL1TGMTClient::bookClone1D: not able to clone histogram "
00310                 << hrefName << std::endl;
00311         me = 0;
00312     }
00313 
00314     return me;
00315 }
00316 
00318 MonitorElement* L1TGMTClient::bookClone1DVB(const std::string& name,
00319         const std::string& title, const std::string& hrefName) {
00320 
00321     MonitorElement* me;
00322 
00323     TH1F* href = get1DHisto(input_dir_ + "/" + hrefName, dbe_);
00324 
00325     if (href) {
00326         LogDebug("TriggerDQM")
00327                 << "\nL1TGMTClient::bookClone1DVB: booking histogram "
00328                 << hrefName << std::endl;
00329         int nbx = href->GetNbinsX();
00330         if (nbx > 99)
00331             nbx = 99;
00332         float xbins[100];
00333         for (int i = 0; i < nbx; i++) {
00334             xbins[i] = href->GetBinLowEdge(i + 1);
00335         }
00336         xbins[nbx] = href->GetXaxis()->GetXmax();
00337 
00338         dbe_->setCurrentFolder(output_dir_);
00339         me = dbe_->book1D(name, title, nbx, xbins);
00340 
00341     } else {
00342         LogDebug("TriggerDQM")
00343                 << "\nL1TGMTClient::bookClone1DVB: not able to clone histogram "
00344                 << hrefName << std::endl;
00345         me = 0;
00346     }
00347 
00348     return me;
00349 }
00350 
00352 MonitorElement* L1TGMTClient::bookClone2D(const std::string& name,
00353         const std::string& title, const std::string& hrefName) {
00354 
00355     MonitorElement* me;
00356 
00357     TH2F* href = get2DHisto(input_dir_ + "/" + hrefName, dbe_);
00358 
00359     if (href) {
00360         LogDebug("TriggerDQM")
00361                 << "\nL1TGMTClient::bookClone2D: booking histogram "
00362                 << hrefName << std::endl;
00363         const unsigned nbx = href->GetNbinsX();
00364         const double xmin = href->GetXaxis()->GetXmin();
00365         const double xmax = href->GetXaxis()->GetXmax();
00366         const unsigned nby = href->GetNbinsY();
00367         const double ymin = href->GetYaxis()->GetXmin();
00368         const double ymax = href->GetYaxis()->GetXmax();
00369         me = dbe_->book2D(name, title, nbx, xmin, xmax, nby, ymin, ymax);
00370     } else {
00371         LogDebug("TriggerDQM")
00372                 << "\nL1TGMTClient::bookClone2D: not able to clone histogram "
00373                 << hrefName << std::endl;
00374         me = 0;
00375     }
00376 
00377     return me;
00378 }
00379 
00381 
00382