CMS 3D CMS Logo

Phase1PixelROCMaps.h
Go to the documentation of this file.
1 #ifndef DQM_TRACKERREMAPPER_PHASE1PIXELROCMAPS_H
2 #define DQM_TRACKERREMAPPER_PHASE1PIXELROCMAPS_H
3 
4 #include "TH1.h"
5 #include "TH2.h"
6 #include "TStyle.h"
7 #include "TCanvas.h"
8 #include "TLine.h"
9 #include "TLatex.h"
10 
11 #include <array>
12 #include <bitset>
13 #include <fmt/printf.h>
14 #include <fstream>
15 #include <boost/tokenizer.hpp>
16 #include <boost/range/adaptor/indexed.hpp>
17 
23 
24 #ifndef PHASE1PIXELMAP_STANDALONE
25 #define LOGDEBUG(x) LogDebug(x)
26 #else
27 #define LOGDEBUG(x) std::cout << x << ": "
28 #endif
29 
30 /*--------------------------------------------------------------------
31 / helper functions to dress plots
32 /--------------------------------------------------------------------*/
33 namespace PixelROCMapHelper {
34  void draw_line(double x1, double x2, double y1, double y2, int width, int style, int color);
35  void dress_plot(
36  TPad*& canv, TH2* h, int lay, int ring, int phase, bool half_shift, bool mark_zero, bool standard_palette);
37 } // namespace PixelROCMapHelper
38 
39 /*--------------------------------------------------------------------
40 / Ancillary struct to contain detector topology coordinates
41 /--------------------------------------------------------------------*/
43  int m_layer;
46  int m_ring;
47  int m_s_blade;
48  int m_s_disk;
49  int m_panel;
51 
52 public:
54  m_layer = -1;
55  m_s_ladder = -1;
56  m_s_module = -1;
57  m_ring = -1;
58  m_s_blade = -1;
59  m_s_disk = -1;
60  m_panel = -1;
61  m_isFlipped = false;
62  }
63 
65  if (this->isBarrel()) {
66  edm::LogPrint("DetCoordinates") << "layer: " << m_layer << " ladder: " << m_s_ladder << " module: " << m_s_module
67  << std::endl;
68  } else {
69  edm::LogPrint("DetCoordinates") << "ring: " << m_ring << " blade: " << m_s_blade << " panel: " << m_panel
70  << " disk: " << m_s_disk << std::endl;
71  }
72  }
73 
74  bool isBarrel() { return (m_layer > 0); }
75  bool isEndcap() { return (m_ring > 0); }
76  bool isFlipped() { return m_isFlipped; }
77 };
78 
79 /*--------------------------------------------------------------------
80 / Ancillary class to build pixel phase-1 tracker maps
81 /--------------------------------------------------------------------*/
83 public:
84  Phase1PixelROCMaps(const char* option, const std::string& zAxisTitle = "")
85  : m_option{option},
86  m_zAxisTitle{zAxisTitle},
88  edm::FileInPath("Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml").fullPath())} {
89  // --------------------- BOOK HISTOGRAMS
90  // barrel
91  for (unsigned int lay = 1; lay <= n_layers; lay++) {
92  int nlad = nlad_list[lay - 1];
93  std::string name = "occ_Layer_" + std::to_string(lay);
94  std::string title = "; Module # ; Ladder #";
95 
96  // if a z-axis title is specified, add the z-axis title
97  if (!m_zAxisTitle.empty()) {
98  title += fmt::sprintf(" ;%s", m_zAxisTitle.c_str());
99  }
100 
101  h_bpix_maps[lay - 1] = std::make_shared<TH2D>(
102  name.c_str(), title.c_str(), 72, -n_layers - 0.5, n_layers + 0.5, (nlad * 4 + 2), -nlad - 0.5, nlad + 0.5);
103  }
104 
105  // endcaps
106  for (unsigned int ring = 1; ring <= n_rings; ring++) {
107  int n = nybins_list[ring - 1];
108  float y = nxbins_list[ring - 1] + 0.5;
109  std::string name = "occ_ring_" + std::to_string(ring);
110  std::string title = "; Disk # ; Blade/Panel #";
111 
112  // if a z-axis title is specified, add the z-axis title
113  if (!m_zAxisTitle.empty()) {
114  title += fmt::sprintf(" ;%s", m_zAxisTitle.c_str());
115  }
116 
117  h_fpix_maps[ring - 1] =
118  std::make_shared<TH2D>(name.c_str(), title.c_str(), 56, -n_rings - 1.5, n_rings + 1.5, n, -y, y);
119  }
120  }
121 
123 
124  // Forward declarations
125  DetCoordinates findDetCoordinates(const uint32_t& t_detid);
126  void fillWholeModule(const uint32_t& detid, double value);
127  void fillSelectedRocs(const uint32_t& detid, const std::bitset<16>& theROCs, double value);
128  void drawBarrelMaps(TCanvas& canvas, const std::string& text = "");
129  void drawForwardMaps(TCanvas& canvas, const std::string& text = "");
130  void drawMaps(TCanvas& canvas, const std::string& text = "");
131 
132  inline std::array<std::shared_ptr<TH2D>, 4> getLayerMaps() { return h_bpix_maps; }
133  inline std::array<std::shared_ptr<TH2D>, 2> getRingMaps() { return h_fpix_maps; }
134 
135 private:
136  Option_t* m_option;
139 
140  // tough luck, we can only do phase-1...
141  static constexpr int numColumns = 416;
142  static constexpr int numRows = 160;
143  static constexpr int n_rings = 2;
144  static constexpr int n_layers = 4;
145 
146  const int nlad_list[n_layers] = {6, 14, 22, 32};
147  const int nybins_list[n_rings] = {92, 140};
148  const int nxbins_list[n_rings] = {11, 17};
149 
150  // maps
151  std::array<std::shared_ptr<TH2D>, n_layers> h_bpix_maps;
152  std::array<std::shared_ptr<TH2D>, n_rings> h_fpix_maps;
153 
154  // options
155  static constexpr const char* kVerbose = "verbose";
156 
157  // Forward declarations of private methods
158  std::vector<std::pair<int, int> > maskedBarrelRocsToBins(DetCoordinates coord);
159  std::vector<std::tuple<int, int, int> > maskedBarrelRocsToBins(DetCoordinates coord, std::bitset<16> myRocs);
160  std::vector<std::pair<int, int> > maskedForwardRocsToBins(DetCoordinates coord);
161  std::vector<std::tuple<int, int, int> > maskedForwardRocsToBins(DetCoordinates coord, std::bitset<16> myRocs);
162 
163 protected:
164  //============================================================================
165  // Taken from pixel naming classes
166  // BmO (-z-x) = 1, BmI (-z+x) = 2 , BpO (+z-x) = 3 , BpI (+z+x) = 4
167  inline int quadrant(const DetId& detid, bool phase_) {
168  if (detid.subdetId() == PixelSubdetector::PixelBarrel) {
169  return PixelBarrelName(detid, &m_trackerTopo, phase_).shell();
170  } else {
171  return PixelEndcapName(detid, &m_trackerTopo, phase_).halfCylinder();
172  }
173  }
174 
175  //============================================================================
176  // Online ladder convention taken from pixel naming class for barrel
177  // Apply sign convention (- sign for BmO and BpO)
178  inline int signed_ladder(const DetId& detid, bool phase_) {
180  return -9999;
181  int signed_ladder = PixelBarrelName(detid, &m_trackerTopo, phase_).ladderName();
182  if (quadrant(detid, phase_) % 2)
183  signed_ladder *= -1;
184  return signed_ladder;
185  }
186 
187  //============================================================================
188  // Online mdoule convention taken from pixel naming class for barrel
189  // Apply sign convention (- sign for BmO and BmI)
190  inline int signed_module(const DetId& detid, bool phase_) {
192  return -9999;
193  int signed_module = PixelBarrelName(detid, &m_trackerTopo, phase_).moduleName();
194  if (quadrant(detid, phase_) < 3)
195  signed_module *= -1;
196  return signed_module;
197  }
198 
199  //============================================================================
200  // Phase 0: Ring was not an existing convention
201  // but the 7 plaquettes were split by HV group
202  // --> Derive Ring 1/2 for them
203  // Panel 1 plq 1-2, Panel 2, plq 1 = Ring 1
204  // Panel 1 plq 3-4, Panel 2, plq 2-3 = Ring 2
205  // Phase 1: Using pixel naming class for endcap
206  inline int ring(const DetId& detid, bool phase_) {
208  return -9999;
209  int ring = -9999;
210  if (phase_ == 0) {
211  ring = 1 + (m_trackerTopo.pxfPanel(detid) + m_trackerTopo.pxfModule(detid) > 3);
212  } else if (phase_ == 1) {
213  ring = PixelEndcapName(detid, &m_trackerTopo, phase_).ringName();
214  }
215  return ring;
216  }
217 
218  //============================================================================
219  // Online blade convention taken from pixel naming class for endcap
220  // Apply sign convention (- sign for BmO and BpO)
221  inline int signed_blade(const DetId& detid, bool phase_) {
223  return -9999;
224  int signed_blade = PixelEndcapName(detid, &m_trackerTopo, phase_).bladeName();
225  if (quadrant(detid, phase_) % 2)
226  signed_blade *= -1;
227  return signed_blade;
228  }
229 
230  //============================================================================
231  inline int signed_blade_panel(const DetId& detid, bool phase_) {
233  return -9999;
234  int signed_blade_panel = signed_blade(detid, phase_) + (m_trackerTopo.pxfPanel(detid) - 1);
235  return signed_blade_panel;
236  }
237 
238  //============================================================================
239  // Online disk convention
240  // Apply sign convention (- sign for BmO and BmI)
241  inline int signed_disk(const DetId& detid, bool phase_) {
243  return -9999;
244  int signed_disk = m_trackerTopo.pxfDisk(DetId(detid));
245  if (quadrant(detid, phase_) < 3)
246  signed_disk *= -1;
247  return signed_disk;
248  }
249 
250  //============================================================================
251  // Determines if the BPix ldder is inner or outer
252  inline bool isBPixOuterLadder(const DetId& detid, bool isPhase0) {
253  bool isOuter = false;
254  int layer = m_trackerTopo.pxbLayer(detid.rawId());
255  bool odd_ladder = m_trackerTopo.pxbLadder(detid.rawId()) % 2;
256  if (isPhase0) {
257  if (layer == 2)
258  isOuter = !odd_ladder;
259  else
260  isOuter = odd_ladder;
261  } else {
262  if (layer == 4)
263  isOuter = odd_ladder;
264  else
265  isOuter = !odd_ladder;
266  }
267  return isOuter;
268  }
269 };
270 
271 #endif
static constexpr const char * kVerbose
void drawMaps(TCanvas &canvas, const std::string &text="")
int signed_blade(const DetId &detid, bool phase_)
unsigned int pxbLayer(const DetId &id) const
int ringName() const
ring Id
int bladeName() const
blade id
int moduleName() const
module id (index in z)
unsigned int pxfModule(const DetId &id) const
DetCoordinates findDetCoordinates(const uint32_t &t_detid)
int signed_module(const DetId &detid, bool phase_)
std::vector< std::pair< int, int > > maskedForwardRocsToBins(DetCoordinates coord)
static constexpr int n_layers
int signed_blade_panel(const DetId &detid, bool phase_)
std::array< std::shared_ptr< TH2D >, n_layers > h_bpix_maps
void drawBarrelMaps(TCanvas &canvas, const std::string &text="")
unsigned int pxbLadder(const DetId &id) const
static std::string to_string(const XMLCh *ch)
void fillWholeModule(const uint32_t &detid, double value)
void draw_line(double x1, double x2, double y1, double y2, int width, int style, int color)
const int nlad_list[n_layers]
constexpr std::array< uint8_t, layerIndexSize< TrackerTraits > > layer
void drawForwardMaps(TCanvas &canvas, const std::string &text="")
bool isBPixOuterLadder(const DetId &detid, bool isPhase0)
int quadrant(const DetId &detid, bool phase_)
static constexpr int n_rings
std::vector< std::pair< int, int > > maskedBarrelRocsToBins(DetCoordinates coord)
const int nxbins_list[n_rings]
Definition: style.py:1
std::array< std::shared_ptr< TH2D >, 4 > getLayerMaps()
Phase1PixelROCMaps(const char *option, const std::string &zAxisTitle="")
const int nybins_list[n_rings]
HalfCylinder halfCylinder() const
unsigned int pxfDisk(const DetId &id) const
int signed_disk(const DetId &detid, bool phase_)
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
Definition: value.py:1
Log< level::Warning, true > LogPrint
int ring(const DetId &detid, bool phase_)
Shell shell() const
unsigned int pxfPanel(const DetId &id) const
Definition: DetId.h:17
static constexpr int numColumns
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
void fillSelectedRocs(const uint32_t &detid, const std::bitset< 16 > &theROCs, double value)
TrackerTopology m_trackerTopo
std::array< std::shared_ptr< TH2D >, n_rings > h_fpix_maps
TrackerTopology fromTrackerParametersXMLFile(const std::string &xmlFileName)
def canvas(sub, attr)
Definition: svgfig.py:482
int signed_ladder(const DetId &detid, bool phase_)
int ladderName() const
ladder id (index in phi)
static constexpr int numRows
void dress_plot(TPad *&canv, TH2 *h, int lay, int ring, int phase, bool half_shift, bool mark_zero, bool standard_palette)
std::array< std::shared_ptr< TH2D >, 2 > getRingMaps()
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4