CMS 3D CMS Logo

SiPixelVCal_PayloadInspector.cc
Go to the documentation of this file.
1 
17 
18 #include <memory>
19 #include <sstream>
20 #include <fmt/printf.h>
21 
22 // include ROOT
23 #include "TH2F.h"
24 #include "TH1F.h"
25 #include "TLegend.h"
26 #include "TCanvas.h"
27 #include "TLine.h"
28 #include "TStyle.h"
29 #include "TLatex.h"
30 #include "TPave.h"
31 #include "TPaveStats.h"
32 #include "TGaxis.h"
33 
34 namespace {
35 
36  using namespace cond::payloadInspector;
37 
38  namespace SiPixelVCalPI {
39  enum type { t_slope = 0, t_offset = 1 };
40  }
41 
42  /************************************************
43  1d histogram of SiPixelVCal of 1 IOV
44  *************************************************/
45  // inherit from one of the predefined plot class: Histogram1D
46  template <SiPixelVCalPI::type myType>
47  class SiPixelVCalValue : public Histogram1D<SiPixelVCal, SINGLE_IOV> {
48  public:
49  SiPixelVCalValue()
50  : Histogram1D<SiPixelVCal, SINGLE_IOV>("SiPixel VCal values",
51  "SiPixel VCal values",
52  100,
53  myType == SiPixelVCalPI::t_slope ? 40. : -700,
54  myType == SiPixelVCalPI::t_slope ? 60. : 0) {}
55 
56  bool fill() override {
57  auto tag = PlotBase::getTag<0>();
58  for (auto const &iov : tag.iovs) {
59  std::shared_ptr<SiPixelVCal> payload = Base::fetchPayload(std::get<1>(iov));
60  if (payload.get()) {
61  auto VCalMap_ = payload->getSlopeAndOffset();
62  for (const auto &element : VCalMap_) {
63  if (myType == SiPixelVCalPI::t_slope) {
64  fillWithValue(element.second.slope);
65  } else if (myType == SiPixelVCalPI::t_offset) {
66  fillWithValue(element.second.offset);
67  }
68  }
69  } // payload
70  } // iovs
71  return true;
72  } // fill
73  };
74 
75  using SiPixelVCalSlopeValue = SiPixelVCalValue<SiPixelVCalPI::t_slope>;
76  using SiPixelVCalOffsetValue = SiPixelVCalValue<SiPixelVCalPI::t_offset>;
77 
78  /************************************************
79  1d histogram of SiPixelVCal of 1 IOV
80  *************************************************/
81  class SiPixelVCalValues : public PlotImage<SiPixelVCal, SINGLE_IOV> {
82  public:
83  SiPixelVCalValues() : PlotImage<SiPixelVCal, SINGLE_IOV>("SiPixelVCal Values") {}
84 
85  bool fill() override {
86  gStyle->SetOptStat("emr");
87 
88  auto tag = PlotBase::getTag<0>();
89  auto iov = tag.iovs.front();
90  std::shared_ptr<SiPixelVCal> payload = fetchPayload(std::get<1>(iov));
91  auto VCalMap_ = payload->getSlopeAndOffset();
92 
93  auto slopes = payload->getAllSlopes();
94  auto offsets = payload->getAllOffsets();
95 
96  /*
97  std::transform(VCalMap_.begin(),
98  VCalMap_.end(),
99  std::inserter(slopes,slopes.end()),
100  [](std::pair<unsigned int, SiPixelVCal::VCal> vcalentry) -> std::pair<uint32_t,float> { return std::make_pair(vcalentry.first,vcalentry.second.slope); });
101 
102  */
103 
104  auto s_extrema = SiPixelPI::findMinMaxInMap(slopes);
105  auto o_extrema = SiPixelPI::findMinMaxInMap(offsets);
106 
107  auto o_range = (o_extrema.second - o_extrema.first) / 10.;
108 
109  TCanvas canvas("Canv", "Canv", 1000, 1000);
110  canvas.Divide(1, 2);
111  canvas.cd();
112  auto h_slope =
113  std::make_shared<TH1F>("slope value",
114  "SiPixel VCal slope value;SiPixel VCal slope value [ADC/VCal units];# modules",
115  50,
116  s_extrema.first * 0.9,
117  s_extrema.second * 1.1);
118 
119  auto h_offset = std::make_shared<TH1F>("offset value",
120  "SiPixel VCal offset value;SiPixel VCal offset value [ADC];# modules",
121  50,
122  o_extrema.first - o_range,
123  o_extrema.second + o_range);
124 
125  for (unsigned int i = 1; i <= 2; i++) {
126  SiPixelPI::adjustCanvasMargins(canvas.cd(i), 0.06, 0.12, 0.08, 0.03);
127  }
128 
129  for (const auto &slope : slopes) {
130  h_slope->Fill(slope.second);
131  }
132 
133  for (const auto &offset : offsets) {
134  h_offset->Fill(offset.second);
135  }
136 
137  canvas.cd(1);
138  adjustHisto(h_slope);
139  canvas.Update();
140 
141  TLegend legend = TLegend(0.40, 0.83, 0.94, 0.93);
142  legend.SetHeader(("Payload hash: #bf{" + (std::get<1>(iov)) + "}").c_str(),
143  "C"); // option "C" allows to center the header
144  legend.AddEntry(h_slope.get(), ("TAG: " + tag.name).c_str(), "F");
145  legend.SetTextSize(0.035);
146  legend.SetLineColor(10);
147  legend.Draw("same");
148 
149  TPaveStats *st = (TPaveStats *)h_slope->FindObject("stats");
150  st->SetTextSize(0.035);
151  SiPixelPI::adjustStats(st, 0.15, 0.83, 0.30, 0.93);
152 
153  auto ltx = TLatex();
154  ltx.SetTextFont(62);
155  ltx.SetTextSize(0.05);
156  ltx.SetTextAlign(11);
157  ltx.DrawLatexNDC(gPad->GetLeftMargin(),
158  1 - gPad->GetTopMargin() + 0.01,
159  ("SiPixel VCal Slope IOV:" + std::to_string(std::get<0>(iov))).c_str());
160 
161  canvas.cd(2);
162  adjustHisto(h_offset);
163  canvas.Update();
164  legend.Draw("same");
165 
166  ltx.DrawLatexNDC(gPad->GetLeftMargin(),
167  1 - gPad->GetTopMargin() + 0.01,
168  ("SiPixel VCal Offset IOV:" + std::to_string(std::get<0>(iov))).c_str());
169 
170  TPaveStats *st2 = (TPaveStats *)h_offset->FindObject("stats");
171  st2->SetTextSize(0.035);
172  SiPixelPI::adjustStats(st2, 0.15, 0.83, 0.30, 0.93);
173 
174  std::string fileName(m_imageFileName);
175  canvas.SaveAs(fileName.c_str());
176 
177  return true;
178  }
179 
180  private:
181  void adjustHisto(const std::shared_ptr<TH1F> &histo) {
182  histo->SetTitle("");
183  histo->GetYaxis()->SetRangeUser(0., histo->GetMaximum() * 1.30);
184  histo->SetFillColor(kRed);
185  histo->SetMarkerStyle(20);
186  histo->SetMarkerSize(1);
187  histo->Draw("bar2");
189  histo->SetStats(true);
190  histo->GetYaxis()->SetTitleOffset(0.9);
191  }
192  };
193 
194  /************************************************
195  1d histogram of SiPixelVCal of 1 IOV per region
196  *************************************************/
197  template <bool isBarrel, SiPixelVCalPI::type myType>
198  class SiPixelVCalValuesPerRegion : public PlotImage<SiPixelVCal, SINGLE_IOV> {
199  public:
200  SiPixelVCalValuesPerRegion() : PlotImage<SiPixelVCal, SINGLE_IOV>("SiPixelVCal Values per region") {}
201 
202  bool fill() override {
203  gStyle->SetOptStat("emr");
204 
205  auto tag = PlotBase::getTag<0>();
206  auto tagname = tag.name;
207  auto iov = tag.iovs.front();
208  std::shared_ptr<SiPixelVCal> payload = fetchPayload(std::get<1>(iov));
210  if (myType == SiPixelVCalPI::t_slope) {
211  Map_ = payload->getAllSlopes();
212  } else {
213  Map_ = payload->getAllOffsets();
214  }
215  auto extrema = SiPixelPI::findMinMaxInMap(Map_);
216  auto range = (extrema.second - extrema.first) / 10.;
217 
218  TCanvas canvas("Canv", "Canv", isBarrel ? 1400 : 1800, 1200);
219  if (Map_.size() > SiPixelPI::phase1size) {
221  std::string fileName(this->m_imageFileName);
222  canvas.SaveAs(fileName.c_str());
223  return false;
224  }
225 
226  canvas.cd();
227 
228  SiPixelPI::PhaseInfo phaseInfo(Map_.size());
229  const char *path_toTopologyXML = phaseInfo.pathToTopoXML();
230 
231  TrackerTopology tTopo =
233 
234  auto myPlots = PixelRegions::PixelRegionContainers(&tTopo, phaseInfo.phase());
235  myPlots.bookAll((myType == SiPixelVCalPI::t_slope) ? "SiPixel VCal slope value" : "SiPixel VCal offset value",
236  (myType == SiPixelVCalPI::t_slope) ? "SiPixel VCal slope value [ADC/VCal units]"
237  : "SiPixel VCal offset value [ADC]",
238  "#modules",
239  50,
240  extrema.first - range,
241  extrema.second + range);
242 
243  canvas.Modified();
244 
245  for (const auto &element : Map_) {
246  myPlots.fill(element.first, element.second);
247  }
248 
249  myPlots.beautify();
250  myPlots.draw(canvas, isBarrel);
251 
252  TLegend legend = TLegend(0.40, 0.88, 0.93, 0.90);
253  legend.SetHeader(("Hash: #bf{" + (std::get<1>(iov)) + "}").c_str(),
254  "C"); // option "C" allows to center the header
255  //legend.AddEntry(h1.get(), ("IOV: " + std::to_string(std::get<0>(iov))).c_str(), "PL");
256  legend.SetTextSize(0.025);
257  legend.SetLineColor(10);
258 
259  unsigned int maxPads = isBarrel ? 4 : 12;
260  for (unsigned int c = 1; c <= maxPads; c++) {
261  canvas.cd(c);
262  SiPixelPI::adjustCanvasMargins(canvas.cd(c), 0.06, 0.12, 0.12, 0.05);
263  legend.Draw("same");
264  canvas.cd(c)->Update();
265  }
266 
267  myPlots.stats();
268 
269  auto ltx = TLatex();
270  ltx.SetTextFont(62);
271  ltx.SetTextSize(0.05);
272  ltx.SetTextAlign(11);
273 
274  for (unsigned int c = 1; c <= maxPads; c++) {
275  auto index = isBarrel ? c - 1 : c + 3;
276 
277  canvas.cd(c);
278  ltx.DrawLatexNDC(gPad->GetLeftMargin(),
279  1 - gPad->GetTopMargin() + 0.01,
280  fmt::sprintf("%s, #color[2]{%s, IOV: %s}",
282  tagname,
283  std::to_string(std::get<0>(iov)))
284  .c_str());
285  }
286 
287  std::string fileName(m_imageFileName);
288  canvas.SaveAs(fileName.c_str());
289 
290  return true;
291  }
292  };
293 
294  using SiPixelVCalSlopeValuesBarrel = SiPixelVCalValuesPerRegion<true, SiPixelVCalPI::t_slope>;
295  using SiPixelVCalSlopeValuesEndcap = SiPixelVCalValuesPerRegion<false, SiPixelVCalPI::t_slope>;
296 
297  using SiPixelVCalOffsetValuesBarrel = SiPixelVCalValuesPerRegion<true, SiPixelVCalPI::t_offset>;
298  using SiPixelVCalOffsetValuesEndcap = SiPixelVCalValuesPerRegion<false, SiPixelVCalPI::t_offset>;
299 
300  /************************************************
301  1d histogram of SiPixelVCal (slope or offset) of 2 IOV per region
302  *************************************************/
303  template <bool isBarrel, SiPixelVCalPI::type myType, IOVMultiplicity nIOVs, int ntags>
304  class SiPixelVCalValuesCompareSubdet : public PlotImage<SiPixelVCal, nIOVs, ntags> {
305  public:
306  SiPixelVCalValuesCompareSubdet()
307  : PlotImage<SiPixelVCal, nIOVs, ntags>(Form("SiPixelVCal Values Comparisons by Subdet %i tags(s)", ntags)) {}
308 
309  bool fill() override {
310  gStyle->SetOptStat("emr");
311 
312  // trick to deal with the multi-ioved tag and two tag case at the same time
313  auto theIOVs = PlotBase::getTag<0>().iovs;
314  auto f_tagname = PlotBase::getTag<0>().name;
315  std::string l_tagname = "";
316  auto firstiov = theIOVs.front();
317  std::tuple<cond::Time_t, cond::Hash> lastiov;
318 
319  // we don't support (yet) comparison with more than 2 tags
320  assert(this->m_plotAnnotations.ntags < 3);
321 
322  if (this->m_plotAnnotations.ntags == 2) {
323  auto tag2iovs = PlotBase::getTag<1>().iovs;
324  l_tagname = PlotBase::getTag<1>().name;
325  lastiov = tag2iovs.front();
326  } else {
327  lastiov = theIOVs.back();
328  }
329 
332 
333  std::shared_ptr<SiPixelVCal> last_payload = this->fetchPayload(std::get<1>(lastiov));
334  if (myType == SiPixelVCalPI::t_slope) {
335  l_Map_ = last_payload->getAllSlopes();
336  } else {
337  l_Map_ = last_payload->getAllOffsets();
338  }
339  auto l_extrema = SiPixelPI::findMinMaxInMap(l_Map_);
340 
341  std::shared_ptr<SiPixelVCal> first_payload = this->fetchPayload(std::get<1>(firstiov));
342  if (myType == SiPixelVCalPI::t_slope) {
343  f_Map_ = first_payload->getAllSlopes();
344  } else {
345  f_Map_ = first_payload->getAllOffsets();
346  }
347 
348  auto f_extrema = SiPixelPI::findMinMaxInMap(f_Map_);
349 
350  auto max = (l_extrema.second > f_extrema.second) ? l_extrema.second : f_extrema.second;
351  auto min = (l_extrema.first < f_extrema.first) ? l_extrema.first : f_extrema.first;
352 
353  std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
354  std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
355 
356  TCanvas canvas("Canv", "Canv", isBarrel ? 1400 : 1800, 1200);
357  if ((f_Map_.size() > SiPixelPI::phase1size) || (l_Map_.size() > SiPixelPI::phase1size)) {
358  SiPixelPI::displayNotSupported(canvas, std::max(f_Map_.size(), l_Map_.size()));
359  std::string fileName(this->m_imageFileName);
360  canvas.SaveAs(fileName.c_str());
361  return false;
362  }
363 
364  canvas.cd();
365 
366  SiPixelPI::PhaseInfo l_phaseInfo(l_Map_.size());
367  SiPixelPI::PhaseInfo f_phaseInfo(f_Map_.size());
368 
369  // deal with last IOV
370 
371  const char *path_toTopologyXML = l_phaseInfo.pathToTopoXML();
372 
373  auto l_tTopo =
375 
376  auto delta = std::abs(max - min) / 10.;
377 
378  auto l_myPlots = PixelRegions::PixelRegionContainers(&l_tTopo, l_phaseInfo.phase());
379  l_myPlots.bookAll(
380  fmt::sprintf("SiPixel VCal %s,last", (myType == SiPixelVCalPI::t_slope ? "slope" : "offset")),
381  fmt::sprintf("SiPixel VCal %s", (myType == SiPixelVCalPI::t_slope ? " slope [ADC/VCal]" : " offset [ADC]")),
382  "#modules",
383  50,
384  min - delta,
385  max + delta);
386 
387  for (const auto &element : l_Map_) {
388  l_myPlots.fill(element.first, element.second);
389  }
390 
391  l_myPlots.beautify();
392  l_myPlots.draw(canvas, isBarrel, "bar2", true);
393 
394  // deal with first IOV
395 
396  path_toTopologyXML = f_phaseInfo.pathToTopoXML();
397 
398  auto f_tTopo =
400 
401  auto f_myPlots = PixelRegions::PixelRegionContainers(&f_tTopo, f_phaseInfo.phase());
402  f_myPlots.bookAll(
403  fmt::sprintf("SiPixel VCal %s,first", (myType == SiPixelVCalPI::t_slope ? "slope" : "offset")),
404  fmt::sprintf("SiPixel VCal %s", (myType == SiPixelVCalPI::t_slope ? " slope [ADC/VCal]" : " offset [ADC]")),
405  "#modules",
406  50,
407  min - delta,
408  max + delta);
409 
410  for (const auto &element : f_Map_) {
411  f_myPlots.fill(element.first, element.second);
412  }
413 
414  f_myPlots.beautify(kAzure, kBlue);
415  f_myPlots.draw(canvas, isBarrel, "HISTsames", true);
416 
417  // rescale the y-axis ranges in order to fit the canvas
418  l_myPlots.rescaleMax(f_myPlots);
419 
420  // done dealing with IOVs
421 
422  auto colorTag = isBarrel ? PixelRegions::L1 : PixelRegions::Rm1l;
423  std::unique_ptr<TLegend> legend;
424  if (this->m_plotAnnotations.ntags == 2) {
425  legend = std::make_unique<TLegend>(0.36, 0.86, 0.94, 0.92);
426  legend->AddEntry(l_myPlots.getHistoFromMap(colorTag).get(), ("#color[2]{" + l_tagname + "}").c_str(), "F");
427  legend->AddEntry(f_myPlots.getHistoFromMap(colorTag).get(), ("#color[4]{" + f_tagname + "}").c_str(), "F");
428  legend->SetTextSize(0.024);
429  } else {
430  legend = std::make_unique<TLegend>(0.58, 0.80, 0.90, 0.92);
431  legend->AddEntry(l_myPlots.getHistoFromMap(colorTag).get(), ("#color[2]{" + lastIOVsince + "}").c_str(), "F");
432  legend->AddEntry(f_myPlots.getHistoFromMap(colorTag).get(), ("#color[4]{" + firstIOVsince + "}").c_str(), "F");
433  legend->SetTextSize(0.040);
434  }
435  legend->SetLineColor(10);
436 
437  unsigned int maxPads = isBarrel ? 4 : 12;
438  for (unsigned int c = 1; c <= maxPads; c++) {
439  canvas.cd(c);
440  SiPixelPI::adjustCanvasMargins(canvas.cd(c), 0.06, 0.12, 0.12, 0.05);
441  legend->Draw("same");
442  canvas.cd(c)->Update();
443  }
444 
445  f_myPlots.stats(0);
446  l_myPlots.stats(1);
447 
448  auto ltx = TLatex();
449  ltx.SetTextFont(62);
450  ltx.SetTextSize(0.05);
451  ltx.SetTextAlign(11);
452 
453  for (unsigned int c = 1; c <= maxPads; c++) {
454  auto index = isBarrel ? c - 1 : c + 3;
455  canvas.cd(c);
456  ltx.DrawLatexNDC(gPad->GetLeftMargin(),
457  1 - gPad->GetTopMargin() + 0.01,
458  (PixelRegions::IDlabels.at(index) + " : #color[4]{" + std::to_string(std::get<0>(firstiov)) +
459  "} vs #color[2]{" + std::to_string(std::get<0>(lastiov)) + "}")
460  .c_str());
461  }
462 
463  std::string fileName(this->m_imageFileName);
464  canvas.SaveAs(fileName.c_str());
465 
466  return true;
467  }
468  };
469 
470  using SiPixelVCalSlopesBarrelCompareSingleTag =
471  SiPixelVCalValuesCompareSubdet<true, SiPixelVCalPI::t_slope, MULTI_IOV, 1>;
472  using SiPixelVCalOffsetsBarrelCompareSingleTag =
473  SiPixelVCalValuesCompareSubdet<true, SiPixelVCalPI::t_offset, MULTI_IOV, 1>;
474 
475  using SiPixelVCalSlopesEndcapCompareSingleTag =
476  SiPixelVCalValuesCompareSubdet<false, SiPixelVCalPI::t_slope, MULTI_IOV, 1>;
477  using SiPixelVCalOffsetsEndcapCompareSingleTag =
478  SiPixelVCalValuesCompareSubdet<false, SiPixelVCalPI::t_offset, MULTI_IOV, 1>;
479 
480  using SiPixelVCalSlopesBarrelCompareTwoTags =
481  SiPixelVCalValuesCompareSubdet<true, SiPixelVCalPI::t_slope, SINGLE_IOV, 2>;
482  using SiPixelVCalOffsetsBarrelCompareTwoTags =
483  SiPixelVCalValuesCompareSubdet<true, SiPixelVCalPI::t_offset, SINGLE_IOV, 2>;
484 
485  using SiPixelVCalSlopesEndcapCompareTwoTags =
486  SiPixelVCalValuesCompareSubdet<false, SiPixelVCalPI::t_slope, SINGLE_IOV, 2>;
487  using SiPixelVCalOffsetsEndcapCompareTwoTags =
488  SiPixelVCalValuesCompareSubdet<false, SiPixelVCalPI::t_offset, SINGLE_IOV, 2>;
489 
490  /************************************************
491  1d histogram of SiPixelVCal of 1 IOV
492  *************************************************/
493 
494  template <SiPixelVCalPI::type myType, IOVMultiplicity nIOVs, int ntags>
495  class SiPixelVCalValueComparisonBase : public PlotImage<SiPixelVCal, nIOVs, ntags> {
496  public:
497  SiPixelVCalValueComparisonBase()
498  : PlotImage<SiPixelVCal, nIOVs, ntags>(Form("SiPixelVCal Synoptic Values Comparison %i tag(s)", ntags)) {}
499 
500  bool fill() override {
501  TH1F::SetDefaultSumw2(true);
502 
503  // trick to deal with the multi-ioved tag and two tag case at the same time
504  auto theIOVs = PlotBase::getTag<0>().iovs;
505  auto f_tagname = PlotBase::getTag<0>().name;
506  std::string l_tagname = "";
507  auto firstiov = theIOVs.front();
508  std::tuple<cond::Time_t, cond::Hash> lastiov;
509 
510  // we don't support (yet) comparison with more than 2 tags
511  assert(this->m_plotAnnotations.ntags < 3);
512 
513  if (this->m_plotAnnotations.ntags == 2) {
514  auto tag2iovs = PlotBase::getTag<1>().iovs;
515  l_tagname = PlotBase::getTag<1>().name;
516  lastiov = tag2iovs.front();
517  } else {
518  lastiov = theIOVs.back();
519  }
520 
523 
524  std::shared_ptr<SiPixelVCal> last_payload = this->fetchPayload(std::get<1>(lastiov));
525  if (myType == SiPixelVCalPI::t_slope) {
526  l_Map_ = last_payload->getAllSlopes();
527  } else {
528  l_Map_ = last_payload->getAllOffsets();
529  }
530  auto l_extrema = SiPixelPI::findMinMaxInMap(l_Map_);
531 
532  std::shared_ptr<SiPixelVCal> first_payload = this->fetchPayload(std::get<1>(firstiov));
533  if (myType == SiPixelVCalPI::t_slope) {
534  f_Map_ = first_payload->getAllSlopes();
535  } else {
536  f_Map_ = first_payload->getAllOffsets();
537  }
538  auto f_extrema = SiPixelPI::findMinMaxInMap(f_Map_);
539 
540  auto max = (l_extrema.second > f_extrema.second) ? l_extrema.second : f_extrema.second;
541  auto min = (l_extrema.first < f_extrema.first) ? l_extrema.first : f_extrema.first;
542  auto delta = std::abs(max - min) / 10.;
543 
544  std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
545  std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
546 
547  TCanvas canvas("Canv", "Canv", 1200, 1000);
548  canvas.cd();
549  auto hfirst = std::unique_ptr<TH1F>(
550  new TH1F("value_first",
551  fmt::sprintf("SiPixel VCal %s value;SiPixel VCal %s ;# modules",
552  (myType == SiPixelVCalPI::t_slope ? "slope" : "offset"),
553  (myType == SiPixelVCalPI::t_slope ? " slope [ADC/VCal]" : " offset [ADC]"))
554  .c_str(),
555  50,
556  min - delta,
557  max + delta));
558  hfirst->SetStats(false);
559 
560  auto hlast = std::unique_ptr<TH1F>(
561  new TH1F("value_last",
562  fmt::sprintf("SiPixel VCal %s value;SiPixel VCal %s ;# modules",
563  (myType == SiPixelVCalPI::t_slope ? "slope" : "offset"),
564  (myType == SiPixelVCalPI::t_slope ? " slope [ADC/VCal]" : " offset [ADC]"))
565  .c_str(),
566  50,
567  min - delta,
568  max + delta));
569  hlast->SetStats(false);
570 
571  SiPixelPI::adjustCanvasMargins(canvas.cd(), 0.06, 0.12, 0.12, 0.05);
572  canvas.Modified();
573 
574  for (const auto &element : f_Map_) {
575  hfirst->Fill(element.second);
576  }
577 
578  for (const auto &element : l_Map_) {
579  hlast->Fill(element.second);
580  }
581 
582  auto extrema = SiPixelPI::getExtrema(hfirst.get(), hlast.get());
583  hfirst->GetYaxis()->SetRangeUser(extrema.first, extrema.second * 1.10);
584 
585  hfirst->SetTitle("");
586  hfirst->SetFillColor(kRed);
587  hfirst->SetBarWidth(0.95);
588  hfirst->Draw("histbar");
589 
590  hlast->SetTitle("");
591  hlast->SetFillColorAlpha(kBlue, 0.20);
592  hlast->SetBarWidth(0.95);
593  hlast->Draw("histbarsame");
594 
595  SiPixelPI::makeNicePlotStyle(hfirst.get());
596  SiPixelPI::makeNicePlotStyle(hlast.get());
597 
598  canvas.Update();
599 
600  TLegend legend = TLegend(0.30, 0.86, 0.95, 0.94);
601  legend.AddEntry(hfirst.get(), ("payload: #color[2]{" + std::get<1>(firstiov) + "}").c_str(), "F");
602  legend.AddEntry(hlast.get(), ("payload: #color[4]{" + std::get<1>(lastiov) + "}").c_str(), "F");
603  legend.SetTextSize(0.025);
604  legend.Draw("same");
605 
606  auto ltx = TLatex();
607  ltx.SetTextFont(62);
608  ltx.SetTextSize(0.040);
609  ltx.SetTextAlign(11);
610  std::string ltxText;
611  if (this->m_plotAnnotations.ntags == 2) {
612  ltxText = fmt::sprintf("#color[2]{%s, %s} vs #color[4]{%s, %s}",
613  f_tagname,
614  std::to_string(std::get<0>(firstiov)),
615  l_tagname,
616  std::to_string(std::get<0>(lastiov)));
617  } else {
618  ltxText = fmt::sprintf("%s IOV: #color[2]{%s} vs IOV: #color[4]{%s}",
619  f_tagname,
620  std::to_string(std::get<0>(firstiov)),
621  std::to_string(std::get<0>(lastiov)));
622  }
623  ltx.DrawLatexNDC(gPad->GetLeftMargin(), 1 - gPad->GetTopMargin() + 0.01, ltxText.c_str());
624 
625  std::string fileName(this->m_imageFileName);
626  canvas.SaveAs(fileName.c_str());
627 
628  return true;
629  }
630  };
631 
632  using SiPixelVCalSlopesComparisonSingleTag = SiPixelVCalValueComparisonBase<SiPixelVCalPI::t_slope, MULTI_IOV, 1>;
633  using SiPixelVCalOffsetsComparisonSingleTag = SiPixelVCalValueComparisonBase<SiPixelVCalPI::t_offset, MULTI_IOV, 1>;
634  using SiPixelVCalSlopesComparisonTwoTags = SiPixelVCalValueComparisonBase<SiPixelVCalPI::t_slope, SINGLE_IOV, 2>;
635  using SiPixelVCalOffsetsComparisonTwoTags = SiPixelVCalValueComparisonBase<SiPixelVCalPI::t_offset, SINGLE_IOV, 2>;
636 
637 } // namespace
638 
640  PAYLOAD_INSPECTOR_CLASS(SiPixelVCalSlopeValue);
641  PAYLOAD_INSPECTOR_CLASS(SiPixelVCalOffsetValue);
642  PAYLOAD_INSPECTOR_CLASS(SiPixelVCalValues);
643  PAYLOAD_INSPECTOR_CLASS(SiPixelVCalSlopeValuesBarrel);
644  PAYLOAD_INSPECTOR_CLASS(SiPixelVCalSlopeValuesEndcap);
645  PAYLOAD_INSPECTOR_CLASS(SiPixelVCalOffsetValuesBarrel);
646  PAYLOAD_INSPECTOR_CLASS(SiPixelVCalOffsetValuesEndcap);
647  PAYLOAD_INSPECTOR_CLASS(SiPixelVCalSlopesBarrelCompareSingleTag);
648  PAYLOAD_INSPECTOR_CLASS(SiPixelVCalOffsetsBarrelCompareSingleTag);
649  PAYLOAD_INSPECTOR_CLASS(SiPixelVCalSlopesEndcapCompareSingleTag);
650  PAYLOAD_INSPECTOR_CLASS(SiPixelVCalOffsetsEndcapCompareSingleTag);
651  PAYLOAD_INSPECTOR_CLASS(SiPixelVCalSlopesBarrelCompareTwoTags);
652  PAYLOAD_INSPECTOR_CLASS(SiPixelVCalOffsetsBarrelCompareTwoTags);
653  PAYLOAD_INSPECTOR_CLASS(SiPixelVCalSlopesEndcapCompareTwoTags);
654  PAYLOAD_INSPECTOR_CLASS(SiPixelVCalOffsetsEndcapCompareTwoTags);
655  PAYLOAD_INSPECTOR_CLASS(SiPixelVCalSlopesComparisonSingleTag);
656  PAYLOAD_INSPECTOR_CLASS(SiPixelVCalOffsetsComparisonSingleTag);
657  PAYLOAD_INSPECTOR_CLASS(SiPixelVCalSlopesComparisonTwoTags);
658  PAYLOAD_INSPECTOR_CLASS(SiPixelVCalOffsetsComparisonTwoTags);
659 }
const std::vector< std::string > IDlabels
static const double slope[3]
std::string to_string(const V &value)
Definition: OMSAccess.h:77
std::pair< float, float > getExtrema(TH1 *h1, TH1 *h2)
assert(be >=bs)
#define PAYLOAD_INSPECTOR_CLASS(CLASS_NAME)
std::pair< T, T > findMinMaxInMap(const std::map< unsigned int, T > &theMap)
static const unsigned int phase1size
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define PAYLOAD_INSPECTOR_MODULE(PAYLOAD_TYPENAME)
void adjustCanvasMargins(TVirtualPad *pad, float top, float bottom, float left, float right)
std::map< uint32_t, float > mapToDetId
Definition: SiPixelVCal.h:14
void adjustStats(TPaveStats *stats, float X1, float Y1, float X2, float Y2)
void bookAll(std::string title_label, std::string x_label, std::string y_label, const int nbins, const float xmin, const float xmax)
void makeNicePlotStyle(TH1 *hist)
void displayNotSupported(TCanvas &canv, const unsigned int size)
TrackerTopology fromTrackerParametersXMLFile(const std::string &xmlFileName)
def canvas(sub, attr)
Definition: svgfig.py:482