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
00024 dbe_ = edm::Service<DQMStore>().operator->();
00025
00026
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
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 if (!eff_eta_dtcsc->getTH1F()->GetSumw2N())
00070 eff_eta_dtcsc->getTH1F()->Sumw2();
00071
00072 }
00073
00074 eff_eta_rpc = bookClone1DVB("eff_eta_rpc", "efficiency RPC vs eta",
00075 "eta_DTCSC_and_RPC");
00076
00077 if (eff_eta_rpc != 0) {
00078 eff_eta_rpc->setAxisTitle("eta", 1);
00079 if (!eff_eta_rpc->getTH1F()->GetSumw2N())
00080 eff_eta_rpc->getTH1F()->Sumw2();
00081
00082 }
00083
00084 eff_phi_dtcsc = bookClone1D("eff_phi_dtcsc", "efficiency DTCSC vs phi",
00085 "phi_DTCSC_and_RPC");
00086
00087 if (eff_phi_dtcsc != 0) {
00088 eff_phi_dtcsc->setAxisTitle("phi (deg)", 1);
00089 if (!eff_phi_dtcsc->getTH1F()->GetSumw2N())
00090 eff_phi_dtcsc->getTH1F()->Sumw2();
00091
00092 }
00093
00094 eff_phi_rpc = bookClone1D("eff_phi_rpc", "efficiency RPC vs phi",
00095 "phi_DTCSC_and_RPC");
00096
00097 if (eff_phi_rpc != 0) {
00098 eff_phi_rpc->setAxisTitle("phi (deg)", 1);
00099 if (!eff_phi_rpc->getTH1F()->GetSumw2N())
00100 eff_phi_rpc->getTH1F()->Sumw2();
00101
00102 }
00103
00104 eff_etaphi_dtcsc = bookClone2D("eff_etaphi_dtcsc",
00105 "efficiency DTCSC vs eta and phi", "etaphi_DTCSC_and_RPC");
00106
00107 if (eff_etaphi_dtcsc != 0) {
00108 eff_etaphi_dtcsc->setAxisTitle("eta", 1);
00109 eff_etaphi_dtcsc->setAxisTitle("phi (deg)", 2);
00110 if (!eff_etaphi_dtcsc->getTH2F()->GetSumw2N())
00111 eff_etaphi_dtcsc->getTH2F()->Sumw2();
00112
00113 }
00114
00115 eff_etaphi_rpc = bookClone2D("eff_etaphi_rpc",
00116 "efficiency RPC vs eta and phi", "etaphi_DTCSC_and_RPC");
00117
00118 if (eff_etaphi_rpc != 0) {
00119 eff_etaphi_rpc->setAxisTitle("eta", 1);
00120 eff_etaphi_rpc->setAxisTitle("phi (deg)", 2);
00121 if (!eff_etaphi_rpc->getTH2F()->GetSumw2N())
00122 eff_etaphi_rpc->getTH2F()->Sumw2();
00123
00124 }
00125 }
00126
00127
00128 void L1TGMTClient::beginLuminosityBlock(const edm::LuminosityBlock& lumiSeg,
00129 const edm::EventSetup& evSetup) {
00130
00131
00132
00133 }
00134
00135
00136 void L1TGMTClient::endLuminosityBlock(const edm::LuminosityBlock& lumiSeg,
00137 const edm::EventSetup& evSetup) {
00138
00139 if (m_runInEndLumi) {
00140
00141 processHistograms();
00142 }
00143
00144 }
00145
00146
00147 void L1TGMTClient::analyze(const edm::Event& iEvent,
00148 const edm::EventSetup& evSetup) {
00149
00150
00151
00152
00153 if (m_runInEventLoop) {
00154
00155 processHistograms();
00156 }
00157
00158 }
00159
00160
00161 void L1TGMTClient::processHistograms() {
00162
00163 LogDebug("TriggerDQM") << "L1TGMTClient: processing..." << std::endl;
00164
00165 makeEfficiency1D(eff_eta_dtcsc, "eta_DTCSC_and_RPC", "eta_RPC_only");
00166 makeEfficiency1D(eff_eta_rpc, "eta_DTCSC_and_RPC", "eta_DTCSC_only");
00167
00168 makeEfficiency1D(eff_phi_dtcsc, "phi_DTCSC_and_RPC", "phi_RPC_only");
00169 makeEfficiency1D(eff_phi_rpc, "phi_DTCSC_and_RPC", "phi_DTCSC_only");
00170
00171 makeEfficiency2D(eff_etaphi_dtcsc, "etaphi_DTCSC_and_RPC",
00172 "etaphi_RPC_only");
00173 makeEfficiency2D(eff_etaphi_rpc, "etaphi_DTCSC_and_RPC",
00174 "etaphi_DTCSC_only");
00175
00176 }
00177
00178
00179 void L1TGMTClient::endRun(const edm::Run& r, const edm::EventSetup& evSetup) {
00180
00181 if (m_runInEndRun) {
00182
00183 processHistograms();
00184 }
00185
00186 }
00187
00188
00189 void L1TGMTClient::endJob() {
00190
00191 if (m_runInEndJob) {
00192
00193 processHistograms();
00194 }
00195
00196 }
00197
00199 void L1TGMTClient::makeRatio1D(MonitorElement* mer, std::string h1Name,
00200 std::string h2Name) {
00201
00202 dbe_->setCurrentFolder(output_dir_);
00203
00204 TH1F* h1 = get1DHisto(input_dir_ + "/" + h1Name, dbe_);
00205 TH1F* h2 = get1DHisto(input_dir_ + "/" + h2Name, dbe_);
00206
00207 if (mer == 0) {
00208 LogDebug("TriggerDQM")
00209 << "\nL1TGMTClient::makeRatio1D: monitoring element zero, not able to retrieve histogram"
00210 << std::endl;
00211 return;
00212 }
00213
00214 TH1F* hr = mer->getTH1F();
00215
00216 if (hr && h1 && h2) {
00217 hr->Divide(h1, h2, 1., 1., " ");
00218 }
00219 }
00220
00222 void L1TGMTClient::makeEfficiency1D(MonitorElement* meeff, std::string heName,
00223 std::string hiName) {
00224
00225 dbe_->setCurrentFolder(output_dir_);
00226
00227 TH1F* he = get1DHisto(input_dir_ + "/" + heName, dbe_);
00228 TH1F* hi = get1DHisto(input_dir_ + "/" + hiName, dbe_);
00229
00230 if (meeff == 0) {
00231 LogDebug("TriggerDQM")
00232 << "L1TGMTClient::makeEfficiency1D: monitoring element zero, not able to retrieve histogram"
00233 << std::endl;
00234 return;
00235 }
00236
00237 TH1F* heff = meeff->getTH1F();
00238
00239 if (heff && he && hi) {
00240 TH1F* hall = (TH1F*) he->Clone("hall");
00241 hall->Add(hi);
00242 heff->Divide(he, hall, 1., 1., "B");
00243 delete hall;
00244 }
00245 }
00246
00248 void L1TGMTClient::makeEfficiency2D(MonitorElement* meeff, std::string heName,
00249 std::string hiName) {
00250
00251 dbe_->setCurrentFolder(output_dir_);
00252
00253 TH2F* he = get2DHisto(input_dir_ + "/" + heName, dbe_);
00254 TH2F* hi = get2DHisto(input_dir_ + "/" + hiName, dbe_);
00255
00256 if (meeff == 0) {
00257 LogDebug("TriggerDQM")
00258 << "\nL1TGMTClient::makeEfficiency2D: monitoring element zero, not able to retrieve histogram"
00259 << std::endl;
00260 return;
00261 }
00262
00263 TH2F* heff = meeff->getTH2F();
00264
00265 if (heff && he && hi) {
00266 TH2F* hall = (TH2F*) he->Clone("hall");
00267 hall->Add(hi);
00268 heff->Divide(he, hall, 1., 1., "B");
00269 delete hall;
00270 }
00271 }
00272
00274 TH1F* L1TGMTClient::get1DHisto(std::string meName, DQMStore* dbi) {
00275
00276 MonitorElement * me_ = dbi->get(meName);
00277
00278 if (!me_) {
00279 LogDebug("TriggerDQM") << "\nL1TGMTClient: " << meName << " NOT FOUND.";
00280 return 0;
00281 }
00282
00283 return me_->getTH1F();
00284 }
00285
00287 TH2F* L1TGMTClient::get2DHisto(std::string meName, DQMStore* dbi) {
00288 MonitorElement * me_ = dbi->get(meName);
00289
00290 if (!me_) {
00291 LogDebug("TriggerDQM") << "\nL1TGMTClient: " << meName << " NOT FOUND.";
00292 return 0;
00293 }
00294 return me_->getTH2F();
00295 }
00296
00298 MonitorElement* L1TGMTClient::bookClone1D(const std::string& name,
00299 const std::string& title, const std::string& hrefName) {
00300
00301 MonitorElement* me;
00302
00303 TH1F* href = get1DHisto(input_dir_ + "/" + hrefName, dbe_);
00304
00305 if (href) {
00306 LogDebug("TriggerDQM")
00307 << "\nL1TGMTClient::bookClone1D: booking histogram "
00308 << hrefName << std::endl;
00309 const unsigned nbx = href->GetNbinsX();
00310 const double xmin = href->GetXaxis()->GetXmin();
00311 const double xmax = href->GetXaxis()->GetXmax();
00312 me = dbe_->book1D(name, title, nbx, xmin, xmax);
00313 } else {
00314 LogDebug("TriggerDQM")
00315 << "\nL1TGMTClient::bookClone1D: not able to clone histogram "
00316 << hrefName << std::endl;
00317 me = 0;
00318 }
00319
00320 return me;
00321 }
00322
00324 MonitorElement* L1TGMTClient::bookClone1DVB(const std::string& name,
00325 const std::string& title, const std::string& hrefName) {
00326
00327 MonitorElement* me;
00328
00329 TH1F* href = get1DHisto(input_dir_ + "/" + hrefName, dbe_);
00330
00331 if (href) {
00332 LogDebug("TriggerDQM")
00333 << "\nL1TGMTClient::bookClone1DVB: booking histogram "
00334 << hrefName << std::endl;
00335 int nbx = href->GetNbinsX();
00336 if (nbx > 99)
00337 nbx = 99;
00338 float xbins[100];
00339 for (int i = 0; i < nbx; i++) {
00340 xbins[i] = href->GetBinLowEdge(i + 1);
00341 }
00342 xbins[nbx] = href->GetXaxis()->GetXmax();
00343
00344 dbe_->setCurrentFolder(output_dir_);
00345 me = dbe_->book1D(name, title, nbx, xbins);
00346
00347 } else {
00348 LogDebug("TriggerDQM")
00349 << "\nL1TGMTClient::bookClone1DVB: not able to clone histogram "
00350 << hrefName << std::endl;
00351 me = 0;
00352 }
00353
00354 return me;
00355 }
00356
00358 MonitorElement* L1TGMTClient::bookClone2D(const std::string& name,
00359 const std::string& title, const std::string& hrefName) {
00360
00361 MonitorElement* me;
00362
00363 TH2F* href = get2DHisto(input_dir_ + "/" + hrefName, dbe_);
00364
00365 if (href) {
00366 LogDebug("TriggerDQM")
00367 << "\nL1TGMTClient::bookClone2D: booking histogram "
00368 << hrefName << std::endl;
00369 const unsigned nbx = href->GetNbinsX();
00370 const double xmin = href->GetXaxis()->GetXmin();
00371 const double xmax = href->GetXaxis()->GetXmax();
00372 const unsigned nby = href->GetNbinsY();
00373 const double ymin = href->GetYaxis()->GetXmin();
00374 const double ymax = href->GetYaxis()->GetXmax();
00375 me = dbe_->book2D(name, title, nbx, xmin, xmax, nby, ymin, ymax);
00376 } else {
00377 LogDebug("TriggerDQM")
00378 << "\nL1TGMTClient::bookClone2D: not able to clone histogram "
00379 << hrefName << std::endl;
00380 me = 0;
00381 }
00382
00383 return me;
00384 }
00385
00387
00388