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