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