CMS 3D CMS Logo

Phase1PixelMaps.cc
Go to the documentation of this file.
1 #include "TH2Poly.h"
2 #include "TGraph.h"
3 #include "TH1.h"
4 #include "TH2.h"
5 #include "TStyle.h"
6 #include "TCanvas.h"
7 
8 #include <fmt/printf.h>
9 #include <fstream>
10 #include <boost/tokenizer.hpp>
11 #include <boost/range/adaptor/indexed.hpp>
12 
17 
18 // set option, but only if not already set
19 //============================================================================
21  if (m_option != nullptr && !m_option[0]) {
22  m_option = option;
23  } else {
24  edm::LogError("Phase1PixelMaps") << "Option has already been set to " << m_option
25  << ". It's not possible to reset it.";
26  }
27 }
28 
29 //============================================================================
30 void Phase1PixelMaps::bookBarrelHistograms(const std::string& currentHistoName, const char* what, const char* zaxis) {
31  std::string histName;
32  std::shared_ptr<TH2Poly> th2p;
33 
34  // check if the passed histogram name already exists, if not store it
35  if (std::find(m_knownNames.begin(), m_knownNames.end(), currentHistoName) == m_knownNames.end()) {
36  m_knownNames.emplace_back(currentHistoName);
37  }
38 
39  for (unsigned i = 0; i < 4; ++i) {
40  histName = "barrel_layer_";
41 
42  th2p = std::make_shared<TH2Poly>(
43  (histName + std::to_string(i + 1)).c_str(), Form("PXBMap of %s - Layer %i", what, i + 1), -15.0, 15.0, 0.0, 5.0);
44 
45  th2p->SetFloat();
46 
47  th2p->GetXaxis()->SetTitle("z [cm]");
48  th2p->GetYaxis()->SetTitle("ladder");
49  th2p->GetZaxis()->SetTitle(zaxis);
50  th2p->GetZaxis()->CenterTitle();
51  th2p->SetStats(false);
52  th2p->SetOption(m_option);
53  pxbTh2PolyBarrel[currentHistoName].push_back(th2p);
54  }
55 
56  th2p = std::make_shared<TH2Poly>("barrel_summary", Form("Barrel Pixel Map of %s", what), -5.0, 5.0, 0.0, 15.0);
57  th2p->SetFloat();
58 
59  th2p->GetXaxis()->SetTitle("");
60  th2p->GetYaxis()->SetTitle("");
61  th2p->GetZaxis()->SetTitle(zaxis);
62  th2p->GetZaxis()->CenterTitle();
63  th2p->SetStats(false);
64  th2p->SetOption(m_option);
65  pxbTh2PolyBarrelSummary[currentHistoName] = th2p;
66 
67  // book the bins
68  bookBarrelBins(currentHistoName);
69 
70  // set the isBooked bit to true;
71  m_isBooked.first = true;
72 }
73 
74 //============================================================================
75 void Phase1PixelMaps::bookForwardHistograms(const std::string& currentHistoName, const char* what, const char* zaxis) {
76  std::string histName;
77  std::shared_ptr<TH2Poly> th2p;
78 
79  // check if the passed histogram name already exists, if not store it
80  if (std::find(m_knownNames.begin(), m_knownNames.end(), currentHistoName) == m_knownNames.end()) {
81  m_knownNames.emplace_back(currentHistoName);
82  }
83 
84  for (unsigned side = 1; side <= 2; ++side) {
85  for (unsigned disk = 1; disk <= 3; ++disk) {
86  histName = "forward_disk_";
87 
88  th2p = std::make_shared<TH2Poly>((histName + std::to_string((side == 1 ? -(int(disk)) : (int)disk))).c_str(),
89  Form("PXFMap of %s - Side %i Disk %i", what, side, disk),
90  -15.0,
91  15.0,
92  -15.0,
93  15.0);
94  th2p->SetFloat();
95  th2p->GetXaxis()->SetTitle("x [cm]");
96  th2p->GetYaxis()->SetTitle("y [cm]");
97  th2p->GetZaxis()->SetTitle(zaxis);
98  th2p->GetZaxis()->CenterTitle();
99  th2p->SetStats(false);
100  th2p->SetOption(m_option);
101  pxfTh2PolyForward[currentHistoName].push_back(th2p);
102  }
103  }
104 
105  th2p = std::make_shared<TH2Poly>("forward_summary", Form("Forward Pixel Map of %s", what), -40.0, 40.0, -20.0, 90.0);
106  th2p->SetFloat();
107 
108  th2p->GetXaxis()->SetTitle("");
109  th2p->GetYaxis()->SetTitle("");
110  th2p->GetZaxis()->SetTitle(zaxis);
111  th2p->GetZaxis()->CenterTitle();
112  th2p->SetStats(false);
113  th2p->SetOption(m_option);
114  pxfTh2PolyForwardSummary[currentHistoName] = th2p;
115 
116  // book the bins
117  bookForwardBins(currentHistoName);
118 
119  m_isBooked.second = true;
120 }
121 
122 //============================================================================
123 void Phase1PixelMaps::bookBarrelBins(const std::string& currentHistoName) {
124  auto theIndexedCorners = this->retrieveCorners(m_cornersBPIX, 4);
125 
126  for (const auto& entry : theIndexedCorners) {
127  auto id = entry.first;
128  auto detid = DetId(id);
129  if (detid.subdetId() != PixelSubdetector::PixelBarrel)
130  continue;
131 
132  int layer = m_trackerTopo.pxbLayer(detid);
133  int ladder = m_trackerTopo.pxbLadder(detid);
134 
135  auto theVectX = entry.second.first;
136  auto theVectY = entry.second.second;
137 
138  float vertX[] = {theVectX[0], theVectX[1], theVectX[2], theVectX[3], theVectX[4]};
139  float vertY[] = {(ladder - 1.0f), (ladder - 1.0f), (float)ladder, (float)ladder, (ladder - 1.0f)};
140 
141  bins[id] = std::make_shared<TGraph>(5, vertX, vertY);
142  bins[id]->SetName(TString::Format("%u", id));
143 
144  // Summary plot
145  for (unsigned k = 0; k < 5; ++k) {
146  vertX[k] += ((layer == 2 || layer == 3) ? 0.0f : -60.0f);
147  vertY[k] += ((layer > 2) ? 30.0f : 0.0f);
148  }
149 
150  binsSummary[id] = std::make_shared<TGraph>(5, vertX, vertY);
151  binsSummary[id]->SetName(TString::Format("%u", id));
152 
153  if (pxbTh2PolyBarrel.find(currentHistoName) != pxbTh2PolyBarrel.end()) {
154  pxbTh2PolyBarrel[currentHistoName][layer - 1]->AddBin(bins[id]->Clone());
155  } else {
156  throw cms::Exception("LogicError") << currentHistoName << " is not found in the Barrel map! Aborting.";
157  }
158 
159  if (pxbTh2PolyBarrelSummary.find(currentHistoName) != pxbTh2PolyBarrelSummary.end()) {
160  pxbTh2PolyBarrelSummary[currentHistoName]->AddBin(binsSummary[id]->Clone());
161  } else {
162  throw cms::Exception("LocalError") << currentHistoName << " is not found in the Barrel Summary map! Aborting.";
163  }
164  }
165 }
166 
167 //============================================================================
168 void Phase1PixelMaps::bookForwardBins(const std::string& currentHistoName) {
169  auto theIndexedCorners = this->retrieveCorners(m_cornersFPIX, 3);
170 
171  for (const auto& entry : theIndexedCorners) {
172  auto id = entry.first;
173  auto detid = DetId(id);
174  if (detid.subdetId() != PixelSubdetector::PixelEndcap)
175  continue;
176 
177  int disk = m_trackerTopo.pxfDisk(detid);
178  int side = m_trackerTopo.pxfSide(detid);
179 
180  unsigned mapIdx = disk + (side - 1) * 3 - 1;
181 
182  auto theVectX = entry.second.first;
183  auto theVectY = entry.second.second;
184 
185  float vertX[] = {theVectX[0], theVectX[1], theVectX[2], theVectX[3]};
186  float vertY[] = {theVectY[0], theVectY[1], theVectY[2], theVectY[3]};
187 
188  bins[id] = std::make_shared<TGraph>(4, vertX, vertY);
189  bins[id]->SetName(TString::Format("%u", id));
190 
191  // Summary plot
192  for (unsigned k = 0; k < 4; ++k) {
193  vertX[k] += (float(side) - 1.5f) * 40.0f;
194  vertY[k] += (disk - 1) * 35.0f;
195  }
196 
197  binsSummary[id] = std::make_shared<TGraph>(4, vertX, vertY);
198  binsSummary[id]->SetName(TString::Format("%u", id));
199 
200  if (pxfTh2PolyForward.find(currentHistoName) != pxfTh2PolyForward.end()) {
201  pxfTh2PolyForward[currentHistoName][mapIdx]->AddBin(bins[id]->Clone());
202  } else {
203  throw cms::Exception("LogicError") << currentHistoName << " is not found in the Forward map! Aborting.";
204  }
205 
206  if (pxfTh2PolyForwardSummary.find(currentHistoName) != pxfTh2PolyForwardSummary.end()) {
207  pxfTh2PolyForwardSummary[currentHistoName]->AddBin(binsSummary[id]->Clone());
208  } else {
209  throw cms::Exception("LogicError") << currentHistoName << " is not found in the Forward Summary map! Aborting.";
210  }
211  }
212 }
213 
214 //============================================================================
215 void Phase1PixelMaps::book(const std::string& currentHistoName, const char* what, const char* zaxis) {
216  bookBarrelHistograms(currentHistoName, what, zaxis);
217  bookForwardHistograms(currentHistoName, what, zaxis);
218  m_isBooked = std::make_pair(true, true);
219 }
220 
221 //============================================================================
222 void Phase1PixelMaps::fill(const std::string& currentHistoName, unsigned int id, double value) {
223  auto detid = DetId(id);
224  if (detid.subdetId() == PixelSubdetector::PixelBarrel) {
225  int layer = m_trackerTopo.pxbLayer(id);
226  if (!m_isBooked.first) {
227  edm::LogError("Phase1PixelMaps") << __func__ << ": trying to fill a histogram not booked";
228  return;
229  }
230 
231  LogDebug("Phase1PixelMaps") << __func__ << " filling barrel with value: " << value << std::endl;
232 
233  pxbTh2PolyBarrel[currentHistoName][layer - 1]->Fill(TString::Format("%u", id), value);
234  pxbTh2PolyBarrelSummary[currentHistoName]->Fill(TString::Format("%u", id), value);
235  } else if (detid.subdetId() == PixelSubdetector::PixelEndcap) {
236  int disk = m_trackerTopo.pxfDisk(id);
237  int side = m_trackerTopo.pxfSide(id);
238  unsigned mapIdx = disk + (side - 1) * 3 - 1;
239  if (!m_isBooked.second) {
240  edm::LogError("Phase1PixelMaps") << __func__ << ": trying to fill a histogram not booked";
241  return;
242  }
243 
244  LogDebug("Phase1PixelMaps") << __func__ << " filling endcaps with value: " << value << std::endl;
245 
246  pxfTh2PolyForward[currentHistoName][mapIdx]->Fill(TString::Format("%u", id), value);
247  pxfTh2PolyForwardSummary[currentHistoName]->Fill(TString::Format("%u", id), value);
248  }
249 }
250 
251 //============================================================================
252 void Phase1PixelMaps::fillBarrelBin(const std::string& currentHistoName, unsigned int id, double value) {
253  auto detid = DetId(id);
254  if (detid.subdetId() != PixelSubdetector::PixelBarrel) {
255  edm::LogError("Phase1PixelMaps") << "fillBarrelBin() The following detid " << id << " is not Pixel Barrel!"
256  << std::endl;
257  return;
258  }
259  if (!m_isBooked.first) {
260  edm::LogError("Phase1PixelMaps") << __func__ << ": trying to fill a histogram not booked";
261  return;
262  }
263  int layer = m_trackerTopo.pxbLayer(id);
264  pxbTh2PolyBarrel[currentHistoName][layer - 1]->Fill(TString::Format("%u", id), value);
265  pxbTh2PolyBarrelSummary[currentHistoName]->Fill(TString::Format("%u", id), value);
266 }
267 
268 //============================================================================
269 void Phase1PixelMaps::fillForwardBin(const std::string& currentHistoName, unsigned int id, double value) {
270  auto detid = DetId(id);
271  if (detid.subdetId() != PixelSubdetector::PixelEndcap) {
272  edm::LogError("Phase1PixelMaps") << "fillForwardBin() The following detid " << id << " is not Pixel Forward!"
273  << std::endl;
274  return;
275  }
276  if (!m_isBooked.second) {
277  edm::LogError("Phase1PixelMaps") << __func__ << ": trying to fill a histogram not booked";
278  return;
279  }
280  int disk = m_trackerTopo.pxfDisk(id);
281  int side = m_trackerTopo.pxfSide(id);
282  unsigned mapIdx = disk + (side - 1) * 3 - 1;
283  pxfTh2PolyForward[currentHistoName][mapIdx]->Fill(TString::Format("%u", id), value);
284  pxfTh2PolyForwardSummary[currentHistoName]->Fill(TString::Format("%u", id), value);
285 }
286 
287 //============================================================================
289  if (!m_isBooked.first && !m_isBooked.second) {
290  edm::LogError("Phase1PixelMaps") << __func__ << ": trying to beautify a histogram not booked";
291  return;
292  }
293 
294  // only if the barrel is booked
295  if (m_isBooked.first) {
296  for (const auto& vec : pxbTh2PolyBarrel) {
297  for (const auto& plot : vec.second) {
298  this->makeNicePlotStyle(plot.get());
299  plot->GetXaxis()->SetTitleOffset(0.9);
300  plot->GetYaxis()->SetTitleOffset(0.9);
301  plot->GetZaxis()->SetTitleOffset(1.2);
302  plot->GetZaxis()->SetTitleSize(0.05);
303  }
304  }
305  }
306 
307  // only if the forwards are booked
308  if (m_isBooked.second) {
309  for (const auto& vec : pxfTh2PolyForward) {
310  for (const auto& plot : vec.second) {
311  this->makeNicePlotStyle(plot.get());
312  plot->GetXaxis()->SetTitleOffset(0.9);
313  plot->GetYaxis()->SetTitleOffset(0.9);
314  plot->GetZaxis()->SetTitleOffset(1.2);
315  plot->GetZaxis()->SetTitleSize(0.05);
316  }
317  }
318  }
319 }
320 
321 //============================================================================
322 void Phase1PixelMaps::setBarrelScale(const std::string& currentHistoName, std::pair<float, float> extrema) {
323  for (auto& histo : pxbTh2PolyBarrel[currentHistoName]) {
324  histo->GetZaxis()->SetRangeUser(extrema.first, extrema.second);
325  }
326 }
327 
328 //============================================================================
329 void Phase1PixelMaps::setForwardScale(const std::string& currentHistoName, std::pair<float, float> extrema) {
330  for (auto& histo : pxfTh2PolyForward[currentHistoName]) {
331  histo->GetZaxis()->SetRangeUser(extrema.first, extrema.second);
332  }
333 }
334 
335 //============================================================================
336 void Phase1PixelMaps::drawBarrelMaps(const std::string& currentHistoName, TCanvas& canvas, const char* drawOption) {
337  canvas.Divide(2, 2);
338  auto found = (std::find(m_knownNames.begin(), m_knownNames.end(), currentHistoName) != m_knownNames.end());
339 
340  if (!m_isBooked.first || !found) {
341  edm::LogError("Phase1PixelMaps") << __func__ << ": trying to draw a histogram not booked";
342  return;
343  }
344 
345  for (int i = 1; i <= 4; i++) {
346  canvas.cd(i);
347  if (strcmp(m_option, "text") == 0) {
348  canvas.cd(i)->SetRightMargin(0.02);
349  pxbTh2PolyBarrel[currentHistoName].at(i - 1)->SetMarkerColor(kRed);
350  } else {
351  if (m_autorescale)
352  rescaleAllBarrel(currentHistoName);
353  adjustCanvasMargins(canvas.cd(i), 0.07, 0.12, 0.10, 0.18);
354  }
355  if (drawOption) {
356  pxbTh2PolyBarrel[currentHistoName].at(i - 1)->Draw("L");
357  pxbTh2PolyBarrel[currentHistoName].at(i - 1)->Draw(fmt::sprintf("%s%ssame", m_option, drawOption).c_str());
358  } else {
359  pxbTh2PolyBarrel[currentHistoName].at(i - 1)->Draw("L");
360  pxbTh2PolyBarrel[currentHistoName].at(i - 1)->Draw(fmt::sprintf("%ssame", m_option).c_str());
361  }
362  }
363 }
364 
365 //============================================================================
366 void Phase1PixelMaps::drawForwardMaps(const std::string& currentHistoName, TCanvas& canvas, const char* drawOption) {
367  canvas.Divide(3, 2);
368  auto found = (std::find(m_knownNames.begin(), m_knownNames.end(), currentHistoName) != m_knownNames.end());
369 
370  if (!m_isBooked.second || !found) {
371  edm::LogError("Phase1PixelMaps") << __func__ << ": trying to draw a histogram not booked";
372  return;
373  }
374 
375  for (int i = 1; i <= 6; i++) {
376  canvas.cd(i);
377  if (strcmp(m_option, "text") == 0) {
378  canvas.cd(i)->SetRightMargin(0.02);
379  pxfTh2PolyForward[currentHistoName].at(i - 1)->SetMarkerColor(kRed);
380  } else {
381  if (m_autorescale)
382  rescaleAllForward(currentHistoName);
383  adjustCanvasMargins(canvas.cd(i), 0.07, 0.12, 0.10, 0.18);
384  }
385  if (drawOption) {
386  pxfTh2PolyForward[currentHistoName].at(i - 1)->Draw("L");
387  pxfTh2PolyForward[currentHistoName].at(i - 1)->Draw(fmt::sprintf("%s%ssame", m_option, drawOption).c_str());
388  } else {
389  pxfTh2PolyForward[currentHistoName].at(i - 1)->Draw("L");
390  pxfTh2PolyForward[currentHistoName].at(i - 1)->Draw(fmt::sprintf("%ssame", m_option).c_str());
391  }
392  }
393 }
394 
395 //============================================================================
396 void Phase1PixelMaps::drawSummaryMaps(const std::string& currentHistoName, TCanvas& canvas, const char* drawOption) {
397  auto found = (std::find(m_knownNames.begin(), m_knownNames.end(), currentHistoName) != m_knownNames.end());
398 
399  if (!m_isBooked.second || !m_isBooked.first || !found) {
400  edm::LogError("Phase1PixelMaps") << __func__ << ": trying to draw a histogram not booked";
401  return;
402  }
403 
404  canvas.Divide(2, 1);
405  canvas.cd(1);
406  std::string temp(m_option); // create a std string
407  auto isText = (temp.find("text") != std::string::npos);
408  adjustCanvasMargins(canvas.cd(1), 0.07, 0.02, 0.01, isText ? 0.05 : 0.15);
409  if (isText) {
410  pxbTh2PolyBarrelSummary[currentHistoName]->SetMarkerColor(kRed);
411  pxbTh2PolyBarrelSummary[currentHistoName]->SetMarkerSize(0.5);
412  }
413  pxbTh2PolyBarrelSummary[currentHistoName]->GetZaxis()->SetTitleOffset(1.4);
414  pxbTh2PolyBarrelSummary[currentHistoName]->Draw("AL");
415  pxbTh2PolyBarrelSummary[currentHistoName]->Draw(fmt::sprintf("%s%ssame", m_option, drawOption).c_str());
416 
417  canvas.cd(2);
418  adjustCanvasMargins(canvas.cd(2), 0.07, 0.02, 0.01, isText ? 0.05 : 0.15);
419  if (isText) {
420  pxfTh2PolyForwardSummary[currentHistoName]->SetMarkerColor(kRed);
421  pxfTh2PolyForwardSummary[currentHistoName]->SetMarkerSize(0.5);
422  }
423  pxfTh2PolyForwardSummary[currentHistoName]->GetZaxis()->SetTitleOffset(1.4);
424  pxfTh2PolyForwardSummary[currentHistoName]->Draw("AL");
425  pxfTh2PolyForwardSummary[currentHistoName]->Draw(fmt::sprintf("%s%ssame", m_option, drawOption).c_str());
426 }
427 
428 //============================================================================
429 const indexedCorners Phase1PixelMaps::retrieveCorners(const std::vector<edm::FileInPath>& cornerFiles,
430  const unsigned int reads) {
431  indexedCorners theOutMap;
432 
433  for (const auto& file : cornerFiles) {
434  auto cornerFileName = file.fullPath();
435  std::ifstream cornerFile(cornerFileName.c_str());
436  if (!cornerFile.good()) {
437  throw cms::Exception("FileError") << "Problem opening corner file: " << cornerFileName;
438  }
440  while (std::getline(cornerFile, line)) {
441  if (!line.empty()) {
442  std::istringstream iss(line);
443  unsigned int id;
445  std::vector<std::string> corners(reads, "");
446  std::vector<float> xP, yP;
447 
448  iss >> id >> name;
449  for (unsigned int i = 0; i < reads; ++i) {
450  iss >> corners.at(i);
451  }
452 
453  LOGDEBUG("Phase1PixelMaps") << id << " : ";
454  for (unsigned int i = 0; i < reads; i++) {
455  // remove the leading and trailing " signs in the corners list
456  (corners[i]).erase(std::remove(corners[i].begin(), corners[i].end(), '"'), corners[i].end());
457  LOGDEBUG("Phase1PixelMaps") << corners.at(i) << " ";
458  typedef boost::tokenizer<boost::char_separator<char>> tokenizer;
459  boost::char_separator<char> sep{","};
460  tokenizer tok{corners.at(i), sep};
461  for (const auto& t : tok | boost::adaptors::indexed(0)) {
462  if (t.index() == 0) {
463  xP.push_back(atof((t.value()).c_str()));
464  } else if (t.index() == 1) {
465  yP.push_back(atof((t.value()).c_str()));
466  } else {
467  edm::LogError("LogicError") << "There should not be any token with index " << t.index() << std::endl;
468  }
469  }
470  }
471  LOGDEBUG("Phase1PixelMaps") << std::endl;
472 
473  xP.push_back(xP.front());
474  yP.push_back(yP.front());
475 
476  for (unsigned int i = 0; i < xP.size(); i++) {
477  LOGDEBUG("Phase1PixelMaps") << "x[" << i << "]=" << xP[i] << " y[" << i << "]" << yP[i] << std::endl;
478  }
479 
480  theOutMap[id] = std::make_pair(xP, yP);
481 
482  } // if line is empty
483  } // loop on lines
484  } // loop on files
485  return theOutMap;
486 }
487 
488 //============================================================================
490  hist->SetStats(kFALSE);
491  hist->SetLineWidth(2);
492  hist->GetXaxis()->CenterTitle(true);
493  hist->GetYaxis()->CenterTitle(true);
494  hist->GetXaxis()->SetTitleFont(42);
495  hist->GetYaxis()->SetTitleFont(42);
496  hist->GetXaxis()->SetTitleSize(0.05);
497  hist->GetYaxis()->SetTitleSize(0.05);
498  hist->GetXaxis()->SetTitleOffset(1.1);
499  hist->GetYaxis()->SetTitleOffset(1.3);
500  hist->GetXaxis()->SetLabelFont(42);
501  hist->GetYaxis()->SetLabelFont(42);
502  hist->GetYaxis()->SetLabelSize(.05);
503  hist->GetXaxis()->SetLabelSize(.05);
504 
505  if (hist->InheritsFrom(TH2::Class())) {
506  hist->GetZaxis()->SetLabelFont(42);
507  hist->GetZaxis()->SetLabelFont(42);
508  hist->GetZaxis()->SetLabelSize(.05);
509  hist->GetZaxis()->SetLabelSize(.05);
510  }
511 }
512 
513 //============================================================================
514 void Phase1PixelMaps::adjustCanvasMargins(TVirtualPad* pad, float top, float bottom, float left, float right) {
515  if (top > 0)
516  pad->SetTopMargin(top);
517  if (bottom > 0)
518  pad->SetBottomMargin(bottom);
519  if (left > 0)
520  pad->SetLeftMargin(left);
521  if (right > 0)
522  pad->SetRightMargin(right);
523 }
524 
525 //============================================================================
526 void Phase1PixelMaps::rescaleAllBarrel(const std::string& currentHistoName) {
527  if (std::find(m_knownNames.begin(), m_knownNames.end(), currentHistoName) == m_knownNames.end()) {
528  edm::LogError("Phase1PixelMaps") << __func__ << ": trying to manipulate a histogram not booked";
529  return;
530  }
531 
532  std::vector<float> maxima;
533  std::transform(pxbTh2PolyBarrel[currentHistoName].begin(),
534  pxbTh2PolyBarrel[currentHistoName].end(),
535  std::back_inserter(maxima),
536  [](std::shared_ptr<TH2Poly> thp) -> float { return thp->GetMaximum(); });
537  std::vector<float> minima;
538  std::transform(pxbTh2PolyBarrel[currentHistoName].begin(),
539  pxbTh2PolyBarrel[currentHistoName].end(),
540  std::back_inserter(minima),
541  [](std::shared_ptr<TH2Poly> thp) -> float { return thp->GetMinimum(); });
542 
543  auto globalMax = *std::max_element(maxima.begin(), maxima.end());
544  auto globalMin = *std::min_element(minima.begin(), minima.end());
545 
546  for (auto& histo : pxbTh2PolyBarrel[currentHistoName]) {
547  histo->GetZaxis()->SetRangeUser(globalMin, globalMax);
548  }
549 }
550 
551 //============================================================================
552 void Phase1PixelMaps::rescaleAllForward(const std::string& currentHistoName) {
553  if (std::find(m_knownNames.begin(), m_knownNames.end(), currentHistoName) == m_knownNames.end()) {
554  edm::LogError("Phase1PixelMaps") << __func__ << ": trying to manipulate a histogram not booked";
555  return;
556  }
557 
558  std::vector<float> maxima;
559  std::transform(pxfTh2PolyForward[currentHistoName].begin(),
560  pxfTh2PolyForward[currentHistoName].end(),
561  std::back_inserter(maxima),
562  [](std::shared_ptr<TH2Poly> thp) -> float { return thp->GetMaximum(); });
563  std::vector<float> minima;
564  std::transform(pxfTh2PolyForward[currentHistoName].begin(),
565  pxfTh2PolyForward[currentHistoName].end(),
566  std::back_inserter(minima),
567  [](std::shared_ptr<TH2Poly> thp) -> float { return thp->GetMinimum(); });
568 
569  auto globalMax = *std::max_element(maxima.begin(), maxima.end());
570  auto globalMin = *std::min_element(minima.begin(), minima.end());
571 
572  for (auto& histo : pxfTh2PolyForward[currentHistoName]) {
573  histo->GetZaxis()->SetRangeUser(globalMin, globalMax);
574  }
575 }
void drawForwardMaps(const std::string &currentHistoName, TCanvas &canvas, const char *drawOption=nullptr)
unsigned int pxbLayer(const DetId &id) const
std::vector< edm::FileInPath > m_cornersFPIX
void setBarrelScale(const std::string &currentHistoName, std::pair< float, float > extrema)
void adjustCanvasMargins(TVirtualPad *pad, float top, float bottom, float left, float right)
void makeNicePlotStyle(TH1 *hist)
Option_t * m_option
void rescaleAllForward(const std::string &currentHistoName)
std::string to_string(const V &value)
Definition: OMSAccess.h:71
unsigned int pxbLadder(const DetId &id) const
Log< level::Error, false > LogError
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
void drawSummaryMaps(const std::string &currentHistoName, TCanvas &canvas, const char *drawOption=nullptr)
constexpr std::array< uint8_t, layerIndexSize > layer
std::map< std::string, std::vector< std::shared_ptr< TH2Poly > > > pxfTh2PolyForward
void fillForwardBin(const std::string &currentHistoName, unsigned int id, double value)
#define LOGDEBUG(x)
std::map< uint32_t, std::shared_ptr< TGraph > > binsSummary
std::vector< edm::FileInPath > m_cornersBPIX
void drawBarrelMaps(const std::string &currentHistoName, TCanvas &canvas, const char *drawOption=nullptr)
std::map< uint32_t, std::shared_ptr< TGraph > > bins
const indexedCorners retrieveCorners(const std::vector< edm::FileInPath > &cornerFiles, const unsigned int reads)
void resetOption(const char *option)
void bookBarrelBins(const std::string &currentHistoName)
TrackerTopology m_trackerTopo
void rescaleAllBarrel(const std::string &currentHistoName)
unsigned int pxfDisk(const DetId &id) const
std::map< unsigned int, std::pair< std::vector< float >, std::vector< float > >> indexedCorners
std::map< std::string, std::vector< std::shared_ptr< TH2Poly > > > pxbTh2PolyBarrel
std::map< std::string, std::shared_ptr< TH2Poly > > pxbTh2PolyBarrelSummary
double f[11][100]
Definition: value.py:1
void bookForwardBins(const std::string &currentHistoName)
void beautifyAllHistograms()
std::pair< bool, bool > m_isBooked
std::map< std::string, std::shared_ptr< TH2Poly > > pxfTh2PolyForwardSummary
__shared__ Hist hist
Definition: DetId.h:17
std::vector< std::string > m_knownNames
unsigned int pxfSide(const DetId &id) const
void bookBarrelHistograms(const std::string &currentHistoName, const char *what, const char *zaxis)
void book(const std::string &currentHistoName, const char *what, const char *zaxis)
def remove(d, key, TELL=False)
Definition: MatrixUtil.py:223
void bookForwardHistograms(const std::string &currentHistoName, const char *what, const char *zaxis)
void fillBarrelBin(const std::string &currentHistoName, unsigned int id, double value)
def canvas(sub, attr)
Definition: svgfig.py:482
void setForwardScale(const std::string &currentHistoName, std::pair< float, float > extrema)
void fill(const std::string &currentHistoName, unsigned int id, double value)
#define LogDebug(id)
unsigned transform(const HcalDetId &id, unsigned transformCode)