CMS 3D CMS Logo

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