CMS 3D CMS Logo

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