File indexing completed on 2024-09-07 04:35:33
0001
0002
0003
0004
0005
0006
0007
0008
0009 #include "CondCore/Utilities/interface/PayloadInspectorModule.h"
0010 #include "CondCore/Utilities/interface/PayloadInspector.h"
0011 #include "CondCore/CondDB/interface/Time.h"
0012 #include "CondFormats/SiPixelObjects/interface/SiPixelVCal.h"
0013 #include "FWCore/ParameterSet/interface/FileInPath.h"
0014 #include "CalibTracker/StandaloneTrackerTopology/interface/StandaloneTrackerTopology.h"
0015 #include "CondCore/SiPixelPlugins/interface/SiPixelPayloadInspectorHelper.h"
0016 #include "CondCore/SiPixelPlugins/interface/PixelRegionContainers.h"
0017
0018 #include <memory>
0019 #include <sstream>
0020 #include <fmt/printf.h>
0021
0022
0023 #include "TH2F.h"
0024 #include "TH1F.h"
0025 #include "TLegend.h"
0026 #include "TCanvas.h"
0027 #include "TLine.h"
0028 #include "TStyle.h"
0029 #include "TLatex.h"
0030 #include "TPave.h"
0031 #include "TPaveStats.h"
0032 #include "TGaxis.h"
0033
0034 namespace {
0035
0036 using namespace cond::payloadInspector;
0037
0038 namespace SiPixelVCalPI {
0039 enum type { t_slope = 0, t_offset = 1 };
0040 }
0041
0042
0043
0044
0045
0046 template <SiPixelVCalPI::type myType>
0047 class SiPixelVCalValue : public Histogram1D<SiPixelVCal, SINGLE_IOV> {
0048 public:
0049 SiPixelVCalValue()
0050 : Histogram1D<SiPixelVCal, SINGLE_IOV>("SiPixel VCal values",
0051 "SiPixel VCal values",
0052 100,
0053 myType == SiPixelVCalPI::t_slope ? 40. : -700,
0054 myType == SiPixelVCalPI::t_slope ? 60. : 0) {}
0055
0056 bool fill() override {
0057 auto tag = PlotBase::getTag<0>();
0058 for (auto const &iov : tag.iovs) {
0059 std::shared_ptr<SiPixelVCal> payload = Base::fetchPayload(std::get<1>(iov));
0060 if (payload.get()) {
0061 auto VCalMap_ = payload->getSlopeAndOffset();
0062 for (const auto &element : VCalMap_) {
0063 if (myType == SiPixelVCalPI::t_slope) {
0064 fillWithValue(element.second.slope);
0065 } else if (myType == SiPixelVCalPI::t_offset) {
0066 fillWithValue(element.second.offset);
0067 }
0068 }
0069 }
0070 }
0071 return true;
0072 }
0073 };
0074
0075 using SiPixelVCalSlopeValue = SiPixelVCalValue<SiPixelVCalPI::t_slope>;
0076 using SiPixelVCalOffsetValue = SiPixelVCalValue<SiPixelVCalPI::t_offset>;
0077
0078
0079
0080
0081 class SiPixelVCalValues : public PlotImage<SiPixelVCal, SINGLE_IOV> {
0082 public:
0083 SiPixelVCalValues() : PlotImage<SiPixelVCal, SINGLE_IOV>("SiPixelVCal Values") {}
0084
0085 bool fill() override {
0086 gStyle->SetOptStat("emr");
0087
0088 auto tag = PlotBase::getTag<0>();
0089 auto iov = tag.iovs.front();
0090 std::shared_ptr<SiPixelVCal> payload = fetchPayload(std::get<1>(iov));
0091 auto VCalMap_ = payload->getSlopeAndOffset();
0092
0093 auto slopes = payload->getAllSlopes();
0094 auto offsets = payload->getAllOffsets();
0095
0096
0097
0098
0099
0100
0101
0102
0103
0104 auto s_extrema = SiPixelPI::findMinMaxInMap(slopes);
0105 auto o_extrema = SiPixelPI::findMinMaxInMap(offsets);
0106
0107 auto o_range = (o_extrema.second - o_extrema.first) / 10.;
0108
0109 TCanvas canvas("Canv", "Canv", 1000, 1000);
0110 canvas.Divide(1, 2);
0111 canvas.cd();
0112 auto h_slope =
0113 std::make_shared<TH1F>("slope value",
0114 "SiPixel VCal slope value;SiPixel VCal slope value [ADC/VCal units];# modules",
0115 50,
0116 s_extrema.first * 0.9,
0117 s_extrema.second * 1.1);
0118
0119 auto h_offset = std::make_shared<TH1F>("offset value",
0120 "SiPixel VCal offset value;SiPixel VCal offset value [ADC];# modules",
0121 50,
0122 o_extrema.first - o_range,
0123 o_extrema.second + o_range);
0124
0125 for (unsigned int i = 1; i <= 2; i++) {
0126 SiPixelPI::adjustCanvasMargins(canvas.cd(i), 0.06, 0.12, 0.08, 0.03);
0127 }
0128
0129 for (const auto &slope : slopes) {
0130 h_slope->Fill(slope.second);
0131 }
0132
0133 for (const auto &offset : offsets) {
0134 h_offset->Fill(offset.second);
0135 }
0136
0137 canvas.cd(1);
0138 adjustHisto(h_slope);
0139 canvas.Update();
0140
0141 TLegend legend = TLegend(0.40, 0.83, 0.94, 0.93);
0142 legend.SetHeader(("Payload hash: #bf{" + (std::get<1>(iov)) + "}").c_str(),
0143 "C");
0144 legend.AddEntry(h_slope.get(), ("TAG: " + tag.name).c_str(), "F");
0145 legend.SetTextSize(0.035);
0146 legend.SetLineColor(10);
0147 legend.Draw("same");
0148
0149 TPaveStats *st = (TPaveStats *)h_slope->FindObject("stats");
0150 st->SetTextSize(0.035);
0151 SiPixelPI::adjustStats(st, 0.15, 0.83, 0.30, 0.93);
0152
0153 auto ltx = TLatex();
0154 ltx.SetTextFont(62);
0155 ltx.SetTextSize(0.05);
0156 ltx.SetTextAlign(11);
0157 ltx.DrawLatexNDC(gPad->GetLeftMargin(),
0158 1 - gPad->GetTopMargin() + 0.01,
0159 ("SiPixel VCal Slope IOV:" + std::to_string(std::get<0>(iov))).c_str());
0160
0161 canvas.cd(2);
0162 adjustHisto(h_offset);
0163 canvas.Update();
0164 legend.Draw("same");
0165
0166 ltx.DrawLatexNDC(gPad->GetLeftMargin(),
0167 1 - gPad->GetTopMargin() + 0.01,
0168 ("SiPixel VCal Offset IOV:" + std::to_string(std::get<0>(iov))).c_str());
0169
0170 TPaveStats *st2 = (TPaveStats *)h_offset->FindObject("stats");
0171 st2->SetTextSize(0.035);
0172 SiPixelPI::adjustStats(st2, 0.15, 0.83, 0.30, 0.93);
0173
0174 std::string fileName(m_imageFileName);
0175 canvas.SaveAs(fileName.c_str());
0176
0177 return true;
0178 }
0179
0180 private:
0181 void adjustHisto(const std::shared_ptr<TH1F> &histo) {
0182 histo->SetTitle("");
0183 histo->GetYaxis()->SetRangeUser(0., histo->GetMaximum() * 1.30);
0184 histo->SetFillColor(kRed);
0185 histo->SetMarkerStyle(20);
0186 histo->SetMarkerSize(1);
0187 histo->Draw("bar2");
0188 SiPixelPI::makeNicePlotStyle(histo.get());
0189 histo->SetStats(true);
0190 histo->GetYaxis()->SetTitleOffset(0.9);
0191 }
0192 };
0193
0194
0195
0196
0197 template <bool isBarrel, SiPixelVCalPI::type myType>
0198 class SiPixelVCalValuesPerRegion : public PlotImage<SiPixelVCal, SINGLE_IOV> {
0199 public:
0200 SiPixelVCalValuesPerRegion() : PlotImage<SiPixelVCal, SINGLE_IOV>("SiPixelVCal Values per region") {}
0201
0202 bool fill() override {
0203 gStyle->SetOptStat("emr");
0204
0205 auto tag = PlotBase::getTag<0>();
0206 auto tagname = tag.name;
0207 auto iov = tag.iovs.front();
0208 std::shared_ptr<SiPixelVCal> payload = fetchPayload(std::get<1>(iov));
0209 SiPixelVCal::mapToDetId Map_;
0210 if (myType == SiPixelVCalPI::t_slope) {
0211 Map_ = payload->getAllSlopes();
0212 } else {
0213 Map_ = payload->getAllOffsets();
0214 }
0215 auto extrema = SiPixelPI::findMinMaxInMap(Map_);
0216 auto range = (extrema.second - extrema.first) / 10.;
0217
0218 TCanvas canvas("Canv", "Canv", isBarrel ? 1400 : 1800, 1200);
0219 if (Map_.size() > SiPixelPI::phase1size) {
0220 SiPixelPI::displayNotSupported(canvas, Map_.size());
0221 std::string fileName(this->m_imageFileName);
0222 canvas.SaveAs(fileName.c_str());
0223 return false;
0224 }
0225
0226 canvas.cd();
0227
0228 SiPixelPI::PhaseInfo phaseInfo(Map_.size());
0229 const char *path_toTopologyXML = phaseInfo.pathToTopoXML();
0230
0231 TrackerTopology tTopo =
0232 StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
0233
0234 auto myPlots = PixelRegions::PixelRegionContainers(&tTopo, phaseInfo.phase());
0235 myPlots.bookAll((myType == SiPixelVCalPI::t_slope) ? "SiPixel VCal slope value" : "SiPixel VCal offset value",
0236 (myType == SiPixelVCalPI::t_slope) ? "SiPixel VCal slope value [ADC/VCal units]"
0237 : "SiPixel VCal offset value [ADC]",
0238 "#modules",
0239 50,
0240 extrema.first - range,
0241 extrema.second + range);
0242
0243 canvas.Modified();
0244
0245 for (const auto &element : Map_) {
0246 myPlots.fill(element.first, element.second);
0247 }
0248
0249 myPlots.beautify();
0250 myPlots.draw(canvas, isBarrel);
0251
0252 TLegend legend = TLegend(0.40, 0.88, 0.93, 0.90);
0253 legend.SetHeader(("Hash: #bf{" + (std::get<1>(iov)) + "}").c_str(),
0254 "C");
0255
0256 legend.SetTextSize(0.025);
0257 legend.SetLineColor(10);
0258
0259 unsigned int maxPads = isBarrel ? 4 : 12;
0260 for (unsigned int c = 1; c <= maxPads; c++) {
0261 canvas.cd(c);
0262 SiPixelPI::adjustCanvasMargins(canvas.cd(c), 0.06, 0.12, 0.12, 0.05);
0263 legend.Draw("same");
0264 canvas.cd(c)->Update();
0265 }
0266
0267 myPlots.stats();
0268
0269 auto ltx = TLatex();
0270 ltx.SetTextFont(62);
0271 ltx.SetTextSize(0.05);
0272 ltx.SetTextAlign(11);
0273
0274 for (unsigned int c = 1; c <= maxPads; c++) {
0275 auto index = isBarrel ? c - 1 : c + 3;
0276
0277 canvas.cd(c);
0278 ltx.DrawLatexNDC(gPad->GetLeftMargin(),
0279 1 - gPad->GetTopMargin() + 0.01,
0280 fmt::sprintf("%s, #color[2]{%s, IOV: %s}",
0281 PixelRegions::IDlabels.at(index),
0282 tagname,
0283 std::to_string(std::get<0>(iov)))
0284 .c_str());
0285 }
0286
0287 std::string fileName(m_imageFileName);
0288 canvas.SaveAs(fileName.c_str());
0289
0290 return true;
0291 }
0292 };
0293
0294 using SiPixelVCalSlopeValuesBarrel = SiPixelVCalValuesPerRegion<true, SiPixelVCalPI::t_slope>;
0295 using SiPixelVCalSlopeValuesEndcap = SiPixelVCalValuesPerRegion<false, SiPixelVCalPI::t_slope>;
0296
0297 using SiPixelVCalOffsetValuesBarrel = SiPixelVCalValuesPerRegion<true, SiPixelVCalPI::t_offset>;
0298 using SiPixelVCalOffsetValuesEndcap = SiPixelVCalValuesPerRegion<false, SiPixelVCalPI::t_offset>;
0299
0300
0301
0302
0303 template <bool isBarrel, SiPixelVCalPI::type myType, IOVMultiplicity nIOVs, int ntags>
0304 class SiPixelVCalValuesCompareSubdet : public PlotImage<SiPixelVCal, nIOVs, ntags> {
0305 public:
0306 SiPixelVCalValuesCompareSubdet()
0307 : PlotImage<SiPixelVCal, nIOVs, ntags>(Form("SiPixelVCal Values Comparisons by Subdet %i tags(s)", ntags)) {}
0308
0309 bool fill() override {
0310 gStyle->SetOptStat("emr");
0311
0312
0313 auto theIOVs = PlotBase::getTag<0>().iovs;
0314 auto f_tagname = PlotBase::getTag<0>().name;
0315 std::string l_tagname = "";
0316 auto firstiov = theIOVs.front();
0317 std::tuple<cond::Time_t, cond::Hash> lastiov;
0318
0319
0320 assert(this->m_plotAnnotations.ntags < 3);
0321
0322 if (this->m_plotAnnotations.ntags == 2) {
0323 auto tag2iovs = PlotBase::getTag<1>().iovs;
0324 l_tagname = PlotBase::getTag<1>().name;
0325 lastiov = tag2iovs.front();
0326 } else {
0327 lastiov = theIOVs.back();
0328 }
0329
0330 SiPixelVCal::mapToDetId l_Map_;
0331 SiPixelVCal::mapToDetId f_Map_;
0332
0333 std::shared_ptr<SiPixelVCal> last_payload = this->fetchPayload(std::get<1>(lastiov));
0334 if (myType == SiPixelVCalPI::t_slope) {
0335 l_Map_ = last_payload->getAllSlopes();
0336 } else {
0337 l_Map_ = last_payload->getAllOffsets();
0338 }
0339 auto l_extrema = SiPixelPI::findMinMaxInMap(l_Map_);
0340
0341 std::shared_ptr<SiPixelVCal> first_payload = this->fetchPayload(std::get<1>(firstiov));
0342 if (myType == SiPixelVCalPI::t_slope) {
0343 f_Map_ = first_payload->getAllSlopes();
0344 } else {
0345 f_Map_ = first_payload->getAllOffsets();
0346 }
0347
0348 auto f_extrema = SiPixelPI::findMinMaxInMap(f_Map_);
0349
0350 auto max = (l_extrema.second > f_extrema.second) ? l_extrema.second : f_extrema.second;
0351 auto min = (l_extrema.first < f_extrema.first) ? l_extrema.first : f_extrema.first;
0352
0353 std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
0354 std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
0355
0356 TCanvas canvas("Canv", "Canv", isBarrel ? 1400 : 1800, 1200);
0357 if ((f_Map_.size() > SiPixelPI::phase1size) || (l_Map_.size() > SiPixelPI::phase1size)) {
0358 SiPixelPI::displayNotSupported(canvas, std::max(f_Map_.size(), l_Map_.size()));
0359 std::string fileName(this->m_imageFileName);
0360 canvas.SaveAs(fileName.c_str());
0361 return false;
0362 }
0363
0364 canvas.cd();
0365
0366 SiPixelPI::PhaseInfo l_phaseInfo(l_Map_.size());
0367 SiPixelPI::PhaseInfo f_phaseInfo(f_Map_.size());
0368
0369
0370
0371 const char *path_toTopologyXML = l_phaseInfo.pathToTopoXML();
0372
0373 auto l_tTopo =
0374 StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
0375
0376 auto delta = std::abs(max - min) / 10.;
0377
0378 auto l_myPlots = PixelRegions::PixelRegionContainers(&l_tTopo, l_phaseInfo.phase());
0379 l_myPlots.bookAll(
0380 fmt::sprintf("SiPixel VCal %s,last", (myType == SiPixelVCalPI::t_slope ? "slope" : "offset")),
0381 fmt::sprintf("SiPixel VCal %s", (myType == SiPixelVCalPI::t_slope ? " slope [ADC/VCal]" : " offset [ADC]")),
0382 "#modules",
0383 50,
0384 min - delta,
0385 max + delta);
0386
0387 for (const auto &element : l_Map_) {
0388 l_myPlots.fill(element.first, element.second);
0389 }
0390
0391 l_myPlots.beautify();
0392 l_myPlots.draw(canvas, isBarrel, "bar2", true);
0393
0394
0395
0396 path_toTopologyXML = f_phaseInfo.pathToTopoXML();
0397
0398 auto f_tTopo =
0399 StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
0400
0401 auto f_myPlots = PixelRegions::PixelRegionContainers(&f_tTopo, f_phaseInfo.phase());
0402 f_myPlots.bookAll(
0403 fmt::sprintf("SiPixel VCal %s,first", (myType == SiPixelVCalPI::t_slope ? "slope" : "offset")),
0404 fmt::sprintf("SiPixel VCal %s", (myType == SiPixelVCalPI::t_slope ? " slope [ADC/VCal]" : " offset [ADC]")),
0405 "#modules",
0406 50,
0407 min - delta,
0408 max + delta);
0409
0410 for (const auto &element : f_Map_) {
0411 f_myPlots.fill(element.first, element.second);
0412 }
0413
0414 f_myPlots.beautify(kAzure, kBlue);
0415 f_myPlots.draw(canvas, isBarrel, "HISTsames", true);
0416
0417
0418 l_myPlots.rescaleMax(f_myPlots);
0419
0420
0421
0422 auto colorTag = isBarrel ? PixelRegions::L1 : PixelRegions::Rm1l;
0423 std::unique_ptr<TLegend> legend;
0424 if (this->m_plotAnnotations.ntags == 2) {
0425 legend = std::make_unique<TLegend>(0.36, 0.86, 0.94, 0.92);
0426 legend->AddEntry(l_myPlots.getHistoFromMap(colorTag).get(), ("#color[2]{" + l_tagname + "}").c_str(), "F");
0427 legend->AddEntry(f_myPlots.getHistoFromMap(colorTag).get(), ("#color[4]{" + f_tagname + "}").c_str(), "F");
0428 legend->SetTextSize(0.024);
0429 } else {
0430 legend = std::make_unique<TLegend>(0.58, 0.80, 0.90, 0.92);
0431 legend->AddEntry(l_myPlots.getHistoFromMap(colorTag).get(), ("#color[2]{" + lastIOVsince + "}").c_str(), "F");
0432 legend->AddEntry(f_myPlots.getHistoFromMap(colorTag).get(), ("#color[4]{" + firstIOVsince + "}").c_str(), "F");
0433 legend->SetTextSize(0.040);
0434 }
0435 legend->SetLineColor(10);
0436
0437 unsigned int maxPads = isBarrel ? 4 : 12;
0438 for (unsigned int c = 1; c <= maxPads; c++) {
0439 canvas.cd(c);
0440 SiPixelPI::adjustCanvasMargins(canvas.cd(c), 0.06, 0.12, 0.12, 0.05);
0441 legend->Draw("same");
0442 canvas.cd(c)->Update();
0443 }
0444
0445 f_myPlots.stats(0);
0446 l_myPlots.stats(1);
0447
0448 auto ltx = TLatex();
0449 ltx.SetTextFont(62);
0450 ltx.SetTextSize(0.05);
0451 ltx.SetTextAlign(11);
0452
0453 for (unsigned int c = 1; c <= maxPads; c++) {
0454 auto index = isBarrel ? c - 1 : c + 3;
0455 canvas.cd(c);
0456 ltx.DrawLatexNDC(gPad->GetLeftMargin(),
0457 1 - gPad->GetTopMargin() + 0.01,
0458 (PixelRegions::IDlabels.at(index) + " : #color[4]{" + std::to_string(std::get<0>(firstiov)) +
0459 "} vs #color[2]{" + std::to_string(std::get<0>(lastiov)) + "}")
0460 .c_str());
0461 }
0462
0463 std::string fileName(this->m_imageFileName);
0464 canvas.SaveAs(fileName.c_str());
0465
0466 return true;
0467 }
0468 };
0469
0470 using SiPixelVCalSlopesBarrelCompareSingleTag =
0471 SiPixelVCalValuesCompareSubdet<true, SiPixelVCalPI::t_slope, MULTI_IOV, 1>;
0472 using SiPixelVCalOffsetsBarrelCompareSingleTag =
0473 SiPixelVCalValuesCompareSubdet<true, SiPixelVCalPI::t_offset, MULTI_IOV, 1>;
0474
0475 using SiPixelVCalSlopesEndcapCompareSingleTag =
0476 SiPixelVCalValuesCompareSubdet<false, SiPixelVCalPI::t_slope, MULTI_IOV, 1>;
0477 using SiPixelVCalOffsetsEndcapCompareSingleTag =
0478 SiPixelVCalValuesCompareSubdet<false, SiPixelVCalPI::t_offset, MULTI_IOV, 1>;
0479
0480 using SiPixelVCalSlopesBarrelCompareTwoTags =
0481 SiPixelVCalValuesCompareSubdet<true, SiPixelVCalPI::t_slope, SINGLE_IOV, 2>;
0482 using SiPixelVCalOffsetsBarrelCompareTwoTags =
0483 SiPixelVCalValuesCompareSubdet<true, SiPixelVCalPI::t_offset, SINGLE_IOV, 2>;
0484
0485 using SiPixelVCalSlopesEndcapCompareTwoTags =
0486 SiPixelVCalValuesCompareSubdet<false, SiPixelVCalPI::t_slope, SINGLE_IOV, 2>;
0487 using SiPixelVCalOffsetsEndcapCompareTwoTags =
0488 SiPixelVCalValuesCompareSubdet<false, SiPixelVCalPI::t_offset, SINGLE_IOV, 2>;
0489
0490
0491
0492
0493
0494 template <SiPixelVCalPI::type myType, IOVMultiplicity nIOVs, int ntags>
0495 class SiPixelVCalValueComparisonBase : public PlotImage<SiPixelVCal, nIOVs, ntags> {
0496 public:
0497 SiPixelVCalValueComparisonBase()
0498 : PlotImage<SiPixelVCal, nIOVs, ntags>(Form("SiPixelVCal Synoptic Values Comparison %i tag(s)", ntags)) {}
0499
0500 bool fill() override {
0501 TH1F::SetDefaultSumw2(true);
0502
0503
0504 auto theIOVs = PlotBase::getTag<0>().iovs;
0505 auto f_tagname = PlotBase::getTag<0>().name;
0506 std::string l_tagname = "";
0507 auto firstiov = theIOVs.front();
0508 std::tuple<cond::Time_t, cond::Hash> lastiov;
0509
0510
0511 assert(this->m_plotAnnotations.ntags < 3);
0512
0513 if (this->m_plotAnnotations.ntags == 2) {
0514 auto tag2iovs = PlotBase::getTag<1>().iovs;
0515 l_tagname = PlotBase::getTag<1>().name;
0516 lastiov = tag2iovs.front();
0517 } else {
0518 lastiov = theIOVs.back();
0519 }
0520
0521 SiPixelVCal::mapToDetId l_Map_;
0522 SiPixelVCal::mapToDetId f_Map_;
0523
0524 std::shared_ptr<SiPixelVCal> last_payload = this->fetchPayload(std::get<1>(lastiov));
0525 if (myType == SiPixelVCalPI::t_slope) {
0526 l_Map_ = last_payload->getAllSlopes();
0527 } else {
0528 l_Map_ = last_payload->getAllOffsets();
0529 }
0530 auto l_extrema = SiPixelPI::findMinMaxInMap(l_Map_);
0531
0532 std::shared_ptr<SiPixelVCal> first_payload = this->fetchPayload(std::get<1>(firstiov));
0533 if (myType == SiPixelVCalPI::t_slope) {
0534 f_Map_ = first_payload->getAllSlopes();
0535 } else {
0536 f_Map_ = first_payload->getAllOffsets();
0537 }
0538 auto f_extrema = SiPixelPI::findMinMaxInMap(f_Map_);
0539
0540 auto max = (l_extrema.second > f_extrema.second) ? l_extrema.second : f_extrema.second;
0541 auto min = (l_extrema.first < f_extrema.first) ? l_extrema.first : f_extrema.first;
0542 auto delta = std::abs(max - min) / 10.;
0543
0544 std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
0545 std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
0546
0547 TCanvas canvas("Canv", "Canv", 1200, 1000);
0548 canvas.cd();
0549 auto hfirst = std::unique_ptr<TH1F>(
0550 new TH1F("value_first",
0551 fmt::sprintf("SiPixel VCal %s value;SiPixel VCal %s ;# modules",
0552 (myType == SiPixelVCalPI::t_slope ? "slope" : "offset"),
0553 (myType == SiPixelVCalPI::t_slope ? " slope [ADC/VCal]" : " offset [ADC]"))
0554 .c_str(),
0555 50,
0556 min - delta,
0557 max + delta));
0558 hfirst->SetStats(false);
0559
0560 auto hlast = std::unique_ptr<TH1F>(
0561 new TH1F("value_last",
0562 fmt::sprintf("SiPixel VCal %s value;SiPixel VCal %s ;# modules",
0563 (myType == SiPixelVCalPI::t_slope ? "slope" : "offset"),
0564 (myType == SiPixelVCalPI::t_slope ? " slope [ADC/VCal]" : " offset [ADC]"))
0565 .c_str(),
0566 50,
0567 min - delta,
0568 max + delta));
0569 hlast->SetStats(false);
0570
0571 SiPixelPI::adjustCanvasMargins(canvas.cd(), 0.06, 0.12, 0.12, 0.05);
0572 canvas.Modified();
0573
0574 for (const auto &element : f_Map_) {
0575 hfirst->Fill(element.second);
0576 }
0577
0578 for (const auto &element : l_Map_) {
0579 hlast->Fill(element.second);
0580 }
0581
0582 auto extrema = SiPixelPI::getExtrema(hfirst.get(), hlast.get());
0583 hfirst->GetYaxis()->SetRangeUser(extrema.first, extrema.second * 1.10);
0584
0585 hfirst->SetTitle("");
0586 hfirst->SetFillColor(kRed);
0587 hfirst->SetBarWidth(0.95);
0588 hfirst->Draw("histbar");
0589
0590 hlast->SetTitle("");
0591 hlast->SetFillColorAlpha(kBlue, 0.20);
0592 hlast->SetBarWidth(0.95);
0593 hlast->Draw("histbarsame");
0594
0595 SiPixelPI::makeNicePlotStyle(hfirst.get());
0596 SiPixelPI::makeNicePlotStyle(hlast.get());
0597
0598 canvas.Update();
0599
0600 TLegend legend = TLegend(0.30, 0.86, 0.95, 0.94);
0601 legend.AddEntry(hfirst.get(), ("payload: #color[2]{" + std::get<1>(firstiov) + "}").c_str(), "F");
0602 legend.AddEntry(hlast.get(), ("payload: #color[4]{" + std::get<1>(lastiov) + "}").c_str(), "F");
0603 legend.SetTextSize(0.025);
0604 legend.Draw("same");
0605
0606 auto ltx = TLatex();
0607 ltx.SetTextFont(62);
0608 ltx.SetTextSize(0.040);
0609 ltx.SetTextAlign(11);
0610 std::string ltxText;
0611 if (this->m_plotAnnotations.ntags == 2) {
0612 ltxText = fmt::sprintf("#color[2]{%s, %s} vs #color[4]{%s, %s}",
0613 f_tagname,
0614 std::to_string(std::get<0>(firstiov)),
0615 l_tagname,
0616 std::to_string(std::get<0>(lastiov)));
0617 } else {
0618 ltxText = fmt::sprintf("%s IOV: #color[2]{%s} vs IOV: #color[4]{%s}",
0619 f_tagname,
0620 std::to_string(std::get<0>(firstiov)),
0621 std::to_string(std::get<0>(lastiov)));
0622 }
0623 ltx.DrawLatexNDC(gPad->GetLeftMargin(), 1 - gPad->GetTopMargin() + 0.01, ltxText.c_str());
0624
0625 std::string fileName(this->m_imageFileName);
0626 canvas.SaveAs(fileName.c_str());
0627
0628 return true;
0629 }
0630 };
0631
0632 using SiPixelVCalSlopesComparisonSingleTag = SiPixelVCalValueComparisonBase<SiPixelVCalPI::t_slope, MULTI_IOV, 1>;
0633 using SiPixelVCalOffsetsComparisonSingleTag = SiPixelVCalValueComparisonBase<SiPixelVCalPI::t_offset, MULTI_IOV, 1>;
0634 using SiPixelVCalSlopesComparisonTwoTags = SiPixelVCalValueComparisonBase<SiPixelVCalPI::t_slope, SINGLE_IOV, 2>;
0635 using SiPixelVCalOffsetsComparisonTwoTags = SiPixelVCalValueComparisonBase<SiPixelVCalPI::t_offset, SINGLE_IOV, 2>;
0636
0637 }
0638
0639 PAYLOAD_INSPECTOR_MODULE(SiPixelVCal) {
0640 PAYLOAD_INSPECTOR_CLASS(SiPixelVCalSlopeValue);
0641 PAYLOAD_INSPECTOR_CLASS(SiPixelVCalOffsetValue);
0642 PAYLOAD_INSPECTOR_CLASS(SiPixelVCalValues);
0643 PAYLOAD_INSPECTOR_CLASS(SiPixelVCalSlopeValuesBarrel);
0644 PAYLOAD_INSPECTOR_CLASS(SiPixelVCalSlopeValuesEndcap);
0645 PAYLOAD_INSPECTOR_CLASS(SiPixelVCalOffsetValuesBarrel);
0646 PAYLOAD_INSPECTOR_CLASS(SiPixelVCalOffsetValuesEndcap);
0647 PAYLOAD_INSPECTOR_CLASS(SiPixelVCalSlopesBarrelCompareSingleTag);
0648 PAYLOAD_INSPECTOR_CLASS(SiPixelVCalOffsetsBarrelCompareSingleTag);
0649 PAYLOAD_INSPECTOR_CLASS(SiPixelVCalSlopesEndcapCompareSingleTag);
0650 PAYLOAD_INSPECTOR_CLASS(SiPixelVCalOffsetsEndcapCompareSingleTag);
0651 PAYLOAD_INSPECTOR_CLASS(SiPixelVCalSlopesBarrelCompareTwoTags);
0652 PAYLOAD_INSPECTOR_CLASS(SiPixelVCalOffsetsBarrelCompareTwoTags);
0653 PAYLOAD_INSPECTOR_CLASS(SiPixelVCalSlopesEndcapCompareTwoTags);
0654 PAYLOAD_INSPECTOR_CLASS(SiPixelVCalOffsetsEndcapCompareTwoTags);
0655 PAYLOAD_INSPECTOR_CLASS(SiPixelVCalSlopesComparisonSingleTag);
0656 PAYLOAD_INSPECTOR_CLASS(SiPixelVCalOffsetsComparisonSingleTag);
0657 PAYLOAD_INSPECTOR_CLASS(SiPixelVCalSlopesComparisonTwoTags);
0658 PAYLOAD_INSPECTOR_CLASS(SiPixelVCalOffsetsComparisonTwoTags);
0659 }