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  class ESChannelStatusDiff : public cond::payloadInspector::PlotImage<ESChannelStatus> {
124  public:
125  ESChannelStatusDiff() : cond::payloadInspector::PlotImage<ESChannelStatus>("ES channel status difference") {
126  setSingleIov(false);
127  }
128  bool fill(const std::vector<std::tuple<cond::Time_t, cond::Hash> >& iovs) 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], irun = 0;
139  int stat[kESChannels];
140  for (auto const& iov : iovs) {
141  std::shared_ptr<ESChannelStatus> payload = fetchPayload(std::get<1>(iov));
142  run[irun] = std::get<0>(iov);
143  // std::cout << " irun " << irun << " IOV " << run[irun] << std::endl;
144  if (payload.get()) {
145  for (int id = 0; id < kESChannels; id++) // looping over all the ES channels
146  if (ESDetId::validHashIndex(id)) {
147  ESDetId myESId = ESDetId::unhashIndex(id);
148  ESChannelStatusCode status_it = (payload->getMap())[myESId];
149  int status = status_it.getStatusCode();
150  if (irun == 0)
151  stat[id] = status;
152  else {
153  int side = myESId.zside(); // -1, 1
154  if (side < 0)
155  side = 1;
156  else
157  side = 0;
158  int plane = myESId.plane() - 1; // 1, 2
159  if (side < 0 || side > 1 || plane < 0 || plane > 1) {
160  std::cout << " channel " << id << " side " << myESId.zside() << " plane " << myESId.plane()
161  << std::endl;
162  return false;
163  }
164  int diff = status - stat[id];
165  if (diff != 0) {
166  if (myESId.strip() == 1) { // we get 32 times the same status, plot it only once!
167  esmap[plane][side]->Fill(myESId.six() - 1, myESId.siy() - 1, diff);
168  escount++;
169  }
170  }
171  } // 2nd IOV
172  } // validHashIndex
173  } // payload
174  irun++;
175  } // loop over IOVs
176 
177  gStyle->SetOptStat(0);
178  gStyle->SetPalette(1);
179  TCanvas canvas("CC map", "CC map", 1680, 1320);
180  TLatex t1;
181  t1.SetNDC();
182  t1.SetTextAlign(26);
183  t1.SetTextSize(0.05);
184  t1.DrawLatex(0.5, 0.96, Form("ES Channel Status, IOV %i - %i", run[1], run[0]));
185  t1.SetTextSize(0.025);
186 
187  float xmi[2] = {0.0, 0.5};
188  float xma[2] = {0.5, 1.0};
189  TPad*** pad = new TPad**[2];
190  for (int plane = 0; plane < 2; plane++) {
191  pad[plane] = new TPad*[2];
192  for (int side = 0; side < 2; side++) {
193  float yma = 0.94 - (0.46 * plane);
194  float ymi = yma - 0.44;
195  pad[plane][side] =
196  new TPad(Form("p_%i_%i", plane, side), Form("p_%i_%i", plane, side), xmi[side], ymi, xma[side], yma);
197  pad[plane][side]->Draw();
198  }
199  }
200 
201  for (int side = 0; side < 2; side++) {
202  for (int plane = 0; plane < 2; plane++) {
203  pad[plane][side]->cd();
204  esmap[plane][side]->Draw("colz1");
205  DrawES(plane, side);
206  }
207  }
208  canvas.cd();
209  t1.SetTextSize(0.025);
210  int Nbdead = escount * 32;
211  // int ipercent = Nbdead * 1000000 / kESChannels;
212  // float percent = (float)ipercent / 1000000.; // keep 2 digits
213  // t1.DrawLatex(0.1, 0.94, Form("Number of dead strips %i (%f)", Nbdead, percent));
214  t1.DrawLatex(0.5, 0.92, Form("Number of different strips %i", Nbdead));
215 
216  std::string ImageName(m_imageFileName);
217  canvas.SaveAs(ImageName.c_str());
218  return true;
219  } // fill method
220  };
221 } // namespace
222 
223 // Register the classes as boost python plugin
225  PAYLOAD_INSPECTOR_CLASS(ESChannelStatusPlot);
226  PAYLOAD_INSPECTOR_CLASS(ESChannelStatusDiff);
227 }
change_name.diff
diff
Definition: change_name.py:13
svgfig.canvas
def canvas(*sub, **attr)
Definition: svgfig.py:482
runGCPTkAlMap.title
string title
Definition: runGCPTkAlMap.py:94
ESChannelStatusCode.h
ESChannelStatusCode
Definition: ESChannelStatusCode.h:9
IX_MAX
Definition: EcalFloatCondObjectContainerUtils.h:20
mps_update.status
status
Definition: mps_update.py:69
DrawES
void DrawES(int plane, int side)
Definition: ESDrawUtils.h:4
IY_MAX
Definition: EcalFloatCondObjectContainerUtils.h:21
PayloadInspector.h
ESDetId::strip
int strip() const
Definition: ESDetId.h:47
gather_cfg.cout
cout
Definition: gather_cfg.py:144
ESDrawUtils.h
ESChannelStatusCode::getStatusCode
uint16_t getStatusCode() const
Definition: ESChannelStatusCode.h:21
PAYLOAD_INSPECTOR_CLASS
#define PAYLOAD_INSPECTOR_CLASS(CLASS_NAME)
Definition: PayloadInspectorModule.h:10
ESDetId
Definition: ESDetId.h:15
ESDetId.h
PayloadInspectorModule.h
RandomServiceHelper.t1
t1
Definition: RandomServiceHelper.py:256
hgcalPlots.stat
stat
Definition: hgcalPlots.py:1119
jets_cff.payload
payload
Definition: jets_cff.py:32
PAYLOAD_INSPECTOR_MODULE
#define PAYLOAD_INSPECTOR_MODULE(PAYLOAD_TYPENAME)
Definition: PayloadInspectorModule.h:8
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
cond
Definition: plugin.cc:23
Time.h
IY_MIN
Definition: EcalFloatCondObjectContainerUtils.h:19
cond::payloadInspector::PlotImage::fetchPayload
std::shared_ptr< PayloadType > fetchPayload(const cond::Hash &payloadHash)
Definition: PayloadInspector.h:863
ESChannelStatus.h
trackerHitRTTI::vector
Definition: trackerHitRTTI.h:21
ESCondObjectContainer
Definition: ESCondObjectContainer.h:11
ESDetId::plane
int plane() const
Definition: ESDetId.h:41
ESDetId::six
int six() const
Definition: ESDetId.h:43
writedatasetfile.run
run
Definition: writedatasetfile.py:27
ESDetId::unhashIndex
static ESDetId unhashIndex(int hi)
get a DetId from a compact index for arrays
Definition: ESDetId.cc:32
triggerObjects_cff.id
id
Definition: triggerObjects_cff.py:31
cond::payloadInspector::PlotImpl::fill
virtual bool fill()=0
cond::payloadInspector::PlotImage
Definition: PayloadInspector.h:852
ESDetId::validHashIndex
static bool validHashIndex(int hi)
Definition: ESDetId.h:59
ESDetId::siy
int siy() const
Definition: ESDetId.h:45
IX_MIN
Definition: EcalFloatCondObjectContainerUtils.h:18
ESDetId::zside
int zside() const
Definition: ESDetId.h:39