CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
L1TGMTClient.cc
Go to the documentation of this file.
2 
8 
9 #include <sstream>
10 
12  parameters_ = ps;
13  initialize();
14 }
15 
17  LogDebug("TriggerDQM") << "[TriggerDQM]: ending... ";
18 }
19 
20 //--------------------------------------------------------
22 
23  // get back-end interface
25 
26  // base folder for the contents of this job
28  "monitorName", "");
29  LogDebug("TriggerDQM") << "Monitor name = " << monitorName_ << std::endl;
30 
32  "");
33  LogDebug("TriggerDQM") << "DQM output dir = " << output_dir_ << std::endl;
34 
36  "");
37  LogDebug("TriggerDQM") << "DQM input dir = " << input_dir_ << std::endl;
38 
40  "runInEventLoop", false);
41  m_runInEndLumi = parameters_.getUntrackedParameter<bool> ("runInEndLumi",
42  false);
43  m_runInEndRun = parameters_.getUntrackedParameter<bool> ("runInEndRun",
44  false);
45  m_runInEndJob = parameters_.getUntrackedParameter<bool> ("runInEndJob",
46  false);
47 
48 }
49 
50 //--------------------------------------------------------
52 
53  LogDebug("TriggerDQM") << "[TriggerDQM]: Begin Job";
54 
55 }
56 
57 //--------------------------------------------------------
58 void L1TGMTClient::beginRun(const edm::Run& r, const edm::EventSetup& evSetup) {
59 
60  // booking histograms in the output_dir_
61 
63 
64  eff_eta_dtcsc = bookClone1DVB("eff_eta_dtcsc", "efficiency DTCSC vs eta",
65  "eta_DTCSC_and_RPC");
66 
67  if (eff_eta_dtcsc != 0) {
68  eff_eta_dtcsc->setAxisTitle("eta", 1);
69  eff_eta_dtcsc->getTH1F()->Sumw2();
70 
71  }
72 
73  eff_eta_rpc = bookClone1DVB("eff_eta_rpc", "efficiency RPC vs eta",
74  "eta_DTCSC_and_RPC");
75 
76  if (eff_eta_rpc != 0) {
77  eff_eta_rpc->setAxisTitle("eta", 1);
78  eff_eta_rpc->getTH1F()->Sumw2();
79 
80  }
81 
82  eff_phi_dtcsc = bookClone1D("eff_phi_dtcsc", "efficiency DTCSC vs phi",
83  "phi_DTCSC_and_RPC");
84 
85  if (eff_phi_dtcsc != 0) {
86  eff_phi_dtcsc->setAxisTitle("phi (deg)", 1);
87  eff_phi_dtcsc->getTH1F()->Sumw2();
88 
89  }
90 
91  eff_phi_rpc = bookClone1D("eff_phi_rpc", "efficiency RPC vs phi",
92  "phi_DTCSC_and_RPC");
93 
94  if (eff_phi_rpc != 0) {
95  eff_phi_rpc->setAxisTitle("phi (deg)", 1);
96  eff_phi_rpc->getTH1F()->Sumw2();
97 
98  }
99 
100  eff_etaphi_dtcsc = bookClone2D("eff_etaphi_dtcsc",
101  "efficiency DTCSC vs eta and phi", "etaphi_DTCSC_and_RPC");
102 
103  if (eff_etaphi_dtcsc != 0) {
104  eff_etaphi_dtcsc->setAxisTitle("eta", 1);
105  eff_etaphi_dtcsc->setAxisTitle("phi (deg)", 2);
106  eff_etaphi_dtcsc->getTH2F()->Sumw2();
107 
108  }
109 
110  eff_etaphi_rpc = bookClone2D("eff_etaphi_rpc",
111  "efficiency RPC vs eta and phi", "etaphi_DTCSC_and_RPC");
112 
113  if (eff_etaphi_rpc != 0) {
114  eff_etaphi_rpc->setAxisTitle("eta", 1);
115  eff_etaphi_rpc->setAxisTitle("phi (deg)", 2);
116  eff_etaphi_rpc->getTH2F()->Sumw2();
117 
118  }
119 }
120 
121 //--------------------------------------------------------
123  const edm::EventSetup& evSetup) {
124 
125  // empty
126 
127 }
128 //--------------------------------------------------------
129 
131  const edm::EventSetup& evSetup) {
132 
133  if (m_runInEndLumi) {
134 
136  }
137 
138 }
139 
140 //--------------------------------------------------------
142  const edm::EventSetup& evSetup) {
143 
144  // there is no loop on events in the offline harvesting step
145  // code here will not be executed offline
146 
147  if (m_runInEventLoop) {
148 
150  }
151 
152 }
153 
154 //--------------------------------------------------------
156 
157  LogDebug("TriggerDQM") << "L1TGMTClient: processing..." << std::endl;
158 
159  makeEfficiency1D(eff_eta_dtcsc, "eta_DTCSC_and_RPC", "eta_RPC_only");
160  makeEfficiency1D(eff_eta_rpc, "eta_DTCSC_and_RPC", "eta_DTCSC_only");
161 
162  makeEfficiency1D(eff_phi_dtcsc, "phi_DTCSC_and_RPC", "phi_RPC_only");
163  makeEfficiency1D(eff_phi_rpc, "phi_DTCSC_and_RPC", "phi_DTCSC_only");
164 
165  makeEfficiency2D(eff_etaphi_dtcsc, "etaphi_DTCSC_and_RPC",
166  "etaphi_RPC_only");
167  makeEfficiency2D(eff_etaphi_rpc, "etaphi_DTCSC_and_RPC",
168  "etaphi_DTCSC_only");
169 
170 }
171 
172 //--------------------------------------------------------
173 void L1TGMTClient::endRun(const edm::Run& r, const edm::EventSetup& evSetup) {
174 
175  if (m_runInEndRun) {
176 
178  }
179 
180 }
181 
182 //--------------------------------------------------------
184 
185  if (m_runInEndJob) {
186 
188  }
189 
190 }
191 
194  std::string h2Name) {
195 
197 
198  TH1F* h1 = get1DHisto(input_dir_ + "/" + h1Name, dbe_);
199  TH1F* h2 = get1DHisto(input_dir_ + "/" + h2Name, dbe_);
200 
201  if (mer == 0) {
202  LogDebug("TriggerDQM")
203  << "\nL1TGMTClient::makeRatio1D: monitoring element zero, not able to retrieve histogram"
204  << std::endl;
205  return;
206  }
207 
208  TH1F* hr = mer->getTH1F();
209 
210  if (hr && h1 && h2) {
211  hr->Divide(h1, h2, 1., 1., " ");
212  }
213 }
214 
217  std::string hiName) {
218 
220 
221  TH1F* he = get1DHisto(input_dir_ + "/" + heName, dbe_);
222  TH1F* hi = get1DHisto(input_dir_ + "/" + hiName, dbe_);
223 
224  if (meeff == 0) {
225  LogDebug("TriggerDQM")
226  << "L1TGMTClient::makeEfficiency1D: monitoring element zero, not able to retrieve histogram"
227  << std::endl;
228  return;
229  }
230 
231  TH1F* heff = meeff->getTH1F();
232 
233  if (heff && he && hi) {
234  TH1F* hall = (TH1F*) he->Clone("hall");
235  hall->Add(hi);
236  heff->Divide(he, hall, 1., 1., "B");
237  delete hall;
238  }
239 }
240 
243  std::string hiName) {
244 
246 
247  TH2F* he = get2DHisto(input_dir_ + "/" + heName, dbe_);
248  TH2F* hi = get2DHisto(input_dir_ + "/" + hiName, dbe_);
249 
250  if (meeff == 0) {
251  LogDebug("TriggerDQM")
252  << "\nL1TGMTClient::makeEfficiency2D: monitoring element zero, not able to retrieve histogram"
253  << std::endl;
254  return;
255  }
256 
257  TH2F* heff = meeff->getTH2F();
258 
259  if (heff && he && hi) {
260  TH2F* hall = (TH2F*) he->Clone("hall");
261  hall->Add(hi);
262  heff->Divide(he, hall, 1., 1., "B");
263  delete hall;
264  }
265 }
266 
269 
270  MonitorElement * me_ = dbi->get(meName);
271 
272  if (!me_) {
273  LogDebug("TriggerDQM") << "\nL1TGMTClient: " << meName << " NOT FOUND.";
274  return 0;
275  }
276 
277  return me_->getTH1F();
278 }
279 
282  MonitorElement * me_ = dbi->get(meName);
283 
284  if (!me_) {
285  LogDebug("TriggerDQM") << "\nL1TGMTClient: " << meName << " NOT FOUND.";
286  return 0;
287  }
288  return me_->getTH2F();
289 }
290 
293  const std::string& title, const std::string& hrefName) {
294 
295  MonitorElement* me;
296 
297  TH1F* href = get1DHisto(input_dir_ + "/" + hrefName, dbe_);
298 
299  if (href) {
300  LogDebug("TriggerDQM")
301  << "\nL1TGMTClient::bookClone1D: booking histogram "
302  << hrefName << std::endl;
303  const unsigned nbx = href->GetNbinsX();
304  const double xmin = href->GetXaxis()->GetXmin();
305  const double xmax = href->GetXaxis()->GetXmax();
306  me = dbe_->book1D(name, title, nbx, xmin, xmax);
307  } else {
308  LogDebug("TriggerDQM")
309  << "\nL1TGMTClient::bookClone1D: not able to clone histogram "
310  << hrefName << std::endl;
311  me = 0;
312  }
313 
314  return me;
315 }
316 
319  const std::string& title, const std::string& hrefName) {
320 
321  MonitorElement* me;
322 
323  TH1F* href = get1DHisto(input_dir_ + "/" + hrefName, dbe_);
324 
325  if (href) {
326  LogDebug("TriggerDQM")
327  << "\nL1TGMTClient::bookClone1DVB: booking histogram "
328  << hrefName << std::endl;
329  int nbx = href->GetNbinsX();
330  if (nbx > 99)
331  nbx = 99;
332  float xbins[100];
333  for (int i = 0; i < nbx; i++) {
334  xbins[i] = href->GetBinLowEdge(i + 1);
335  }
336  xbins[nbx] = href->GetXaxis()->GetXmax();
337 
339  me = dbe_->book1D(name, title, nbx, xbins);
340 
341  } else {
342  LogDebug("TriggerDQM")
343  << "\nL1TGMTClient::bookClone1DVB: not able to clone histogram "
344  << hrefName << std::endl;
345  me = 0;
346  }
347 
348  return me;
349 }
350 
353  const std::string& title, const std::string& hrefName) {
354 
355  MonitorElement* me;
356 
357  TH2F* href = get2DHisto(input_dir_ + "/" + hrefName, dbe_);
358 
359  if (href) {
360  LogDebug("TriggerDQM")
361  << "\nL1TGMTClient::bookClone2D: booking histogram "
362  << hrefName << std::endl;
363  const unsigned nbx = href->GetNbinsX();
364  const double xmin = href->GetXaxis()->GetXmin();
365  const double xmax = href->GetXaxis()->GetXmax();
366  const unsigned nby = href->GetNbinsY();
367  const double ymin = href->GetYaxis()->GetXmin();
368  const double ymax = href->GetYaxis()->GetXmax();
369  me = dbe_->book2D(name, title, nbx, xmin, xmax, nby, ymin, ymax);
370  } else {
371  LogDebug("TriggerDQM")
372  << "\nL1TGMTClient::bookClone2D: not able to clone histogram "
373  << hrefName << std::endl;
374  me = 0;
375  }
376 
377  return me;
378 }
379 
381 
382 
#define LogDebug(id)
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
std::string monitorName_
Definition: L1TGMTClient.h:68
void endJob()
Endjob.
void endLuminosityBlock(const edm::LuminosityBlock &, const edm::EventSetup &)
DQM Client Diagnostic.
MonitorElement * bookClone1DVB(const std::string &name, const std::string &title, const std::string &hrefName)
void beginRun(const edm::Run &, const edm::EventSetup &)
BeginRun.
Definition: L1TGMTClient.cc:58
const double xbins[]
void beginLuminosityBlock(const edm::LuminosityBlock &, const edm::EventSetup &)
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:954
MonitorElement * eff_eta_rpc
Definition: L1TGMTClient.h:79
bool m_runInEndJob
Definition: L1TGMTClient.h:75
MonitorElement * eff_etaphi_dtcsc
Definition: L1TGMTClient.h:82
DQMStore * dbe_
Definition: L1TGMTClient.h:67
virtual ~L1TGMTClient()
Destructor.
Definition: L1TGMTClient.cc:16
void initialize()
Definition: L1TGMTClient.cc:21
std::string input_dir_
Definition: L1TGMTClient.h:69
int iEvent
Definition: GenABIO.cc:230
bool m_runInEventLoop
Definition: L1TGMTClient.h:72
void analyze(const edm::Event &, const edm::EventSetup &)
Fake Analyze.
susybsm::HSCParticleRef hr
Definition: classes.h:26
bool m_runInEndLumi
Definition: L1TGMTClient.h:73
MonitorElement * eff_etaphi_rpc
Definition: L1TGMTClient.h:83
void endRun(const edm::Run &, const edm::EventSetup &)
EndRun.
MonitorElement * get(const std::string &path) const
get ME from full pathname (e.g. &quot;my/long/dir/my_histo&quot;)
Definition: DQMStore.cc:1708
MonitorElement * eff_phi_rpc
Definition: L1TGMTClient.h:81
MonitorElement * eff_eta_dtcsc
Definition: L1TGMTClient.h:78
MonitorElement * eff_phi_dtcsc
Definition: L1TGMTClient.h:80
void makeEfficiency1D(MonitorElement *meeff, std::string heName, std::string hiName)
TH1F * getTH1F(void) const
bool m_runInEndRun
Definition: L1TGMTClient.h:74
L1TGMTClient(const edm::ParameterSet &)
Constructor.
Definition: L1TGMTClient.cc:11
MonitorElement * bookClone1D(const std::string &name, const std::string &title, const std::string &hrefName)
TH1F * get1DHisto(std::string meName, DQMStore *dbi)
void makeRatio1D(MonitorElement *mer, std::string h1Name, std::string h2Name)
void beginJob()
BeginJob.
Definition: L1TGMTClient.cc:51
MonitorElement * bookClone2D(const std::string &name, const std::string &title, const std::string &hrefName)
std::string output_dir_
Definition: L1TGMTClient.h:70
void makeEfficiency2D(MonitorElement *meeff, std::string heName, std::string hiName)
TH2F * get2DHisto(std::string meName, DQMStore *dbi)
void processHistograms()
TH2F * getTH2F(void) const
MonitorElement * book2D(const char *name, const char *title, int nchX, double lowX, double highX, int nchY, double lowY, double highY)
Book 2D histogram.
Definition: DQMStore.cc:1082
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:667
edm::ParameterSet parameters_
Definition: L1TGMTClient.h:66
Definition: Run.h:41