CMS 3D CMS Logo

PixelRegionContainers.h
Go to the documentation of this file.
1 #ifndef CONDCORE_SIPIXELPLUGINS_PIXELREGIONCONTAINERS_H
2 #define CONDCORE_SIPIXELPLUGINS_PIXELREGIONCONTAINERS_H
3 
11 #include <boost/range/adaptor/indexed.hpp>
12 #include <cstdlib>
13 #include "TH1.h"
14 #include "TLatex.h"
15 
16 namespace PixelRegions {
17 
18  inline std::string itoa(int i) {
19  char temp[20];
20  sprintf(temp, "%d", i);
21  return ((std::string)temp);
22  }
23 
24  // "PixelId"
25  // BPix: 1000*(subdetId=1) + 100*(layer=1,2,3,4)
26  // FPix: 1000*(subdetId=2) + 100*(side=1,2) + 10*(disk=1,2,3) + 1*(ring=1,2)
27  // Phase-2 FPix: 1000*(subdetId=3) + 10*(disk=1...12)
28 
29  // clang-format off
30  enum PixelId {
31  // BPix
32  L1 = 1100, L2 = 1200, L3 = 1300, L4 = 1400, // common to phase1 and phase2
33  // FPix minus
34  Rm1l = 2111, Rm1u = 2112, Rm2l = 2121, Rm2u = 2122, Rm3l = 2131, Rm3u = 2132, // phase1-only
35  // FPix plus
36  Rp1l = 2211, Rp1u = 2212, Rp2l = 2221, Rp2u = 2222, Rp3l = 2231, Rp3u = 2232, // phase1-only
37  // Phase-2 endcaps
38  Ph2EmR1 = 3101, Ph2EmR2 = 3102, Ph2EmR3 = 3103, Ph2EmR4 = 3104, // phase-2 EPix
39  Ph2EpR1 = 3201, Ph2EpR2 = 3202, Ph2EpR3 = 3203, Ph2EpR4 = 3204,
40  Ph2FmR1 = 4101, Ph2FmR2 = 4102, Ph2FmR3 = 4103, Ph2FmR4 = 4104, Ph2FmR5 = 4105, //phase-2 FPix
41  Ph2FpR1 = 4201, Ph2FpR2 = 4202, Ph2FpR3 = 4203, Ph2FpR4 = 4204, Ph2FpR5 = 4205,
42  End = 99999
43  };
44 
45  const std::vector<PixelId> PixelIDs = {
46  // BPix
48  // Phase-1 FPix minus
50  // Phase-1 FPix plus
52  // Phase-2 endcaps
58  };
59 
60  // clang-format on
61 
62  const std::vector<std::string> IDlabels = {"Barrel Pixel L1", //1
63  "Barrel Pixel L2", //2
64  "Barrel Pixel L3", //3
65  "Barrel Pixel L4", //4
66  "FPIX(-) Disk 1 inner ring", //5
67  "FPIX(-) Disk 1 outer ring", //6
68  "FPIX(-) Disk 2 inner ring", //7
69  "FPIX(-) Disk 2 outer ring", //8
70  "FPIX(-) Disk 3 inner ring", //9
71  "FPIX(-) Disk 3 outer ring", //10
72  "FPIX(+) Disk 1 inner ring", //11
73  "FPIX(+) Disk 1 outer ring", //12
74  "FPIX(+) Disk 2 inner ring", //13
75  "FPIX(+) Disk 2 outer ring", //14
76  "FPIX(+) Disk 3 inner ring", //15
77  "FPIX(+) Disk 3 outer ring", //16
78  "EPix(-) Ring 1", //17
79  "EPix(-) Ring 2", //18
80  "EPix(-) Ring 3", //19
81  "EPix(-) Ring 4", //20
82  "EPix(+) Ring 1", //21
83  "EPix(+) Ring 2", //22
84  "EPix(+) Ring 3", //23
85  "EPix(+) Ring 4", //24
86  "FPix(-) Ring 1", //25
87  "FPix(-) Ring 2", //26
88  "FPix(-) Ring 3", //27
89  "FPix(-) Ring 4", //28
90  "FPix(-) Ring 5", //29
91  "FPix(+) Ring 1", //30
92  "FPix(+) Ring 2", //31
93  "FPix(+) Ring 3", //32
94  "FPix(+) Ring 4", //33
95  "FPix(+) Ring 5", //34
96  "END"}; //35
97 
98  //============================================================================
99  [[maybe_unused]] static const std::vector<std::string> getIDLabels(const SiPixelPI::phase& ph, bool isBarrel) {
100  std::vector<std::string> out;
101  for (const auto& label : IDlabels) {
102  if (isBarrel) {
103  if (label.find("Barrel") != std::string::npos) {
104  out.push_back(label);
105  }
106  } else {
107  if (ph == SiPixelPI::phase::two) {
108  if (label.find("Ring") != std::string::npos) {
109  out.push_back(label);
110  }
111  } else {
112  if (label.find("ring") != std::string::npos) {
113  out.push_back(label);
114  }
115  }
116  }
117  }
118  return out;
119  }
120 
121  //============================================================================
122  static const PixelId calculateBPixID(const unsigned int layer) {
123  // BPix: 1000*(subdetId=1) + 100*(layer=1,2,3,4)
124  PixelId bpixLayer = static_cast<PixelId>(1000 + 100 * layer);
125  return bpixLayer;
126  }
127 
128  //============================================================================
130  const unsigned int side,
131  const unsigned int disk,
132  const unsigned int ring) {
133  // FPix: 1000*(subdetId=2) + 100*(side=1,2) + 10*(disk=1,2,3) + 1*(ring=1,2)
134  using namespace SiPixelPI;
135  unsigned int prefix(2000);
136  unsigned int disk_(0); // if that's phase-2 set the disk n. to zero
137  if (ph > phase::one) {
138  if (disk < 9) {
139  prefix += 1000; // if that's EPix id starts with 3k
140  } else {
141  prefix += 2000; // if that's FPix id starts with 4k
142  }
143  } else {
144  disk_ = disk; // if that's not phase-2 set the disk to real disk n.
145  }
146  PixelId fpixRing = static_cast<PixelId>(prefix + 100 * side + 10 * disk_ + ring);
147  return fpixRing;
148  }
149 
150  //============================================================================
151  static const PixelId detIdToPixelId(const unsigned int detid,
152  const TrackerTopology* trackTopo,
153  const SiPixelPI::phase& ph) {
154  using namespace SiPixelPI;
155  DetId detId = DetId(detid);
156  unsigned int subid = detId.subdetId();
157  unsigned int pixid = 0;
158  if (subid == PixelSubdetector::PixelBarrel) {
159  int layer = trackTopo->pxbLayer(detId); // 1, 2, 3, 4
160  pixid = calculateBPixID(layer);
161  } else if (subid == PixelSubdetector::PixelEndcap) {
162  int side = trackTopo->pxfSide(detId); // 1 (-z), 2 for (+z)
163  int disk = trackTopo->pxfDisk(detId); // 1, 2, 3
164  int ring(0);
165  switch (ph) {
166  case phase::zero:
167  ring = SiPixelPI::ring(detid, *trackTopo, false); // 1 (lower), 2 (upper)
168  break;
169  case phase::one:
170  ring = SiPixelPI::ring(detid, *trackTopo, true); // 1 (lower), 2 (upper)
171  break;
172  case phase::two:
173  // I know, that's funny but go look at:
174  // https://github.com/cms-sw/cmssw/blob/master/Geometry/TrackerNumberingBuilder/README.md
175  ring = trackTopo->pxfBlade(detId);
176  break;
177  default:
178  throw cms::Exception("LogicalError") << " there is not such phase as " << ph;
179  }
180  pixid = calculateFPixID(ph, side, disk, ring);
181  }
182  PixelId pixID = static_cast<PixelId>(pixid);
183  return pixID;
184  }
185 
186  //============================================================================
187  [[maybe_unused]] static const std::vector<uint32_t> attachedDets(const PixelRegions::PixelId theId,
188  const TrackerTopology* trackTopo,
189  const SiPixelPI::phase& ph) {
190  using namespace SiPixelPI;
191  std::vector<uint32_t> out = {};
192  edm::FileInPath m_fp;
193 
194  switch (ph) {
195  case phase::zero:
196  // phase-1 skimmed geometry from release
197  m_fp = edm::FileInPath("CalibTracker/SiPixelESProducers/data/PixelSkimmedGeometry.txt");
198  break;
199  case phase::one:
200  // phase-0 skimmed geometry from release
201  m_fp = edm::FileInPath("SLHCUpgradeSimulations/Geometry/data/PhaseI/PixelSkimmedGeometry_phase1.txt");
202  break;
203  case phase::two:
204  m_fp = edm::FileInPath("SLHCUpgradeSimulations/Geometry/data/PhaseII/Tilted/PixelSkimmedGeometryT14.txt");
205  break;
206  default:
207  throw cms::Exception("LogicalError") << " there is not such phase as " << ph;
208  }
209 
210  SiPixelDetInfoFileReader pxlreader(m_fp.fullPath());
211  const std::vector<uint32_t>& pxldetids = pxlreader.getAllDetIds();
212  for (const auto& d : pxldetids) {
213  auto ID = detIdToPixelId(d, trackTopo, ph);
214  if (ID == theId) {
215  out.push_back(d);
216  }
217  }
218 
219  // in case no DetIds are assigned, fill with UINT32_MAX
220  // in order to be able to tell a part the default case
221  // in SiPixelGainCalibHelper::fillTheHitsto
222  if (out.empty()) {
223  out.push_back(0xFFFFFFFF);
224  }
225 
226  COUT << "ID:" << theId << " ";
227  for (const auto& entry : out) {
228  COUT << entry << ",";
229  }
230  COUT << std::endl;
231 
232  return out;
233  }
234 
235  /*--------------------------------------------------------------------
236  / Ancillary class to build pixel phase0/phase1 plots per region
237  /--------------------------------------------------------------------*/
239  public:
241  : m_trackerTopo(t_topo), m_Phase(ph) {
242  // set log scale by default to false
243  m_isLog = false;
244  if (m_trackerTopo) {
245  m_isTrackerTopologySet = true;
246  } else {
247  m_isTrackerTopologySet = false;
248  }
249  }
250 
252 
253  //============================================================================
254  void setTheTopo(const TrackerTopology* t_topo) {
255  m_trackerTopo = t_topo;
256  m_isTrackerTopologySet = true;
257  }
258 
259  //============================================================================
261 
262  //============================================================================
263  void bookAll(std::string title_label,
264  std::string x_label,
265  std::string y_label,
266  const int nbins,
267  const float xmin,
268  const float xmax) {
269  using namespace SiPixelPI;
270  for (const auto& pixelId : PixelIDs | boost::adaptors::indexed(0)) {
271  if (m_Phase == phase::two && // if that's phase-2
272  pixelId.value() > PixelId::L4 && // if it's end-cap
273  pixelId.value() < PixelId::Ph2EmR1 // it's a phase-1 ring
274  ) {
275  continue;
276  }
277 
278  m_theMap[pixelId.value()] = std::make_shared<TH1F>((title_label + itoa(pixelId.value())).c_str(),
279  Form("%s %s;%s;%s",
280  (IDlabels.at(pixelId.index())).c_str(),
281  title_label.c_str(),
282  x_label.c_str(),
283  y_label.c_str()),
284  nbins,
285  xmin,
286  xmax);
287  }
288  }
289 
290  //============================================================================
291  void fill(const unsigned int detid, const float value) {
292  // first check that the topology is set
294 
295  // convert from detid to pixelid
297  if (m_theMap.find(myId) != m_theMap.end()) {
298  m_theMap[myId]->Fill(value);
299  } else {
300  edm::LogError("PixelRegionContainers")
301  << detid << " :=> " << myId << " is not a recongnized PixelId enumerator! \n"
302  << m_trackerTopo->print(detid);
303  }
304  }
305 
306  //============================================================================
307  void draw(TCanvas& canv, bool isBarrel, const char* option = "bar2", bool isPhase1Comparison = false) {
308  using namespace SiPixelPI;
309  if (isBarrel) {
310  if (canv.GetListOfPrimitives()->GetSize() == 0) {
311  // divide only if it was not already divided
312  canv.Divide(2, 2);
313  }
314  for (int j = 1; j <= 4; j++) {
315  if (!m_isLog) {
316  canv.cd(j);
317  } else {
318  canv.cd(j)->SetLogy();
319  }
320  if ((j == 4) && (m_Phase < 1) && !isPhase1Comparison) {
321  m_theMap.at(PixelIDs[j - 1])->Draw("AXIS");
322  TLatex t2;
323  t2.SetTextAlign(22);
324  t2.SetTextSize(0.1);
325  t2.SetTextAngle(45);
326  t2.SetTextFont(61);
327  t2.SetTextColor(kBlack);
328  t2.DrawLatexNDC(0.5, 0.5, "Not in Phase-0!");
329  } else {
330  m_theMap.at(PixelIDs[j - 1])->Draw(option);
331  }
332  }
333  } else { // forward
334  // pattern where to insert the plots in the canvas
335  // clang-format off
336  const std::array<int, 18> phase2Pattern = {{1, 2, 3, 4,
337  6, 7, 8, 9,
338  11, 12, 13, 14, 15,
339  16, 17, 18, 19, 20}};
340  // clang-format on
341 
342  if (canv.GetListOfPrimitives()->GetSize() == 0) {
343  canv.Divide(m_Phase == phase::two ? 5 : 4, m_Phase == phase::two ? 4 : 3);
344  }
345  const int maxPads = (m_Phase == phase::two) ? 18 : 12;
346  for (int j = 1; j <= maxPads; j++) {
347  unsigned int mapIndex = m_Phase == 2 ? j + 15 : j + 3;
348  if (!m_isLog) {
349  canv.cd((m_Phase == phase::two) ? phase2Pattern[j - 1] : j);
350  //canv.cd(j);
351  } else {
352  canv.cd(j)->SetLogy();
353  }
354  if ((j % 6 == 5 || j % 6 == 0) && (m_Phase < 1) && !isPhase1Comparison) {
355  m_theMap.at(PixelIDs[mapIndex])->Draw("AXIS");
356  TLatex t2;
357  t2.SetTextAlign(22);
358  t2.SetTextSize(0.1);
359  t2.SetTextAngle(45);
360  t2.SetTextFont(61);
361  t2.SetTextColor(kBlack);
362  t2.DrawLatexNDC(0.5, 0.5, "Not in Phase-0!");
363  } else {
364  m_theMap.at(PixelIDs[mapIndex])->Draw(option);
365  }
366  }
367  }
368  }
369 
370  //============================================================================
371  void beautify(const int linecolor = kBlack, const int fillcolor = kRed) {
372  for (const auto& plot : m_theMap) {
373  plot.second->SetTitle("");
374  if (!m_isLog) {
375  plot.second->GetYaxis()->SetRangeUser(0., plot.second->GetMaximum() * 1.30);
376  } else {
377  plot.second->GetYaxis()->SetRangeUser(0.1, plot.second->GetMaximum() * 100.);
378  }
379  plot.second->SetLineColor(linecolor);
380  if (fillcolor > 0) {
381  plot.second->SetFillColor(fillcolor);
382  }
383  plot.second->SetMarkerStyle(20);
384  plot.second->SetMarkerSize(1);
385  SiPixelPI::makeNicePlotStyle(plot.second.get());
386  plot.second->SetStats(true);
387  }
388  }
389 
390  //============================================================================
391  void setLogScale() { m_isLog = true; }
392 
393  //============================================================================
394  void stats(int index = 0) {
395  for (const auto& plot : m_theMap) {
396  TPaveStats* st = (TPaveStats*)plot.second->FindObject("stats");
397  if (st) {
398  st->SetTextSize(0.03);
399  st->SetLineColor(10);
400  if (plot.second->GetFillColor() != 0) {
401  st->SetTextColor(plot.second->GetFillColor());
402  } else {
403  st->SetTextColor(plot.second->GetLineColor());
404  }
405  SiPixelPI::adjustStats(st, 0.13, 0.85 - index * 0.08, 0.36, 0.93 - index * 0.08);
406  }
407  }
408  }
409 
410  //============================================================================
411  std::shared_ptr<TH1F>& getHistoFromMap(const PixelRegions::PixelId& theId) {
412  auto it = m_theMap.find(theId);
413  if (it != m_theMap.end()) {
414  return it->second;
415  } else {
416  throw cms::Exception("LogicError")
417  << "PixelRegionContainer::getHistoFromMap(): No histogram is available for PixelId: " << theId << "\n";
418  }
419  }
420 
421  //============================================================================
422  void rescaleMax(PixelRegionContainers& the2ndContainer) {
423  for (const auto& plot : m_theMap) {
424  auto thePixId = plot.first;
425  auto extrema = SiPixelPI::getExtrema((plot.second).get(), the2ndContainer.getHistoFromMap(thePixId).get());
426  plot.second->GetYaxis()->SetRangeUser(extrema.first, extrema.second);
427  the2ndContainer.getHistoFromMap(thePixId)->GetYaxis()->SetRangeUser(extrema.first, extrema.second);
428  }
429  }
430 
431  private:
435  std::map<PixelId, std::shared_ptr<TH1F>> m_theMap;
436  int m_nbins;
437  float m_xim, m_xmax;
438  bool m_isLog;
439  };
440 } // namespace PixelRegions
441 #endif
std::map< PixelId, std::shared_ptr< TH1F > > m_theMap
unsigned int pxbLayer(const DetId &id) const
std::string fullPath() const
Definition: FileInPath.cc:161
unsigned int pxfBlade(const DetId &id) const
const std::vector< std::string > IDlabels
void rescaleMax(PixelRegionContainers &the2ndContainer)
std::string print(DetId detid) const
uint32_t ID
Definition: Definitions.h:24
void beautify(const int linecolor=kBlack, const int fillcolor=kRed)
Log< level::Error, false > LogError
std::pair< float, float > getExtrema(TH1 *h1, TH1 *h2)
assert(be >=bs)
char const * label
std::string itoa(int i)
std::shared_ptr< TH1F > & getHistoFromMap(const PixelRegions::PixelId &theId)
void draw(TCanvas &canv, bool isBarrel, const char *option="bar2", bool isPhase1Comparison=false)
void setTheTopo(const TrackerTopology *t_topo)
unsigned int pxfDisk(const DetId &id) const
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
const std::vector< PixelId > PixelIDs
Definition: value.py:1
#define COUT
void fill(const unsigned int detid, const float value)
d
Definition: ztail.py:151
void adjustStats(TPaveStats *stats, float X1, float Y1, float X2, float Y2)
void bookAll(std::string title_label, std::string x_label, std::string y_label, const int nbins, const float xmin, const float xmax)
int ring(const DetId &detid, const TrackerTopology &tTopo_, bool phase_)
Definition: DetId.h:17
unsigned int pxfSide(const DetId &id) const
void makeNicePlotStyle(TH1 *hist)
static const std::vector< std::string > getIDLabels(const SiPixelPI::phase &ph, bool isBarrel)
const std::vector< uint32_t > & getAllDetIds() const
static const PixelId calculateFPixID(const SiPixelPI::phase &ph, const unsigned int side, const unsigned int disk, const unsigned int ring)
static const PixelId detIdToPixelId(const unsigned int detid, const TrackerTopology *trackTopo, const SiPixelPI::phase &ph)
static const char disk_[]
static const std::vector< uint32_t > attachedDets(const PixelRegions::PixelId theId, const TrackerTopology *trackTopo, const SiPixelPI::phase &ph)
PixelRegionContainers(const TrackerTopology *t_topo, const SiPixelPI::phase &ph)
static const PixelId calculateBPixID(const unsigned int layer)