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  // base folder for the contents of this job
25  "monitorName", "");
26  LogDebug("TriggerDQM") << "Monitor name = " << monitorName_ << std::endl;
27 
29  "");
30  LogDebug("TriggerDQM") << "DQM output dir = " << output_dir_ << std::endl;
31 
33  "");
34  LogDebug("TriggerDQM") << "DQM input dir = " << input_dir_ << std::endl;
35 
37  "runInEventLoop", false);
38  m_runInEndLumi = parameters_.getUntrackedParameter<bool> ("runInEndLumi",
39  false);
40  m_runInEndRun = parameters_.getUntrackedParameter<bool> ("runInEndRun",
41  false);
42  m_runInEndJob = parameters_.getUntrackedParameter<bool> ("runInEndJob",
43  false);
44 
45 }
46 
47 //--------------------------------------------------------
49 
50  // booking histograms in the output_dir_
51 
53 
54  eff_eta_dtcsc = bookClone1DVB(ibooker, igetter, "eff_eta_dtcsc", "efficiency DTCSC vs eta",
55  "eta_DTCSC_and_RPC");
56 
57  if (eff_eta_dtcsc != 0) {
58  eff_eta_dtcsc->setAxisTitle("eta", 1);
59  if(eff_eta_dtcsc->getTH1F()->GetSumw2N() == 0) eff_eta_dtcsc->getTH1F()->Sumw2();
60 
61  }
62 
63  eff_eta_rpc = bookClone1DVB(ibooker, igetter, "eff_eta_rpc", "efficiency RPC vs eta",
64  "eta_DTCSC_and_RPC");
65 
66  if (eff_eta_rpc != 0) {
67  eff_eta_rpc->setAxisTitle("eta", 1);
68  if(eff_eta_rpc->getTH1F()->GetSumw2N() == 0) eff_eta_rpc->getTH1F()->Sumw2();
69 
70  }
71 
72  eff_phi_dtcsc = bookClone1D(ibooker, igetter, "eff_phi_dtcsc", "efficiency DTCSC vs phi",
73  "phi_DTCSC_and_RPC");
74 
75  if (eff_phi_dtcsc != 0) {
76  eff_phi_dtcsc->setAxisTitle("phi (deg)", 1);
77  if(eff_phi_dtcsc->getTH1F()->GetSumw2N() == 0) eff_phi_dtcsc->getTH1F()->Sumw2();
78 
79  }
80 
81  eff_phi_rpc = bookClone1D(ibooker, igetter, "eff_phi_rpc", "efficiency RPC vs phi",
82  "phi_DTCSC_and_RPC");
83 
84  if (eff_phi_rpc != 0) {
85  eff_phi_rpc->setAxisTitle("phi (deg)", 1);
86  if(eff_phi_rpc->getTH1F()->GetSumw2N() == 0) eff_phi_rpc->getTH1F()->Sumw2();
87 
88  }
89 
90  eff_etaphi_dtcsc = bookClone2D(ibooker, igetter, "eff_etaphi_dtcsc",
91  "efficiency DTCSC vs eta and phi", "etaphi_DTCSC_and_RPC");
92 
93  if (eff_etaphi_dtcsc != 0) {
94  eff_etaphi_dtcsc->setAxisTitle("eta", 1);
95  eff_etaphi_dtcsc->setAxisTitle("phi (deg)", 2);
96  if(eff_etaphi_dtcsc->getTH2F()->GetSumw2N() == 0) eff_etaphi_dtcsc->getTH2F()->Sumw2();
97 
98  }
99 
100  eff_etaphi_rpc = bookClone2D(ibooker, igetter, "eff_etaphi_rpc",
101  "efficiency RPC vs eta and phi", "etaphi_DTCSC_and_RPC");
102 
103  if (eff_etaphi_rpc != 0) {
104  eff_etaphi_rpc->setAxisTitle("eta", 1);
105  eff_etaphi_rpc->setAxisTitle("phi (deg)", 2);
106  if(eff_etaphi_rpc->getTH2F()->GetSumw2N() == 0) eff_etaphi_rpc->getTH2F()->Sumw2();
107 
108  }
109 
110  processHistograms(ibooker, igetter);
111 }
112 
113 //--------------------------------------------------------
114 
116  }
117 
118 //--------------------------------------------------------
120 
121  LogDebug("TriggerDQM") << "L1TGMTClient: processing..." << std::endl;
122 
123  makeEfficiency1D(ibooker, igetter, eff_eta_dtcsc, "eta_DTCSC_and_RPC", "eta_RPC_only");
124  makeEfficiency1D(ibooker, igetter, eff_eta_rpc, "eta_DTCSC_and_RPC", "eta_DTCSC_only");
125 
126  makeEfficiency1D(ibooker, igetter, eff_phi_dtcsc, "phi_DTCSC_and_RPC", "phi_RPC_only");
127  makeEfficiency1D(ibooker, igetter, eff_phi_rpc, "phi_DTCSC_and_RPC", "phi_DTCSC_only");
128 
129  makeEfficiency2D(ibooker, igetter, eff_etaphi_dtcsc, "etaphi_DTCSC_and_RPC", "etaphi_RPC_only");
130  makeEfficiency2D(ibooker, igetter, eff_etaphi_rpc, "etaphi_DTCSC_and_RPC", "etaphi_DTCSC_only");
131 }
132 
135 
136  igetter.setCurrentFolder(output_dir_);
137 
138  TH1F* h1 = get1DHisto(input_dir_ + "/" + h1Name, igetter);
139  TH1F* h2 = get1DHisto(input_dir_ + "/" + h2Name, igetter);
140 
141  if (mer == 0) {
142  LogDebug("TriggerDQM")
143  << "\nL1TGMTClient::makeRatio1D: monitoring element zero, not able to retrieve histogram"
144  << std::endl;
145  return;
146  }
147 
148  TH1F* hr = mer->getTH1F();
149 
150  if (hr && h1 && h2) {
151  hr->Divide(h1, h2, 1., 1., " ");
152  }
153 }
154 
157 
158  igetter.setCurrentFolder(output_dir_);
159 
160  TH1F* he = get1DHisto(input_dir_ + "/" + heName, igetter);
161  TH1F* hi = get1DHisto(input_dir_ + "/" + hiName, igetter);
162 
163  if (meeff == 0) {
164  LogDebug("TriggerDQM")
165  << "L1TGMTClient::makeEfficiency1D: monitoring element zero, not able to retrieve histogram"
166  << std::endl;
167  return;
168  }
169 
170  TH1F* heff = meeff->getTH1F();
171 
172  if (heff && he && hi) {
173  TH1F* hall = (TH1F*) he->Clone("hall");
174  hall->Add(hi);
175  heff->Divide(he, hall, 1., 1., "B");
176  delete hall;
177  }
178 }
179 
182 
183  igetter.setCurrentFolder(output_dir_);
184 
185  TH2F* he = get2DHisto(input_dir_ + "/" + heName, igetter);
186  TH2F* hi = get2DHisto(input_dir_ + "/" + hiName, igetter);
187 
188  if (meeff == 0) {
189  LogDebug("TriggerDQM")
190  << "\nL1TGMTClient::makeEfficiency2D: monitoring element zero, not able to retrieve histogram"
191  << std::endl;
192  return;
193  }
194 
195  TH2F* heff = meeff->getTH2F();
196 
197  if (heff && he && hi) {
198  TH2F* hall = (TH2F*) he->Clone("hall");
199  hall->Add(hi);
200  heff->Divide(he, hall, 1., 1., "B");
201  delete hall;
202  }
203 }
204 
207 
208  MonitorElement * me_ = igetter.get(meName);
209 
210  if (!me_) {
211  LogDebug("TriggerDQM") << "\nL1TGMTClient: " << meName << " NOT FOUND.";
212  return 0;
213  }
214 
215  return me_->getTH1F();
216 }
217 
220  MonitorElement * me_ = igetter.get(meName);
221 
222  if (!me_) {
223  LogDebug("TriggerDQM") << "\nL1TGMTClient: " << meName << " NOT FOUND.";
224  return 0;
225  }
226  return me_->getTH2F();
227 }
228 
231 
232  MonitorElement* me;
233 
234  TH1F* href = get1DHisto(input_dir_ + "/" + hrefName, igetter);
235 
236  if (href) {
237  LogDebug("TriggerDQM")
238  << "\nL1TGMTClient::bookClone1D: booking histogram "
239  << hrefName << std::endl;
240  const unsigned nbx = href->GetNbinsX();
241  const double xmin = href->GetXaxis()->GetXmin();
242  const double xmax = href->GetXaxis()->GetXmax();
243  ibooker.setCurrentFolder(output_dir_);
244  me = ibooker.book1D(name, title, nbx, xmin, xmax);
245  } else {
246  LogDebug("TriggerDQM")
247  << "\nL1TGMTClient::bookClone1D: not able to clone histogram "
248  << hrefName << std::endl;
249  me = 0;
250  }
251 
252  return me;
253 }
254 
257 
258  MonitorElement* me;
259 
260  TH1F* href = get1DHisto(input_dir_ + "/" + hrefName, igetter);
261 
262  if (href) {
263  LogDebug("TriggerDQM")
264  << "\nL1TGMTClient::bookClone1DVB: booking histogram "
265  << hrefName << std::endl;
266  int nbx = href->GetNbinsX();
267  if (nbx > 99)
268  nbx = 99;
269  float xbins[100];
270  for (int i = 0; i < nbx; i++) {
271  xbins[i] = href->GetBinLowEdge(i + 1);
272  }
273  xbins[nbx] = href->GetXaxis()->GetXmax();
274 
275  ibooker.setCurrentFolder(output_dir_);
276  me = ibooker.book1D(name, title, nbx, xbins);
277 
278  } else {
279  LogDebug("TriggerDQM")
280  << "\nL1TGMTClient::bookClone1DVB: not able to clone histogram "
281  << hrefName << std::endl;
282  me = 0;
283  }
284 
285  return me;
286 }
287 
290  const std::string& title, const std::string& hrefName) {
291 
292  MonitorElement* me;
293 
294  TH2F* href = get2DHisto(input_dir_ + "/" + hrefName, igetter);
295 
296  if (href) {
297  LogDebug("TriggerDQM")
298  << "\nL1TGMTClient::bookClone2D: booking histogram "
299  << hrefName << std::endl;
300  const unsigned nbx = href->GetNbinsX();
301  const double xmin = href->GetXaxis()->GetXmin();
302  const double xmax = href->GetXaxis()->GetXmax();
303  const unsigned nby = href->GetNbinsY();
304  const double ymin = href->GetYaxis()->GetXmin();
305  const double ymax = href->GetYaxis()->GetXmax();
306  ibooker.setCurrentFolder(output_dir_);
307  me = ibooker.book2D(name, title, nbx, xmin, xmax, nby, ymin, ymax);
308  } else {
309  LogDebug("TriggerDQM")
310  << "\nL1TGMTClient::bookClone2D: not able to clone histogram "
311  << hrefName << std::endl;
312  me = 0;
313  }
314 
315  return me;
316 }
317 
319 
320 
#define LogDebug(id)
void makeEfficiency2D(DQMStore::IBooker &ibooker, DQMStore::IGetter &igetter, MonitorElement *meeff, std::string heName, std::string hiName)
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
void processHistograms(DQMStore::IBooker &ibooker, DQMStore::IGetter &igetter)
std::string monitorName_
Definition: L1TGMTClient.h:43
const double xbins[]
TH2F * get2DHisto(std::string meName, DQMStore::IGetter &igetter)
MonitorElement * eff_eta_rpc
Definition: L1TGMTClient.h:54
void makeEfficiency1D(DQMStore::IBooker &ibooker, DQMStore::IGetter &igetter, MonitorElement *meeff, std::string heName, std::string hiName)
bool m_runInEndJob
Definition: L1TGMTClient.h:50
MonitorElement * eff_etaphi_dtcsc
Definition: L1TGMTClient.h:57
MonitorElement * bookClone1DVB(DQMStore::IBooker &ibooker, DQMStore::IGetter &igetter, const std::string &name, const std::string &title, const std::string &hrefName)
MonitorElement * get(const std::string &path)
Definition: DQMStore.cc:301
virtual ~L1TGMTClient()
Destructor.
Definition: L1TGMTClient.cc:16
void initialize()
Definition: L1TGMTClient.cc:21
MonitorElement * bookClone1D(DQMStore::IBooker &ibooker, DQMStore::IGetter &igetter, const std::string &name, const std::string &title, const std::string &hrefName)
std::string input_dir_
Definition: L1TGMTClient.h:44
bool m_runInEventLoop
Definition: L1TGMTClient.h:47
susybsm::HSCParticleRef hr
Definition: classes.h:26
bool m_runInEndLumi
Definition: L1TGMTClient.h:48
MonitorElement * bookClone2D(DQMStore::IBooker &ibooker, DQMStore::IGetter &igetter, const std::string &name, const std::string &title, const std::string &hrefName)
MonitorElement * eff_etaphi_rpc
Definition: L1TGMTClient.h:58
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
MonitorElement * eff_phi_rpc
Definition: L1TGMTClient.h:56
MonitorElement * eff_eta_dtcsc
Definition: L1TGMTClient.h:53
MonitorElement * eff_phi_dtcsc
Definition: L1TGMTClient.h:55
void makeRatio1D(DQMStore::IBooker &ibooker, DQMStore::IGetter &igetter, MonitorElement *mer, std::string h1Name, std::string h2Name)
virtual void dqmEndJob(DQMStore::IBooker &ibooker, DQMStore::IGetter &igetter) override
Definition: L1TGMTClient.cc:48
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:273
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:133
TH1F * getTH1F(void) const
bool m_runInEndRun
Definition: L1TGMTClient.h:49
L1TGMTClient(const edm::ParameterSet &)
Constructor.
Definition: L1TGMTClient.cc:11
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:329
virtual void dqmEndLuminosityBlock(DQMStore::IBooker &ibooker, DQMStore::IGetter &igetter, const edm::LuminosityBlock &lumiSeg, const edm::EventSetup &evSetup)
std::string output_dir_
Definition: L1TGMTClient.h:45
TH2F * getTH2F(void) const
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
TH1F * get1DHisto(std::string meName, DQMStore::IGetter &igetter)
edm::ParameterSet parameters_
Definition: L1TGMTClient.h:42