CMS 3D CMS Logo

JetResolution_PayloadInspector.cc
Go to the documentation of this file.
2 
6 
7 // the data format of the condition to be inspected
11 
12 #include <memory>
13 #include <sstream>
14 #include <fstream>
15 #include <iostream>
16 
17 // include ROOT
18 #include "TH2F.h"
19 #include "TLegend.h"
20 #include "TCanvas.h"
21 #include "TLine.h"
22 #include "TStyle.h"
23 #include "TLatex.h"
24 #include "TPave.h"
25 #include "TPaveStats.h"
26 
27 #define MIN_ETA -5.05
28 #define MAX_ETA 5.05
29 #define NBIN_ETA 51
30 #define MIN_PT -5.
31 #define MAX_PT 3005.
32 #define NBIN_PT 301
33 
34 namespace JME {
35 
36  using namespace cond::payloadInspector;
37 
38  enum index { NORM = 0, DOWN = 1, UP = 2 };
39 
40  /*******************************************************
41  *
42  * 1d histogram of JetResolution of 1 IOV
43  *
44  *******************************************************/
45 
46  // inherit from one of the predefined plot class: Histogram1D
47 
48  class JetResolutionVsEta : public cond::payloadInspector::Histogram1D<JetResolutionObject, SINGLE_IOV> {
49  public:
51  : cond::payloadInspector::Histogram1D<JetResolutionObject, SINGLE_IOV>(
52  "Jet Resolution", "#eta", NBIN_ETA, MIN_ETA, MAX_ETA, "Resolution") {
55  }
56 
57  bool fill() override {
58  double par_Pt = 100.;
59  double par_Eta = 1.;
60  double par_Rho = 20.;
61 
62  // Default values will be used if no input parameters
64  auto ip = paramValues.find("Jet_Pt");
65  if (ip != paramValues.end()) {
66  par_Pt = std::stod(ip->second);
67  edm::LogPrint("JER_PI") << "Jet Pt: " << par_Pt;
68  }
69  ip = paramValues.find("Jet_Rho");
70  if (ip != paramValues.end()) {
71  par_Rho = std::stod(ip->second);
72  edm::LogPrint("JER_PI") << "Jet Rho: " << par_Rho;
73  }
74 
75  auto tag = PlotBase::getTag<0>();
76  for (auto const& iov : tag.iovs) {
77  std::shared_ptr<JetResolutionObject> payload = Base::fetchPayload(std::get<1>(iov));
78  if (payload.get()) {
79  if (!payload->getRecords().empty() && // No formula for SF
80  payload->getDefinition().getFormulaString().compare("") == 0)
81  return false;
82 
83  for (size_t idx = 0; idx <= NBIN_ETA; idx++) {
84  par_Eta = (idx + 0.5) * (MAX_ETA - MIN_ETA) / NBIN_ETA + MIN_ETA;
85 
86  JetParameters j_param;
87  j_param.setJetPt(par_Pt);
88  j_param.setRho(par_Rho);
89  j_param.setJetEta(par_Eta);
90 
91  if (payload->getRecord(j_param) == nullptr) {
92  continue;
93  }
94 
95  const JetResolutionObject::Record record = *(payload->getRecord(j_param));
96  float res = payload->evaluateFormula(record, j_param);
97  fillWithBinAndValue(idx, res);
98  } // x-axis
99  } else
100  return false;
101  }
102  return true;
103  }
104  }; // class
105 
106  class JetResolutionVsPt : public cond::payloadInspector::Histogram1D<JetResolutionObject, SINGLE_IOV> {
107  public:
109  : cond::payloadInspector::Histogram1D<JetResolutionObject, SINGLE_IOV>(
110  "Jet Energy Resolution", "p_T", NBIN_PT, MIN_PT, MAX_PT, "Resolution") {
113  }
114 
115  bool fill() override {
116  double par_Pt = 100.;
117  double par_Eta = 1.;
118  double par_Rho = 20.;
119 
121  auto ip = paramValues.find("Jet_Eta");
122  if (ip != paramValues.end()) {
123  par_Eta = std::stod(ip->second);
124  edm::LogPrint("JER_PI") << "Jet Eta: " << par_Eta;
125  }
126  ip = paramValues.find("Jet_Rho");
127  if (ip != paramValues.end()) {
128  par_Rho = std::stod(ip->second);
129  edm::LogPrint("JER_PI") << "Jet Rho: " << par_Rho;
130  }
131 
132  auto tag = PlotBase::getTag<0>();
133  for (auto const& iov : tag.iovs) {
134  std::shared_ptr<JetResolutionObject> payload = Base::fetchPayload(std::get<1>(iov));
135  if (payload.get()) {
136  if (!payload->getRecords().empty() && // No formula for SF
137  payload->getDefinition().getFormulaString().compare("") == 0)
138  return false;
139 
140  for (size_t idx = 0; idx <= NBIN_PT; idx++) {
141  par_Pt = (idx + 0.5) * (MAX_PT - MIN_PT) / NBIN_PT + MIN_PT;
142 
143  JetParameters j_param;
144  j_param.setJetEta(par_Eta);
145  j_param.setRho(par_Rho);
146  j_param.setJetPt(par_Pt);
147 
148  if (payload->getRecord(j_param) == nullptr) {
149  continue;
150  }
151 
152  const JetResolutionObject::Record record = *(payload->getRecord(j_param));
153  float res = payload->evaluateFormula(record, j_param);
154  fillWithBinAndValue(idx + 1, res);
155  } // x-axis
156  } else
157  return false;
158  }
159  return true;
160  } // fill
161  };
162 
163  class JetResolutionSummary : public cond::payloadInspector::PlotImage<JetResolutionObject, SINGLE_IOV> {
164  public:
166  : cond::payloadInspector::PlotImage<JetResolutionObject, SINGLE_IOV>("Jet Resolution Summary") {
170  }
171 
172  bool fill() override {
173  double par_Pt = 100.;
174  double par_Eta = 1.;
175  double par_Rho = 20.;
176 
178  auto ip = paramValues.find("Jet_Pt");
179  if (ip != paramValues.end()) {
180  par_Pt = std::stod(ip->second);
181  }
182  ip = paramValues.find("Jet_Eta");
183  if (ip != paramValues.end()) {
184  par_Eta = std::stod(ip->second);
185  }
186  ip = paramValues.find("Jet_Rho");
187  if (ip != paramValues.end()) {
188  par_Rho = std::stod(ip->second);
189  }
190 
191  TH1D* resol_eta = new TH1D("Jet Resolution vs #eta", "", NBIN_ETA, MIN_ETA, MAX_ETA);
192  TH1D* resol_pt = new TH1D("Jet Resolution vs p_T", "", NBIN_PT, MIN_PT, MAX_PT);
193  TLegend* leg_eta = new TLegend(0.26, 0.73, 0.935, 0.90);
194  TLegend* leg_pt = new TLegend(0.26, 0.73, 0.935, 0.90);
195 
196  leg_eta->SetBorderSize(0);
197  leg_eta->SetLineStyle(0);
198  leg_eta->SetFillStyle(0);
199 
200  leg_eta->SetTextFont(42);
201  leg_pt->SetBorderSize(0);
202  leg_pt->SetLineStyle(0);
203  leg_pt->SetFillStyle(0);
204 
205  auto tag = PlotBase::getTag<0>();
206  auto iov = tag.iovs.front();
207  std::shared_ptr<JetResolutionObject> payload = fetchPayload(std::get<1>(iov));
208  unsigned int run = std::get<0>(iov);
209  std::string tagname = tag.name;
210  std::stringstream ss_tagname(tag.name);
211  std::string stmp;
212 
213  std::string tag_ver;
214  std::string tag_res;
215  std::string tag_jet;
216 
217  getline(ss_tagname, stmp, '_'); // drop first
218  getline(ss_tagname, stmp, '_'); // year
219  tag_ver = stmp;
220  getline(ss_tagname, stmp, '_'); // ver
221  tag_ver += '_' + stmp;
222  getline(ss_tagname, stmp, '_'); // cmssw
223  tag_ver += '_' + stmp;
224  getline(ss_tagname, stmp, '_'); // data/mc
225  tag_ver += '_' + stmp;
226  getline(ss_tagname, stmp, '_'); // bin
227  tag_res = stmp;
228  getline(ss_tagname, stmp, '_'); // jet algorithm
229  tag_jet = stmp;
230 
231  if (payload.get()) {
232  if (!payload->getRecords().empty() && // No formula for SF
233  payload->getDefinition().getFormulaString().compare("") == 0)
234  return false;
235 
236  for (size_t idx = 0; idx <= NBIN_ETA; idx++) {
237  double x_axis = (idx + 0.5) * (MAX_ETA - MIN_ETA) / NBIN_ETA + MIN_ETA;
238 
239  JetParameters j_param;
240  j_param.setJetPt(par_Pt);
241  j_param.setRho(par_Rho);
242  j_param.setJetEta(x_axis);
243 
244  if (payload->getRecord(j_param) == nullptr) {
245  continue;
246  }
247 
248  const JetResolutionObject::Record record = *(payload->getRecord(j_param));
249  float res = payload->evaluateFormula(record, j_param);
250  resol_eta->SetBinContent(idx + 1, res);
251  } // x-axis
252 
253  for (size_t idx = 0; idx <= NBIN_PT; idx++) {
254  double x_axis = (idx + 0.5) * (MAX_PT - MIN_PT) / NBIN_PT + MIN_PT;
255 
256  JetParameters j_param;
257  j_param.setJetEta(par_Eta);
258  j_param.setRho(par_Rho);
259  j_param.setJetPt(x_axis);
260 
261  if (payload->getRecord(j_param) == nullptr) {
262  continue;
263  }
264 
265  const JetResolutionObject::Record record = *(payload->getRecord(j_param));
266  float res = payload->evaluateFormula(record, j_param);
267  resol_pt->SetBinContent(idx + 1, res);
268  } // x-axis
269 
270  gStyle->SetOptStat(0);
271  gStyle->SetLabelFont(42, "XYZ");
272  gStyle->SetLabelSize(0.05, "XYZ");
273  gStyle->SetFrameLineWidth(3);
274 
275  std::string title = Form("Summary Run %i", run);
276  TCanvas canvas("Jet Resolution Summary", title.c_str(), 800, 1200);
277  canvas.Divide(1, 2);
278 
279  canvas.cd(1);
280  resol_eta->SetTitle("Jet Resolution vs. #eta");
281  resol_eta->SetXTitle("#eta");
282  resol_eta->SetYTitle("Resolution");
283  resol_eta->SetLineWidth(3);
284  resol_eta->SetMaximum(resol_eta->GetMaximum() * 1.25);
285  resol_eta->Draw("");
286 
287  leg_eta->AddEntry(resol_eta, (tag_ver + '_' + tag_jet).c_str(), "l");
288  leg_eta->AddEntry((TObject*)nullptr, Form("JetPt=%.2f; JetRho=%.2f", par_Pt, par_Rho), "");
289  leg_eta->Draw();
290 
291  canvas.cd(2);
292  resol_pt->SetTitle("Jet Resolution vs. p_{T}");
293  resol_pt->SetXTitle("p_{T} [GeV]");
294  resol_pt->SetYTitle("Resolution");
295  resol_pt->SetLineWidth(3);
296  resol_pt->Draw("][");
297 
298  leg_pt->AddEntry(resol_pt, (tag_ver + '_' + tag_jet).c_str(), "l");
299  leg_pt->AddEntry((TObject*)nullptr, Form("JetEta=%.2f; JetRho=%.2f", par_Eta, par_Rho), "");
300  leg_pt->Draw();
301 
302  canvas.SaveAs(m_imageFileName.c_str());
303 
304  return true;
305  } else // no payload.get()
306  return false;
307  } // fill
308 
309  }; // class
310 
311  class JetResolutionCompare : public cond::payloadInspector::PlotImage<JetResolutionObject, SINGLE_IOV, 2> {
312  public:
314  : cond::payloadInspector::PlotImage<JetResolutionObject, SINGLE_IOV, 2>("Jet Resolution Compare") {
318  }
319 
320  bool fill() override {
321  double par_Pt = 100.;
322  double par_Eta = 1.;
323  double par_Rho = 20.;
324 
326  auto ip = paramValues.find("Jet_Pt");
327  if (ip != paramValues.end()) {
328  par_Pt = std::stod(ip->second);
329  }
330  ip = paramValues.find("Jet_Eta");
331  if (ip != paramValues.end()) {
332  par_Eta = std::stod(ip->second);
333  }
334  ip = paramValues.find("Jet_Rho");
335  if (ip != paramValues.end()) {
336  par_Rho = std::stod(ip->second);
337  }
338 
339  TH1D* resol_eta1 = new TH1D("Jet Resolution vs #eta one", "", NBIN_ETA, MIN_ETA, MAX_ETA);
340  TH1D* resol_eta2 = new TH1D("Jet Resolution vs #eta two", "", NBIN_ETA, MIN_ETA, MAX_ETA);
341  TH1D* resol_pt1 = new TH1D("Jet Resolution vs p_T one", "", NBIN_PT, MIN_PT, MAX_PT);
342  TH1D* resol_pt2 = new TH1D("Jet Resolution vs p_T two", "", NBIN_PT, MIN_PT, MAX_PT);
343  TLegend* leg_eta = new TLegend(0.26, 0.73, 0.935, 0.90);
344  TLegend* leg_pt = new TLegend(0.26, 0.73, 0.935, 0.90);
345 
346  leg_eta->SetBorderSize(0);
347  leg_eta->SetLineStyle(0);
348  leg_eta->SetFillStyle(0);
349 
350  leg_eta->SetTextFont(42);
351  leg_pt->SetBorderSize(0);
352  leg_pt->SetLineStyle(0);
353  leg_pt->SetFillStyle(0);
354 
355  auto tag1 = PlotBase::getTag<0>();
356  auto iov1 = tag1.iovs.front();
357  std::shared_ptr<JetResolutionObject> payload1 = fetchPayload(std::get<1>(iov1));
358  auto tag2 = PlotBase::getTag<1>();
359  auto iov2 = tag2.iovs.front();
360  std::shared_ptr<JetResolutionObject> payload2 = fetchPayload(std::get<1>(iov2));
361 
362  std::stringstream ss_tagname1(tag1.name);
363  std::stringstream ss_tagname2(tag2.name);
364  std::string stmp;
365 
366  std::string tag_ver1;
367  std::string tag_ver2;
368 
369  getline(ss_tagname1, stmp, '_'); // drop first
370  getline(ss_tagname1, stmp); // year
371  tag_ver1 = stmp;
372 
373  getline(ss_tagname2, stmp, '_'); // drop first
374  getline(ss_tagname2, stmp); // year
375  tag_ver2 = stmp;
376 
377  if (payload1.get() && payload2.get()) {
378  if (!payload1->getRecords().empty() && // No formula for SF
379  payload1->getDefinition().getFormulaString().compare("") == 0)
380  return false;
381 
382  for (size_t idx = 0; idx <= NBIN_ETA; idx++) {
383  double x_axis = (idx + 0.5) * (MAX_ETA - MIN_ETA) / NBIN_ETA + MIN_ETA;
384 
385  JetParameters j_param;
386  j_param.setJetPt(par_Pt);
387  j_param.setRho(par_Rho);
388  j_param.setJetEta(x_axis);
389 
390  if (payload1->getRecord(j_param) == nullptr || payload2->getRecord(j_param) == nullptr) {
391  continue;
392  }
393 
394  const JetResolutionObject::Record record1 = *(payload1->getRecord(j_param));
395  const JetResolutionObject::Record record2 = *(payload2->getRecord(j_param));
396  float res1 = payload1->evaluateFormula(record1, j_param);
397  resol_eta1->SetBinContent(idx + 1, res1);
398  float res2 = payload2->evaluateFormula(record2, j_param);
399  resol_eta2->SetBinContent(idx + 1, res2);
400  } // x-axis
401 
402  for (size_t idx = 0; idx <= NBIN_PT; idx++) {
403  double x_axis = (idx + 0.5) * (MAX_PT - MIN_PT) / NBIN_PT + MIN_PT;
404 
405  JetParameters j_param;
406  j_param.setJetEta(par_Eta);
407  j_param.setRho(par_Rho);
408  j_param.setJetPt(x_axis);
409 
410  if (payload1->getRecord(j_param) == nullptr || payload2->getRecord(j_param) == nullptr) {
411  continue;
412  }
413 
414  const JetResolutionObject::Record record1 = *(payload1->getRecord(j_param));
415  const JetResolutionObject::Record record2 = *(payload2->getRecord(j_param));
416  float res1 = payload1->evaluateFormula(record1, j_param);
417  resol_pt1->SetBinContent(idx + 1, res1);
418  float res2 = payload2->evaluateFormula(record2, j_param);
419  resol_pt2->SetBinContent(idx + 1, res2);
420  } // x-axis
421 
422  gStyle->SetOptStat(0);
423  gStyle->SetLabelFont(42, "XYZ");
424  gStyle->SetLabelSize(0.05, "XYZ");
425  gStyle->SetFrameLineWidth(3);
426 
427  std::string title = Form("Comparison between %s and %s", tag1.name.c_str(), tag2.name.c_str());
428  TCanvas canvas("Jet Resolution Comparison", title.c_str(), 800, 1200);
429  canvas.Divide(1, 2);
430 
431  canvas.cd(1);
432  resol_eta1->SetTitle("Jet Resolution Comparison vs. #eta");
433  resol_eta1->SetXTitle("#eta");
434  resol_eta1->SetYTitle("Resolution");
435  resol_eta1->SetLineWidth(3);
436  resol_eta1->SetMaximum(resol_eta1->GetMaximum() * 1.25);
437  resol_eta1->Draw("][");
438 
439  resol_eta2->SetLineColor(2);
440  resol_eta2->SetLineWidth(3);
441  resol_eta2->SetLineStyle(2);
442  resol_eta2->Draw("][same");
443 
444  leg_eta->AddEntry(resol_eta1, tag_ver1.c_str(), "l");
445  leg_eta->AddEntry(resol_eta2, tag_ver2.c_str(), "l");
446  leg_eta->AddEntry((TObject*)nullptr, Form("JetPt=%.2f; JetRho=%.2f", par_Pt, par_Rho), "");
447  leg_eta->Draw();
448 
449  canvas.cd(2);
450  resol_pt1->SetTitle("Jet Resolution Comparison vs. p_{T}");
451  resol_pt1->SetXTitle("p_{T} [GeV]");
452  resol_pt1->SetYTitle("Resolution");
453  resol_pt1->SetLineWidth(3);
454  resol_pt1->Draw("][");
455 
456  resol_pt2->SetLineColor(2);
457  resol_pt2->SetLineWidth(3);
458  resol_pt2->SetLineStyle(2);
459  resol_pt2->Draw("][same");
460 
461  leg_pt->AddEntry(resol_pt1, tag_ver1.c_str(), "l");
462  leg_pt->AddEntry(resol_pt2, tag_ver2.c_str(), "l");
463  leg_pt->AddEntry((TObject*)nullptr, Form("JetEta=%.2f; JetRho=%.2f", par_Eta, par_Rho), "");
464  leg_pt->Draw();
465 
466  canvas.SaveAs(m_imageFileName.c_str());
467 
468  return true;
469  } else // no payload.get()
470  return false;
471  } // fill
472 
473  }; // class
474 
475  template <index ii>
476  class JetScaleFactorVsEta : public cond::payloadInspector::Histogram1D<JetResolutionObject, SINGLE_IOV> {
477  public:
479  : cond::payloadInspector::Histogram1D<JetResolutionObject, SINGLE_IOV>(
480  "Jet Energy Scale Factor", "#eta", NBIN_ETA, MIN_ETA, MAX_ETA, "Scale Factor") {
483  }
484 
485  bool fill() override {
486  double par_Pt = 100.;
487  double par_Eta = 1.;
488 
490  auto ip = paramValues.find("Jet_Pt");
491  if (ip != paramValues.end()) {
492  par_Pt = std::stod(ip->second);
493  edm::LogPrint("JER_PI") << "Jet Pt: " << par_Pt;
494  }
495  ip = paramValues.find("Jet_Eta");
496  if (ip != paramValues.end()) {
497  par_Eta = std::stod(ip->second);
498  edm::LogPrint("JER_PI") << "Jet Eta: " << par_Eta;
499  }
500 
501  auto tag = PlotBase::getTag<0>();
502  for (auto const& iov : tag.iovs) {
503  std::shared_ptr<JetResolutionObject> payload = Base::fetchPayload(std::get<1>(iov));
504  if (payload.get()) {
505  if (!payload->getRecords().empty() && // No formula for SF
506  payload->getDefinition().getFormulaString().compare("") != 0)
507  return false;
508 
509  for (size_t idx = 0; idx <= NBIN_ETA; idx++) {
510  par_Eta = (idx + 0.5) * (MAX_ETA - MIN_ETA) / NBIN_ETA + MIN_ETA;
511 
512  JetParameters j_param;
513  j_param.setJetPt(par_Pt);
514  j_param.setJetEta(par_Eta);
515 
516  if (payload->getRecord(j_param) == nullptr) {
517  continue;
518  }
519 
520  const JetResolutionObject::Record record = *(payload->getRecord(j_param));
521  if (record.getParametersValues().size() == 3) { // norm, down, up
522  double sf = record.getParametersValues()[ii];
523  fillWithBinAndValue(idx + 1, sf);
524  }
525  } // x-axis
526  } else
527  return false;
528  } // for
529  return true;
530  } // fill
531  }; // class
532 
536 
537  template <index ii>
538  class JetScaleFactorVsPt : public cond::payloadInspector::Histogram1D<JetResolutionObject, SINGLE_IOV> {
539  public:
541  : cond::payloadInspector::Histogram1D<JetResolutionObject, SINGLE_IOV>(
542  "Jet Energy Scale Factor", "p_T", NBIN_PT, MIN_PT, MAX_PT, "Scale Factor") {
545  }
546 
547  bool fill() override {
548  double par_Pt = 100.;
549  double par_Eta = 1.;
550 
552  auto ip = paramValues.find("Jet_Pt");
553  if (ip != paramValues.end()) {
554  par_Pt = std::stod(ip->second);
555  edm::LogPrint("JER_PI") << "Jet Pt: " << par_Pt;
556  }
557  ip = paramValues.find("Jet_Eta");
558  if (ip != paramValues.end()) {
559  par_Eta = std::stod(ip->second);
560  edm::LogPrint("JER_PI") << "Jet Eta: " << par_Eta;
561  }
562 
563  auto tag = PlotBase::getTag<0>();
564  for (auto const& iov : tag.iovs) {
565  std::shared_ptr<JetResolutionObject> payload = Base::fetchPayload(std::get<1>(iov));
566  if (payload.get()) {
567  if (!payload->getRecords().empty() && // No formula for SF
568  payload->getDefinition().getFormulaString().compare("") != 0)
569  return false;
570 
571  for (size_t idx = 0; idx <= NBIN_PT; idx++) {
572  par_Pt = (idx + 0.5) * (MAX_PT - MIN_PT) / NBIN_PT + MIN_PT;
573 
574  JetParameters j_param;
575  j_param.setJetEta(par_Eta);
576  j_param.setJetPt(par_Pt);
577 
578  if (payload->getRecord(j_param) == nullptr) {
579  continue;
580  }
581 
582  const JetResolutionObject::Record record = *(payload->getRecord(j_param));
583  if (record.getParametersValues().size() == 3) { // norm, down, up
584  double sf = record.getParametersValues()[ii];
585  fillWithBinAndValue(idx + 1, sf);
586  }
587  } // x-axis
588  } else
589  return false;
590  } // for
591  return true;
592  } // fill
593  }; // class
594 
598 
599  class JetScaleFactorSummary : public cond::payloadInspector::PlotImage<JetResolutionObject, SINGLE_IOV> {
600  public:
602  : cond::payloadInspector::PlotImage<JetResolutionObject, SINGLE_IOV>("Jet ScaleFactor Summary") {
605  }
606 
607  bool fill() override {
608  double par_Pt = 100.;
609  double par_Eta = 1.;
610 
612  auto ip = paramValues.find("Jet_Pt");
613  if (ip != paramValues.end()) {
614  par_Pt = std::stod(ip->second);
615  }
616  ip = paramValues.find("Jet_Eta");
617  if (ip != paramValues.end()) {
618  par_Eta = std::stod(ip->second);
619  }
620 
621  TH1D* sf_eta_norm = new TH1D("Jet SF vs #eta NORM", "", NBIN_ETA, MIN_ETA, MAX_ETA);
622  TH1D* sf_eta_down = new TH1D("Jet SF vs #eta DOWN", "", NBIN_ETA, MIN_ETA, MAX_ETA);
623  TH1D* sf_eta_up = new TH1D("Jet SF vs #eta UP", "", NBIN_ETA, MIN_ETA, MAX_ETA);
624  TH1D* sf_pt_norm = new TH1D("Jet SF vs p_T NORM", "", NBIN_PT, MIN_PT, MAX_PT);
625  TH1D* sf_pt_down = new TH1D("Jet SF vs p_T DOWN", "", NBIN_PT, MIN_PT, MAX_PT);
626  TH1D* sf_pt_up = new TH1D("Jet SF vs p_T UP", "", NBIN_PT, MIN_PT, MAX_PT);
627 
628  TLegend* leg_eta = new TLegend(0.26, 0.73, 0.935, 0.90);
629  TLegend* leg_pt = new TLegend(0.26, 0.73, 0.935, 0.90);
630 
631  leg_eta->SetBorderSize(0);
632  leg_eta->SetLineStyle(0);
633  leg_eta->SetFillStyle(0);
634 
635  leg_eta->SetTextFont(42);
636  leg_pt->SetBorderSize(0);
637  leg_pt->SetLineStyle(0);
638  leg_pt->SetFillStyle(0);
639 
640  auto tag = PlotBase::getTag<0>();
641  auto iov = tag.iovs.front();
642  std::shared_ptr<JetResolutionObject> payload = fetchPayload(std::get<1>(iov));
643  unsigned int run = std::get<0>(iov);
644  std::string tagname = tag.name;
645  std::stringstream ss_tagname(tag.name);
646  std::string stmp;
647 
648  std::string tag_ver;
649  std::string tag_res;
650  std::string tag_jet;
651 
652  getline(ss_tagname, stmp, '_'); // drop first
653  getline(ss_tagname, stmp, '_'); // year
654  tag_ver = stmp;
655  getline(ss_tagname, stmp, '_'); // ver
656  tag_ver += '_' + stmp;
657  getline(ss_tagname, stmp, '_'); // cmssw
658  tag_ver += '_' + stmp;
659  getline(ss_tagname, stmp, '_'); // data/mc
660  tag_ver += '_' + stmp;
661  getline(ss_tagname, stmp, '_'); // bin
662  tag_res = stmp;
663  getline(ss_tagname, stmp, '_'); // jet algorithm
664  tag_jet = stmp;
665 
666  bool is_2bin = false;
667 
668  if (payload.get()) {
669  if (!payload->getRecords().empty() && // No formula for SF
670  payload->getDefinition().getFormulaString().compare("") != 0)
671  return false;
672 
673  is_2bin = false;
674  for (size_t idx = 0; idx <= NBIN_ETA; idx++) {
675  double x_axis = (idx + 0.5) * (MAX_ETA - MIN_ETA) / NBIN_ETA + MIN_ETA;
676 
677  JetParameters j_param;
678  j_param.setJetPt(par_Pt);
679  j_param.setJetEta(x_axis);
680 
681  if (payload->getRecord(j_param) == nullptr) {
682  continue;
683  }
684 
685  const JetResolutionObject::Record record = *(payload->getRecord(j_param));
686  if (record.getBinsRange().size() > 1)
687  is_2bin = true;
688 
689  if (record.getParametersValues().size() == 3) { // norm, down, up
690  sf_eta_norm->SetBinContent(idx + 1, record.getParametersValues()[0]);
691  sf_eta_down->SetBinContent(idx + 1, record.getParametersValues()[1]);
692  sf_eta_up->SetBinContent(idx + 1, record.getParametersValues()[2]);
693  }
694  } // x-axis
695 
696  for (size_t idx = 0; idx <= NBIN_PT; idx++) {
697  double x_axis = (idx + 0.5) * (MAX_PT - MIN_PT) / NBIN_PT + MIN_PT;
698 
699  JetParameters j_param;
700  j_param.setJetEta(par_Eta);
701  j_param.setJetPt(x_axis);
702 
703  if (payload->getRecord(j_param) == nullptr) {
704  continue;
705  }
706 
707  const JetResolutionObject::Record record = *(payload->getRecord(j_param));
708  if (record.getParametersValues().size() == 3) { // norm, down, up
709  sf_pt_norm->SetBinContent(idx + 1, record.getParametersValues()[0]);
710  sf_pt_down->SetBinContent(idx + 1, record.getParametersValues()[1]);
711  sf_pt_up->SetBinContent(idx + 1, record.getParametersValues()[2]);
712  }
713  } // x-axis
714 
715  gStyle->SetOptStat(0);
716  gStyle->SetLabelFont(42, "XYZ");
717  gStyle->SetLabelSize(0.05, "XYZ");
718  gStyle->SetFrameLineWidth(3);
719 
720  std::string title = Form("Summary Run %i", run);
721  TCanvas canvas("Jet ScaleFactor Summary", title.c_str(), 800, 1200);
722  canvas.Divide(1, 2);
723 
724  canvas.cd(1);
725  sf_eta_up->SetTitle("ScaleFactor vs. #eta");
726  sf_eta_up->SetXTitle("#eta");
727  sf_eta_up->SetYTitle("Scale Factor");
728  sf_eta_up->SetLineStyle(7);
729  sf_eta_up->SetLineWidth(3);
730  sf_eta_up->SetFillColorAlpha(kGray, 0.5);
731  sf_eta_up->SetMaximum(sf_eta_up->GetMaximum() * 1.25);
732  sf_eta_up->SetMinimum(0.);
733  sf_eta_up->Draw("][");
734 
735  sf_eta_down->SetLineStyle(7);
736  sf_eta_down->SetLineWidth(3);
737  sf_eta_down->SetFillColorAlpha(kWhite, 1);
738  sf_eta_down->Draw("][ same");
739 
740  sf_eta_norm->SetLineStyle(1);
741  sf_eta_norm->SetLineWidth(5);
742  sf_eta_norm->SetFillColor(0);
743  sf_eta_norm->Draw("][ same");
744  sf_eta_norm->Draw("axis same");
745 
746  leg_eta->AddEntry(sf_eta_norm, (tag_ver + '_' + tag_jet).c_str(), "l");
747  leg_eta->AddEntry((TObject*)nullptr, Form("JetPt=%.2f", par_Pt), "");
748  leg_eta->Draw();
749 
750  if (is_2bin == true) {
751  canvas.cd(2);
752  sf_pt_up->SetTitle("ScaleFactor vs. p_{T}");
753  sf_pt_up->SetXTitle("p_{T} [GeV]");
754  sf_pt_up->SetYTitle("Scale Factor");
755  sf_pt_up->SetLineStyle(7);
756  sf_pt_up->SetLineWidth(3);
757  sf_pt_up->SetFillColorAlpha(kGray, 0.5);
758  sf_pt_up->SetMaximum(sf_pt_up->GetMaximum() * 1.25);
759  sf_pt_up->SetMinimum(0.);
760  sf_pt_up->Draw("][");
761 
762  sf_pt_down->SetLineStyle(7);
763  sf_pt_down->SetLineWidth(3);
764  sf_pt_down->SetFillColorAlpha(kWhite, 1);
765  sf_pt_down->Draw("][ same");
766 
767  sf_pt_norm->SetLineStyle(1);
768  sf_pt_norm->SetLineWidth(5);
769  sf_pt_norm->SetFillColor(0);
770  sf_pt_norm->Draw("][ same");
771  sf_pt_norm->Draw("axis same");
772 
773  leg_pt->AddEntry(sf_pt_norm, (tag_ver + '_' + tag_jet).c_str(), "l");
774  leg_pt->AddEntry((TObject*)nullptr, Form("JetEta=%.2f", par_Eta), "");
775  leg_pt->Draw();
776  }
777 
778  canvas.SaveAs(m_imageFileName.c_str());
779 
780  return true;
781  } else // no payload.get()
782  return false;
783  } // fill
784 
785  }; // class
786 
787  class JetScaleFactorCompare : public cond::payloadInspector::PlotImage<JetResolutionObject, SINGLE_IOV, 2> {
788  public:
790  : cond::payloadInspector::PlotImage<JetResolutionObject, SINGLE_IOV, 2>("Jet ScaleFactor Summary") {
793  }
794 
795  bool fill() override {
796  double par_Pt = 100.;
797  double par_Eta = 1.;
798 
800  auto ip = paramValues.find("Jet_Pt");
801  if (ip != paramValues.end()) {
802  par_Pt = std::stod(ip->second);
803  }
804  ip = paramValues.find("Jet_Eta");
805  if (ip != paramValues.end()) {
806  par_Eta = std::stod(ip->second);
807  }
808 
809  TH1D* sf_eta_one = new TH1D("Jet SF vs #eta one", "", NBIN_ETA, MIN_ETA, MAX_ETA);
810  TH1D* sf_eta_two = new TH1D("Jet SF vs #eta two", "", NBIN_ETA, MIN_ETA, MAX_ETA);
811  TH1D* sf_pt_one = new TH1D("Jet SF vs p_T one", "", NBIN_PT, MIN_PT, MAX_PT);
812  TH1D* sf_pt_two = new TH1D("Jet SF vs p_T two", "", NBIN_PT, MIN_PT, MAX_PT);
813 
814  TLegend* leg_eta = new TLegend(0.26, 0.73, 0.935, 0.90);
815  TLegend* leg_pt = new TLegend(0.26, 0.73, 0.935, 0.90);
816 
817  leg_eta->SetBorderSize(0);
818  leg_eta->SetLineStyle(0);
819  leg_eta->SetFillStyle(0);
820 
821  leg_eta->SetTextFont(42);
822  leg_pt->SetBorderSize(0);
823  leg_pt->SetLineStyle(0);
824  leg_pt->SetFillStyle(0);
825 
826  auto tag1 = PlotBase::getTag<0>();
827  auto iov1 = tag1.iovs.front();
828  std::shared_ptr<JetResolutionObject> payload1 = fetchPayload(std::get<1>(iov1));
829 
830  auto tag2 = PlotBase::getTag<1>();
831  auto iov2 = tag2.iovs.front();
832  std::shared_ptr<JetResolutionObject> payload2 = fetchPayload(std::get<1>(iov2));
833 
834  std::stringstream ss_tagname1(tag1.name);
835  std::stringstream ss_tagname2(tag2.name);
836  std::string stmp;
837 
838  std::string tag_ver1;
839  std::string tag_ver2;
840 
841  getline(ss_tagname1, stmp, '_'); // drop first
842  getline(ss_tagname1, stmp); // year
843  tag_ver1 = stmp;
844 
845  getline(ss_tagname2, stmp, '_'); // drop first
846  getline(ss_tagname2, stmp); // year
847  tag_ver2 = stmp;
848 
849  bool is_2bin = false;
850 
851  if (payload1.get() && payload2.get()) {
852  if (!payload1->getRecords().empty() && // No formula for SF
853  payload1->getDefinition().getFormulaString().compare("") != 0)
854  return false;
855 
856  is_2bin = false;
857  for (size_t idx = 0; idx <= NBIN_ETA; idx++) {
858  double x_axis = (idx + 0.5) * (MAX_ETA - MIN_ETA) / NBIN_ETA + MIN_ETA;
859 
860  JetParameters j_param;
861  j_param.setJetPt(par_Pt);
862  j_param.setJetEta(x_axis);
863 
864  const JetResolutionObject::Record* rcd_p1 = payload1->getRecord(j_param);
865  const JetResolutionObject::Record* rcd_p2 = payload2->getRecord(j_param);
866  if (rcd_p1 == nullptr || rcd_p2 == nullptr) {
867  continue;
868  }
869 
870  const JetResolutionObject::Record record1 = *(rcd_p1);
871  const JetResolutionObject::Record record2 = *(rcd_p2);
872 
873  if (record1.getBinsRange().size() > 1 && record2.getBinsRange().size() > 1)
874  is_2bin = true;
875 
876  if (record1.getParametersValues().size() == 3) { // norm, down, up
877  sf_eta_one->SetBinContent(idx + 1, record1.getParametersValues()[0]);
878  }
879  if (record2.getParametersValues().size() == 3) { // norm, down, up
880  sf_eta_two->SetBinContent(idx + 1, record2.getParametersValues()[0]);
881  }
882  } // x-axis
883 
884  for (size_t idx = 0; idx <= NBIN_PT; idx++) {
885  double x_axis = (idx + 0.5) * (MAX_PT - MIN_PT) / NBIN_PT + MIN_PT;
886 
887  JetParameters j_param;
888  j_param.setJetEta(par_Eta);
889  j_param.setJetPt(x_axis);
890 
891  const JetResolutionObject::Record* rcd_p1 = payload1->getRecord(j_param);
892  const JetResolutionObject::Record* rcd_p2 = payload2->getRecord(j_param);
893  if (rcd_p1 == nullptr || rcd_p2 == nullptr) {
894  continue;
895  }
896 
897  const JetResolutionObject::Record record1 = *(rcd_p1);
898  const JetResolutionObject::Record record2 = *(rcd_p2);
899 
900  if (record1.getParametersValues().size() == 3) { // norm, down, up
901  sf_pt_one->SetBinContent(idx + 1, record1.getParametersValues()[0]);
902  }
903  if (record2.getParametersValues().size() == 3) { // norm, down, up
904  sf_pt_two->SetBinContent(idx + 1, record2.getParametersValues()[0]);
905  }
906  } // x-axis
907 
908  gStyle->SetOptStat(0);
909  gStyle->SetLabelFont(42, "XYZ");
910  gStyle->SetLabelSize(0.05, "XYZ");
911  gStyle->SetFrameLineWidth(3);
912 
913  std::string title = Form("Comparison");
914  TCanvas canvas("Jet ScaleFactor", title.c_str(), 800, 1200);
915  canvas.Divide(1, 2);
916 
917  canvas.cd(1);
918  sf_eta_one->SetTitle("Jet ScaleFactor Comparison vs. #eta");
919  sf_eta_one->SetXTitle("#eta");
920  sf_eta_one->SetYTitle("Scale Factor");
921  sf_eta_one->SetLineWidth(3);
922  sf_eta_one->SetMaximum(sf_eta_one->GetMaximum() * 1.25);
923  sf_eta_one->SetMinimum(0.);
924  sf_eta_one->Draw("][");
925 
926  sf_eta_two->SetLineColor(2);
927  sf_eta_two->SetLineWidth(3);
928  sf_eta_two->SetLineStyle(2);
929  sf_eta_two->Draw("][ same");
930 
931  leg_eta->AddEntry(sf_eta_one, tag_ver1.c_str(), "l");
932  leg_eta->AddEntry(sf_eta_two, tag_ver2.c_str(), "l");
933  leg_eta->AddEntry((TObject*)nullptr, Form("JetPt=%.2f", par_Pt), "");
934  leg_eta->Draw();
935 
936  if (is_2bin == true) {
937  canvas.cd(2);
938  sf_pt_one->SetTitle("Jet ScaleFactor Comparison vs. p_{T}");
939  sf_pt_one->SetXTitle("p_{T} [GeV]");
940  sf_pt_one->SetYTitle("Scale Factor");
941  sf_pt_one->SetLineWidth(3);
942  sf_pt_one->SetMaximum(sf_pt_one->GetMaximum() * 1.25);
943  sf_pt_one->SetMinimum(0.);
944  sf_pt_one->Draw("][");
945 
946  sf_pt_two->SetLineColor(2);
947  sf_pt_two->SetLineWidth(3);
948  sf_pt_two->SetLineStyle(2);
949  sf_pt_two->Draw("][ same");
950 
951  leg_pt->AddEntry(sf_pt_one, tag_ver1.c_str(), "l");
952  leg_pt->AddEntry(sf_pt_two, tag_ver2.c_str(), "l");
953  leg_pt->AddEntry((TObject*)nullptr, Form("JetEta=%.2f", par_Eta), "");
954  leg_pt->Draw();
955  }
956 
957  canvas.SaveAs(m_imageFileName.c_str());
958 
959  return true;
960  } else // no payload.get()
961  return false;
962  } // fill
963 
964  }; // class
965 
966  // Register the classes as boost python plugin
976 
981  }
982 
983 } // namespace JME
const std::vector< Range > & getBinsRange() const
JetParameters & setJetEta(float eta)
JetParameters & setRho(float rho)
JetScaleFactorVsEta< NORM > JetScaleFactorVsEtaNORM
JetScaleFactorVsEta< UP > JetScaleFactorVsEtaUP
JetScaleFactorVsPt< NORM > JetScaleFactorVsPtNORM
Definition: Electron.h:6
#define PAYLOAD_INSPECTOR_CLASS(CLASS_NAME)
JetScaleFactorVsPt< DOWN > JetScaleFactorVsPtDOWN
#define PAYLOAD_INSPECTOR_MODULE(PAYLOAD_TYPENAME)
Log< level::Warning, true > LogPrint
const std::vector< float > & getParametersValues() const
ii
Definition: cuy.py:589
const std::map< std::string, std::string > & inputParamValues() const
void addInputParam(const std::string &paramName)
JetScaleFactorVsPt< UP > JetScaleFactorVsPtUP
def canvas(sub, attr)
Definition: svgfig.py:482
JetParameters & setJetPt(float pt)
JetScaleFactorVsEta< DOWN > JetScaleFactorVsEtaDOWN