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 
10 #include <boost/range/adaptor/indexed.hpp>
11 #include "TH1.h"
12 #include <cstdlib>
13 
14 namespace PixelRegions {
15 
17  char temp[20];
18  sprintf(temp, "%d", i);
19  return ((std::string)temp);
20  }
21 
22  // "PixelId"
23  // BPix: 1000*(subdetId=1) + 100*(layer=1,2,3,4)
24  // FPix: 1000*(subdetId=2) + 100*(side=1,2) + 10*(disk=1,2,3) + 1*(ring=1,2)
25 
26  enum PixelId {
27  L1 = 1100, // BPix
28  L2 = 1200,
29  L3 = 1300,
30  L4 = 1400,
31  Rm1l = 2111, // FPix minus
32  Rm1u = 2112,
33  Rm2l = 2121,
34  Rm2u = 2122,
35  Rm3l = 2131,
36  Rm3u = 2132,
37  Rp1l = 2211, // FPix plus
38  Rp1u = 2212,
39  Rp2l = 2221,
40  Rp2u = 2222,
41  Rp3l = 2231,
42  Rp3u = 2232,
43  End = 99999
44  };
45 
46  const std::vector<PixelId> PixelIDs = {PixelId::L1, // BPix
50  PixelId::Rm1l, // FPix minus
56  PixelId::Rp1l, // FPix plus
62  PixelId::End};
63 
64  const std::vector<std::string> IDlabels = {"Barrel Pixel L1",
65  "Barrel Pixel L2",
66  "Barrel Pixel L3",
67  "Barrel Pixel L4",
68  "FPIX(-) Disk 1 inner ring",
69  "FPIX(-) Disk 1 outer ring",
70  "FPIX(-) Disk 2 inner ring",
71  "FPIX(-) Disk 2 outer ring",
72  "FPIX(-) Disk 3 inner ring",
73  "FPIX(-) Disk 3 outer ring",
74  "FPIX(+) Disk 1 inner ring",
75  "FPIX(+) Disk 1 outer ring",
76  "FPIX(+) Disk 2 inner ring",
77  "FPIX(+) Disk 2 outer ring",
78  "FPIX(+) Disk 3 inner ring",
79  "FPIX(+) Disk 3 outer ring",
80  "END"};
81 
82  //============================================================================
83  static const PixelId calculateBPixID(const unsigned int layer) {
84  // BPix: 1000*(subdetId=1) + 100*(layer=1,2,3,4)
85  PixelId bpixLayer = static_cast<PixelId>(1000 + 100 * layer);
86  return bpixLayer;
87  }
88 
89  //============================================================================
90  static const PixelId calculateFPixID(const unsigned int side, const unsigned int disk, const unsigned int ring) {
91  // FPix: 1000*(subdetId=2) + 100*(side=1,2) + 10*(disk=1,2,3) + 1*(ring=1,2)
92  PixelId fpixRing = static_cast<PixelId>(2000 + 100 * side + 10 * disk + ring);
93  return fpixRing;
94  }
95 
96  //============================================================================
97  static const PixelId detIdToPixelId(const unsigned int detid, const TrackerTopology* trackTopo, const bool phase1) {
98  DetId detId = DetId(detid);
99  unsigned int subid = detId.subdetId();
100  unsigned int pixid = 0;
101  if (subid == PixelSubdetector::PixelBarrel) {
102  PixelBarrelName bpix(detId, trackTopo, phase1);
103  int layer = bpix.layerName();
104  pixid = calculateBPixID(layer);
105  } else if (subid == PixelSubdetector::PixelEndcap) {
106  PixelEndcapName fpix(detId, trackTopo, phase1);
107  int side = trackTopo->pxfSide(detId); // 1 (-z), 2 for (+z)
108  int disk = fpix.diskName(); //trackTopo->pxfDisk(detId); // 1, 2, 3
109  // This works only in case of the Phase-1 detector
110  //int ring = fpix.ringName(); // 1 (lower), 2 (upper)
111  int ring = SiPixelPI::ring(detid, *trackTopo, phase1); // 1 (lower), 2 (upper)
112  pixid = calculateFPixID(side, disk, ring);
113  }
114  PixelId pixID = static_cast<PixelId>(pixid);
115  return pixID;
116  }
117 
118  //============================================================================
119  const std::vector<uint32_t> attachedDets(const PixelRegions::PixelId theId,
120  const TrackerTopology* trackTopo,
121  const bool phase1) {
122  std::vector<uint32_t> out = {};
123  edm::FileInPath m_fp;
124  if (phase1) {
125  // phase-1 skimmed geometry from release
126  m_fp = edm::FileInPath("SLHCUpgradeSimulations/Geometry/data/PhaseI/PixelSkimmedGeometry_phase1.txt");
127  } else {
128  // phase-0 skimmed geometry from release
129  m_fp = edm::FileInPath("CalibTracker/SiPixelESProducers/data/PixelSkimmedGeometry.txt");
130  }
131 
132  SiPixelDetInfoFileReader pxlreader(m_fp.fullPath());
133  const std::vector<uint32_t>& pxldetids = pxlreader.getAllDetIds();
134  for (const auto& d : pxldetids) {
135  auto ID = detIdToPixelId(d, trackTopo, phase1);
136  if (ID == theId) {
137  out.push_back(d);
138  }
139  }
140 
141  // in case no DetIds are assigned, fill with UINT32_MAX
142  // in order to be able to tell a part the default case
143  // in SiPixelGainCalibHelper::fillTheHitsto
144  if (out.empty()) {
145  out.push_back(0xFFFFFFFF);
146  }
147 
148  COUT << "ID:" << theId << " ";
149  for (const auto& entry : out) {
150  COUT << entry << ",";
151  }
152  COUT << std::endl;
153 
154  return out;
155  }
156 
157  /*--------------------------------------------------------------------
158  / Ancillary class to build pixel phase0/phase1 plots per region
159  /--------------------------------------------------------------------*/
161  public:
162  PixelRegionContainers(const TrackerTopology* t_topo, const bool isPhase1)
163  : m_trackerTopo(t_topo), m_isPhase1(isPhase1) {
164  // set log scale by default to false
165  m_isLog = false;
166  }
167 
169 
170  //============================================================================
171  void bookAll(std::string title_label,
172  std::string x_label,
173  std::string y_label,
174  const int nbins,
175  const float xmin,
176  const float xmax) {
177  for (const auto& pixelId : PixelIDs | boost::adaptors::indexed(0)) {
178  m_theMap[pixelId.value()] = std::make_shared<TH1F>((title_label + itoa(pixelId.value())).c_str(),
179  Form("%s %s;%s;%s",
180  (IDlabels.at(pixelId.index())).c_str(),
181  title_label.c_str(),
182  x_label.c_str(),
183  y_label.c_str()),
184  nbins,
185  xmin,
186  xmax);
187  }
188  }
189 
190  //============================================================================
191  void fill(const unsigned int detid, const float value) {
192  // convert from detid to pixelid
194  if (m_theMap.find(myId) != m_theMap.end()) {
195  m_theMap[myId]->Fill(value);
196  } else {
197  edm::LogError("PixelRegionContainers")
198  << detid << " :=> " << myId << " is not a recongnized PixelId enumerator!" << std::endl;
199  }
200  }
201 
202  //============================================================================
203  void draw(TCanvas& canv, bool isBarrel, const char* option = "bar2", bool isPhase1Comparison = false) {
204  if (isBarrel) {
205  for (int j = 1; j <= 4; j++) {
206  if (!m_isLog) {
207  canv.cd(j);
208  } else {
209  canv.cd(j)->SetLogy();
210  }
211  if ((j == 4) && !m_isPhase1 && !isPhase1Comparison) {
212  m_theMap.at(PixelIDs[j - 1])->Draw("AXIS");
213  TLatex t2;
214  t2.SetTextAlign(22);
215  t2.SetTextSize(0.1);
216  t2.SetTextAngle(45);
217  t2.SetTextFont(61);
218  t2.SetTextColor(kBlack);
219  t2.DrawLatexNDC(0.5, 0.5, "Not in Phase-0!");
220  } else {
221  m_theMap.at(PixelIDs[j - 1])->Draw(option);
222  }
223  }
224  } else { // forward
225  for (int j = 1; j <= 12; j++) {
226  if (!m_isLog) {
227  canv.cd(j);
228  } else {
229  canv.cd(j)->SetLogy();
230  }
231  if ((j % 6 == 5 || j % 6 == 0) && !m_isPhase1 && !isPhase1Comparison) {
232  m_theMap.at(PixelIDs[j + 3])->Draw("AXIS");
233  TLatex t2;
234  t2.SetTextAlign(22);
235  t2.SetTextSize(0.1);
236  t2.SetTextAngle(45);
237  t2.SetTextFont(61);
238  t2.SetTextColor(kBlack);
239  t2.DrawLatexNDC(0.5, 0.5, "Not in Phase-0!");
240  } else {
241  m_theMap.at(PixelIDs[j + 3])->Draw(option);
242  }
243  }
244  }
245  }
246 
247  //============================================================================
248  void beautify(const int linecolor = kBlack, const int fillcolor = kRed) {
249  for (const auto& plot : m_theMap) {
250  plot.second->SetTitle("");
251  if (!m_isLog) {
252  plot.second->GetYaxis()->SetRangeUser(0., plot.second->GetMaximum() * 1.30);
253  } else {
254  plot.second->GetYaxis()->SetRangeUser(0.1, plot.second->GetMaximum() * 100.);
255  }
256  plot.second->SetLineColor(linecolor);
257  if (fillcolor > 0) {
258  plot.second->SetFillColor(fillcolor);
259  }
260  plot.second->SetMarkerStyle(20);
261  plot.second->SetMarkerSize(1);
262  SiPixelPI::makeNicePlotStyle(plot.second.get());
263  plot.second->SetStats(true);
264  }
265  }
266 
267  //============================================================================
268  void setLogScale() { m_isLog = true; }
269 
270  //============================================================================
271  void stats(int index = 0) {
272  for (const auto& plot : m_theMap) {
273  TPaveStats* st = (TPaveStats*)plot.second->FindObject("stats");
274  if (st) {
275  st->SetTextSize(0.03);
276  st->SetLineColor(10);
277  if (plot.second->GetFillColor() != 0) {
278  st->SetTextColor(plot.second->GetFillColor());
279  } else {
280  st->SetTextColor(plot.second->GetLineColor());
281  }
282  SiPixelPI::adjustStats(st, 0.13, 0.85 - index * 0.08, 0.36, 0.93 - index * 0.08);
283  }
284  }
285  }
286 
287  //============================================================================
288  std::shared_ptr<TH1F>& getHistoFromMap(const PixelRegions::PixelId& theId) {
289  auto it = m_theMap.find(theId);
290  if (it != m_theMap.end()) {
291  return it->second;
292  } else {
293  throw cms::Exception("PixelRegionContainer") << "No histogram is available for PixelId" << theId << "\n";
294  }
295  }
296 
297  //============================================================================
298  void rescaleMax(PixelRegionContainers& the2ndContainer) {
299  for (const auto& plot : m_theMap) {
300  auto thePixId = plot.first;
301  auto extrema = SiPixelPI::getExtrema((plot.second).get(), the2ndContainer.getHistoFromMap(thePixId).get());
302  plot.second->GetYaxis()->SetRangeUser(extrema.first, extrema.second);
303  the2ndContainer.getHistoFromMap(thePixId)->GetYaxis()->SetRangeUser(extrema.first, extrema.second);
304  }
305  }
306 
307  private:
310 
311  std::map<PixelId, std::shared_ptr<TH1F>> m_theMap;
312  int m_nbins;
313  float m_xim, m_xmax;
314  bool m_isLog;
315  };
316 } // namespace PixelRegions
317 #endif
PixelRegions::L3
Definition: PixelRegionContainers.h:29
SiPixelDetInfoFileReader::getAllDetIds
const std::vector< uint32_t > & getAllDetIds() const
Definition: SiPixelDetInfoFileReader.cc:85
PixelRegions::L1
Definition: PixelRegionContainers.h:27
RandomServiceHelper.t2
t2
Definition: RandomServiceHelper.py:257
mps_fire.i
i
Definition: mps_fire.py:355
PixelRegions::L2
Definition: PixelRegionContainers.h:28
PixelSubdetector.h
PixelBarrelName.h
PixelSubdetector::PixelEndcap
Definition: PixelSubdetector.h:11
PixelSubdetector::PixelBarrel
Definition: PixelSubdetector.h:11
PixelRegions::PixelIDs
const std::vector< PixelId > PixelIDs
Definition: PixelRegionContainers.h:46
TrackerTopology::pxfSide
unsigned int pxfSide(const DetId &id) const
Definition: TrackerTopology.h:192
SiPixelPI::ring
int ring(const DetId &detid, const TrackerTopology &tTopo_, bool phase_)
Definition: SiPixelPayloadInspectorHelper.h:93
PixelRegions::PixelRegionContainers::m_trackerTopo
const TrackerTopology * m_trackerTopo
Definition: PixelRegionContainers.h:308
PixelRegions::PixelRegionContainers::draw
void draw(TCanvas &canv, bool isBarrel, const char *option="bar2", bool isPhase1Comparison=false)
Definition: PixelRegionContainers.h:203
mps_splice.entry
entry
Definition: mps_splice.py:68
TrackerTopology
Definition: TrackerTopology.h:16
PixelRegions::Rm1l
Definition: PixelRegionContainers.h:31
pixelClusterTagInfos_cfi.isPhase1
isPhase1
Definition: pixelClusterTagInfos_cfi.py:7
PixelRegions::Rp1l
Definition: PixelRegionContainers.h:37
PixelRegions::Rp3u
Definition: PixelRegionContainers.h:42
PixelBarrelName
Definition: PixelBarrelName.h:16
PixelRegions::PixelRegionContainers::m_nbins
int m_nbins
Definition: PixelRegionContainers.h:312
PixelRegions::PixelRegionContainers::stats
void stats(int index=0)
Definition: PixelRegionContainers.h:271
PixelRegions::PixelId
PixelId
Definition: PixelRegionContainers.h:26
PixelBarrelName::layerName
int layerName() const
layer id
Definition: PixelBarrelName.h:43
groupFilesInBlocks.temp
list temp
Definition: groupFilesInBlocks.py:142
plotFactory.plot
plot
Definition: plotFactory.py:109
fileinputsource_cfi.option
option
Definition: fileinputsource_cfi.py:87
FileInPath.h
DetId
Definition: DetId.h:17
SiPixelDetInfoFileReader
Definition: SiPixelDetInfoFileReader.h:28
edm::FileInPath
Definition: FileInPath.h:64
TrackerTopology.h
PixelRegions::PixelRegionContainers::fill
void fill(const unsigned int detid, const float value)
Definition: PixelRegionContainers.h:191
PixelRegions::PixelRegionContainers::m_isPhase1
bool m_isPhase1
Definition: PixelRegionContainers.h:309
PixelRegions::Rp1u
Definition: PixelRegionContainers.h:38
PixelEndcapName
Definition: PixelEndcapName.h:16
PixelRegions::calculateFPixID
static const PixelId calculateFPixID(const unsigned int side, const unsigned int disk, const unsigned int ring)
Definition: PixelRegionContainers.h:90
PixelRegions::PixelRegionContainers::~PixelRegionContainers
~PixelRegionContainers()
Definition: PixelRegionContainers.h:168
PixelRegions::Rp2u
Definition: PixelRegionContainers.h:40
PixelRegions::PixelRegionContainers
Definition: PixelRegionContainers.h:160
PixelRegions::PixelRegionContainers::rescaleMax
void rescaleMax(PixelRegionContainers &the2ndContainer)
Definition: PixelRegionContainers.h:298
LaserClient_cfi.nbins
nbins
Definition: LaserClient_cfi.py:51
PixelRegions::PixelRegionContainers::PixelRegionContainers
PixelRegionContainers(const TrackerTopology *t_topo, const bool isPhase1)
Definition: PixelRegionContainers.h:162
PixelRegions::Rm2l
Definition: PixelRegionContainers.h:33
DetId::subdetId
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector's numbering enum)
Definition: DetId.h:48
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
PixelRegions::PixelRegionContainers::getHistoFromMap
std::shared_ptr< TH1F > & getHistoFromMap(const PixelRegions::PixelId &theId)
Definition: PixelRegionContainers.h:288
PixelRegions::Rp3l
Definition: PixelRegionContainers.h:41
PixelPluginsPhase0_cfi.isBarrel
isBarrel
Definition: PixelPluginsPhase0_cfi.py:17
PixelRegions::Rm2u
Definition: PixelRegionContainers.h:34
edm::LogError
Definition: MessageLogger.h:183
align::ID
uint32_t ID
Definition: Definitions.h:24
PixelRegions::PixelRegionContainers::m_theMap
std::map< PixelId, std::shared_ptr< TH1F > > m_theMap
Definition: PixelRegionContainers.h:311
PixelRegions::IDlabels
const std::vector< std::string > IDlabels
Definition: PixelRegionContainers.h:64
PixelRegions::Rm3u
Definition: PixelRegionContainers.h:36
PixelRegions::L4
Definition: PixelRegionContainers.h:30
PixelRegions::PixelRegionContainers::setLogScale
void setLogScale()
Definition: PixelRegionContainers.h:268
value
Definition: value.py:1
PixelRegions::Rm1u
Definition: PixelRegionContainers.h:32
PixelRegions::calculateBPixID
static const PixelId calculateBPixID(const unsigned int layer)
Definition: PixelRegionContainers.h:83
SiPixelPI::adjustStats
void adjustStats(TPaveStats *stats, float X1, float Y1, float X2, float Y2)
Definition: SiPixelPayloadInspectorHelper.h:449
SiPixelDetInfoFileReader.h
PixelRegions::Rm3l
Definition: PixelRegionContainers.h:35
COUT
#define COUT
Definition: PVValidationHelpers.h:13
PixelRegions::PixelRegionContainers::m_xmax
float m_xmax
Definition: PixelRegionContainers.h:313
PixelEndcapName::diskName
int diskName() const
disk id
Definition: PixelEndcapName.h:45
PixelRegions::detIdToPixelId
static const PixelId detIdToPixelId(const unsigned int detid, const TrackerTopology *trackTopo, const bool phase1)
Definition: PixelRegionContainers.h:97
TrackerOfflineValidation_Dqm_cff.xmax
xmax
Definition: TrackerOfflineValidation_Dqm_cff.py:11
Exception
Definition: hltDiff.cc:246
PixelRegions::PixelRegionContainers::m_xim
float m_xim
Definition: PixelRegionContainers.h:313
relativeConstraints.ring
ring
Definition: relativeConstraints.py:68
SiPixelPI::makeNicePlotStyle
void makeNicePlotStyle(TH1 *hist)
Definition: SiPixelPayloadInspectorHelper.h:475
PixelRegions::PixelRegionContainers::m_isLog
bool m_isLog
Definition: PixelRegionContainers.h:314
PixelRegions::End
Definition: PixelRegionContainers.h:43
SiPixelPI::getExtrema
std::pair< float, float > getExtrema(TH1 *h1, TH1 *h2)
Definition: SiPixelPayloadInspectorHelper.h:459
AlignmentPI::index
index
Definition: AlignmentPayloadInspectorHelper.h:46
MillePedeFileConverter_cfg.out
out
Definition: MillePedeFileConverter_cfg.py:31
PixelRegions::itoa
std::string itoa(int i)
Definition: PixelRegionContainers.h:16
ztail.d
d
Definition: ztail.py:151
PixelRegions::attachedDets
const std::vector< uint32_t > attachedDets(const PixelRegions::PixelId theId, const TrackerTopology *trackTopo, const bool phase1)
Definition: PixelRegionContainers.h:119
TrackerOfflineValidation_Dqm_cff.xmin
xmin
Definition: TrackerOfflineValidation_Dqm_cff.py:10
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
PixelRegions::Rp2l
Definition: PixelRegionContainers.h:39
PixelRegions
Definition: PixelRegionContainers.h:14
PixelRegions::PixelRegionContainers::bookAll
void bookAll(std::string title_label, std::string x_label, std::string y_label, const int nbins, const float xmin, const float xmax)
Definition: PixelRegionContainers.h:171
edm::FileInPath::fullPath
std::string fullPath() const
Definition: FileInPath.cc:163
PixelRegions::PixelRegionContainers::beautify
void beautify(const int linecolor=kBlack, const int fillcolor=kRed)
Definition: PixelRegionContainers.h:248
PixelEndcapName.h