CMS 3D CMS Logo

Phase1PixelSummaryMap.h
Go to the documentation of this file.
1 #ifndef DQM_TRACKERREMAPPER_PHASE1PIXELSUMMARYMAP_H
2 #define DQM_TRACKERREMAPPER_PHASE1PIXELSUMMARYMAP_H
3 
4 #include "TArrow.h"
5 #include "TCanvas.h"
6 #include "TGraph.h"
7 #include "TH1.h"
8 #include "TH2.h"
9 #include "TH2Poly.h"
10 #include "TLatex.h"
11 #include "TStyle.h"
12 
13 #include <array>
14 #include <fmt/printf.h>
15 #include <fstream>
16 #include <memory>
17 #include <boost/tokenizer.hpp>
18 #include <boost/range/adaptor/indexed.hpp>
19 
23 
24 #ifndef PH1PSUMMARYMAP_STANDALONE
25 #define LOGDEBUG(x) LogDebug(x)
26 #define LOGINFO(x) edm::LogInfo(x)
27 #define LOGPRINT(x) edm::LogPrint(x)
28 #else
29 #define LOGDEBUG(x) std::cout << x << " Debug : "
30 #define LOGINFO(x) std::cout << x << " Info : "
31 #define LOGPRINT(x) std::cout << x << " : "
32 #endif
33 
34 using indexedCorners = std::map<unsigned int, std::pair<std::vector<float>, std::vector<float>>>;
35 
37  //============================================================================
38  // utility to tokenize std::string
39  //============================================================================
40  inline std::vector<std::string> tokenize(std::string line, char delimiter) {
41  // Vector of string to save tokens
42  std::vector<std::string> tokens;
43  std::stringstream check1(line);
44  std::string intermediate;
45 
46  // Tokenizing w.r.t. delimiter
47  while (getline(check1, intermediate, delimiter)) {
48  tokens.push_back(intermediate);
49  }
50  return tokens;
51  }
52 } // namespace Ph1PMapSummaryHelper
53 
54 /*--------------------------------------------------------------------
55 / Ancillary class to build pixel phase-1 tracker maps
56 /--------------------------------------------------------------------*/
58 public:
60  : m_option{option},
61  m_title{title},
62  m_zAxisTitle{zAxisTitle},
64  edm::FileInPath("Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml").fullPath())} {
65  // store the file in path for the corners (BPIX)
66  for (unsigned int i = 1; i <= 4; i++) {
67  m_cornersBPIX.push_back(edm::FileInPath(Form("DQM/SiStripMonitorClient/data/Geometry/vertices_barrel_%i", i)));
68  }
69 
70  // store the file in path for the corners (BPIX)
71  for (int j : {-3, -2, -1, 1, 2, 3}) {
72  m_cornersFPIX.push_back(edm::FileInPath(Form("DQM/SiStripMonitorClient/data/Geometry/vertices_forward_%i", j)));
73  }
74  }
75 
76  ~Phase1PixelSummaryMap() = default;
77 
78  void resetOption(const char* option);
79  void createTrackerBaseMap();
80  void printTrackerMap(TCanvas& canvas, const float topMargin = 0.02, int index = 0);
81  bool fillTrackerMap(unsigned int id, double value);
82  void setZAxisRange(const double min, const double max);
83  const std::pair<float, float> getZAxisRange() const;
84 
85 protected:
86  void addNamedBins(edm::FileInPath geoFile, int tX, int tY, int sX, int sY, bool applyModuleRotation = false);
87 
88 private:
89  Option_t* m_option;
92 
94  std::shared_ptr<TH2Poly> m_BaseTrackerMap;
95  std::map<uint32_t, std::shared_ptr<TGraph>> bins;
96 
97  std::vector<edm::FileInPath> m_cornersBPIX;
98  std::vector<edm::FileInPath> m_cornersFPIX;
99 
100  static const unsigned int maxPxBarrel = 4;
101  static const unsigned int maxPxForward = 3;
102  const std::array<int, maxPxBarrel> barrelLadderShift = {{0, 14, 44, 90}};
103  const std::array<int, maxPxForward> forwardDiskXShift = {{25, 75, 125}};
104 
105  const int forwardDiskYShift = 45; //# to make +DISK on top in the 'strip-like' layout
106 
107  const int plotWidth = 3000;
108  const int plotHeight = 2000;
109 
111 };
112 
113 #endif
const std::string m_zAxisTitle
void resetOption(const char *option)
static const unsigned int maxPxBarrel
const std::array< int, maxPxBarrel > barrelLadderShift
std::map< uint32_t, std::shared_ptr< TGraph > > bins
const std::array< int, maxPxForward > forwardDiskXShift
std::vector< std::string > tokenize(std::string line, char delimiter)
const std::pair< float, float > getZAxisRange() const
~Phase1PixelSummaryMap()=default
static const unsigned int maxPxForward
std::shared_ptr< TH2Poly > m_BaseTrackerMap
std::vector< edm::FileInPath > m_cornersBPIX
std::map< unsigned int, std::pair< std::vector< float >, std::vector< float > >> indexedCorners
void setZAxisRange(const double min, const double max)
Definition: value.py:1
void addNamedBins(edm::FileInPath geoFile, int tX, int tY, int sX, int sY, bool applyModuleRotation=false)
std::vector< edm::FileInPath > m_cornersFPIX
void printTrackerMap(TCanvas &canvas, const float topMargin=0.02, int index=0)
bool fillTrackerMap(unsigned int id, double value)
TrackerTopology fromTrackerParametersXMLFile(const std::string &xmlFileName)
def canvas(sub, attr)
Definition: svgfig.py:482
Phase1PixelSummaryMap(const char *option, std::string title, std::string zAxisTitle)