CMS 3D CMS Logo

SiPixelPayloadInspectorHelper.h
Go to the documentation of this file.
1 #ifndef CONDCORE_SIPIXELPLUGINS_SIPIXELPAYLOADINSPECTORHELPER_H
2 #define CONDCORE_SIPIXELPLUGINS_SIPIXELPAYLOADINSPECTORHELPER_H
3 
4 #include <vector>
5 #include <numeric>
6 #include <fstream> // std::ifstream
7 #include <string>
8 #include <bitset>
9 
10 #include "TGraph.h"
11 #include "TH1.h"
12 #include "TH2.h"
13 #include "TLatex.h"
14 #include "TLine.h"
15 #include "TPave.h"
16 #include "TPaveStats.h"
17 #include "TPaveText.h"
18 #include "TStyle.h"
19 #include "TCanvas.h"
20 
28 
29 //#define MMDEBUG
30 #ifdef MMDEBUG
31 #include <iostream>
32 #define COUT std::cout << "MM "
33 #else
34 #define COUT edm::LogVerbatim("")
35 #endif
36 
37 namespace SiPixelPI {
38 
39  enum phase { zero = 0, one = 1, two = 2, undefined = 999 };
40 
41  // size of the phase-0 pixel detID list
42  static const unsigned int phase0size = 1440;
43  static const unsigned int phase1size = 1856;
44  static const unsigned int phase2size = 3892;
45  static const unsigned int mismatched = 9999;
46 
47  //============================================================================
48  // struct to store info useful to construct topology based on the detid list
49  struct PhaseInfo {
50  // construct with det size
51  PhaseInfo(unsigned int size) : m_detsize(size) {}
52  // construct passing the phase
53  PhaseInfo(const phase& thePhase) {
54  switch (thePhase) {
55  case phase::zero:
57  break;
58  case phase::one:
60  break;
61  case phase::two:
63  break;
64  default:
65  m_detsize = 99999;
66  edm::LogError("PhaseInfo") << "undefined phase: " << thePhase;
67  }
68  }
69  virtual ~PhaseInfo() { edm::LogInfo("PhaseInfo") << "PhaseInfo::~PhaseInfo()\n"; }
70  const SiPixelPI::phase phase() const {
71  if (m_detsize == phase0size)
72  return phase::zero;
73  else if (m_detsize == phase1size)
74  return phase::one;
75  else if (m_detsize > phase1size)
76  return phase::two;
77  else {
78  throw cms::Exception("LogicError") << "this detId list size: " << m_detsize << "should not exist!";
79  }
80  }
81 
82  const char* pathToTopoXML() {
83  if (m_detsize == phase0size)
84  return "Geometry/TrackerCommonData/data/trackerParameters.xml";
85  else if (m_detsize == phase1size)
86  return "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
87  else if (m_detsize > phase1size)
88  return "Geometry/TrackerCommonData/data/PhaseII/trackerParameters.xml";
89  else {
90  throw cms::Exception("LogicError") << "this detId list size: " << m_detsize << "should not exist!";
91  }
92  }
93 
94  const void print(std::stringstream& ss) {
95  ss << "---------------------------------------------------------------\n"
96  << " PhaseInfo Data \n\n"
97  << " Phase : " << phase() << "\n"
98  << " DetSide: " << m_detsize << "\n"
99  << " pathToXML: " << pathToTopoXML() << "\n"
100  << "-------------------------------------------------------------\n\n";
101  }
102 
103  const bool isPhase1Comparison(const PhaseInfo& theOtherPhase) const {
104  if (phase() == phase::one || theOtherPhase.phase() == phase::one)
105  return true;
106  else
107  return false;
108  }
109 
110  const bool isComparedWithPhase2(const PhaseInfo& theOtherPhase) const {
111  if ((phase() == phase::two && theOtherPhase.phase() != phase::two) ||
112  (phase() != phase::two && theOtherPhase.phase() == phase::two)) {
113  return true;
114  } else {
115  return false;
116  }
117  }
118 
119  private:
120  size_t m_detsize;
121  };
122 
123  //============================================================================
124  // add ostream for PhaseInfo
125  inline std::ostream& operator<<(std::ostream& os, PhaseInfo phInfo) {
126  std::stringstream ss;
127  phInfo.print(ss);
128  os << ss.str();
129  return os;
130  }
131 
132  //============================================================================
133  inline std::pair<unsigned int, unsigned int> unpack(cond::Time_t since) {
134  auto kLowMask = 0XFFFFFFFF;
135  auto run = (since >> 32);
136  auto lumi = (since & kLowMask);
137  return std::make_pair(run, lumi);
138  }
139 
140  //============================================================================
141  // Taken from pixel naming classes
142  // BmO (-z-x) = 1, BmI (-z+x) = 2 , BpO (+z-x) = 3 , BpI (+z+x) = 4
143  inline int quadrant(const DetId& detid, const TrackerTopology* tTopo_, bool phase_) {
144  if (detid.subdetId() == PixelSubdetector::PixelBarrel) {
145  return PixelBarrelName(detid, tTopo_, phase_).shell();
146  } else {
147  return PixelEndcapName(detid, tTopo_, phase_).halfCylinder();
148  }
149  }
150 
151  //============================================================================
152  // Online ladder convention taken from pixel naming class for barrel
153  // Apply sign convention (- sign for BmO and BpO)
154  inline int signed_ladder(const DetId& detid, const TrackerTopology& tTopo_, bool phase_) {
156  return -9999;
157  int signed_ladder = PixelBarrelName(detid, &tTopo_, phase_).ladderName();
158  if (quadrant(detid, &tTopo_, phase_) % 2)
159  signed_ladder *= -1;
160  return signed_ladder;
161  }
162 
163  //============================================================================
164  // Online mdoule convention taken from pixel naming class for barrel
165  // Apply sign convention (- sign for BmO and BmI)
166  inline int signed_module(const DetId& detid, const TrackerTopology& tTopo_, bool phase_) {
168  return -9999;
169  int signed_module = PixelBarrelName(detid, &tTopo_, phase_).moduleName();
170  if (quadrant(detid, &tTopo_, phase_) < 3)
171  signed_module *= -1;
172  return signed_module;
173  }
174 
175  //============================================================================
176  // Phase 0: Ring was not an existing convention
177  // but the 7 plaquettes were split by HV group
178  // --> Derive Ring 1/2 for them
179  // Panel 1 plq 1-2, Panel 2, plq 1 = Ring 1
180  // Panel 1 plq 3-4, Panel 2, plq 2-3 = Ring 2
181  // Phase 1: Using pixel naming class for endcap
182  inline int ring(const DetId& detid, const TrackerTopology& tTopo_, bool phase_) {
184  return -9999;
185  int ring = -9999;
186  if (phase_ == 0) {
187  ring = 1 + (tTopo_.pxfPanel(detid) + tTopo_.pxfModule(detid) > 3);
188  } else if (phase_ == 1) {
189  ring = PixelEndcapName(detid, &tTopo_, phase_).ringName();
190  }
191  return ring;
192  }
193 
194  //============================================================================
195  // Online blade convention taken from pixel naming class for endcap
196  // Apply sign convention (- sign for BmO and BpO)
197  inline int signed_blade(const DetId& detid, const TrackerTopology& tTopo_, bool phase_) {
199  return -9999;
200  int signed_blade = PixelEndcapName(detid, &tTopo_, phase_).bladeName();
201  if (quadrant(detid, &tTopo_, phase_) % 2)
202  signed_blade *= -1;
203  return signed_blade;
204  }
205 
206  //============================================================================
207  inline int signed_blade_panel(const DetId& detid, const TrackerTopology& tTopo_, bool phase_) {
209  return -9999;
210  int signed_blade_panel = signed_blade(detid, tTopo_, phase_) + (tTopo_.pxfPanel(detid) - 1);
211  return signed_blade_panel;
212  }
213 
214  //============================================================================
215  // Online disk convention
216  // Apply sign convention (- sign for BmO and BmI)
217  inline int signed_disk(const DetId& detid, const TrackerTopology& tTopo_, bool phase_) {
219  return -9999;
220  int signed_disk = tTopo_.pxfDisk(DetId(detid));
221  if (quadrant(detid, &tTopo_, phase_) < 3)
222  signed_disk *= -1;
223  return signed_disk;
224  }
225 
226  //============================================================================
227  inline void draw_line(double x1, double x2, double y1, double y2, int width = 2, int style = 1, int color = 1) {
228  TLine* l = new TLine(x1, y1, x2, y2);
229  l->SetBit(kCanDelete);
230  l->SetLineWidth(width);
231  l->SetLineStyle(style);
232  l->SetLineColor(color);
233  l->Draw();
234  }
235 
236  //============================================================================
237  inline void dress_occup_plot(TCanvas& canv,
238  TH2* h,
239  int lay,
240  int ring = 0,
241  int phase = 0,
242  bool half_shift = true,
243  bool mark_zero = true,
244  bool standard_palette = true) {
245  std::string s_title;
246 
247  if (lay > 0) {
248  canv.cd(lay);
249  s_title = "Barrel Pixel Layer " + std::to_string(lay);
250  } else {
251  canv.cd(ring);
252  if (ring > 4) {
253  ring = ring - 4;
254  }
255  s_title = "Forward Pixel Ring " + std::to_string(ring);
256  }
257 
258  gStyle->SetPadRightMargin(0.125);
259 
260  if (standard_palette) {
261  gStyle->SetPalette(1);
262  } else {
263  // this is the fine gradient palette
264  const Int_t NRGBs = 5;
265  const Int_t NCont = 255;
266 
267  Double_t stops[NRGBs] = {0.00, 0.34, 0.61, 0.84, 1.00};
268  Double_t red[NRGBs] = {0.00, 0.00, 0.87, 1.00, 0.51};
269  Double_t green[NRGBs] = {0.00, 0.81, 1.00, 0.20, 0.00};
270  Double_t blue[NRGBs] = {0.51, 1.00, 0.12, 0.00, 0.00};
271  TColor::CreateGradientColorTable(NRGBs, stops, red, green, blue, NCont);
272  gStyle->SetNumberContours(NCont);
273  }
274 
275  h->SetMarkerSize(0.7);
276  h->Draw("colz1");
277 
278  auto ltx = TLatex();
279  ltx.SetTextFont(62);
280  ltx.SetTextColor(1);
281  ltx.SetTextSize(0.06);
282  ltx.SetTextAlign(31);
283  ltx.DrawLatexNDC(1 - gPad->GetRightMargin(), 1 - gPad->GetTopMargin() + 0.01, (s_title).c_str());
284 
285  // Draw Lines around modules
286  if (lay > 0) {
287  std::vector<std::vector<int>> nladder = {{10, 16, 22}, {6, 14, 22, 32}};
288  int nlad = nladder[phase][lay - 1];
289  for (int xsign = -1; xsign <= 1; xsign += 2)
290  for (int ysign = -1; ysign <= 1; ysign += 2) {
291  float xlow = xsign * (half_shift * 0.5);
292  float xhigh = xsign * (half_shift * 0.5 + 4);
293  float ylow = ysign * (half_shift * 0.5 + (phase == 0) * 0.5);
294  float yhigh = ysign * (half_shift * 0.5 - (phase == 0) * 0.5 + nlad);
295  // Outside box
296  draw_line(xlow, xhigh, ylow, ylow, 1); // bottom
297  draw_line(xlow, xhigh, yhigh, yhigh, 1); // top
298  draw_line(xlow, xlow, ylow, yhigh, 1); // left
299  draw_line(xhigh, xhigh, ylow, yhigh, 1); // right
300  // Inner Horizontal lines
301  for (int lad = 1; lad < nlad; ++lad) {
302  float y = ysign * (lad + half_shift * 0.5);
303  draw_line(xlow, xhigh, y, y, 1);
304  }
305  for (int lad = 1; lad <= nlad; ++lad)
306  if (!(phase == 0 && (lad == 1 || lad == nlad))) {
307  float y = ysign * (lad + half_shift * 0.5 - 0.5);
308  draw_line(xlow, xhigh, y, y, 1, 3);
309  }
310  // Inner Vertical lines
311  for (int mod = 1; mod < 4; ++mod) {
312  float x = xsign * (mod + half_shift * 0.5);
313  draw_line(x, x, ylow, yhigh, 1);
314  }
315  // Make a BOX around ROC 0
316  // Phase 0 - ladder +1 is always non-flipped
317  // Phase 1 - ladder +1 is always flipped
318  if (mark_zero) {
319  for (int mod = 1; mod <= 4; ++mod)
320  for (int lad = 1; lad <= nlad; ++lad) {
321  bool flipped = ysign == 1 ? lad % 2 == 0 : lad % 2 == 1;
322  if (phase == 1)
323  flipped = !flipped;
324  int roc0_orientation = flipped ? -1 : 1;
325  if (xsign == -1)
326  roc0_orientation *= -1;
327  if (ysign == -1)
328  roc0_orientation *= -1;
329  float x1 = xsign * (mod + half_shift * 0.5);
330  float x2 = xsign * (mod + half_shift * 0.5 - 1. / 8);
331  float y1 = ysign * (lad + half_shift * 0.5 - 0.5);
332  float y2 = ysign * (lad + half_shift * 0.5 - 0.5 + roc0_orientation * 1. / 2);
333  if (!(phase == 0 && (lad == 1 || lad == nlad) && xsign == -1)) {
334  if (lay == 1 && xsign <= -1) {
335  float x1 = xsign * ((mod - 1) + half_shift * 0.5);
336  float x2 = xsign * ((mod - 1) + half_shift * 0.5 + 1. / 8);
337  float y1 = ysign * (lad + half_shift * 0.5 - 0.5 + roc0_orientation);
338  float y2 = ysign * (lad + half_shift * 0.5 - 0.5 + roc0_orientation * 3. / 2);
339  draw_line(x1, x2, y1, y1, 1);
340  draw_line(x2, x2, y1, y2, 1);
341  } else {
342  draw_line(x1, x2, y1, y1, 1);
343  //draw_line(x1, x2, y2, y2, 1);
344  //draw_line(x1, x1, y1, y2, 1);
345  draw_line(x2, x2, y1, y2, 1);
346  }
347  }
348  }
349  }
350  }
351  } else {
352  // FPIX
353  for (int dsk = 1, ndsk = 2 + (phase == 1); dsk <= ndsk; ++dsk) {
354  for (int xsign = -1; xsign <= 1; xsign += 2)
355  for (int ysign = -1; ysign <= 1; ysign += 2) {
356  if (phase == 0) {
357  int first_roc = 3, nbin = 16;
358  for (int bld = 1, nbld = 12; bld <= nbld; ++bld) {
359  // Horizontal lines
360  for (int plq = 1, nplq = 7; plq <= nplq; ++plq) {
361  float xlow =
362  xsign * (half_shift * 0.5 + dsk - 1 + (first_roc - 3 + 2 * plq + (plq == 1)) / (float)nbin);
363  float xhigh =
364  xsign * (half_shift * 0.5 + dsk - 1 + (first_roc - 3 + 2 * (plq + 1) - (plq == 7)) / (float)nbin);
365  float ylow = ysign * (half_shift * 0.5 + (bld - 0.5) - (2 + plq / 2) * 0.1);
366  float yhigh = ysign * (half_shift * 0.5 + (bld - 0.5) + (2 + plq / 2) * 0.1);
367  draw_line(xlow, xhigh, ylow, ylow, 1); // bottom
368  draw_line(xlow, xhigh, yhigh, yhigh, 1); // top
369  }
370  // Vertical lines
371  for (int plq = 1, nplq = 7 + 1; plq <= nplq; ++plq) {
372  float x = xsign * (half_shift * 0.5 + dsk - 1 +
373  (first_roc - 3 + 2 * plq + (plq == 1) - (plq == 8)) / (float)nbin);
374  float ylow = ysign * (half_shift * 0.5 + (bld - 0.5) - (2 + (plq - (plq == 8)) / 2) * 0.1);
375  float yhigh = ysign * (half_shift * 0.5 + (bld - 0.5) + (2 + (plq - (plq == 8)) / 2) * 0.1);
376  draw_line(x, x, ylow, yhigh, 1);
377  }
378  // Panel 2 has dashed mid-plane
379  for (int plq = 2, nplq = 6; plq <= nplq; ++plq)
380  if (plq % 2 == 0) {
381  float x = xsign * (half_shift * 0.5 + dsk - 1 +
382  (first_roc - 3 + 2 * plq + (plq == 1) - (plq == 8) + 1) / (float)nbin);
383  float ylow = ysign * (half_shift * 0.5 + (bld - 0.5) - (2 + (plq - (plq == 8)) / 2) * 0.1);
384  float yhigh = ysign * (half_shift * 0.5 + (bld - 0.5) + (2 + (plq - (plq == 8)) / 2) * 0.1);
385  draw_line(x, x, ylow, yhigh, 1, 2);
386  }
387  // Make a BOX around ROC 0
388  for (int plq = 1, nplq = 7; plq <= nplq; ++plq) {
389  float x1 =
390  xsign * (half_shift * 0.5 + dsk - 1 + (first_roc - 3 + 2 * plq + (plq == 1)) / (float)nbin);
391  float x2 =
392  xsign * (half_shift * 0.5 + dsk - 1 + (first_roc - 3 + 2 * plq + (plq == 1) + 1) / (float)nbin);
393  int sign = xsign * ysign * ((plq % 2) ? 1 : -1);
394  float y1 = ysign * (half_shift * 0.5 + (bld - 0.5) + sign * (2 + plq / 2) * 0.1);
395  float y2 = ysign * (half_shift * 0.5 + (bld - 0.5) + sign * (plq / 2) * 0.1);
396  //draw_line(x1, x2, y1, y1, 1);
397  draw_line(x1, x2, y2, y2, 1);
398  //draw_line(x1, x1, y1, y2, 1);
399  draw_line(x2, x2, y1, y2, 1);
400  }
401  }
402  } else if (phase == 1) {
403  if (ring == 0) { // both
404  for (int ring = 1; ring <= 2; ++ring)
405  for (int bld = 1, nbld = 5 + ring * 6; bld <= nbld; ++bld) {
406  float scale = (ring == 1) ? 1.5 : 1;
407  Color_t p1_color = 1, p2_color = 1;
408  // Horizontal lines
409  // Panel 2 has dashed mid-plane
410  float x1 = xsign * (half_shift * 0.5 + dsk - 1 + (ring - 1) * 0.5);
411  float x2 = xsign * (half_shift * 0.5 + dsk - 1 + ring * 0.5);
412  int sign = ysign;
413  float y1 = ysign * (half_shift * 0.5 - 0.5 + scale * bld + sign * 0.5);
414  //float yp1_mid = ysign * (half_shift*0.5 - 0.5 + scale*bld + sign*0.25);
415  float y2 = ysign * (half_shift * 0.5 - 0.5 + scale * bld);
416  float yp2_mid = ysign * (half_shift * 0.5 - 0.5 + scale * bld - sign * 0.25);
417  float y3 = ysign * (half_shift * 0.5 - 0.5 + scale * bld - sign * 0.5);
418  draw_line(x1, x2, y1, y1, 1, 1, p1_color);
419  //draw_line(x1, x2, yp1_mid, yp1_mid, 1, 3);
420  draw_line(x1, x2, y2, y2, 1, 1, p1_color);
421  draw_line(x1, x2, yp2_mid, yp2_mid, 1, 2);
422  draw_line(x1, x2, y3, y3, 1, 1, p2_color);
423  // Vertical lines
424  float x = xsign * (half_shift * 0.5 + dsk - 1 + (ring - 1) * 0.5);
425  draw_line(x, x, y1, y2, 1, 1, p1_color);
426  draw_line(x, x, y2, y3, 1, 1, p2_color);
427  if (ring == 2) {
428  //draw_line(x, x, y2, y3, 1, 1, p1_color);
429  x = xsign * (half_shift * 0.5 + dsk);
430  draw_line(x, x, y1, y2, 1, 1, p1_color);
431  draw_line(x, x, y2, y3, 1, 1, p2_color);
432  }
433  // Make a BOX around ROC 0
434  x1 = xsign * (half_shift * 0.5 + dsk - 1 + ring * 0.5 - 1 / 16.);
435  x2 = xsign * (half_shift * 0.5 + dsk - 1 + ring * 0.5);
436  float y1_p1 = ysign * (half_shift * 0.5 - 0.5 + scale * bld + sign * 0.25);
437  float y2_p1 = ysign * (half_shift * 0.5 - 0.5 + scale * bld + sign * 0.25 + xsign * ysign * 0.25);
438  draw_line(x1, x2, y1_p1, y1_p1, 1);
439  //draw_line(x1, x2, y2_p1, y2_p1, 1);
440  draw_line(x1, x1, y1_p1, y2_p1, 1);
441  //draw_line(x2, x2, y1_p1, y2_p1, 1);
442  float y1_p2 = ysign * (half_shift * 0.5 - 0.5 + scale * bld - sign * 0.25);
443  float y2_p2 = ysign * (half_shift * 0.5 - 0.5 + scale * bld - sign * 0.25 - xsign * ysign * 0.25);
444  draw_line(x1, x2, y1_p2, y1_p2, 1);
445  //draw_line(x1, x2, y2_p2, y2_p2, 1);
446  draw_line(x1, x1, y1_p2, y2_p2, 1);
447  //draw_line(x2, x2, y1_p2, y2_p2, 1);
448  }
449  } else { // only one ring, 1 or 2
450  for (int bld = 1, nbld = 5 + ring * 6; bld <= nbld; ++bld) {
451  Color_t p1_color = 1, p2_color = 1;
452  // Horizontal lines
453  // Panel 2 has dashed mid-plane
454  float x1 = xsign * (half_shift * 0.5 + dsk - 1);
455  float x2 = xsign * (half_shift * 0.5 + dsk);
456  int sign = ysign;
457  float y1 = ysign * (half_shift * 0.5 - 0.5 + bld + sign * 0.5);
458  //float yp1_mid = ysign * (half_shift*0.5 - 0.5 + bld + sign*0.25);
459  float y2 = ysign * (half_shift * 0.5 - 0.5 + bld);
460  float yp2_mid = ysign * (half_shift * 0.5 - 0.5 + bld - sign * 0.25);
461  float y3 = ysign * (half_shift * 0.5 - 0.5 + bld - sign * 0.5);
462  draw_line(x1, x2, y1, y1, 1, 1, p1_color);
463  //draw_line(x1, x2, yp1_mid, yp1_mid, 1, 3);
464  draw_line(x1, x2, y2, y2, 1, 1, p1_color);
465  draw_line(x1, x2, yp2_mid, yp2_mid, 1, 2);
466  draw_line(x1, x2, y3, y3, 1, 1, p2_color);
467  // Vertical lines
468  float x = xsign * (half_shift * 0.5 + dsk - 1);
469  draw_line(x, x, y1, y2, 1, 1, p1_color);
470  draw_line(x, x, y2, y3, 1, 1, p2_color);
471  if (ring == 2) {
472  //draw_line(x, x, y2, y3, 1, 1, p1_color);
473  x = xsign * (half_shift * 0.5 + dsk);
474  draw_line(x, x, y1, y2, 1, 1, p1_color);
475  draw_line(x, x, y2, y3, 1, 1, p2_color);
476  }
477  // Make a BOX around ROC 0
478  x1 = xsign * (half_shift * 0.5 + dsk - 1 / 8.);
479  x2 = xsign * (half_shift * 0.5 + dsk);
480  float y1_p1 = ysign * (half_shift * 0.5 - 0.5 + bld + sign * 0.25);
481  float y2_p1 = ysign * (half_shift * 0.5 - 0.5 + bld + sign * 0.25 + xsign * ysign * 0.25);
482  draw_line(x1, x2, y1_p1, y1_p1, 1);
483  //draw_line(x1, x2, y2_p1, y2_p1, 1);
484  draw_line(x1, x1, y1_p1, y2_p1, 1);
485  //draw_line(x2, x2, y1_p1, y2_p1, 1);
486  float y1_p2 = ysign * (half_shift * 0.5 - 0.5 + bld - sign * 0.25);
487  float y2_p2 = ysign * (half_shift * 0.5 - 0.5 + bld - sign * 0.25 - xsign * ysign * 0.25);
488  draw_line(x1, x2, y1_p2, y1_p2, 1);
489  //draw_line(x1, x2, y2_p2, y2_p2, 1);
490  draw_line(x1, x1, y1_p2, y2_p2, 1);
491  //draw_line(x2, x2, y1_p2, y2_p2, 1);
492  }
493  }
494  }
495  }
496  }
497  // Special shifted "rebin" for Phase 0
498  // Y axis should always have at least half-roc granularity because
499  // there are half-ROC size shifts implemented in the coordinates
500  // To remove this and show full ROC granularity
501  // We merge bin contents in each pair of bins corresponding to one ROC
502  // TODO: make sure this works for Profiles
503  if (phase == 0 && h->GetNbinsY() == 250 && h->GetNbinsX() == 80) {
504  int nentries = h->GetEntries();
505  for (int binx = 1; binx <= 80; ++binx) {
506  double sum = 0;
507  for (int biny = 1; biny <= 250; ++biny) {
508  bool odd_nrocy = (binx - 1 < 40) != (((binx - 1) / 4) % 2);
509  if (biny % 2 == odd_nrocy)
510  sum += h->GetBinContent(binx, biny);
511  else {
512  sum += h->GetBinContent(binx, biny);
513  if (sum) {
514  h->SetBinContent(binx, biny, sum);
515  h->SetBinContent(binx, biny - 1, sum);
516  }
517  sum = 0;
518  }
519  }
520  }
521  h->SetEntries(nentries);
522  }
523  }
524  }
525 
526  /*--------------------------------------------------------------------*/
527  inline void adjustCanvasMargins(TVirtualPad* pad, float top, float bottom, float left, float right)
528  /*--------------------------------------------------------------------*/
529  {
530  if (top > 0)
531  pad->SetTopMargin(top);
532  if (bottom > 0)
533  pad->SetBottomMargin(bottom);
534  if (left > 0)
535  pad->SetLeftMargin(left);
536  if (right > 0)
537  pad->SetRightMargin(right);
538  }
539 
540  /*--------------------------------------------------------------------*/
541  inline void adjustStats(TPaveStats* stats, float X1, float Y1, float X2, float Y2)
542  /*--------------------------------------------------------------------*/
543  {
544  stats->SetX1NDC(X1); //new x start position
545  stats->SetY1NDC(Y1); //new y start position
546  stats->SetX2NDC(X2); //new x end position
547  stats->SetY2NDC(Y2); //new y end position
548  }
549 
550  /*--------------------------------------------------------------------*/
551  inline std::pair<float, float> getExtrema(TH1* h1, TH1* h2)
552  /*--------------------------------------------------------------------*/
553  {
554  float theMax(-9999.);
555  float theMin(9999.);
556  theMax = h1->GetMaximum() > h2->GetMaximum() ? h1->GetMaximum() : h2->GetMaximum();
557  theMin = h1->GetMinimum() < h2->GetMaximum() ? h1->GetMinimum() : h2->GetMinimum();
558 
559  float add_min = theMin > 0. ? -0.05 : 0.05;
560  float add_max = theMax > 0. ? 0.05 : -0.05;
561 
562  auto result = std::make_pair(theMin * (1 + add_min), theMax * (1 + add_max));
563  return result;
564  }
565 
566  /*--------------------------------------------------------------------*/
567  inline void makeNicePlotStyle(TH1* hist)
568  /*--------------------------------------------------------------------*/
569  {
570  hist->SetStats(kFALSE);
571  hist->SetLineWidth(2);
572  hist->GetXaxis()->CenterTitle(true);
573  hist->GetYaxis()->CenterTitle(true);
574  hist->GetXaxis()->SetTitleFont(42);
575  hist->GetYaxis()->SetTitleFont(42);
576  hist->GetXaxis()->SetTitleSize(0.05);
577  hist->GetYaxis()->SetTitleSize(0.05);
578  hist->GetXaxis()->SetTitleOffset(1.1);
579  hist->GetYaxis()->SetTitleOffset(1.3);
580  hist->GetXaxis()->SetLabelFont(42);
581  hist->GetYaxis()->SetLabelFont(42);
582  hist->GetYaxis()->SetLabelSize(.05);
583  hist->GetXaxis()->SetLabelSize(.05);
584 
585  if (hist->InheritsFrom(TH2::Class())) {
586  hist->GetZaxis()->SetLabelFont(42);
587  hist->GetZaxis()->SetLabelFont(42);
588  hist->GetZaxis()->SetLabelSize(.05);
589  hist->GetZaxis()->SetLabelSize(.05);
590  }
591  }
592 
593  enum DetType { t_barrel = 0, t_forward = 1, t_all = 2 };
594  const std::array<std::string, 3> DetNames = {{"Barrel", "End Caps", "Whole"}};
595 
596  enum regions {
597  BPixL1o, //0 Barrel Pixel Layer 1 outer
598  BPixL1i, //1 Barrel Pixel Layer 1 inner
599  BPixL2o, //2 Barrel Pixel Layer 2 outer
600  BPixL2i, //3 Barrel Pixel Layer 2 inner
601  BPixL3o, //4 Barrel Pixel Layer 3 outer
602  BPixL3i, //5 Barrel Pixel Layer 3 inner
603  BPixL4o, //6 Barrel Pixel Layer 4 outer
604  BPixL4i, //7 Barrel Pixel Layer 4 inner
605  FPixmL1, //8 Forward Pixel Minus side Disk 1
606  FPixmL2, //9 Forward Pixel Minus side Disk 2
607  FPixmL3, //10 Forward Pixel Minus side Disk 3
608  FPixpL1, //11 Forward Pixel Plus side Disk 1
609  FPixpL2, //12 Forward Pixel Plus side Disk 2
610  FPixpL3, //13 Forward Pixel Plus side Disk 3
611  NUM_OF_REGIONS //14 -- default
612  };
613 
614  /*--------------------------------------------------------------------*/
616  /*--------------------------------------------------------------------*/
617  {
618  switch (e) {
619  case SiPixelPI::BPixL1o:
620  return "BPix L1/o";
621  case SiPixelPI::BPixL1i:
622  return "BPix L1/i";
623  case SiPixelPI::BPixL2o:
624  return "BPix L2/o";
625  case SiPixelPI::BPixL2i:
626  return "BPix L2/i";
627  case SiPixelPI::BPixL3o:
628  return "BPix L3/o";
629  case SiPixelPI::BPixL3i:
630  return "BPix L3/i";
631  case SiPixelPI::BPixL4o:
632  return "BPix L4/o";
633  case SiPixelPI::BPixL4i:
634  return "BPix L4/i";
635  case SiPixelPI::FPixmL1:
636  return "FPix- D1";
637  case SiPixelPI::FPixmL2:
638  return "FPix- D2";
639  case SiPixelPI::FPixmL3:
640  return "FPix- D3";
641  case SiPixelPI::FPixpL1:
642  return "FPix+ D1";
643  case SiPixelPI::FPixpL2:
644  return "FPix+ D2";
645  case SiPixelPI::FPixpL3:
646  return "FPix+ D3";
647  default:
648  edm::LogWarning("LogicError") << "Unknown partition: " << e;
649  return "";
650  }
651  }
652 
653  /*--------------------------------------------------------------------*/
654  inline bool isBPixOuterLadder(const DetId& detid, const TrackerTopology& tTopo, bool isPhase0)
655  /*--------------------------------------------------------------------*/
656  {
657  // Using TrackerTopology
658  // Ladders have a staggered structure
659  // Non-flipped ladders are on the outer radius
660  // Phase 0: Outer ladders are odd for layer 1,3 and even for layer 2
661  // Phase 1: Outer ladders are odd for layer 1,2,3 and even for layer 4
662  bool isOuter = false;
663  int layer = tTopo.pxbLayer(detid.rawId());
664  bool odd_ladder = tTopo.pxbLadder(detid.rawId()) % 2;
665  if (isPhase0) {
666  if (layer == 2)
667  isOuter = !odd_ladder;
668  else
669  isOuter = odd_ladder;
670  } else {
671  if (layer == 4)
672  isOuter = !odd_ladder;
673  else
674  isOuter = odd_ladder;
675  }
676  return isOuter;
677  }
678 
679  // ancillary struct to manage the topology
680  // info in a more compact way
681 
682  struct topolInfo {
683  private:
684  uint32_t m_rawid;
686  int m_layer;
687  int m_side;
688  int m_ring;
691 
692  public:
693  void init();
694  void fillGeometryInfo(const DetId& detId, const TrackerTopology& tTopo, const SiPixelPI::phase& ph);
696  bool sanityCheck();
697  int subDetId() { return m_subdetid; }
698  int layer() { return m_layer; }
699  int side() { return m_side; }
700  int ring() { return m_ring; }
701  bool isInternal() { return m_isInternal; }
702 
703  void printAll(std::stringstream& ss) const;
704  virtual ~topolInfo() {}
705  };
706 
707  /*--------------------------------------------------------------------*/
708  inline void topolInfo::printAll(std::stringstream& ss) const
709  /*--------------------------------------------------------------------*/
710  {
711  ss << " detId: " << m_rawid << " subdetid: " << m_subdetid << " layer: " << m_layer << " side: " << m_side
712  << " ring: " << m_ring << " isInternal: " << m_isInternal;
713  }
714 
715  /*--------------------------------------------------------------------*/
716  inline void topolInfo::init()
717  /*--------------------------------------------------------------------*/
718  {
720  m_rawid = 0;
721  m_subdetid = -1;
722  m_layer = -1;
723  m_side = -1;
724  m_ring = -1;
725  m_isInternal = false;
726  };
727 
728  /*--------------------------------------------------------------------*/
730  /*--------------------------------------------------------------------*/
731  {
732  if (m_layer == 0 || (m_subdetid == 1 && m_layer > 4) || (m_subdetid == 2 && m_layer > 3)) {
733  return false;
734  } else {
735  return true;
736  }
737  }
738  /*--------------------------------------------------------------------*/
739  inline void topolInfo::fillGeometryInfo(const DetId& detId, const TrackerTopology& tTopo, const SiPixelPI::phase& ph)
740  /*--------------------------------------------------------------------*/
741  {
742  // set the phase
743  m_Phase = ph;
744  unsigned int subdetId = static_cast<unsigned int>(detId.subdetId());
745 
746  m_rawid = detId.rawId();
747  m_subdetid = subdetId;
748  if (subdetId == PixelSubdetector::PixelBarrel) {
749  m_layer = tTopo.pxbLayer(detId.rawId());
751  } else if (subdetId == PixelSubdetector::PixelEndcap) {
752  m_layer = tTopo.pxfDisk(detId.rawId());
753  m_side = tTopo.pxfSide(detId.rawId());
754  } else
755  edm::LogWarning("LogicError") << "Unknown subdetid: " << subdetId;
756  }
757 
758  // ------------ method to assign a partition based on the topology struct info ---------------
759 
760  /*--------------------------------------------------------------------*/
762  /*--------------------------------------------------------------------*/
763  {
765 
766  if (m_Phase == SiPixelPI::undefined) {
767  throw cms::Exception("LogicError") << "Cannot call filterThePartition BEFORE filling the geometry info!";
768  }
769 
770  // BPix
771  if (m_subdetid == 1) {
772  switch (m_layer) {
773  case 1:
775  break;
776  case 2:
778  break;
779  case 3:
781  break;
782  case 4:
784  break;
785  default:
786  edm::LogWarning("LogicError") << "Unknow BPix layer: " << m_layer;
787  break;
788  }
789  // FPix
790  } else if (m_subdetid == 2) {
791  switch (m_layer) {
792  case 1:
794  break;
795  case 2:
797  break;
798  case 3:
800  break;
801  default:
803  // warning message only if the phase2 is < 2
804  edm::LogWarning("LogicError") << "Unknow FPix disk: " << m_layer;
805  }
806  break;
807  }
808  }
809  return ret;
810  }
811 
812  /*--------------------------------------------------------------------*/
813  inline void displayNotSupported(TCanvas& canv, const unsigned int size)
814  /*--------------------------------------------------------------------*/
815  {
816  std::string phase = (size < SiPixelPI::phase1size) ? "Phase-0" : "Phase-2";
817  canv.cd();
818  TLatex t2;
819  t2.SetTextAlign(21);
820  t2.SetTextSize(0.1);
821  t2.SetTextAngle(45);
822  t2.SetTextColor(kRed);
823  if (size != SiPixelPI::mismatched) {
824  t2.DrawLatexNDC(0.6, 0.50, Form("%s NOT SUPPORTED!", phase.c_str()));
825  } else {
826  t2.DrawLatexNDC(0.6, 0.50, "MISMATCHED PAYLOAD SIZE!");
827  }
828  }
829 
830  /*--------------------------------------------------------------------*/
831  template <typename T>
832  std::pair<T, T> findMinMaxInMap(const std::map<unsigned int, T>& theMap)
833  /*--------------------------------------------------------------------*/
834  {
835  using pairtype = std::pair<unsigned int, T>;
836  auto max = *std::max_element(
837  theMap.begin(), theMap.end(), [](const pairtype& p1, const pairtype& p2) { return p1.second < p2.second; });
838  auto min = *std::min_element(
839  theMap.begin(), theMap.end(), [](const pairtype& p1, const pairtype& p2) { return p1.second < p2.second; });
840  return std::make_pair(min.second, max.second);
841  }
842 
843  /*--------------------------------------------------------------------*/
844  inline bool checkAnswerOK(std::string& answer, bool& result)
845  /*--------------------------------------------------------------------*/
846  {
847  std::transform(answer.begin(), answer.end(), answer.begin(), [](unsigned char x) { return ::tolower(x); });
848 
849  bool answer_valid = (answer == "y") || (answer == "n") || (answer == "yes") || (answer == "no") ||
850  (answer == "true") || (answer == "false") || (answer == "1") || (answer == "0");
851 
852  result = answer_valid && (answer[0] == 'y' || answer[0] == 't' || answer[0] == '1');
853  return answer_valid;
854  }
855 }; // namespace SiPixelPI
856 #endif
unsigned int pxbLayer(const DetId &id) const
int ringName() const
ring Id
const void print(std::stringstream &ss)
PhaseInfo(const phase &thePhase)
int signed_blade_panel(const DetId &detid, const TrackerTopology &tTopo_, bool phase_)
ret
prodAgent to be discontinued
static const unsigned int phase0size
int bladeName() const
blade id
int moduleName() const
module id (index in z)
unsigned int pxfModule(const DetId &id) const
int quadrant(const DetId &detid, const TrackerTopology *tTopo_, bool phase_)
const Time_t kLowMask(0xFFFFFFFF)
void fillGeometryInfo(const DetId &detId, const TrackerTopology &tTopo, const SiPixelPI::phase &ph)
unsigned int pxbLadder(const DetId &id) const
Log< level::Error, false > LogError
std::pair< float, float > getExtrema(TH1 *h1, TH1 *h2)
bool checkAnswerOK(std::string &answer, bool &result)
static std::string to_string(const XMLCh *ch)
const SiPixelPI::phase phase() const
int signed_ladder(const DetId &detid, const TrackerTopology &tTopo_, bool phase_)
unsigned long long Time_t
Definition: Time.h:14
void printAll(std::stringstream &ss) const
std::pair< T, T > findMinMaxInMap(const std::map< unsigned int, T > &theMap)
int signed_disk(const DetId &detid, const TrackerTopology &tTopo_, bool phase_)
std::pair< unsigned int, unsigned int > unpack(cond::Time_t since)
Definition: style.py:1
int signed_blade(const DetId &detid, const TrackerTopology &tTopo_, bool phase_)
HalfCylinder halfCylinder() const
static const unsigned int phase1size
unsigned int pxfDisk(const DetId &id) const
SiPixelPI::regions filterThePartition()
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
void draw_line(double x1, double x2, double y1, double y2, int width=2, int style=1, int color=1)
void adjustCanvasMargins(TVirtualPad *pad, float top, float bottom, float left, float right)
std::ostream & operator<<(std::ostream &os, PhaseInfo phInfo)
void adjustStats(TPaveStats *stats, float X1, float Y1, float X2, float Y2)
Shell shell() const
unsigned int pxfPanel(const DetId &id) const
Log< level::Info, false > LogInfo
int ring(const DetId &detid, const TrackerTopology &tTopo_, bool phase_)
Definition: DetId.h:17
const std::array< std::string, 3 > DetNames
unsigned int pxfSide(const DetId &id) const
void makeNicePlotStyle(TH1 *hist)
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
int signed_module(const DetId &detid, const TrackerTopology &tTopo_, bool phase_)
void displayNotSupported(TCanvas &canv, const unsigned int size)
const bool isComparedWithPhase2(const PhaseInfo &theOtherPhase) const
bool isBPixOuterLadder(const DetId &detid, const TrackerTopology &tTopo, bool isPhase0)
float x
int ladderName() const
ladder id (index in phi)
std::string getStringFromRegionEnum(SiPixelPI::regions e)
void dress_occup_plot(TCanvas &canv, TH2 *h, int lay, int ring=0, int phase=0, bool half_shift=true, bool mark_zero=true, bool standard_palette=true)
Log< level::Warning, false > LogWarning
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
T mod(const T &a, const T &b)
Definition: ecalDccMap.h:4
static const unsigned int mismatched
static const unsigned int phase2size
const bool isPhase1Comparison(const PhaseInfo &theOtherPhase) const
unsigned transform(const HcalDetId &id, unsigned transformCode)