CMS 3D CMS Logo

EcalIntercalibConstants_PayloadInspector.cc
Go to the documentation of this file.
8 // the data format of the condition to be inspected
10 
11 #include "TH2F.h"
12 #include "TCanvas.h"
13 #include "TStyle.h"
14 #include "TLine.h"
15 #include "TLatex.h"
16 
17 #include <memory>
18 #include <sstream>
19 #include <array>
20 
21 namespace {
22  enum { kEBChannels = 61200, kEEChannels = 14648 };
23  enum {
24  MIN_IETA = 1,
25  MIN_IPHI = 1,
26  MAX_IETA = 85,
27  MAX_IPHI = 360,
28  EBhistEtaMax = 171
29  }; // barrel lower and upper bounds on eta and phi
30  enum {
31  IX_MIN = 1,
32  IY_MIN = 1,
33  IX_MAX = 100,
34  IY_MAX = 100,
35  EEhistXMax = 220
36  }; // endcaps lower and upper bounds on x and y
37 
38  /**************************************************************
39  2d histogram of ECAL barrel Intercalib Constants of 1 IOV
40  ***************************************************************/
41 
42  // inherit from one of the predefined plot class: Histogram2D
43  class EcalIntercalibConstantsEBMap : public cond::payloadInspector::Histogram2D<EcalIntercalibConstants> {
44  public:
45  EcalIntercalibConstantsEBMap()
46  : cond::payloadInspector::Histogram2D<EcalIntercalibConstants>("ECAL Barrel Intercalib Constants - map ",
47  "iphi",
48  MAX_IPHI,
49  MIN_IPHI,
50  MAX_IPHI + 1,
51  "ieta",
53  -MAX_IETA,
54  MAX_IETA + 1) {
55  Base::setSingleIov(true);
56  }
57 
58  // Histogram2D::fill (virtual) needs be overridden - the implementation should use fillWithValue
59  bool fill() override {
60  auto tag = PlotBase::getTag<0>();
61  for (auto const& iov : tag.iovs) {
62  std::shared_ptr<EcalIntercalibConstants> payload = Base::fetchPayload(std::get<1>(iov));
63  if (payload.get()) {
64  // looping over the EB channels, via the dense-index, mapped into EBDetId's
65  if (payload->barrelItems().empty())
66  return false;
67  // set to 1 for ieta 0 (no crystal)
68  for (int iphi = MIN_IPHI; iphi < MAX_IPHI + 1; iphi++)
69  fillWithValue(iphi, 0, 1);
70 
71  for (int cellid = EBDetId::MIN_HASH; cellid < EBDetId::kSizeForDenseIndexing; ++cellid) {
72  uint32_t rawid = EBDetId::unhashIndex(cellid);
73 
74  // check the existence of ECAL Intercalib Constants, for a given ECAL barrel channel
76  if (value_ptr == payload->end())
77  continue; // cell absent from payload
78  float weight = (float)(*value_ptr);
79 
80  // fill the Histogram2D here
81  fillWithValue((EBDetId(rawid)).iphi(), (EBDetId(rawid)).ieta(), weight);
82  } // loop over cellid
83  } // if payload.get()
84  } // loop over IOV's (1 in this case)
85 
86  return true;
87  } // fill method
88  };
89 
90  class EcalIntercalibConstantsEEMap : public cond::payloadInspector::Histogram2D<EcalIntercalibConstants> {
91  private:
92  int EEhistSplit = 20;
93 
94  public:
95  EcalIntercalibConstantsEEMap()
96  : cond::payloadInspector::Histogram2D<EcalIntercalibConstants>("ECAL Endcap Intercalib Constants - map ",
97  "ix",
98  EEhistXMax,
99  IX_MIN,
100  EEhistXMax + 1,
101  "iy",
102  IY_MAX,
103  IY_MIN,
104  IY_MAX + 1) {
105  Base::setSingleIov(true);
106  }
107 
108  bool fill() override {
109  auto tag = PlotBase::getTag<0>();
110  for (auto const& iov : tag.iovs) {
111  std::shared_ptr<EcalIntercalibConstants> payload = Base::fetchPayload(std::get<1>(iov));
112  if (payload.get()) {
113  if (payload->endcapItems().empty())
114  return false;
115 
116  // set to 0 everywhwere
117  for (int ix = IX_MIN; ix < EEhistXMax + 1; ix++)
118  for (int iy = IY_MIN; iy < IY_MAX + 1; iy++)
119  fillWithValue(ix, iy, 0);
120 
121  for (int cellid = 0; cellid < EEDetId::kSizeForDenseIndexing; ++cellid) { // loop on EE cells
122  if (EEDetId::validHashIndex(cellid)) {
123  uint32_t rawid = EEDetId::unhashIndex(cellid);
124  EcalFloatCondObjectContainer::const_iterator value_ptr = payload->find(rawid);
125  if (value_ptr == payload->end())
126  continue; // cell absent from payload
127  float weight = (float)(*value_ptr);
128  EEDetId myEEId(rawid);
129  if (myEEId.zside() == -1)
130  fillWithValue(myEEId.ix(), myEEId.iy(), weight);
131  else
132  fillWithValue(myEEId.ix() + IX_MAX + EEhistSplit, myEEId.iy(), weight);
133  } // validDetId
134  } // loop over cellid
135  } // payload
136  } // loop over IOV's (1 in this case)
137  return true;
138  } // fill method
139  };
140 
141  /*************************************************
142  2d plot of ECAL IntercalibConstants of 1 IOV
143  *************************************************/
144  class EcalIntercalibConstantsPlot : public cond::payloadInspector::PlotImage<EcalIntercalibConstants> {
145  public:
146  EcalIntercalibConstantsPlot()
147  : cond::payloadInspector::PlotImage<EcalIntercalibConstants>("ECAL Intercalib Constants - map ") {
148  setSingleIov(true);
149  }
150 
151  bool fill(const std::vector<std::tuple<cond::Time_t, cond::Hash> >& iovs) override {
152  TH2F* barrel = new TH2F("EB", "mean EB", MAX_IPHI, 0, MAX_IPHI, 2 * MAX_IETA, -MAX_IETA, MAX_IETA);
153  TH2F* endc_p = new TH2F("EE+", "mean EE+", IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
154  TH2F* endc_m = new TH2F("EE-", "mean EE-", IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
155 
156  auto iov = iovs.front();
157  std::shared_ptr<EcalIntercalibConstants> payload = fetchPayload(std::get<1>(iov));
158  unsigned int run = std::get<0>(iov);
159 
160  if (payload.get()) {
161  if (payload->barrelItems().empty())
162  return false;
163 
164  fillEBMap_SingleIOV<EcalIntercalibConstants>(payload, barrel);
165 
166  if (payload->endcapItems().empty())
167  return false;
168 
169  fillEEMap_SingleIOV<EcalIntercalibConstants>(payload, endc_m, endc_p);
170 
171  } // payload
172 
173  gStyle->SetPalette(1);
174  gStyle->SetOptStat(0);
175  TCanvas canvas("CC map", "CC map", 1600, 450);
176  TLatex t1;
177  t1.SetNDC();
178  t1.SetTextAlign(26);
179  t1.SetTextSize(0.05);
180  t1.DrawLatex(0.5, 0.96, Form("Ecal IntercalibConstants, IOV %i", run));
181 
182  float xmi[3] = {0.0, 0.24, 0.76};
183  float xma[3] = {0.24, 0.76, 1.00};
184  std::array<std::unique_ptr<TPad>, 3> pad;
185  for (int obj = 0; obj < 3; obj++) {
186  pad[obj] = std::make_unique<TPad>(Form("p_%i", obj), Form("p_%i", obj), xmi[obj], 0.0, xma[obj], 0.94);
187  pad[obj]->Draw();
188  }
189  // EcalDrawMaps ICMap;
190  pad[0]->cd();
191  // ICMap.DrawEE(endc_m, 0., 2.);
192  DrawEE(endc_m, 0., 2.5);
193  pad[1]->cd();
194  // ICMap.DrawEB(barrel, 0., 2.);
195  DrawEB(barrel, 0., 2.5);
196  pad[2]->cd();
197  // ICMap.DrawEE(endc_p, 0., 2.);
198  DrawEE(endc_p, 0., 2.5);
199 
200  std::string ImageName(m_imageFileName);
201  canvas.SaveAs(ImageName.c_str());
202  return true;
203  } // fill method
204  };
205 
206  /**************************************************************************
207  2d plot of ECAL IntercalibConstants difference or ratio between 2 IOVs
208  ***************************************************************************/
209  template <cond::payloadInspector::IOVMultiplicity nIOVs, int ntags, int method>
210  class EcalIntercalibConstantsBase : public cond::payloadInspector::PlotImage<EcalIntercalibConstants, nIOVs, ntags> {
211  public:
212  EcalIntercalibConstantsBase()
213  : cond::payloadInspector::PlotImage<EcalIntercalibConstants, nIOVs, ntags>(
214  "ECAL Intercalib Constants comparison") {}
215 
216  bool fill() override {
217  TH2F* barrel = new TH2F("EB", "mean EB", MAX_IPHI, 0, MAX_IPHI, 2 * MAX_IETA, -MAX_IETA, MAX_IETA);
218  TH2F* endc_p = new TH2F("EE+", "mean EE+", IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
219  TH2F* endc_m = new TH2F("EE-", "mean EE-", IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
220  float pEBmin, pEEmin, pEBmax, pEEmax;
221  pEBmin = 10.;
222  pEEmin = 10.;
223  pEBmax = -10.;
224  pEEmax = -10.;
225 
226  unsigned int run[2];
227  float pEB[kEBChannels], pEE[kEEChannels];
228  std::string l_tagname[2];
229  auto iovs = cond::payloadInspector::PlotBase::getTag<0>().iovs;
230  l_tagname[0] = cond::payloadInspector::PlotBase::getTag<0>().name;
231  auto firstiov = iovs.front();
232  run[0] = std::get<0>(firstiov);
233  std::tuple<cond::Time_t, cond::Hash> lastiov;
234  if (ntags == 2) {
235  auto tag2iovs = cond::payloadInspector::PlotBase::getTag<1>().iovs;
236  l_tagname[1] = cond::payloadInspector::PlotBase::getTag<1>().name;
237  lastiov = tag2iovs.front();
238  } else {
239  lastiov = iovs.back();
240  l_tagname[1] = l_tagname[0];
241  }
242  run[1] = std::get<0>(lastiov);
243  for (int irun = 0; irun < nIOVs; irun++) {
244  std::shared_ptr<EcalIntercalibConstants> payload;
245  if (irun == 0) {
246  payload = this->fetchPayload(std::get<1>(firstiov));
247  } else {
248  payload = this->fetchPayload(std::get<1>(lastiov));
249  }
250  if (payload.get()) {
251  if (payload->barrelItems().empty())
252  return false;
253  fillEBMap_TwoIOVs<EcalIntercalibConstants>(payload, barrel, irun, pEB, pEBmin, pEBmax, method);
254 
255  if (payload->endcapItems().empty())
256  return false;
257  fillEEMap_TwoIOVs<EcalIntercalibConstants>(payload, endc_m, endc_p, irun, pEE, pEEmin, pEEmax, method);
258  } // payload
259  } // loop over IOVs
260 
261  gStyle->SetPalette(1);
262  gStyle->SetOptStat(0);
263  TCanvas canvas("CC map", "CC map", 1600, 450);
264  TLatex t1;
265  t1.SetNDC();
266  t1.SetTextAlign(26);
267  int len = l_tagname[0].length() + l_tagname[1].length();
268  std::string dr[2] = {"-", "/"};
269  if (ntags == 2) {
270  if (len < 170) {
271  t1.SetTextSize(0.05);
272  t1.DrawLatex(0.5,
273  0.96,
274  Form("%s IOV %i %s %s IOV %i",
275  l_tagname[1].c_str(),
276  run[1],
277  dr[method].c_str(),
278  l_tagname[0].c_str(),
279  run[0]));
280  } else {
281  t1.SetTextSize(0.05);
282  t1.DrawLatex(0.5, 0.96, Form("Ecal IntercalibConstants, IOV %i %s %i", run[1], dr[method].c_str(), run[0]));
283  }
284  } else {
285  t1.SetTextSize(0.05);
286  t1.DrawLatex(0.5, 0.96, Form("%s, IOV %i %s %i", l_tagname[0].c_str(), run[1], dr[method].c_str(), run[0]));
287  }
288  float xmi[3] = {0.0, 0.24, 0.76};
289  float xma[3] = {0.24, 0.76, 1.00};
290  std::array<std::unique_ptr<TPad>, 3> pad;
291  for (int obj = 0; obj < 3; obj++) {
292  pad[obj] = std::make_unique<TPad>(Form("p_%i", obj), Form("p_%i", obj), xmi[obj], 0.0, xma[obj], 0.94);
293  pad[obj]->Draw();
294  }
295  pad[0]->cd();
296  DrawEE(endc_m, pEEmin, pEEmax);
297  pad[1]->cd();
298  DrawEB(barrel, pEBmin, pEBmax);
299  pad[2]->cd();
300  DrawEE(endc_p, pEEmin, pEEmax);
301 
302  std::string ImageName(this->m_imageFileName);
303  canvas.SaveAs(ImageName.c_str());
304  return true;
305  } // fill method
306  }; // class EcalIntercalibConstantsDiffBase
307  using EcalIntercalibConstantsDiffOneTag = EcalIntercalibConstantsBase<cond::payloadInspector::SINGLE_IOV, 1, 0>;
308  using EcalIntercalibConstantsDiffTwoTags = EcalIntercalibConstantsBase<cond::payloadInspector::SINGLE_IOV, 2, 0>;
309  using EcalIntercalibConstantsRatioOneTag = EcalIntercalibConstantsBase<cond::payloadInspector::SINGLE_IOV, 1, 1>;
310  using EcalIntercalibConstantsRatioTwoTags = EcalIntercalibConstantsBase<cond::payloadInspector::SINGLE_IOV, 2, 1>;
311 
312  /*********************************************************
313  2d plot of Ecal Intercalib Constants Summary of 1 IOV
314  *********************************************************/
315  class EcalIntercalibConstantsSummaryPlot : public cond::payloadInspector::PlotImage<EcalIntercalibConstants> {
316  public:
317  EcalIntercalibConstantsSummaryPlot()
318  : cond::payloadInspector::PlotImage<EcalIntercalibConstants>("Ecal Intercalib Constants Summary - map ") {
319  setSingleIov(true);
320  }
321 
322  bool fill(const std::vector<std::tuple<cond::Time_t, cond::Hash> >& iovs) override {
323  auto iov = iovs.front();
324  std::shared_ptr<EcalIntercalibConstants> payload = fetchPayload(std::get<1>(iov));
325  unsigned int run = std::get<0>(iov);
326  TH2F* align;
327  int NbRows;
328 
329  if (payload.get()) {
330  NbRows = 2;
331  align = new TH2F("", "", 0, 0, 0, 0, 0, 0);
332 
333  float mean_x_EB = 0.0f;
334  float mean_x_EE = 0.0f;
335 
336  float rms_EB = 0.0f;
337  float rms_EE = 0.0f;
338 
339  int num_x_EB = 0;
340  int num_x_EE = 0;
341 
342  payload->summary(mean_x_EB, rms_EB, num_x_EB, mean_x_EE, rms_EE, num_x_EE);
344  align, "Ecal Intercalib Constants", mean_x_EB, rms_EB, num_x_EB, mean_x_EE, rms_EE, num_x_EE);
345  } else
346  return false;
347 
348  gStyle->SetPalette(1);
349  gStyle->SetOptStat(0);
350  TCanvas canvas("CC map", "CC map", 1000, 1000);
351  TLatex t1;
352  t1.SetNDC();
353  t1.SetTextAlign(26);
354  t1.SetTextSize(0.04);
355  t1.SetTextColor(2);
356  t1.DrawLatex(0.5, 0.96, Form("Ecal Intercalib Constants Summary, IOV %i", run));
357 
358  TPad pad("pad", "pad", 0.0, 0.0, 1.0, 0.94);
359  pad.Draw();
360  pad.cd();
361  align->Draw("TEXT");
362  TLine* l = new TLine;
363  l->SetLineWidth(1);
364 
365  drawTable(NbRows, 4);
366 
367  std::string ImageName(m_imageFileName);
368  canvas.SaveAs(ImageName.c_str());
369 
370  return true;
371  }
372  };
373 
374 } // namespace
375 
376 // Register the classes as boost python plugin
378  PAYLOAD_INSPECTOR_CLASS(EcalIntercalibConstantsEBMap);
379  PAYLOAD_INSPECTOR_CLASS(EcalIntercalibConstantsEEMap);
380  PAYLOAD_INSPECTOR_CLASS(EcalIntercalibConstantsPlot);
381  PAYLOAD_INSPECTOR_CLASS(EcalIntercalibConstantsDiffOneTag);
382  PAYLOAD_INSPECTOR_CLASS(EcalIntercalibConstantsDiffTwoTags);
383  PAYLOAD_INSPECTOR_CLASS(EcalIntercalibConstantsRatioOneTag);
384  PAYLOAD_INSPECTOR_CLASS(EcalIntercalibConstantsRatioTwoTags);
385  PAYLOAD_INSPECTOR_CLASS(EcalIntercalibConstantsSummaryPlot);
386 }
void fillTableWithSummary(TH2F *&align, std::string title, const float &mean_x_EB, const float &rms_EB, const int &num_x_EB, const float &mean_x_EE, const float &rms_EE, const int &num_x_EE)
void fillWithValue(float xvalue, float yvalue, float weight=1)
Definition: weight.py:1
static EEDetId unhashIndex(int hi)
Definition: EEDetId.cc:65
static bool validHashIndex(int i)
Definition: EEDetId.h:239
void DrawEE(TH2F *endc, float min, float max)
Definition: EcalDrawUtils.h:29
void DrawEB(TH2F *ebmap, float min, float max)
Definition: EcalDrawUtils.h:4
#define PAYLOAD_INSPECTOR_CLASS(CLASS_NAME)
static const int MIN_HASH
Definition: EBDetId.h:149
#define PAYLOAD_INSPECTOR_MODULE(PAYLOAD_TYPENAME)
ALPAKA_FN_ACC ALPAKA_FN_INLINE uint32_t ix(uint32_t id)
std::vector< Item >::const_iterator const_iterator
static EBDetId unhashIndex(int hi)
get a DetId from a compact index for arrays
Definition: EBDetId.h:110
def canvas(sub, attr)
Definition: svgfig.py:482
std::shared_ptr< PayloadType > fetchPayload(const cond::Hash &payloadHash)
ALPAKA_FN_ACC ALPAKA_FN_INLINE uint32_t iy(uint32_t id)
void drawTable(int nbRows, int nbColumns)
Definition: EcalDrawUtils.h:91
std::shared_ptr< PayloadType > fetchPayload(const cond::Hash &payloadHash)