CMS 3D CMS Logo

ESChannelStatus_PayloadInspector.cc
Go to the documentation of this file.
4 
5 // the data format of the condition to be inspected
10 
11 #include <memory>
12 #include <sstream>
13 
14 #include "TStyle.h"
15 #include "TH2F.h"
16 #include "TCanvas.h"
17 #include "TLine.h"
18 #include "TLatex.h"
19 
20 namespace {
21  enum { kESChannels = 137216 };
22  enum { IX_MIN = 1, IY_MIN = 1, IX_MAX = 40, IY_MAX = 40 }; // endcaps lower and upper bounds on x and y
23 
24  /********************************************
25  2d plot of ES channel status of 1 IOV
26  *********************************************/
27  class ESChannelStatusPlot : public cond::payloadInspector::PlotImage<ESChannelStatus> {
28  public:
29  ESChannelStatusPlot() : cond::payloadInspector::PlotImage<ESChannelStatus>("ES channel status") {
30  setSingleIov(true);
31  }
32  bool fill(const std::vector<std::tuple<cond::Time_t, cond::Hash> >& iovs) override {
33  TH2F*** esmap = new TH2F**[2];
34  std::string title[2][2] = {{"ES+F", "ES-F"}, {"ES+R", "ES-R"}};
35  for (int plane = 0; plane < 2; plane++) {
36  esmap[plane] = new TH2F*[2];
37  for (int side = 0; side < 2; side++)
38  esmap[plane][side] = new TH2F(
39  Form("esmap%i%i", plane, side), title[plane][side].c_str(), IX_MAX, 0, IX_MAX, IY_MAX, 0, IY_MAX);
40  }
41  Int_t escount = 0;
42  unsigned int run = 0;
43  auto iov = iovs.front();
44  std::shared_ptr<ESChannelStatus> payload = fetchPayload(std::get<1>(iov));
45  run = std::get<0>(iov);
46  if (payload.get()) {
47  // looping over all the ES channels
48  for (int id = 0; id < kESChannels; id++)
49  if (ESDetId::validHashIndex(id)) {
50  ESDetId myESId = ESDetId::unhashIndex(id);
51  int side = myESId.zside(); // -1, 1
52  if (side < 0)
53  side = 1;
54  else
55  side = 0;
56  int plane = myESId.plane() - 1; // 1, 2
57  if (side < 0 || side > 1 || plane < 0 || plane > 1) {
58  std::cout << " channel " << id << " side " << myESId.zside() << " plane " << myESId.plane() << std::endl;
59  return false;
60  }
61  ESChannelStatusCode status_it = (payload->getMap())[myESId];
62  int status = status_it.getStatusCode();
63  if (status != 0) {
64  if (myESId.strip() == 1) { // we get 32 times the same status, plot it only once!
65  esmap[plane][side]->Fill(myESId.six() - 1, myESId.siy() - 1, status);
66  // std::cout << " channel " << id << " side " << myESId.zside() << " plane " << myESId.plane()
67  // << " x " << myESId.six() << " y " << myESId.siy() << " strip " << myESId.strip()
68  // << " status " << status << std::endl;
69  escount++;
70  }
71  }
72  } // validHashIndex
73  } // payload
74 
75  gStyle->SetOptStat(0);
76  gStyle->SetPalette(1);
77  TCanvas canvas("CC map", "CC map", 1680, 1320);
78  TLatex t1;
79  t1.SetNDC();
80  t1.SetTextAlign(26);
81  t1.SetTextSize(0.05);
82  t1.DrawLatex(0.5, 0.96, Form("ES Channel Status, IOV %i", run));
83  t1.SetTextSize(0.025);
84 
85  float xmi[2] = {0.0, 0.5};
86  float xma[2] = {0.5, 1.0};
87  TPad*** pad = new TPad**[2];
88  for (int plane = 0; plane < 2; plane++) {
89  pad[plane] = new TPad*[2];
90  for (int side = 0; side < 2; side++) {
91  float yma = 0.94 - (0.46 * plane);
92  float ymi = yma - 0.44;
93  pad[plane][side] =
94  new TPad(Form("p_%i_%i", plane, side), Form("p_%i_%i", plane, side), xmi[side], ymi, xma[side], yma);
95  pad[plane][side]->Draw();
96  }
97  }
98 
99  for (int side = 0; side < 2; side++) {
100  for (int plane = 0; plane < 2; plane++) {
101  pad[plane][side]->cd();
102  esmap[plane][side]->Draw("colz1");
103  DrawES(plane, side);
104  }
105  }
106  canvas.cd();
107  t1.SetTextSize(0.025);
108  int Nbdead = escount * 32;
109  // int ipercent = Nbdead * 1000000 / kESChannels;
110  // float percent = (float)ipercent / 1000000.; // keep 2 digits
111  // t1.DrawLatex(0.1, 0.94, Form("Number of dead strips %i (%f)", Nbdead, percent));
112  t1.DrawLatex(0.5, 0.92, Form("Number of dead strips %i", Nbdead));
113 
114  std::string ImageName(m_imageFileName);
115  canvas.SaveAs(ImageName.c_str());
116  return true;
117  } // fill method
118  };
119 
120  /*************************************************************
121  2d plot of ES channel status difference between 2 IOVs
122  **************************************************************/
123  template <cond::payloadInspector::IOVMultiplicity nIOVs, int ntags>
124  class ESChannelStatusDiffBase : public cond::payloadInspector::PlotImage<ESChannelStatus, nIOVs, ntags> {
125  public:
126  ESChannelStatusDiffBase()
127  : cond::payloadInspector::PlotImage<ESChannelStatus, nIOVs, ntags>("ES channel status difference") {}
128  bool fill() override {
129  TH2F*** esmap = new TH2F**[2];
130  std::string title[2][2] = {{"ES+F", "ES-F"}, {"ES+R", "ES-R"}};
131  for (int plane = 0; plane < 2; plane++) {
132  esmap[plane] = new TH2F*[2];
133  for (int side = 0; side < 2; side++)
134  esmap[plane][side] = new TH2F(
135  Form("esmap%i%i", plane, side), title[plane][side].c_str(), IX_MAX, 0, IX_MAX, IY_MAX, 0, IY_MAX);
136  }
137  Int_t escount = 0;
138  unsigned int run[2] = {0, 0};
139  int stat[kESChannels];
140  std::string l_tagname[2];
141  // std::cout << " running with " << nIOVs << " IOVs and " << ntags << " tags " << std::endl;
142  auto iovs = cond::payloadInspector::PlotBase::getTag<0>().iovs;
143  l_tagname[0] = cond::payloadInspector::PlotBase::getTag<0>().name;
144  auto firstiov = iovs.front();
145  run[0] = std::get<0>(firstiov);
146  std::tuple<cond::Time_t, cond::Hash> lastiov;
147  if (ntags == 2) {
148  auto tag2iovs = cond::payloadInspector::PlotBase::getTag<1>().iovs;
149  l_tagname[1] = cond::payloadInspector::PlotBase::getTag<1>().name;
150  lastiov = tag2iovs.front();
151  } else {
152  lastiov = iovs.back();
153  l_tagname[1] = l_tagname[0];
154  }
155  run[1] = std::get<0>(lastiov);
156  for (int irun = 0; irun < nIOVs; irun++) {
157  std::shared_ptr<ESChannelStatus> payload;
158  if (irun == 0) {
159  payload = this->fetchPayload(std::get<1>(firstiov));
160  } else {
161  payload = this->fetchPayload(std::get<1>(lastiov));
162  }
163  // std::cout << " irun " << irun << " tag " << l_tagname[irun] << " IOV " << run[irun] << std ::endl;
164  if (payload.get()) {
165  for (int id = 0; id < kESChannels; id++) // looping over all the ES channels
166  if (ESDetId::validHashIndex(id)) {
167  ESDetId myESId = ESDetId::unhashIndex(id);
168  ESChannelStatusCode status_it = (payload->getMap())[myESId];
169  int status = status_it.getStatusCode();
170  if (irun == 0)
171  stat[id] = status;
172  else {
173  int side = myESId.zside(); // -1, 1
174  if (side < 0)
175  side = 1;
176  else
177  side = 0;
178  int plane = myESId.plane() - 1; // 1, 2
179  if (side < 0 || side > 1 || plane < 0 || plane > 1) {
180  std::cout << " channel " << id << " side " << myESId.zside() << " plane " << myESId.plane()
181  << std::endl;
182  return false;
183  }
184  int diff = status - stat[id];
185  if (diff != 0) {
186  if (myESId.strip() == 1) { // we get 32 times the same status, plot it only once!
187  esmap[plane][side]->Fill(myESId.six() - 1, myESId.siy() - 1, diff);
188  escount++;
189  }
190  }
191  } // 2nd IOV
192  } // validHashIndex
193  } // payload
194  } // loop over IOVs
195 
196  gStyle->SetOptStat(0);
197  gStyle->SetPalette(1);
198  TCanvas canvas("CC map", "CC map", 1680, 1320);
199  TLatex t1;
200  t1.SetNDC();
201  t1.SetTextAlign(26);
202  int len = l_tagname[0].length() + l_tagname[1].length();
203  if (ntags == 2) {
204  if (len < 60) {
205  t1.SetTextSize(0.03);
206  t1.DrawLatex(
207  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]));
208  } else {
209  t1.SetTextSize(0.05);
210  t1.DrawLatex(0.5, 0.96, Form("ES Channel Status, IOV %i - %i", run[1], run[0]));
211  }
212  } else {
213  t1.SetTextSize(0.05);
214  t1.DrawLatex(0.5, 0.96, Form("%s IOV %i - %i", l_tagname[0].c_str(), run[1], run[0]));
215  }
216 
217  float xmi[2] = {0.0, 0.5};
218  float xma[2] = {0.5, 1.0};
219  TPad*** pad = new TPad**[2];
220  for (int plane = 0; plane < 2; plane++) {
221  pad[plane] = new TPad*[2];
222  for (int side = 0; side < 2; side++) {
223  float yma = 0.94 - (0.46 * plane);
224  float ymi = yma - 0.44;
225  pad[plane][side] =
226  new TPad(Form("p_%i_%i", plane, side), Form("p_%i_%i", plane, side), xmi[side], ymi, xma[side], yma);
227  pad[plane][side]->Draw();
228  }
229  }
230 
231  for (int side = 0; side < 2; side++) {
232  for (int plane = 0; plane < 2; plane++) {
233  pad[plane][side]->cd();
234  esmap[plane][side]->Draw("colz1");
235  DrawES(plane, side);
236  }
237  }
238  canvas.cd();
239  t1.SetTextSize(0.025);
240  int Nbdead = escount * 32;
241  // int ipercent = Nbdead * 1000000 / kESChannels;
242  // float percent = (float)ipercent / 1000000.; // keep 2 digits
243  // t1.DrawLatex(0.1, 0.94, Form("Number of dead strips %i (%f)", Nbdead, percent));
244  t1.DrawLatex(0.5, 0.92, Form("Number of different strips %i", Nbdead));
245 
246  std::string ImageName(this->m_imageFileName);
247  canvas.SaveAs(ImageName.c_str());
248  return true;
249  } // fill method
250  };
251  using ESChannelStatusDiffOneTag = ESChannelStatusDiffBase<cond::payloadInspector::SINGLE_IOV, 1>;
252  using ESChannelStatusDiffTwoTags = ESChannelStatusDiffBase<cond::payloadInspector::SINGLE_IOV, 2>;
253 } // namespace
254 
255 // Register the classes as boost python plugin
257  PAYLOAD_INSPECTOR_CLASS(ESChannelStatusPlot);
258  PAYLOAD_INSPECTOR_CLASS(ESChannelStatusDiffOneTag);
259  PAYLOAD_INSPECTOR_CLASS(ESChannelStatusDiffTwoTags);
260 }
int zside() const
Definition: ESDetId.h:39
#define PAYLOAD_INSPECTOR_CLASS(CLASS_NAME)
int plane() const
Definition: ESDetId.h:41
#define PAYLOAD_INSPECTOR_MODULE(PAYLOAD_TYPENAME)
int siy() const
Definition: ESDetId.h:45
int six() const
Definition: ESDetId.h:43
void DrawES(int plane, int side)
Definition: ESDrawUtils.h:4
int strip() const
Definition: ESDetId.h:47
def canvas(sub, attr)
Definition: svgfig.py:482
static ESDetId unhashIndex(int hi)
get a DetId from a compact index for arrays
Definition: ESDetId.cc:32
uint16_t getStatusCode() const
static bool validHashIndex(int hi)
Definition: ESDetId.h:59
std::shared_ptr< PayloadType > fetchPayload(const cond::Hash &payloadHash)