CMS 3D CMS Logo

List of all members | Public Member Functions | Protected Member Functions | Private Attributes | Static Private Attributes
Phase1PixelSummaryMap Class Reference

#include <Phase1PixelSummaryMap.h>

Public Member Functions

void createTrackerBaseMap ()
 
bool fillTrackerMap (unsigned int id, double value)
 
const std::pair< float, float > getZAxisRange () const
 
 Phase1PixelSummaryMap (const char *option, std::string title, std::string zAxisTitle)
 
void printTrackerMap (TCanvas &canvas, const float topMargin=0.02, int index=0)
 
void resetOption (const char *option)
 
void setZAxisRange (const double min, const double max)
 
 ~Phase1PixelSummaryMap ()=default
 

Protected Member Functions

void addNamedBins (edm::FileInPath geoFile, int tX, int tY, int sX, int sY, bool applyModuleRotation=false)
 

Private Attributes

TArrow arrow
 
const std::array< int, maxPxBarrelbarrelLadderShift = {{0, 14, 44, 90}}
 
std::map< uint32_t, std::shared_ptr< TGraph > > bins
 
const std::array< int, maxPxForwardforwardDiskXShift = {{25, 75, 125}}
 
const int forwardDiskYShift = 45
 
std::shared_ptr< TH2Poly > m_BaseTrackerMap
 
std::vector< edm::FileInPathm_cornersBPIX
 
std::vector< edm::FileInPathm_cornersFPIX
 
Option_t * m_option
 
const std::string m_title
 
TrackerTopology m_trackerTopo
 
const std::string m_zAxisTitle
 
TArrow phiArrow
 
const int plotHeight = 2000
 
const int plotWidth = 3000
 
TArrow xArrow
 
TArrow yArrow
 

Static Private Attributes

static const unsigned int maxPxBarrel = 4
 
static const unsigned int maxPxForward = 3
 

Detailed Description

Definition at line 57 of file Phase1PixelSummaryMap.h.

Constructor & Destructor Documentation

◆ Phase1PixelSummaryMap()

Phase1PixelSummaryMap::Phase1PixelSummaryMap ( const char *  option,
std::string  title,
std::string  zAxisTitle 
)
inline

Definition at line 59 of file Phase1PixelSummaryMap.h.

References fileinputsource_cfi::option.

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  }
const std::string m_zAxisTitle
std::vector< edm::FileInPath > m_cornersBPIX
std::vector< edm::FileInPath > m_cornersFPIX
TrackerTopology fromTrackerParametersXMLFile(const std::string &xmlFileName)

◆ ~Phase1PixelSummaryMap()

Phase1PixelSummaryMap::~Phase1PixelSummaryMap ( )
default

Member Function Documentation

◆ addNamedBins()

void Phase1PixelSummaryMap::addNamedBins ( edm::FileInPath  geoFile,
int  tX,
int  tY,
int  sX,
int  sY,
bool  applyModuleRotation = false 
)
protected

Definition at line 180 of file Phase1PixelSummaryMap.cc.

References bins, hcalRecHitTable_cff::detId, Exception, edm::FileInPath::fullPath(), mps_fire::i, mps_splice::line, LOGDEBUG, m_BaseTrackerMap, N, AlCaHLTBitMon_ParallelJobs::p, submitPVResolutionJobs::q, AlCaHLTBitMon_QueryRunRegistry::string, Ph1PMapSummaryHelper::tokenize(), and geometryCSVtoXML::xy.

Referenced by createTrackerBaseMap().

181  {
182  auto cornerFileName = geoFile.fullPath();
183  std::ifstream cornerFile(cornerFileName.c_str());
184  if (!cornerFile.good()) {
185  throw cms::Exception("FileError") << "Problem opening corner file: " << cornerFileName;
186  }
188  while (std::getline(cornerFile, line)) {
189  if (!line.empty()) {
190  std::istringstream iss(line);
191 
192  auto tokens = Ph1PMapSummaryHelper::tokenize(line, '"');
193  // Printing the token vector
194  for (unsigned int i = 0; i < tokens.size(); i++)
195  LOGDEBUG("Phase1PixelSummaryMap") << tokens[i] << '\n';
196 
197  auto detInfo = Ph1PMapSummaryHelper::tokenize(tokens[0], ' ');
198  unsigned int detId = stoi(detInfo[0]);
199  std::string detIdName = detInfo[1];
200  auto xy = Ph1PMapSummaryHelper::tokenize(tokens[1], ' ');
201  unsigned int verNum = 1;
202  std::vector<float> xP, yP;
203  for (const auto& coord : xy) {
204  auto coordSpl = Ph1PMapSummaryHelper::tokenize(coord, ',');
205  if (applyModuleRotation) {
206  xP.push_back(-(std::stof(coordSpl[0]) * sX + tX));
207  yP.push_back(((std::stof(coordSpl[1]) * sY + tY)));
208  } else {
209  xP.push_back(std::stof(coordSpl[0]) * sX + tX);
210  yP.push_back(std::stof(coordSpl[1]) * sY + tY);
211  }
212  verNum++;
213  }
214  //close the polygon
215  xP.push_back(xP[0]);
216  yP.push_back(yP[0]);
217 
218  LOGDEBUG("Phase1PixelSummaryMap") << detId << "[";
219  for (const auto& p : xP) {
220  LOGDEBUG("Phase1PixelSummaryMap") << p << ",";
221  }
222  LOGDEBUG("Phase1PixelSummaryMap") << "] [ ";
223  for (const auto& q : yP) {
224  LOGDEBUG("Phase1PixelSummaryMap") << q << ",";
225  }
226  LOGDEBUG("Phase1PixelSummaryMap") << "]" << std::endl;
227 
228  const unsigned int N = verNum;
229  if (applyModuleRotation) {
230  bins[detId] = std::make_shared<TGraph>(N, &yP[0], &xP[0]);
231  } else {
232  bins[detId] = std::make_shared<TGraph>(N, &xP[0], &yP[0]);
233  //bins[detId] = std::make_shared<TGraph>(N, &yP[0], &xP[0]); // rotation by 90 deg (so that it had the same layout as for the strips)
234  }
235 
236  bins[detId]->SetName(detInfo[0].c_str());
237  m_BaseTrackerMap->AddBin(bins[detId]->Clone());
238  }
239  }
240  return;
241 }
std::map< uint32_t, std::shared_ptr< TGraph > > bins
std::vector< std::string > tokenize(std::string line, char delimiter)
std::shared_ptr< TH2Poly > m_BaseTrackerMap
#define LOGDEBUG(x)
#define N
Definition: blowfish.cc:9
const std::string & fullPath() const
Definition: FileInPath.cc:144

◆ createTrackerBaseMap()

void Phase1PixelSummaryMap::createTrackerBaseMap ( )

Definition at line 33 of file Phase1PixelSummaryMap.cc.

References addNamedBins(), barrelLadderShift, forwardDiskXShift, forwardDiskYShift, mps_fire::i, dqmiolumiharvest::j, isotrackApplyRegressor::k, LOGINFO, m_BaseTrackerMap, m_cornersBPIX, m_cornersFPIX, m_zAxisTitle, maxPxBarrel, and maxPxForward.

Referenced by templateHelper::SiPixelFullPixelIDMap< PayloadType, StoreType, TransientType >::fill(), and main().

33  {
34  m_BaseTrackerMap = std::make_shared<TH2Poly>("Summary", "", -10, 160, -70, 70);
35  m_BaseTrackerMap->SetFloat(true);
36  m_BaseTrackerMap->GetXaxis()->SetTitle("");
37  m_BaseTrackerMap->GetYaxis()->SetTitle("");
38  m_BaseTrackerMap->GetZaxis()->SetTitle(m_zAxisTitle.c_str());
39  m_BaseTrackerMap->GetZaxis()->CenterTitle();
40  m_BaseTrackerMap->GetZaxis()->SetTitleOffset(1.2);
41  m_BaseTrackerMap->SetOption("COLZ L");
42  m_BaseTrackerMap->SetStats(false);
43 
44  //BARREL FIRST
45  for (unsigned int i = 0; i < maxPxBarrel; i++) {
46  LOGINFO("Phase1PixelSummaryMap") << "barrel, shift: " << i << " corner: " << i << std::endl;
47  LOGINFO("Phase1PixelSummaryMap") << "translate x: " << 0 << std::endl;
48  LOGINFO("Phase1PixelSummaryMap") << "translate y: " << barrelLadderShift[i] << std::endl;
49 
50  int currBarrelTranslateX = 0;
51  int currBarrelTranslateY = barrelLadderShift[i];
52  addNamedBins(m_cornersBPIX[i], currBarrelTranslateX, currBarrelTranslateY, 1, 1, true);
53  }
54 
55  //MINUS FORWARD
56  for (int j : {-3, -2, -1}) {
57  LOGINFO("Phase1PixelSummaryMap") << "negative fwd, shift: " << -j - 1 << " corner: " << maxPxForward + j
58  << std::endl;
59  LOGINFO("Phase1PixelSummaryMap") << "translate x: " << forwardDiskXShift[-j - 1] << std::endl;
60  LOGINFO("Phase1PixelSummaryMap") << "translate y: " << -forwardDiskYShift << std::endl;
61 
62  int currForwardTranslateX = forwardDiskXShift[-j - 1];
63  int currForwardTranslateY = -forwardDiskYShift;
64  addNamedBins(m_cornersFPIX[maxPxForward + j], currForwardTranslateX, currForwardTranslateY, 1, 1);
65  }
66 
67  //PLUS FORWARD
68  for (int k : {1, 2, 3}) {
69  LOGINFO("Phase1PixelSummaryMap") << "positive fwd, shift: " << k << " corner: " << maxPxForward + k - 1
70  << std::endl;
71  LOGINFO("Phase1PixelSummaryMap") << "translate x: " << forwardDiskXShift[k - 1] << std::endl;
72  LOGINFO("Phase1PixelSummaryMap") << "translate y: " << forwardDiskYShift << std::endl;
73 
74  int currForwardTranslateX = forwardDiskXShift[k - 1];
75  int currForwardTranslateY = forwardDiskYShift;
76  addNamedBins(m_cornersFPIX[maxPxForward + k - 1], currForwardTranslateX, currForwardTranslateY, 1, 1);
77  }
78 
79  edm::LogPrint("Phase1PixelSummaryMap") << "Base Tracker Map: constructed" << std::endl;
80  return;
81 }
const std::string m_zAxisTitle
static const unsigned int maxPxBarrel
const std::array< int, maxPxBarrel > barrelLadderShift
const std::array< int, maxPxForward > forwardDiskXShift
#define LOGINFO(x)
static const unsigned int maxPxForward
std::shared_ptr< TH2Poly > m_BaseTrackerMap
std::vector< edm::FileInPath > m_cornersBPIX
Log< level::Warning, true > LogPrint
void addNamedBins(edm::FileInPath geoFile, int tX, int tY, int sX, int sY, bool applyModuleRotation=false)
std::vector< edm::FileInPath > m_cornersFPIX

◆ fillTrackerMap()

bool Phase1PixelSummaryMap::fillTrackerMap ( unsigned int  id,
double  value 
)

Definition at line 157 of file Phase1PixelSummaryMap.cc.

References ALCARECOPPSCalTrackBasedSel_cff::detid, m_BaseTrackerMap, PixelSubdetector::PixelBarrel, and PixelSubdetector::PixelEndcap.

Referenced by templateHelper::SiPixelFullPixelIDMap< PayloadType, StoreType, TransientType >::fill(), and main().

157  {
158  auto detid = DetId(id);
159  if (detid.subdetId() != PixelSubdetector::PixelBarrel && detid.subdetId() != PixelSubdetector::PixelEndcap) {
160  edm::LogError("Phase1PixelSummaryMap")
161  << __func__ << " The following detid " << id << " is not Pixel!" << std::endl;
162  return false;
163  } else {
164  m_BaseTrackerMap->Fill(TString::Format("%u", id), value);
165  return true;
166  }
167 }
Log< level::Error, false > LogError
std::shared_ptr< TH2Poly > m_BaseTrackerMap
Definition: value.py:1
Definition: DetId.h:17

◆ getZAxisRange()

const std::pair< float, float > Phase1PixelSummaryMap::getZAxisRange ( ) const

Definition at line 170 of file Phase1PixelSummaryMap.cc.

References m_BaseTrackerMap.

170  {
171  return std::make_pair(m_BaseTrackerMap->GetMinimum(), m_BaseTrackerMap->GetMaximum());
172 }
std::shared_ptr< TH2Poly > m_BaseTrackerMap

◆ printTrackerMap()

void Phase1PixelSummaryMap::printTrackerMap ( TCanvas &  canvas,
const float  topMargin = 0.02,
int  index = 0 
)

Definition at line 84 of file Phase1PixelSummaryMap.cc.

References arrow, svgfig::canvas(), m_BaseTrackerMap, m_title, phiArrow, xArrow, and yArrow.

Referenced by templateHelper::SiPixelFullPixelIDMap< PayloadType, StoreType, TransientType >::fill(), and main().

84  {
85  //canvas = TCanvas("c1","c1",plotWidth,plotHeight);
86  if (index != 0)
87  canvas.cd(index);
88  else
89  canvas.cd();
90 
91  if (index == 0) {
92  canvas.SetTopMargin(topMargin);
93  canvas.SetBottomMargin(0.02);
94  canvas.SetLeftMargin(0.02);
95  canvas.SetRightMargin(0.14);
96  } else {
97  m_BaseTrackerMap->GetZaxis()->SetTitleOffset(1.5);
98  canvas.cd(index)->SetTopMargin(topMargin);
99  canvas.cd(index)->SetBottomMargin(0.02);
100  canvas.cd(index)->SetLeftMargin(0.02);
101  canvas.cd(index)->SetRightMargin(0.14);
102  }
103 
104  m_BaseTrackerMap->Draw("AL");
105  m_BaseTrackerMap->Draw("AC COLZ0 L SAME");
106 
107  //### z arrow
108  arrow = TArrow(0.05, 27.0, 0.05, -30.0, 0.02, "|>");
109  arrow.SetLineWidth(4);
110  arrow.Draw();
111  arrow.SetAngle(30);
112  //### phi arrow
113  phiArrow = TArrow(0.0, 27.0, 30.0, 27.0, 0.02, "|>");
114  phiArrow.SetLineWidth(4);
115  phiArrow.Draw();
116  phiArrow.SetAngle(30);
117  //### x arrow
118  xArrow = TArrow(25.0, 44.5, 50.0, 44.5, 0.02, "|>");
119  xArrow.SetLineWidth(4);
120  xArrow.Draw();
121  xArrow.SetAngle(30);
122  //### y arrow
123  yArrow = TArrow(25.0, 44.5, 25.0, 69.5, 0.02, "|>");
124  yArrow.SetLineWidth(4);
125  yArrow.Draw();
126  yArrow.SetAngle(30);
127 
128  //###################################################
129  //# add some captions
130  auto txt = TLatex();
131  txt.SetNDC();
132  txt.SetTextFont(1);
133  txt.SetTextColor(1);
134  txt.SetTextAlign(22);
135  txt.SetTextAngle(0);
136 
137  //# draw new-style title
138  txt.SetTextSize(0.03);
139  txt.DrawLatex(0.5, ((index == 0) ? 0.95 : 0.93), (fmt::sprintf("Pixel Tracker Map: %s", m_title)).c_str());
140  txt.SetTextSize(0.03);
141 
142  txt.DrawLatex(0.55, 0.125, "-DISK");
143  txt.DrawLatex(0.55, 0.875, "+DISK");
144 
145  txt.DrawLatex(0.08, 0.28, "+z");
146  txt.DrawLatex(0.25, 0.70, "+phi");
147  txt.DrawLatex((index == 0) ? 0.31 : 0.33, 0.78, "+x");
148  txt.DrawLatex((index == 0) ? 0.21 : 0.22, ((index == 0) ? 0.96 : 0.94), "+y");
149 
150  txt.SetTextAngle(90);
151  txt.DrawLatex(0.04, 0.5, "BARREL");
152 
153  edm::LogPrint("Phase1PixelSummaryMap") << "Base Tracker Map: printed" << std::endl;
154 }
std::shared_ptr< TH2Poly > m_BaseTrackerMap
Log< level::Warning, true > LogPrint
def canvas(sub, attr)
Definition: svgfig.py:482

◆ resetOption()

void Phase1PixelSummaryMap::resetOption ( const char *  option)

Definition at line 23 of file Phase1PixelSummaryMap.cc.

References m_option, and fileinputsource_cfi::option.

23  {
24  if (m_option != nullptr && !m_option[0]) {
25  m_option = option;
26  } else {
27  edm::LogError("Phase1PixelSummaryMap")
28  << "Option has already been set to " << m_option << ". It's not possible to reset it.";
29  }
30 }
Log< level::Error, false > LogError

◆ setZAxisRange()

void Phase1PixelSummaryMap::setZAxisRange ( const double  min,
const double  max 
)

Definition at line 175 of file Phase1PixelSummaryMap.cc.

References m_BaseTrackerMap, WZElectronSkims53X_cff::max, and SiStripPI::min.

175  {
176  m_BaseTrackerMap->GetZaxis()->SetRangeUser(min, max);
177 }
std::shared_ptr< TH2Poly > m_BaseTrackerMap

Member Data Documentation

◆ arrow

TArrow Phase1PixelSummaryMap::arrow
private

Definition at line 110 of file Phase1PixelSummaryMap.h.

Referenced by printTrackerMap().

◆ barrelLadderShift

const std::array<int, maxPxBarrel> Phase1PixelSummaryMap::barrelLadderShift = {{0, 14, 44, 90}}
private

Definition at line 102 of file Phase1PixelSummaryMap.h.

Referenced by createTrackerBaseMap().

◆ bins

std::map<uint32_t, std::shared_ptr<TGraph> > Phase1PixelSummaryMap::bins
private

Definition at line 95 of file Phase1PixelSummaryMap.h.

Referenced by addNamedBins().

◆ forwardDiskXShift

const std::array<int, maxPxForward> Phase1PixelSummaryMap::forwardDiskXShift = {{25, 75, 125}}
private

Definition at line 103 of file Phase1PixelSummaryMap.h.

Referenced by createTrackerBaseMap().

◆ forwardDiskYShift

const int Phase1PixelSummaryMap::forwardDiskYShift = 45
private

Definition at line 105 of file Phase1PixelSummaryMap.h.

Referenced by createTrackerBaseMap().

◆ m_BaseTrackerMap

std::shared_ptr<TH2Poly> Phase1PixelSummaryMap::m_BaseTrackerMap
private

◆ m_cornersBPIX

std::vector<edm::FileInPath> Phase1PixelSummaryMap::m_cornersBPIX
private

Definition at line 97 of file Phase1PixelSummaryMap.h.

Referenced by createTrackerBaseMap().

◆ m_cornersFPIX

std::vector<edm::FileInPath> Phase1PixelSummaryMap::m_cornersFPIX
private

Definition at line 98 of file Phase1PixelSummaryMap.h.

Referenced by createTrackerBaseMap().

◆ m_option

Option_t* Phase1PixelSummaryMap::m_option
private

Definition at line 89 of file Phase1PixelSummaryMap.h.

Referenced by resetOption().

◆ m_title

const std::string Phase1PixelSummaryMap::m_title
private

Definition at line 90 of file Phase1PixelSummaryMap.h.

Referenced by printTrackerMap().

◆ m_trackerTopo

TrackerTopology Phase1PixelSummaryMap::m_trackerTopo
private

Definition at line 93 of file Phase1PixelSummaryMap.h.

◆ m_zAxisTitle

const std::string Phase1PixelSummaryMap::m_zAxisTitle
private

Definition at line 91 of file Phase1PixelSummaryMap.h.

Referenced by createTrackerBaseMap().

◆ maxPxBarrel

const unsigned int Phase1PixelSummaryMap::maxPxBarrel = 4
staticprivate

Definition at line 100 of file Phase1PixelSummaryMap.h.

Referenced by createTrackerBaseMap().

◆ maxPxForward

const unsigned int Phase1PixelSummaryMap::maxPxForward = 3
staticprivate

Definition at line 101 of file Phase1PixelSummaryMap.h.

Referenced by createTrackerBaseMap().

◆ phiArrow

TArrow Phase1PixelSummaryMap::phiArrow
private

Definition at line 110 of file Phase1PixelSummaryMap.h.

Referenced by printTrackerMap().

◆ plotHeight

const int Phase1PixelSummaryMap::plotHeight = 2000
private

Definition at line 108 of file Phase1PixelSummaryMap.h.

◆ plotWidth

const int Phase1PixelSummaryMap::plotWidth = 3000
private

Definition at line 107 of file Phase1PixelSummaryMap.h.

◆ xArrow

TArrow Phase1PixelSummaryMap::xArrow
private

Definition at line 110 of file Phase1PixelSummaryMap.h.

Referenced by printTrackerMap().

◆ yArrow

TArrow Phase1PixelSummaryMap::yArrow
private

Definition at line 110 of file Phase1PixelSummaryMap.h.

Referenced by printTrackerMap().