CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
EcalTimeOffsetConstant_PayloadInspector.cc
Go to the documentation of this file.
7 
8 // the data format of the condition to be inspected
10 
11 #include "TH2F.h" // a 2-D histogram with four bytes per cell (float)
12 #include "TCanvas.h"
13 #include "TLine.h"
14 #include "TStyle.h"
15 #include "TLatex.h" //write mathematical equations.
16 #include "TPave.h"
17 #include "TPaveStats.h"
18 #include <string>
19 #include <fstream>
20 
21 namespace {
22 
23  /*******************************************************
24  2d plot of Ecal Time Offset Constant of 1 IOV
25  *******************************************************/
26  class EcalTimeOffsetConstantPlot : public cond::payloadInspector::PlotImage<EcalTimeOffsetConstant> {
27  public:
28  EcalTimeOffsetConstantPlot()
29  : cond::payloadInspector::PlotImage<EcalTimeOffsetConstant>("Ecal Time Offset Constant - map ") {
30  setSingleIov(true);
31  }
32 
33  bool fill(const std::vector<std::tuple<cond::Time_t, cond::Hash> >& iovs) override {
34  auto iov = iovs.front();
35  std::shared_ptr<EcalTimeOffsetConstant> payload = fetchPayload(std::get<1>(iov));
36  unsigned int run = std::get<0>(iov);
37  TH2F* align;
38  int NbRows;
39 
40  if (payload.get()) {
41  NbRows = 1;
42  align = new TH2F("Time Offset Constant [ns]", "EB EE", 2, 0, 2, NbRows, 0, NbRows);
43  EcalTimeOffsetConstant it = (*payload);
44 
45  double row = NbRows - 0.5;
46 
47  align->Fill(0.5, row, it.getEBValue());
48  align->Fill(1.5, row, it.getEEValue());
49  } else
50  return false;
51 
52  gStyle->SetPalette(1);
53  gStyle->SetOptStat(0);
54  TCanvas canvas("CC map", "CC map", 1000, 1000);
55  TLatex t1;
56  t1.SetNDC();
57  t1.SetTextAlign(26);
58  t1.SetTextSize(0.05);
59  t1.SetTextColor(2);
60  t1.DrawLatex(0.5, 0.96, Form("Ecal Time Offset Constant, IOV %i", run));
61 
62  TPad* pad = new TPad("pad", "pad", 0.0, 0.0, 1.0, 0.94);
63  pad->Draw();
64  pad->cd();
65  align->Draw("TEXT");
66 
67  drawTable(NbRows, 2);
68 
69  align->GetXaxis()->SetTickLength(0.);
70  align->GetXaxis()->SetLabelSize(0.);
71  align->GetYaxis()->SetTickLength(0.);
72  align->GetYaxis()->SetLabelSize(0.);
73 
74  std::string ImageName(m_imageFileName);
75  canvas.SaveAs(ImageName.c_str());
76 
77  return true;
78  }
79  };
80 
81  /*****************************************************************
82  2d plot of Ecal Time Offset Constant difference between 2 IOVs
83  *****************************************************************/
84  template <cond::payloadInspector::IOVMultiplicity nIOVs, int ntags>
85  class EcalTimeOffsetConstantDiffBase
86  : public cond::payloadInspector::PlotImage<EcalTimeOffsetConstant, nIOVs, ntags> {
87  public:
88  EcalTimeOffsetConstantDiffBase()
89  : cond::payloadInspector::PlotImage<EcalTimeOffsetConstant, nIOVs, ntags>(
90  "Ecal Time Offset Constant difference") {}
91 
92  bool fill() override {
93  unsigned int run[2], NbRows = 0;
94  float val[2] = {};
95  TH2F* align = new TH2F("", "", 1, 0., 1., 1, 0., 1.); // pseudo creation
96  std::string l_tagname[2];
97  auto iovs = cond::payloadInspector::PlotBase::getTag<0>().iovs;
98  l_tagname[0] = cond::payloadInspector::PlotBase::getTag<0>().name;
99  auto firstiov = iovs.front();
100  run[0] = std::get<0>(firstiov);
101  std::tuple<cond::Time_t, cond::Hash> lastiov;
102  if (ntags == 2) {
103  auto tag2iovs = cond::payloadInspector::PlotBase::getTag<1>().iovs;
104  l_tagname[1] = cond::payloadInspector::PlotBase::getTag<1>().name;
105  lastiov = tag2iovs.front();
106  } else {
107  lastiov = iovs.back();
108  l_tagname[1] = l_tagname[0];
109  }
110  run[1] = std::get<0>(lastiov);
111  for (int irun = 0; irun < nIOVs; irun++) {
112  std::shared_ptr<EcalTimeOffsetConstant> payload;
113  if (irun == 0) {
114  payload = this->fetchPayload(std::get<1>(firstiov));
115  } else {
116  payload = this->fetchPayload(std::get<1>(lastiov));
117  }
118  if (payload.get()) {
119  NbRows = 1;
120 
121  if (irun == 1)
122  align = new TH2F("Ecal Time Offset Constant [ns]", "EB EE", 2, 0, 2, NbRows, 0, NbRows);
123 
124  EcalTimeOffsetConstant it = (*payload);
125 
126  if (irun == 0) {
127  val[0] = it.getEBValue();
128  val[1] = it.getEEValue();
129 
130  } else {
131  double row = NbRows - 0.5;
132  align->Fill(0.5, row, it.getEBValue() - val[0]);
133  align->Fill(1.5, row, it.getEEValue() - val[1]);
134  }
135 
136  } // if payload.get()
137  else
138  return false;
139  } // loop over IOVs
140 
141  gStyle->SetPalette(1);
142  gStyle->SetOptStat(0);
143  TCanvas canvas("CC map", "CC map", 1000, 1000);
144  TLatex t1;
145  t1.SetNDC();
146  t1.SetTextAlign(26);
147  t1.SetTextColor(2);
148  int len = l_tagname[0].length() + l_tagname[1].length();
149  if (ntags == 2) {
150  if (len < 80) {
151  t1.SetTextSize(0.02);
152  t1.DrawLatex(0.5, 0.96, Form("%s %i - %s %i", l_tagname[1].c_str(), run[1], l_tagname[0].c_str(), run[0]));
153  } else {
154  t1.SetTextSize(0.03);
155  t1.DrawLatex(0.5, 0.96, Form("Ecal Time Offset Constant, IOV %i - %i", run[1], run[0]));
156  }
157  } else {
158  t1.SetTextSize(0.03);
159  t1.DrawLatex(0.5, 0.96, Form("%s, IOV %i - %i", l_tagname[0].c_str(), run[1], run[0]));
160  }
161 
162  TPad* pad = new TPad("pad", "pad", 0.0, 0.0, 1.0, 0.94);
163  pad->Draw();
164  pad->cd();
165  align->Draw("TEXT");
166 
167  drawTable(NbRows, 2);
168 
169  align->GetXaxis()->SetTickLength(0.);
170  align->GetXaxis()->SetLabelSize(0.);
171  align->GetYaxis()->SetTickLength(0.);
172  align->GetYaxis()->SetLabelSize(0.);
173 
174  std::string ImageName(this->m_imageFileName);
175  canvas.SaveAs(ImageName.c_str());
176  return true;
177  }
178  }; // class EcalTimeOffsetConstantDiffBase
179  using EcalTimeOffsetConstantDiffOneTag = EcalTimeOffsetConstantDiffBase<cond::payloadInspector::SINGLE_IOV, 1>;
180  using EcalTimeOffsetConstantDiffTwoTags = EcalTimeOffsetConstantDiffBase<cond::payloadInspector::SINGLE_IOV, 2>;
181 
182 } // namespace
183 
184 // Register the classes as boost python plugin
186  PAYLOAD_INSPECTOR_CLASS(EcalTimeOffsetConstantPlot);
187  PAYLOAD_INSPECTOR_CLASS(EcalTimeOffsetConstantDiffOneTag);
188  PAYLOAD_INSPECTOR_CLASS(EcalTimeOffsetConstantDiffTwoTags);
189 }
def canvas
Definition: svgfig.py:482
#define PAYLOAD_INSPECTOR_CLASS(CLASS_NAME)
#define PAYLOAD_INSPECTOR_MODULE(PAYLOAD_TYPENAME)
void drawTable(int nbRows, int nbColumns)
Definition: EcalDrawUtils.h:91
std::shared_ptr< PayloadType > fetchPayload(const cond::Hash &payloadHash)