CMS 3D CMS Logo

EcalDQMTowerStatus_PayloadInspector.cc
Go to the documentation of this file.
8 
9 // the data format of the condition to be inspected
11 
12 #include "TH2F.h"
13 #include "TCanvas.h"
14 #include "TStyle.h"
15 #include "TLine.h"
16 #include "TLatex.h"
17 
18 #include <string>
19 
20 namespace {
21 
22  enum { kEBTotalTowers = 2448, kEETotalTowers = 1584 };
23  enum { MIN_IETA = 1, MIN_IPHI = 1, MAX_IETA = 17, MAX_IPHI = 72 }; // barrel lower and upper bounds on eta and phi
24  enum { IX_MIN = 1, IY_MIN = 1, IX_MAX = 20, IY_MAX = 20 }; // endcaps lower and upper bounds on x and y
25  /***********************************************
26  2d plot of EcalDQMTowerStatus of 1 IOV
27 ************************************************/
28  class EcalDQMTowerStatusPlot : public cond::payloadInspector::PlotImage<EcalDQMTowerStatus> {
29  public:
30  EcalDQMTowerStatusPlot() : cond::payloadInspector::PlotImage<EcalDQMTowerStatus>("EcalDQMTowerStatus - map ") {
31  setSingleIov(true);
32  }
33 
34  bool fill(const std::vector<std::tuple<cond::Time_t, cond::Hash> >& iovs) override {
35  TH2F* barrel = new TH2F("EB", "EB DQM Tower Status", 72, 0, 72, 34, -17, 17);
36  TH2F* endc_p = new TH2F("EE+", "EE+ DQM Tower Status", 22, 0, 22, 22, 0, 22);
37  TH2F* endc_m = new TH2F("EE-", "EE- DQM Tower Status", 22, 0, 22, 22, 0, 22);
38 
39  auto iov = iovs.front();
40  std::shared_ptr<EcalDQMTowerStatus> payload = fetchPayload(std::get<1>(iov));
41  unsigned int run = std::get<0>(iov);
42  double maxEB = 0, maxEE = 0;
43 
44  if (payload.get()) {
45  for (uint cellid = 0; cellid < EcalTrigTowerDetId::kEBTotalTowers; ++cellid) {
46  if (payload->barrelItems().empty())
47  break;
49  if ((*payload).find(rawid) == (*payload).end())
50  continue;
51 
52  int ieta = rawid.ieta();
53  if (ieta > 0)
54  ieta--; // 1 to 17
55  int iphi = rawid.iphi() - 1; // 0 to 71
56  barrel->Fill(iphi, ieta, (*payload)[rawid].getStatusCode());
57 
58  if (maxEB < (*payload)[rawid].getStatusCode())
59  maxEB = (*payload)[rawid].getStatusCode();
60  }
61 
62  if (payload->endcapItems().empty())
63  return false;
64 
65  for (uint cellid = 0; cellid < EcalTrigTowerDetId::kEETotalTowers; ++cellid) {
66  if (EcalScDetId::validHashIndex(cellid)) {
67  EcalScDetId rawid = EcalScDetId::unhashIndex(cellid);
68  if ((*payload).find(rawid) == (*payload).end())
69  continue;
70  int ix = rawid.ix(); // 1 to 20
71  int iy = rawid.iy(); // 1 to 20
72  int side = rawid.zside();
73  if (side == -1)
74  endc_m->Fill(ix, iy, (*payload)[rawid].getStatusCode());
75  else
76  endc_p->Fill(ix, iy, (*payload)[rawid].getStatusCode());
77 
78  if (maxEE < (*payload)[rawid].getStatusCode())
79  maxEE = (*payload)[rawid].getStatusCode();
80  }
81  }
82 
83  } // payload
84 
85  TCanvas canvas("CC map", "CC map", 800, 800);
86  TLatex t1;
87  t1.SetNDC();
88  t1.SetTextAlign(26);
89  t1.SetTextSize(0.05);
90  t1.DrawLatex(0.5, 0.96, Form("Ecal DQM Tower Status, IOV %i", run));
91 
92  //TPad* padb = new TPad("padb","padb", 0., 0.55, 1., 1.);
93  TPad* padb = new TPad("padb", "padb", 0., 0.45, 1., 0.9);
94  padb->Draw();
95  TPad* padem = new TPad("padem", "padem", 0., 0., 0.45, 0.45);
96  padem->Draw();
97  TPad* padep = new TPad("padep", "padep", 0.55, 0., 1., 0.45);
98  padep->Draw();
99 
100  TLine* l = new TLine(0., 0., 0., 0.);
101  l->SetLineWidth(1);
102  padb->cd();
103  barrel->SetStats(false);
104  barrel->SetMaximum(maxEB);
105  barrel->SetMinimum(0);
106  barrel->Draw("colz");
107  //barrel->Draw("col");
108  for (int i = 0; i < 17; i++) {
109  Double_t x = 4. + (i * 4);
110  l = new TLine(x, -17., x, 17.);
111  l->Draw();
112  }
113 
114  l = new TLine(0., 0., 72., 0.);
115  l->Draw();
116 
117  padem->cd();
118  DrawEE_Tower(endc_m, l, 0, maxEE);
119 
120  padep->cd();
121  DrawEE_Tower(endc_p, l, 0, maxEE);
122 
123  std::string ImageName(m_imageFileName);
124  canvas.SaveAs(ImageName.c_str());
125  return true;
126  } // fill method
127  };
128 
129  /***************************************************
130  2d plot of EcalDQMTowerStatus Diff between 2 IOV
131  ****************************************************/
132  template <cond::payloadInspector::IOVMultiplicity nIOVs, int ntags>
133  class EcalDQMTowerStatusDiffBase : public cond::payloadInspector::PlotImage<EcalDQMTowerStatus, nIOVs, ntags> {
134  public:
135  EcalDQMTowerStatusDiffBase()
136  : cond::payloadInspector::PlotImage<EcalDQMTowerStatus, nIOVs, ntags>("EcalDQMTowerStatusDiff - map ") {}
137  bool fill() override {
138  TH2F* barrel = new TH2F("EB", "EB DQM Tower Status", 72, 0, 72, 34, -17, 17);
139  TH2F* endc_p = new TH2F("EE+", "EE+ DQM Tower Status", 22, 0, 22, 22, 0, 22);
140  TH2F* endc_m = new TH2F("EE-", "EE- DQM Tower Status", 22, 0, 22, 22, 0, 22);
141 
142  unsigned int run[2];
143  float pEB[kEBTotalTowers], pEE[kEETotalTowers];
144  std::string l_tagname[2];
145 
146  auto iovs = cond::payloadInspector::PlotBase::getTag<0>().iovs;
147  l_tagname[0] = cond::payloadInspector::PlotBase::getTag<0>().name;
148  auto firstiov = iovs.front();
149  run[0] = std::get<0>(firstiov);
150  std::tuple<cond::Time_t, cond::Hash> lastiov;
151  if (ntags == 2) {
152  auto tag2iovs = cond::payloadInspector::PlotBase::getTag<1>().iovs;
153  l_tagname[1] = cond::payloadInspector::PlotBase::getTag<1>().name;
154  lastiov = tag2iovs.front();
155  } else {
156  lastiov = iovs.back();
157  l_tagname[1] = l_tagname[0];
158  }
159  run[1] = std::get<0>(lastiov);
160  for (int irun = 0; irun < nIOVs; irun++) {
161  std::shared_ptr<EcalDQMTowerStatus> payload;
162  if (irun == 0) {
163  payload = this->fetchPayload(std::get<1>(firstiov));
164  } else {
165  payload = this->fetchPayload(std::get<1>(lastiov));
166  }
167  if (payload.get()) {
168  for (uint cellid = 0; cellid < EcalTrigTowerDetId::kEBTotalTowers; ++cellid) {
169  if (payload->barrelItems().empty())
170  break;
171 
173  if ((*payload).find(rawid) == (*payload).end())
174  continue;
175 
176  float weight = (*payload)[rawid].getStatusCode();
177 
178  if (irun == 0) {
179  pEB[cellid] = weight;
180  } else {
181  int ieta = rawid.ieta();
182  if (ieta > 0)
183  ieta--; // 1 to 17
184  int iphi = rawid.iphi() - 1; // 0 to 71
185  unsigned int new_status = (*payload)[rawid].getStatusCode();
186  if (new_status != pEB[cellid]) {
187  int tmp3 = 0;
188 
189  if (new_status > pEB[cellid])
190  tmp3 = 1;
191  else
192  tmp3 = -1;
193 
194  barrel->Fill(iphi, ieta, 0.05 + 0.95 * (tmp3 > 0));
195  }
196  }
197  }
198 
199  if (payload->endcapItems().empty())
200  return false;
201 
202  for (uint cellid = 0; cellid < EcalTrigTowerDetId::kEETotalTowers; ++cellid) {
203  if (EcalScDetId::validHashIndex(cellid)) {
204  EcalScDetId rawid = EcalScDetId::unhashIndex(cellid);
205  if ((*payload).find(rawid) == (*payload).end())
206  continue;
207 
208  float weight = (*payload)[rawid].getStatusCode();
209 
210  if (irun == 0) {
211  pEE[cellid] = weight;
212  } else {
213  int ix = rawid.ix(); // 1 to 20
214  int iy = rawid.iy(); // 1 to 20
215  int side = rawid.zside();
216 
217  unsigned int new_status = (*payload)[rawid].getStatusCode();
218  if (new_status != pEE[cellid]) {
219  int tmp3 = 0;
220 
221  if (new_status > pEE[cellid])
222  tmp3 = 1;
223  else
224  tmp3 = -1;
225 
226  if (side == -1)
227  endc_m->Fill(ix, iy, 0.05 + 0.95 * (tmp3 > 0));
228  else
229  endc_p->Fill(ix, iy, 0.05 + 0.95 * (tmp3 > 0));
230  }
231  }
232  }
233  }
234  } // payload
235  }
236 
237  TCanvas canvas("CC map", "CC map", 800, 800);
238  TLatex t1;
239  t1.SetNDC();
240  t1.SetTextAlign(26);
241  int len = l_tagname[0].length() + l_tagname[1].length();
242  if (ntags == 2 && len < 58) {
243  t1.SetTextSize(0.025);
244  t1.DrawLatex(
245  0.5, 0.96, Form("%s IOV %i - %s IOV %i", l_tagname[1].c_str(), run[1], l_tagname[0].c_str(), run[0]));
246  } else {
247  t1.SetTextSize(0.04);
248  t1.DrawLatex(0.5, 0.96, Form("Ecal DQM Tower Status (Diff), IOV %i vs %i", run[0], run[1]));
249  }
250  TPad* padb = new TPad("padb", "padb", 0., 0.45, 1., 0.9);
251  padb->Draw();
252  TPad* padem = new TPad("padem", "padem", 0., 0., 0.45, 0.45);
253  padem->Draw();
254  TPad* padep = new TPad("padep", "padep", 0.55, 0., 1., 0.45);
255  padep->Draw();
256 
257  TLine* l = new TLine(0., 0., 0., 0.);
258  l->SetLineWidth(1);
259  padb->cd();
260  barrel->SetStats(false);
261  barrel->SetMaximum(1.15);
262  barrel->SetMinimum(0);
263  barrel->Draw("colz");
264  //barrel->Draw("col");
265  for (int i = 0; i < 17; i++) {
266  Double_t x = 4. + (i * 4);
267  l = new TLine(x, -17., x, 17.);
268  l->Draw();
269  }
270 
271  l = new TLine(0., 0., 72., 0.);
272  l->Draw();
273 
274  padem->cd();
275  DrawEE_Tower(endc_m, l, 0, 1.15);
276 
277  padep->cd();
278  DrawEE_Tower(endc_p, l, 0, 1.15);
279 
280  std::string ImageName(this->m_imageFileName);
281  canvas.SaveAs(ImageName.c_str());
282  return true;
283  } // fill method
284  }; // class EcalDQMTowerStatusDiffBase
285  using EcalDQMTowerStatusDiffOneTag = EcalDQMTowerStatusDiffBase<cond::payloadInspector::SINGLE_IOV, 1>;
286  using EcalDQMTowerStatusDiffTwoTags = EcalDQMTowerStatusDiffBase<cond::payloadInspector::SINGLE_IOV, 2>;
287 
288 } // namespace
289 
290 // Register the classes as boost python plugin
292  PAYLOAD_INSPECTOR_CLASS(EcalDQMTowerStatusPlot);
293  PAYLOAD_INSPECTOR_CLASS(EcalDQMTowerStatusDiffOneTag);
294  PAYLOAD_INSPECTOR_CLASS(EcalDQMTowerStatusDiffTwoTags);
295 }
static EcalTrigTowerDetId detIdFromDenseIndex(uint32_t di)
int zside() const
Definition: EcalScDetId.h:64
int ieta() const
get the tower ieta
static EcalScDetId unhashIndex(int hi)
Definition: EcalScDetId.h:117
Definition: weight.py:1
#define PAYLOAD_INSPECTOR_CLASS(CLASS_NAME)
void DrawEE_Tower(TH2F *endc, TLine *l, double minScale, double maxScale)
Definition: EcalDrawUtils.h:65
#define PAYLOAD_INSPECTOR_MODULE(PAYLOAD_TYPENAME)
int iy() const
Definition: EcalScDetId.h:76
ALPAKA_FN_ACC ALPAKA_FN_INLINE uint32_t ix(uint32_t id)
static bool validHashIndex(int hi)
Definition: EcalScDetId.h:139
def canvas(sub, attr)
Definition: svgfig.py:482
ALPAKA_FN_ACC ALPAKA_FN_INLINE uint32_t iy(uint32_t id)
int iphi() const
get the tower iphi
std::shared_ptr< PayloadType > fetchPayload(const cond::Hash &payloadHash)
int ix() const
Definition: EcalScDetId.h:70