File indexing completed on 2024-09-07 04:35:32
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/SiPixelLorentzAngle.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 #include "DQM/TrackerRemapper/interface/Phase1PixelROCMaps.h"
0018 #include "DQM/TrackerRemapper/interface/Phase1PixelSummaryMap.h"
0019
0020 #include <memory>
0021 #include <sstream>
0022
0023
0024 #include "TH2F.h"
0025 #include "TH1F.h"
0026 #include "TLegend.h"
0027 #include "TCanvas.h"
0028 #include "TLine.h"
0029 #include "TStyle.h"
0030 #include "TLatex.h"
0031 #include "TPave.h"
0032 #include "TPaveStats.h"
0033 #include "TGaxis.h"
0034
0035 namespace {
0036
0037 using namespace cond::payloadInspector;
0038
0039
0040
0041
0042
0043
0044 class SiPixelLorentzAngleValue : public Histogram1D<SiPixelLorentzAngle, SINGLE_IOV> {
0045 public:
0046 SiPixelLorentzAngleValue()
0047 : Histogram1D<SiPixelLorentzAngle, SINGLE_IOV>(
0048 "SiPixel LorentzAngle values", "SiPixel LorentzAngle values", 100, 0.0, 0.1) {}
0049
0050 bool fill() override {
0051 auto tag = PlotBase::getTag<0>();
0052 for (auto const &iov : tag.iovs) {
0053 std::shared_ptr<SiPixelLorentzAngle> payload = Base::fetchPayload(std::get<1>(iov));
0054 if (payload.get()) {
0055 std::map<uint32_t, float> LAMap_ = payload->getLorentzAngles();
0056
0057 for (const auto &element : LAMap_) {
0058 fillWithValue(element.second);
0059 }
0060 }
0061 }
0062 return true;
0063 }
0064 };
0065
0066
0067
0068
0069 class SiPixelLorentzAngleValues : public PlotImage<SiPixelLorentzAngle, SINGLE_IOV> {
0070 public:
0071 SiPixelLorentzAngleValues() : PlotImage<SiPixelLorentzAngle, SINGLE_IOV>("SiPixelLorentzAngle Values") {}
0072
0073 bool fill() override {
0074 gStyle->SetOptStat("emr");
0075
0076 auto tag = PlotBase::getTag<0>();
0077 auto tagname = PlotBase::getTag<0>().name;
0078 auto iov = tag.iovs.front();
0079 std::string IOVsince = std::to_string(std::get<0>(iov));
0080 std::shared_ptr<SiPixelLorentzAngle> payload = fetchPayload(std::get<1>(iov));
0081 std::map<uint32_t, float> LAMap_ = payload->getLorentzAngles();
0082 auto extrema = SiPixelPI::findMinMaxInMap(LAMap_);
0083
0084 TCanvas canvas("Canv", "Canv", 1200, 1000);
0085 canvas.cd();
0086 auto h1 = std::make_unique<TH1F>("value",
0087 "SiPixel LA value;SiPixel LorentzAngle #mu_{H}(tan#theta_{L}/B) [1/T];# modules",
0088 50,
0089 (extrema.first * 0.9),
0090 (extrema.second * 1.1));
0091
0092 SiPixelPI::adjustCanvasMargins(canvas.cd(), 0.06, 0.12, 0.12, 0.05);
0093 canvas.Modified();
0094
0095 for (const auto &element : LAMap_) {
0096 h1->Fill(element.second);
0097 }
0098
0099 h1->SetTitle("");
0100 h1->GetYaxis()->SetRangeUser(0., h1->GetMaximum() * 1.30);
0101 h1->SetFillColor(kRed);
0102 h1->SetMarkerStyle(20);
0103 h1->SetMarkerSize(1);
0104 h1->Draw("bar2");
0105
0106 SiPixelPI::makeNicePlotStyle(h1.get());
0107 h1->SetStats(true);
0108
0109 canvas.Update();
0110
0111 TLegend legend = TLegend(0.35, 0.88, 0.94, 0.93);
0112 legend.SetHeader(("#splitline{Payload hash: #bf{" + (std::get<1>(iov)) + "}}{Tag,IOV: #bf{#color[2]{" + tagname +
0113 "}," + IOVsince + "}}")
0114 .c_str(),
0115 "C");
0116 legend.SetTextSize(0.023);
0117 legend.SetLineColor(10);
0118 legend.Draw("same");
0119
0120 TPaveStats *st = (TPaveStats *)h1->FindObject("stats");
0121 st->SetTextSize(0.03);
0122 SiPixelPI::adjustStats(st, 0.15, 0.83, 0.34, 0.93);
0123
0124 auto ltx = TLatex();
0125 ltx.SetTextFont(62);
0126
0127 ltx.SetTextSize(0.05);
0128 ltx.SetTextAlign(11);
0129 ltx.DrawLatexNDC(gPad->GetLeftMargin(), 1 - gPad->GetTopMargin() + 0.01, "SiPixel Lorentz Angle Value");
0130
0131 std::string fileName(m_imageFileName);
0132 canvas.SaveAs(fileName.c_str());
0133
0134 return true;
0135 }
0136 };
0137
0138
0139
0140
0141 template <bool isBarrel>
0142 class SiPixelLorentzAngleValuesPerRegion : public PlotImage<SiPixelLorentzAngle, SINGLE_IOV> {
0143 public:
0144 SiPixelLorentzAngleValuesPerRegion()
0145 : PlotImage<SiPixelLorentzAngle, SINGLE_IOV>("SiPixelLorentzAngle Values per region") {}
0146
0147 bool fill() override {
0148 gStyle->SetOptStat("emr");
0149
0150 auto tag = PlotBase::getTag<0>();
0151 auto iov = tag.iovs.front();
0152 std::shared_ptr<SiPixelLorentzAngle> payload = fetchPayload(std::get<1>(iov));
0153 std::map<uint32_t, float> LAMap_ = payload->getLorentzAngles();
0154 auto extrema = SiPixelPI::findMinMaxInMap(LAMap_);
0155
0156 TCanvas canvas("Canv", "Canv", isBarrel ? 1400 : 1800, 1200);
0157 canvas.cd();
0158
0159 SiPixelPI::PhaseInfo phaseInfo(LAMap_.size());
0160 const char *path_toTopologyXML = phaseInfo.pathToTopoXML();
0161
0162 TrackerTopology tTopo =
0163 StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
0164
0165 auto myPlots = PixelRegions::PixelRegionContainers(&tTopo, phaseInfo.phase());
0166 myPlots.bookAll("SiPixel LA",
0167 "SiPixel LorentzAngle #mu_{H}(tan#theta_{L}/B) [1/T]",
0168 "#modules",
0169 50,
0170 extrema.first * 0.9,
0171 extrema.second * 1.1);
0172
0173 canvas.Modified();
0174
0175 for (const auto &element : LAMap_) {
0176 myPlots.fill(element.first, element.second);
0177 }
0178
0179 myPlots.beautify();
0180 myPlots.draw(canvas, isBarrel);
0181
0182 TLegend legend = TLegend(0.40, 0.88, 0.93, 0.90);
0183 legend.SetHeader(("Hash: #bf{" + (std::get<1>(iov)) + "}").c_str(),
0184 "C");
0185
0186 legend.SetTextSize(0.025);
0187 legend.SetLineColor(10);
0188
0189 unsigned int maxPads = canvas.GetListOfPrimitives()->GetSize();
0190 for (unsigned int c = 1; c <= maxPads; c++) {
0191 if (phaseInfo.phase() == SiPixelPI::phase::two && (c == 5 || c == 10))
0192 continue;
0193 canvas.cd(c);
0194 SiPixelPI::adjustCanvasMargins(canvas.cd(c), 0.06, 0.12, 0.12, 0.05);
0195 legend.Draw("same");
0196 canvas.cd(c)->Update();
0197 }
0198
0199 myPlots.stats();
0200
0201 auto ltx = TLatex();
0202 ltx.SetTextFont(62);
0203 ltx.SetTextSize(0.05);
0204 ltx.SetTextAlign(11);
0205
0206 int index = 0;
0207 for (unsigned int c = 1; c <= maxPads; c++) {
0208 if (phaseInfo.phase() == SiPixelPI::phase::two && (c == 5 || c == 10))
0209 continue;
0210 canvas.cd(c);
0211 ltx.DrawLatexNDC(gPad->GetLeftMargin(),
0212 1 - gPad->GetTopMargin() + 0.01,
0213 (PixelRegions::getIDLabels(phaseInfo.phase(), isBarrel)[index] +
0214 ", IOV:" + std::to_string(std::get<0>(iov)))
0215 .c_str());
0216
0217 index++;
0218 }
0219
0220 std::string fileName(m_imageFileName);
0221 canvas.SaveAs(fileName.c_str());
0222
0223 return true;
0224 }
0225 };
0226
0227 using SiPixelLorentzAngleValuesBarrel = SiPixelLorentzAngleValuesPerRegion<true>;
0228 using SiPixelLorentzAngleValuesEndcap = SiPixelLorentzAngleValuesPerRegion<false>;
0229
0230
0231
0232
0233 template <bool isBarrel, IOVMultiplicity nIOVs, int ntags>
0234 class SiPixelLorentzAngleValuesComparisonPerRegion : public PlotImage<SiPixelLorentzAngle, nIOVs, ntags> {
0235 public:
0236 SiPixelLorentzAngleValuesComparisonPerRegion()
0237 : PlotImage<SiPixelLorentzAngle, nIOVs, ntags>(
0238 Form("SiPixelLorentzAngle Values Comparisons per region %i tag(s)", ntags)) {}
0239
0240 bool fill() override {
0241 gStyle->SetOptStat("emr");
0242
0243
0244 auto theIOVs = PlotBase::getTag<0>().iovs;
0245 auto f_tagname = PlotBase::getTag<0>().name;
0246 std::string l_tagname = "";
0247 auto firstiov = theIOVs.front();
0248 std::tuple<cond::Time_t, cond::Hash> lastiov;
0249
0250
0251 assert(this->m_plotAnnotations.ntags < 3);
0252
0253 if (this->m_plotAnnotations.ntags == 2) {
0254 auto tag2iovs = PlotBase::getTag<1>().iovs;
0255 l_tagname = PlotBase::getTag<1>().name;
0256 lastiov = tag2iovs.front();
0257 } else {
0258 lastiov = theIOVs.back();
0259 }
0260
0261 std::shared_ptr<SiPixelLorentzAngle> last_payload = this->fetchPayload(std::get<1>(lastiov));
0262 std::map<uint32_t, float> l_LAMap_ = last_payload->getLorentzAngles();
0263 auto l_extrema = SiPixelPI::findMinMaxInMap(l_LAMap_);
0264
0265 std::shared_ptr<SiPixelLorentzAngle> first_payload = this->fetchPayload(std::get<1>(firstiov));
0266 std::map<uint32_t, float> f_LAMap_ = first_payload->getLorentzAngles();
0267 auto f_extrema = SiPixelPI::findMinMaxInMap(f_LAMap_);
0268
0269 auto max = (l_extrema.second > f_extrema.second) ? l_extrema.second : f_extrema.second;
0270 auto min = (l_extrema.first < f_extrema.first) ? l_extrema.first : f_extrema.first;
0271
0272 std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
0273 std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
0274
0275 TCanvas canvas("Canv", "Canv", isBarrel ? 1400 : 1800, 1200);
0276 canvas.cd();
0277
0278 SiPixelPI::PhaseInfo l_phaseInfo(l_LAMap_.size());
0279 SiPixelPI::PhaseInfo f_phaseInfo(f_LAMap_.size());
0280
0281 if (l_phaseInfo.isComparedWithPhase2(f_phaseInfo)) {
0282 SiPixelPI::displayNotSupported(canvas, std::max(f_LAMap_.size(), l_LAMap_.size()));
0283 std::string fileName(this->m_imageFileName);
0284 canvas.SaveAs(fileName.c_str());
0285 return false;
0286 }
0287
0288
0289 const char *path_toTopologyXML = l_phaseInfo.pathToTopoXML();
0290
0291 edm::LogPrint("SiPixelLorentzAngleValuesComparisonPerRegion") << l_phaseInfo;
0292
0293 auto l_tTopo =
0294 StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
0295
0296 auto l_myPlots = PixelRegions::PixelRegionContainers(&l_tTopo, l_phaseInfo.phase());
0297 l_myPlots.bookAll("SiPixel LA,last",
0298 "SiPixel LorentzAngle #mu_{H}(tan#theta_{L}/B) [1/T]",
0299 "#modules",
0300 50,
0301 min * 0.9,
0302 max * 1.1);
0303
0304 for (const auto &element : l_LAMap_) {
0305 l_myPlots.fill(element.first, element.second);
0306 }
0307
0308 l_myPlots.beautify();
0309 l_myPlots.draw(canvas, isBarrel, "bar2", f_phaseInfo.isPhase1Comparison(l_phaseInfo));
0310
0311
0312 path_toTopologyXML = f_phaseInfo.pathToTopoXML();
0313
0314 auto f_tTopo =
0315 StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
0316
0317 auto f_myPlots = PixelRegions::PixelRegionContainers(&f_tTopo, f_phaseInfo.phase());
0318 f_myPlots.bookAll("SiPixel LA,first",
0319 "SiPixel LorentzAngle #mu_{H}(tan#theta_{L}/B) [1/T]",
0320 "#modules",
0321 50,
0322 min * 0.9,
0323 max * 1.1);
0324
0325 for (const auto &element : f_LAMap_) {
0326 f_myPlots.fill(element.first, element.second);
0327 }
0328
0329 f_myPlots.beautify(kAzure, kBlue);
0330 f_myPlots.draw(canvas, isBarrel, "HISTsames", f_phaseInfo.isPhase1Comparison(l_phaseInfo));
0331
0332
0333 l_myPlots.rescaleMax(f_myPlots);
0334
0335
0336
0337 auto colorTag = PixelRegions::L1;
0338 std::unique_ptr<TLegend> legend;
0339 if (this->m_plotAnnotations.ntags == 2) {
0340 legend = std::make_unique<TLegend>(0.36, 0.86, 0.94, 0.92);
0341 legend->AddEntry(l_myPlots.getHistoFromMap(colorTag).get(), ("#color[2]{" + l_tagname + "}").c_str(), "F");
0342 legend->AddEntry(f_myPlots.getHistoFromMap(colorTag).get(), ("#color[4]{" + f_tagname + "}").c_str(), "F");
0343 legend->SetTextSize(0.024);
0344 } else {
0345 legend = std::make_unique<TLegend>(0.58, 0.80, 0.90, 0.92);
0346 legend->AddEntry(l_myPlots.getHistoFromMap(colorTag).get(), ("#color[2]{" + lastIOVsince + "}").c_str(), "F");
0347 legend->AddEntry(f_myPlots.getHistoFromMap(colorTag).get(), ("#color[4]{" + firstIOVsince + "}").c_str(), "F");
0348 legend->SetTextSize(0.040);
0349 }
0350 legend->SetLineColor(10);
0351
0352 unsigned int maxPads = canvas.GetListOfPrimitives()->GetSize();
0353 for (unsigned int c = 1; c <= maxPads; c++) {
0354 if (l_phaseInfo.phase() == SiPixelPI::phase::two && (c == 5 || c == 10))
0355 continue;
0356 canvas.cd(c);
0357 SiPixelPI::adjustCanvasMargins(canvas.cd(c), 0.06, 0.12, 0.12, 0.05);
0358 legend->Draw("same");
0359 canvas.cd(c)->Update();
0360 }
0361
0362 f_myPlots.stats(0);
0363 l_myPlots.stats(1);
0364
0365 auto ltx = TLatex();
0366 ltx.SetTextFont(62);
0367 ltx.SetTextSize(0.05);
0368 ltx.SetTextAlign(11);
0369
0370 int index = 0;
0371 for (unsigned int c = 1; c <= maxPads; c++) {
0372 if (l_phaseInfo.phase() == SiPixelPI::phase::two && (c == 5 || c == 10))
0373 continue;
0374 canvas.cd(c);
0375
0376 COUT << "c:" << c << " index:" << index << " : "
0377 << PixelRegions::getIDLabels(l_phaseInfo.phase(), isBarrel)[index] << "\n";
0378
0379 ltx.DrawLatexNDC(
0380 gPad->GetLeftMargin(),
0381 1 - gPad->GetTopMargin() + 0.01,
0382 (PixelRegions::getIDLabels(l_phaseInfo.phase(), isBarrel)[index] + " : #color[4]{" +
0383 std::to_string(std::get<0>(firstiov)) + "} vs #color[2]{" + std::to_string(std::get<0>(lastiov)) + "}")
0384 .c_str());
0385 index++;
0386 }
0387
0388 std::string fileName(this->m_imageFileName);
0389 canvas.SaveAs(fileName.c_str());
0390
0391 #ifdef MMDEBUG
0392 canvas.SaveAs("DEBUG.root");
0393 #endif
0394
0395 return true;
0396 }
0397 };
0398
0399 using SiPixelLorentzAngleValuesBarrelCompareSingleTag =
0400 SiPixelLorentzAngleValuesComparisonPerRegion<true, MULTI_IOV, 1>;
0401 using SiPixelLorentzAngleValuesEndcapCompareSingleTag =
0402 SiPixelLorentzAngleValuesComparisonPerRegion<false, MULTI_IOV, 1>;
0403
0404 using SiPixelLorentzAngleValuesBarrelCompareTwoTags =
0405 SiPixelLorentzAngleValuesComparisonPerRegion<true, SINGLE_IOV, 2>;
0406 using SiPixelLorentzAngleValuesEndcapCompareTwoTags =
0407 SiPixelLorentzAngleValuesComparisonPerRegion<false, SINGLE_IOV, 2>;
0408
0409
0410
0411
0412 template <IOVMultiplicity nIOVs, int ntags>
0413 class SiPixelLorentzAngleValueComparisonBase : public PlotImage<SiPixelLorentzAngle, nIOVs, ntags> {
0414 public:
0415 SiPixelLorentzAngleValueComparisonBase()
0416 : PlotImage<SiPixelLorentzAngle, nIOVs, ntags>(Form("SiPixelLorentzAngle Values Comparison %i tag(s)", ntags)) {
0417 }
0418
0419 bool fill() override {
0420 TH1F::SetDefaultSumw2(true);
0421
0422
0423 auto theIOVs = PlotBase::getTag<0>().iovs;
0424 auto f_tagname = PlotBase::getTag<0>().name;
0425 std::string l_tagname = "";
0426 auto firstiov = theIOVs.front();
0427 std::tuple<cond::Time_t, cond::Hash> lastiov;
0428
0429
0430 assert(this->m_plotAnnotations.ntags < 3);
0431
0432 if (this->m_plotAnnotations.ntags == 2) {
0433 auto tag2iovs = PlotBase::getTag<1>().iovs;
0434 l_tagname = PlotBase::getTag<1>().name;
0435 lastiov = tag2iovs.front();
0436 } else {
0437 lastiov = theIOVs.back();
0438 }
0439
0440 std::shared_ptr<SiPixelLorentzAngle> last_payload = this->fetchPayload(std::get<1>(lastiov));
0441 std::map<uint32_t, float> l_LAMap_ = last_payload->getLorentzAngles();
0442 auto l_extrema = SiPixelPI::findMinMaxInMap(l_LAMap_);
0443
0444 std::shared_ptr<SiPixelLorentzAngle> first_payload = this->fetchPayload(std::get<1>(firstiov));
0445 std::map<uint32_t, float> f_LAMap_ = first_payload->getLorentzAngles();
0446 auto f_extrema = SiPixelPI::findMinMaxInMap(f_LAMap_);
0447
0448 auto max = (l_extrema.second > f_extrema.second) ? l_extrema.second : f_extrema.second;
0449 auto min = (l_extrema.first < f_extrema.first) ? l_extrema.first : f_extrema.first;
0450
0451 std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
0452 std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
0453
0454 TCanvas canvas("Canv", "Canv", 1200, 1000);
0455 canvas.cd();
0456 auto hfirst =
0457 std::make_unique<TH1F>("value_first",
0458 "SiPixel LA value;SiPixel LorentzAngle #mu_{H}(tan#theta_{L}/B) [1/T];# modules",
0459 50,
0460 min * 0.9,
0461 max * 1.1);
0462 hfirst->SetStats(false);
0463
0464 auto hlast =
0465 std::make_unique<TH1F>("value_last",
0466 "SiPixel LA value;SiPixel LorentzAngle #mu_{H}(tan#theta_{L}/B) [1/T];# modules",
0467 50,
0468 min * 0.9,
0469 max * 1.1);
0470 hlast->SetStats(false);
0471
0472 SiPixelPI::adjustCanvasMargins(canvas.cd(), 0.06, 0.12, 0.12, 0.05);
0473 canvas.Modified();
0474
0475 for (const auto &element : f_LAMap_) {
0476 hfirst->Fill(element.second);
0477 }
0478
0479 for (const auto &element : l_LAMap_) {
0480 hlast->Fill(element.second);
0481 }
0482
0483 auto extrema = SiPixelPI::getExtrema(hfirst.get(), hlast.get());
0484 hfirst->GetYaxis()->SetRangeUser(extrema.first, extrema.second * 1.10);
0485
0486 hfirst->SetTitle("");
0487 hfirst->SetFillColor(kRed);
0488 hfirst->SetBarWidth(0.95);
0489 hfirst->Draw("histbar");
0490
0491 hlast->SetTitle("");
0492 hlast->SetFillColorAlpha(kBlue, 0.20);
0493 hlast->SetBarWidth(0.95);
0494 hlast->Draw("histbarsame");
0495
0496 SiPixelPI::makeNicePlotStyle(hfirst.get());
0497 SiPixelPI::makeNicePlotStyle(hlast.get());
0498
0499 canvas.Update();
0500
0501 TLegend legend = TLegend(0.30, 0.86, 0.95, 0.94);
0502
0503
0504
0505 legend.AddEntry(hfirst.get(), ("payload: #color[2]{" + std::get<1>(firstiov) + "}").c_str(), "F");
0506 legend.AddEntry(hlast.get(), ("payload: #color[4]{" + std::get<1>(lastiov) + "}").c_str(), "F");
0507 legend.SetTextSize(0.025);
0508 legend.Draw("same");
0509
0510 auto ltx = TLatex();
0511 ltx.SetTextFont(62);
0512
0513 ltx.SetTextSize(0.047);
0514 ltx.SetTextAlign(11);
0515 std::string ltxText;
0516 if (this->m_plotAnnotations.ntags == 2) {
0517 ltxText = fmt::sprintf("#color[2]{%s, %s} vs #color[4]{%s, %s}",
0518 f_tagname,
0519 std::to_string(std::get<0>(firstiov)),
0520 l_tagname,
0521 std::to_string(std::get<0>(lastiov)));
0522 } else {
0523 ltxText = fmt::sprintf("%s IOV: #color[2]{%s} vs IOV: #color[4]{%s}",
0524 f_tagname,
0525 std::to_string(std::get<0>(firstiov)),
0526 std::to_string(std::get<0>(lastiov)));
0527 }
0528 ltx.DrawLatexNDC(gPad->GetLeftMargin(), 1 - gPad->GetTopMargin() + 0.01, ltxText.c_str());
0529
0530 std::string fileName(this->m_imageFileName);
0531 canvas.SaveAs(fileName.c_str());
0532
0533 return true;
0534 }
0535 };
0536
0537 using SiPixelLorentzAngleValueComparisonSingleTag = SiPixelLorentzAngleValueComparisonBase<MULTI_IOV, 1>;
0538 using SiPixelLorentzAngleValueComparisonTwoTags = SiPixelLorentzAngleValueComparisonBase<SINGLE_IOV, 2>;
0539
0540
0541
0542
0543 template <IOVMultiplicity nIOVs, int ntags>
0544 class SiPixelLorentzAngleByRegionComparisonBase : public PlotImage<SiPixelLorentzAngle, nIOVs, ntags> {
0545 public:
0546 SiPixelLorentzAngleByRegionComparisonBase()
0547 : PlotImage<SiPixelLorentzAngle, nIOVs, ntags>(
0548 Form("SiPixelLorentzAngle Comparison by Region %i tag(s)", ntags)) {}
0549
0550 bool fill() override {
0551 gStyle->SetPaintTextFormat(".3f");
0552
0553
0554 auto theIOVs = PlotBase::getTag<0>().iovs;
0555 auto f_tagname = PlotBase::getTag<0>().name;
0556 std::string l_tagname = "";
0557 auto firstiov = theIOVs.front();
0558 std::tuple<cond::Time_t, cond::Hash> lastiov;
0559
0560
0561 assert(this->m_plotAnnotations.ntags < 3);
0562
0563 if (this->m_plotAnnotations.ntags == 2) {
0564 auto tag2iovs = PlotBase::getTag<1>().iovs;
0565 l_tagname = PlotBase::getTag<1>().name;
0566 lastiov = tag2iovs.front();
0567 } else {
0568 lastiov = theIOVs.back();
0569 }
0570
0571 std::shared_ptr<SiPixelLorentzAngle> last_payload = this->fetchPayload(std::get<1>(lastiov));
0572 std::map<uint32_t, float> l_LAMap_ = last_payload->getLorentzAngles();
0573 std::shared_ptr<SiPixelLorentzAngle> first_payload = this->fetchPayload(std::get<1>(firstiov));
0574 std::map<uint32_t, float> f_LAMap_ = first_payload->getLorentzAngles();
0575
0576 std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
0577 std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
0578
0579 SiPixelPI::PhaseInfo l_phaseInfo(l_LAMap_.size());
0580 SiPixelPI::PhaseInfo f_phaseInfo(f_LAMap_.size());
0581
0582 TCanvas canvas("Comparison", "Comparison", 1600, 800);
0583 std::map<SiPixelPI::regions, std::shared_ptr<TH1F>> FirstLA_spectraByRegion;
0584 std::map<SiPixelPI::regions, std::shared_ptr<TH1F>> LastLA_spectraByRegion;
0585 std::shared_ptr<TH1F> summaryFirst;
0586 std::shared_ptr<TH1F> summaryLast;
0587
0588
0589 for (int r = SiPixelPI::BPixL1o; r != SiPixelPI::NUM_OF_REGIONS; r++) {
0590 SiPixelPI::regions part = static_cast<SiPixelPI::regions>(r);
0591 std::string s_part = SiPixelPI::getStringFromRegionEnum(part);
0592
0593 FirstLA_spectraByRegion[part] = std::make_shared<TH1F>(Form("hfirstLA_%s", s_part.c_str()),
0594 Form(";%s #mu_{H} [1/T];n. of modules", s_part.c_str()),
0595 1000,
0596 0.,
0597 1000.);
0598 LastLA_spectraByRegion[part] = std::make_shared<TH1F>(Form("hlastLA_%s", s_part.c_str()),
0599 Form(";%s #mu_{H} [1/T];n. of modules", s_part.c_str()),
0600 1000,
0601 0.,
0602 1000.);
0603 }
0604
0605 summaryFirst = std::make_shared<TH1F>("first Summary",
0606 "Summary for #LT tan#theta_{L}/B #GT;;average LA #LT #mu_{H} #GT [1/T]",
0607 FirstLA_spectraByRegion.size(),
0608 0,
0609 FirstLA_spectraByRegion.size());
0610 summaryLast = std::make_shared<TH1F>("last Summary",
0611 "Summary for #LT tan#theta_{L}/B #GT;;average LA #LT #mu_{H} #GT [1/T]",
0612 LastLA_spectraByRegion.size(),
0613 0,
0614 LastLA_spectraByRegion.size());
0615
0616
0617 const char *path_toTopologyXML = f_phaseInfo.pathToTopoXML();
0618
0619 auto f_tTopo =
0620 StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
0621
0622
0623
0624
0625 for (const auto &it : f_LAMap_) {
0626 if (DetId(it.first).det() != DetId::Tracker) {
0627 edm::LogWarning("SiPixelLorentzAngle_PayloadInspector")
0628 << "Encountered invalid Tracker DetId:" << it.first << " - terminating ";
0629 return false;
0630 }
0631
0632 SiPixelPI::topolInfo t_info_fromXML;
0633 t_info_fromXML.init();
0634 DetId detid(it.first);
0635 t_info_fromXML.fillGeometryInfo(detid, f_tTopo, f_phaseInfo.phase());
0636
0637 SiPixelPI::regions thePart = t_info_fromXML.filterThePartition();
0638 if (thePart != SiPixelPI::NUM_OF_REGIONS) {
0639 FirstLA_spectraByRegion[thePart]->Fill(it.second);
0640 }
0641 }
0642
0643
0644 path_toTopologyXML = l_phaseInfo.pathToTopoXML();
0645
0646 auto l_tTopo =
0647 StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
0648
0649
0650
0651
0652 for (const auto &it : l_LAMap_) {
0653 if (DetId(it.first).det() != DetId::Tracker) {
0654 edm::LogWarning("SiPixelLorentzAngle_PayloadInspector")
0655 << "Encountered invalid Tracker DetId:" << it.first << " - terminating ";
0656 return false;
0657 }
0658
0659 SiPixelPI::topolInfo t_info_fromXML;
0660 t_info_fromXML.init();
0661 DetId detid(it.first);
0662 t_info_fromXML.fillGeometryInfo(detid, l_tTopo, l_phaseInfo.phase());
0663
0664 SiPixelPI::regions thePart = t_info_fromXML.filterThePartition();
0665 if (thePart != SiPixelPI::NUM_OF_REGIONS) {
0666 LastLA_spectraByRegion[thePart]->Fill(it.second);
0667 }
0668 }
0669
0670
0671 int bin = 1;
0672 for (int r = SiPixelPI::BPixL1o; r != SiPixelPI::NUM_OF_REGIONS; r++) {
0673 SiPixelPI::regions part = static_cast<SiPixelPI::regions>(r);
0674
0675 summaryFirst->GetXaxis()->SetBinLabel(bin, SiPixelPI::getStringFromRegionEnum(part).c_str());
0676
0677 float f_mean =
0678 FirstLA_spectraByRegion[part]->GetMean() > 10.e-6 ? FirstLA_spectraByRegion[part]->GetMean() : 0.;
0679 summaryFirst->SetBinContent(bin, f_mean);
0680
0681
0682 summaryLast->GetXaxis()->SetBinLabel(bin, SiPixelPI::getStringFromRegionEnum(part).c_str());
0683
0684 float l_mean = LastLA_spectraByRegion[part]->GetMean() > 10.e-6 ? LastLA_spectraByRegion[part]->GetMean() : 0.;
0685 summaryLast->SetBinContent(bin, l_mean);
0686
0687 bin++;
0688 }
0689
0690 SiPixelPI::makeNicePlotStyle(summaryFirst.get());
0691 summaryFirst->SetMarkerColor(kRed);
0692 summaryFirst->GetXaxis()->LabelsOption("v");
0693 summaryFirst->GetXaxis()->SetLabelSize(0.05);
0694 summaryFirst->GetYaxis()->SetTitleOffset(0.9);
0695
0696 SiPixelPI::makeNicePlotStyle(summaryLast.get());
0697 summaryLast->SetMarkerColor(kBlue);
0698 summaryLast->GetYaxis()->SetTitleOffset(0.9);
0699 summaryLast->GetXaxis()->LabelsOption("v");
0700 summaryLast->GetXaxis()->SetLabelSize(0.05);
0701
0702 canvas.cd()->SetGridy();
0703
0704 canvas.SetBottomMargin(0.18);
0705 canvas.SetLeftMargin(0.11);
0706 canvas.SetRightMargin(0.02);
0707 canvas.Modified();
0708
0709 summaryFirst->SetFillColor(kRed);
0710 summaryLast->SetFillColor(kBlue);
0711
0712 summaryFirst->SetBarWidth(0.45);
0713 summaryFirst->SetBarOffset(0.1);
0714
0715 summaryLast->SetBarWidth(0.4);
0716 summaryLast->SetBarOffset(0.55);
0717
0718 summaryLast->SetMarkerSize(1.5);
0719 summaryFirst->SetMarkerSize(1.5);
0720
0721 float max = (summaryFirst->GetMaximum() > summaryLast->GetMaximum()) ? summaryFirst->GetMaximum()
0722 : summaryLast->GetMaximum();
0723
0724 summaryFirst->GetYaxis()->SetRangeUser(0., std::max(0., max * 1.40));
0725
0726 summaryFirst->Draw("b text0");
0727 summaryLast->Draw("b text0 same");
0728
0729 TLegend legend = TLegend(0.52, 0.80, 0.98, 0.9);
0730 legend.SetHeader("#mu_{H} value comparison", "C");
0731 std::string l_tagOrHash, f_tagOrHash;
0732 if (this->m_plotAnnotations.ntags == 2) {
0733 l_tagOrHash = l_tagname;
0734 f_tagOrHash = f_tagname;
0735 } else {
0736 l_tagOrHash = std::get<1>(lastiov);
0737 f_tagOrHash = std::get<1>(firstiov);
0738 }
0739
0740 legend.AddEntry(
0741 summaryLast.get(),
0742 ("IOV: #scale[1.2]{" + std::to_string(std::get<0>(lastiov)) + "} | #color[4]{" + l_tagOrHash + "}").c_str(),
0743 "F");
0744 legend.AddEntry(
0745 summaryFirst.get(),
0746 ("IOV: #scale[1.2]{" + std::to_string(std::get<0>(firstiov)) + "} | #color[2]{" + f_tagOrHash + "}").c_str(),
0747 "F");
0748
0749 legend.SetTextSize(0.025);
0750 legend.Draw("same");
0751
0752 std::string fileName(this->m_imageFileName);
0753 canvas.SaveAs(fileName.c_str());
0754 return true;
0755 }
0756 };
0757
0758 using SiPixelLorentzAngleByRegionComparisonSingleTag = SiPixelLorentzAngleByRegionComparisonBase<MULTI_IOV, 1>;
0759 using SiPixelLorentzAngleByRegionComparisonTwoTags = SiPixelLorentzAngleByRegionComparisonBase<SINGLE_IOV, 2>;
0760
0761
0762
0763
0764
0765 template <SiPixelPI::DetType myType>
0766 class SiPixelLorentzAngleMap : public PlotImage<SiPixelLorentzAngle, SINGLE_IOV> {
0767 public:
0768 SiPixelLorentzAngleMap()
0769 : PlotImage<SiPixelLorentzAngle, SINGLE_IOV>("SiPixelLorentzAngle Pixel Map"),
0770 m_trackerTopo{StandaloneTrackerTopology::fromTrackerParametersXMLFile(
0771 edm::FileInPath("Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml").fullPath())} {}
0772
0773 bool fill() override {
0774 auto tag = PlotBase::getTag<0>();
0775 auto iov = tag.iovs.front();
0776 auto tagname = tag.name;
0777 std::shared_ptr<SiPixelLorentzAngle> payload = fetchPayload(std::get<1>(iov));
0778
0779 Phase1PixelROCMaps thePixLAMap("");
0780
0781 std::map<uint32_t, float> LAMap_ = payload->getLorentzAngles();
0782 if (LAMap_.size() != SiPixelPI::phase1size) {
0783 edm::LogError("SiPixelLorentzAngle_PayloadInspector")
0784 << "SiPixelLorentzAngle maps are not supported for non-Phase1 Pixel geometries !";
0785 TCanvas canvas("Canv", "Canv", 1200, 1000);
0786 SiPixelPI::displayNotSupported(canvas, LAMap_.size());
0787 std::string fileName(m_imageFileName);
0788 canvas.SaveAs(fileName.c_str());
0789 return false;
0790 }
0791
0792
0793 std::array<double, n_layers> b_minima = {{999., 999., 999., 999.}};
0794 std::array<double, n_rings> f_minima = {{999., 999.}};
0795
0796 for (const auto &element : LAMap_) {
0797 int subid = DetId(element.first).subdetId();
0798 if (subid == PixelSubdetector::PixelBarrel) {
0799 auto layer = m_trackerTopo.pxbLayer(DetId(element.first));
0800 if (element.second < b_minima.at(layer - 1)) {
0801 b_minima.at(layer - 1) = element.second;
0802 }
0803 } else if (subid == PixelSubdetector::PixelEndcap) {
0804 auto ring = SiPixelPI::ring(DetId(element.first), m_trackerTopo, true);
0805 if (element.second < f_minima.at(ring - 1)) {
0806 f_minima.at(ring - 1) = element.second;
0807 }
0808 }
0809 thePixLAMap.fillWholeModule(element.first, element.second);
0810 }
0811
0812 gStyle->SetOptStat(0);
0813
0814 TCanvas canvas("Summary", "Summary", 1200, k_height[myType]);
0815 canvas.cd();
0816
0817 auto unpacked = SiPixelPI::unpack(std::get<0>(iov));
0818
0819 std::string IOVstring = (unpacked.first == 0)
0820 ? std::to_string(unpacked.second)
0821 : (std::to_string(unpacked.first) + "," + std::to_string(unpacked.second));
0822
0823 const auto headerText = fmt::sprintf("#color[4]{%s}, IOV: #color[4]{%s}", tagname, IOVstring);
0824
0825 switch (myType) {
0826 case SiPixelPI::t_barrel:
0827 thePixLAMap.drawBarrelMaps(canvas, headerText);
0828 break;
0829 case SiPixelPI::t_forward:
0830 thePixLAMap.drawForwardMaps(canvas, headerText);
0831 break;
0832 case SiPixelPI::t_all:
0833 thePixLAMap.drawMaps(canvas, headerText);
0834 break;
0835 default:
0836 throw cms::Exception("SiPixelLorentzAngleMap") << "\nERROR: unrecognized Pixel Detector part " << std::endl;
0837 }
0838
0839 if (myType == SiPixelPI::t_barrel || myType == SiPixelPI::t_all) {
0840
0841 for (unsigned int lay = 1; lay <= n_layers; lay++) {
0842 auto h_bpix_LA = thePixLAMap.getLayerMaps();
0843
0844 COUT << " layer:" << lay << " max:" << h_bpix_LA[lay - 1]->GetMaximum() << " min: " << b_minima.at(lay - 1)
0845 << std::endl;
0846
0847 h_bpix_LA[lay - 1]->GetZaxis()->SetRangeUser(b_minima.at(lay - 1) - 0.001,
0848 h_bpix_LA[lay - 1]->GetMaximum() + 0.001);
0849 }
0850 }
0851
0852 if (myType == SiPixelPI::t_forward || myType == SiPixelPI::t_all) {
0853
0854 for (unsigned int ring = 1; ring <= n_rings; ring++) {
0855 auto h_fpix_LA = thePixLAMap.getRingMaps();
0856
0857 COUT << " ringer:" << ring << " max:" << h_fpix_LA[ring - 1]->GetMaximum()
0858 << " min: " << f_minima.at(ring - 1) << std::endl;
0859
0860 h_fpix_LA[ring - 1]->GetZaxis()->SetRangeUser(f_minima.at(ring - 1) - 0.001,
0861 h_fpix_LA[ring - 1]->GetMaximum() + 0.001);
0862 }
0863 }
0864
0865 std::string fileName(m_imageFileName);
0866 canvas.SaveAs(fileName.c_str());
0867 #ifdef MMDEBUG
0868 canvas.SaveAs("outPixLA.root");
0869 #endif
0870
0871 return true;
0872 }
0873
0874 private:
0875 TrackerTopology m_trackerTopo;
0876 static constexpr std::array<int, 3> k_height = {{1200, 600, 1600}};
0877 static constexpr int n_layers = 4;
0878 static constexpr int n_rings = 2;
0879 };
0880
0881 using SiPixelBPixLorentzAngleMap = SiPixelLorentzAngleMap<SiPixelPI::t_barrel>;
0882 using SiPixelFPixLorentzAngleMap = SiPixelLorentzAngleMap<SiPixelPI::t_forward>;
0883 using SiPixelFullLorentzAngleMapByROC = SiPixelLorentzAngleMap<SiPixelPI::t_all>;
0884
0885
0886
0887
0888 class SiPixelLorentzAngleFullPixelMap : public PlotImage<SiPixelLorentzAngle, SINGLE_IOV> {
0889 public:
0890 SiPixelLorentzAngleFullPixelMap() : PlotImage<SiPixelLorentzAngle, SINGLE_IOV>("SiPixelLorentzAngle Map") {
0891 label_ = "SiPixelLorentzAngleFullPixelMap";
0892 payloadString = "Lorentz Angle";
0893 }
0894
0895 bool fill() override {
0896 gStyle->SetPalette(1);
0897 auto tag = PlotBase::getTag<0>();
0898 auto iov = tag.iovs.front();
0899 std::shared_ptr<SiPixelLorentzAngle> payload = this->fetchPayload(std::get<1>(iov));
0900
0901 if (payload.get()) {
0902 Phase1PixelSummaryMap fullMap(
0903 "", fmt::sprintf("%s", payloadString), fmt::sprintf("%s #mu_{H} [1/T]", payloadString));
0904 fullMap.createTrackerBaseMap();
0905
0906 std::map<uint32_t, float> LAMap_ = payload->getLorentzAngles();
0907
0908 if (LAMap_.size() == SiPixelPI::phase0size || LAMap_.size() > SiPixelPI::phase1size) {
0909 edm::LogError(label_) << "There are " << LAMap_.size()
0910 << " DetIds in this payload. SiPixelLorentzAngleFullPixelMap maps are not supported "
0911 "for non-Phase1 Pixel geometries !";
0912 TCanvas canvas("Canv", "Canv", 1200, 1000);
0913 SiPixelPI::displayNotSupported(canvas, LAMap_.size());
0914 std::string fileName(this->m_imageFileName);
0915 canvas.SaveAs(fileName.c_str());
0916 return false;
0917 } else {
0918 if (LAMap_.size() < SiPixelPI::phase1size) {
0919 edm::LogWarning(label_) << "\n ********* WARNING! ********* \n There are " << LAMap_.size()
0920 << " DetIds in this payload !"
0921 << "\n **************************** \n";
0922 }
0923 }
0924
0925 for (const auto &entry : LAMap_) {
0926 fullMap.fillTrackerMap(entry.first, entry.second);
0927 }
0928
0929 TCanvas canvas("Canv", "Canv", 3000, 2000);
0930 fullMap.printTrackerMap(canvas);
0931
0932 auto ltx = TLatex();
0933 ltx.SetTextFont(62);
0934 ltx.SetTextSize(0.025);
0935 ltx.SetTextAlign(11);
0936 ltx.DrawLatexNDC(
0937 gPad->GetLeftMargin() + 0.01,
0938 gPad->GetBottomMargin() + 0.01,
0939 ("#color[4]{" + tag.name + "}, IOV: #color[4]{" + std::to_string(std::get<0>(iov)) + "}").c_str());
0940
0941 std::string fileName(this->m_imageFileName);
0942 canvas.SaveAs(fileName.c_str());
0943 }
0944 return true;
0945 }
0946
0947 protected:
0948 std::string payloadString;
0949 std::string label_;
0950 };
0951
0952
0953
0954
0955 template <IOVMultiplicity nIOVs, int ntags>
0956 class SiPixelLorentzAngleFullMapCompareBase : public PlotImage<SiPixelLorentzAngle, nIOVs, ntags> {
0957 public:
0958 SiPixelLorentzAngleFullMapCompareBase() : PlotImage<SiPixelLorentzAngle, nIOVs, ntags>("SiPixelLorentzAngle Map") {
0959 label_ = "SiPixelLorentzAngleFullPixelMap";
0960 payloadString = "Lorentz Angle";
0961 }
0962
0963 bool fill() override {
0964 gStyle->SetPalette(kBlueRedYellow);
0965
0966
0967 auto theIOVs = PlotBase::getTag<0>().iovs;
0968 auto f_tagname = PlotBase::getTag<0>().name;
0969 std::string l_tagname = "";
0970 auto firstiov = theIOVs.front();
0971 std::tuple<cond::Time_t, cond::Hash> lastiov;
0972
0973
0974 assert(this->m_plotAnnotations.ntags < 3);
0975
0976 if (this->m_plotAnnotations.ntags == 2) {
0977 auto tag2iovs = PlotBase::getTag<1>().iovs;
0978 l_tagname = PlotBase::getTag<1>().name;
0979 lastiov = tag2iovs.front();
0980 } else {
0981 lastiov = theIOVs.back();
0982 }
0983
0984 std::shared_ptr<SiPixelLorentzAngle> last_payload = this->fetchPayload(std::get<1>(lastiov));
0985 std::shared_ptr<SiPixelLorentzAngle> first_payload = this->fetchPayload(std::get<1>(firstiov));
0986
0987 if (first_payload.get() && last_payload.get()) {
0988
0989 Phase1PixelSummaryMap fullMap("",
0990 fmt::sprintf("%s Diff", payloadString),
0991 fmt::sprintf("%s difference #Delta#mu_{H}/#mu_{H} [%%]", payloadString));
0992 fullMap.createTrackerBaseMap();
0993
0994 std::map<uint32_t, float> l_LAMap_ = last_payload->getLorentzAngles();
0995 std::map<uint32_t, float> f_LAMap_ = first_payload->getLorentzAngles();
0996
0997 std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
0998 std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
0999
1000 if (l_LAMap_.size() != f_LAMap_.size()) {
1001 edm::LogError(label_) << "There are " << l_LAMap_.size() << "dets in payload" << std::get<1>(lastiov)
1002 << "and " << f_LAMap_.size() << " dets in payload" << std::get<1>(firstiov)
1003 << "display is not possible!";
1004
1005 TCanvas canvas("Canv", "Canv", 1200, 1000);
1006 SiPixelPI::displayNotSupported(canvas, SiPixelPI::mismatched);
1007 std::string fileName(this->m_imageFileName);
1008 canvas.SaveAs(fileName.c_str());
1009 return false;
1010 }
1011
1012 if (l_LAMap_.size() == SiPixelPI::phase0size || l_LAMap_.size() > SiPixelPI::phase1size) {
1013 edm::LogError(label_) << "There are " << l_LAMap_.size()
1014 << " DetIds in this payload. SiPixelLorentzAngleFullPixelMap maps are not supported "
1015 "for non-Phase1 Pixel geometries !";
1016 TCanvas canvas("Canv", "Canv", 1200, 1000);
1017 SiPixelPI::displayNotSupported(canvas, l_LAMap_.size());
1018 std::string fileName(this->m_imageFileName);
1019 canvas.SaveAs(fileName.c_str());
1020 return false;
1021 } else {
1022 if (l_LAMap_.size() < SiPixelPI::phase1size) {
1023 edm::LogWarning(label_) << "\n ********* WARNING! ********* \n There are " << l_LAMap_.size()
1024 << " DetIds in this payload !"
1025 << "\n **************************** \n";
1026 }
1027 }
1028
1029
1030 for (const auto &[id, value] : l_LAMap_) {
1031 assert(value != 0.);
1032 const auto &diff = (value - f_LAMap_[id]) * 100.f / value;
1033 fullMap.fillTrackerMap(id, diff);
1034 }
1035
1036
1037 TCanvas canvas("Canv", "Canv", 3000, 2000);
1038 fullMap.printTrackerMap(canvas, 0.03);
1039
1040
1041 auto ltx = TLatex();
1042 ltx.SetTextFont(62);
1043 ltx.SetTextSize(0.025);
1044 ltx.SetTextAlign(11);
1045 std::string ltxText;
1046 if (this->m_plotAnnotations.ntags == 2) {
1047 ltxText = fmt::sprintf("#color[2]{%s, %s} vs #color[4]{%s, %s}",
1048 f_tagname,
1049 std::to_string(std::get<0>(firstiov)),
1050 l_tagname,
1051 std::to_string(std::get<0>(lastiov)));
1052 } else {
1053 ltxText = fmt::sprintf("%s IOV: #color[2]{%s} vs IOV: #color[4]{%s}",
1054 f_tagname,
1055 std::to_string(std::get<0>(firstiov)),
1056 std::to_string(std::get<0>(lastiov)));
1057 }
1058 ltx.DrawLatexNDC(gPad->GetLeftMargin(), 1 - gPad->GetTopMargin() + 0.01, ltxText.c_str());
1059
1060 std::string fileName(this->m_imageFileName);
1061 canvas.SaveAs(fileName.c_str());
1062 }
1063 return true;
1064 }
1065
1066 protected:
1067 std::string payloadString;
1068 std::string label_;
1069 };
1070
1071 using SiPixelLorentzAngleFullMapCompareSingleTag = SiPixelLorentzAngleFullMapCompareBase<MULTI_IOV, 1>;
1072 using SiPixelLorentzAngleFullMapCompareTwoTags = SiPixelLorentzAngleFullMapCompareBase<SINGLE_IOV, 2>;
1073
1074 }
1075
1076 PAYLOAD_INSPECTOR_MODULE(SiPixelLorentzAngle) {
1077 PAYLOAD_INSPECTOR_CLASS(SiPixelLorentzAngleValue);
1078 PAYLOAD_INSPECTOR_CLASS(SiPixelLorentzAngleValues);
1079 PAYLOAD_INSPECTOR_CLASS(SiPixelLorentzAngleValuesBarrel);
1080 PAYLOAD_INSPECTOR_CLASS(SiPixelLorentzAngleValuesEndcap);
1081 PAYLOAD_INSPECTOR_CLASS(SiPixelLorentzAngleValuesBarrelCompareSingleTag);
1082 PAYLOAD_INSPECTOR_CLASS(SiPixelLorentzAngleValuesEndcapCompareSingleTag);
1083 PAYLOAD_INSPECTOR_CLASS(SiPixelLorentzAngleValuesBarrelCompareTwoTags);
1084 PAYLOAD_INSPECTOR_CLASS(SiPixelLorentzAngleValuesEndcapCompareTwoTags);
1085 PAYLOAD_INSPECTOR_CLASS(SiPixelLorentzAngleValueComparisonSingleTag);
1086 PAYLOAD_INSPECTOR_CLASS(SiPixelLorentzAngleValueComparisonTwoTags);
1087 PAYLOAD_INSPECTOR_CLASS(SiPixelLorentzAngleByRegionComparisonSingleTag);
1088 PAYLOAD_INSPECTOR_CLASS(SiPixelLorentzAngleByRegionComparisonTwoTags);
1089 PAYLOAD_INSPECTOR_CLASS(SiPixelBPixLorentzAngleMap);
1090 PAYLOAD_INSPECTOR_CLASS(SiPixelFPixLorentzAngleMap);
1091 PAYLOAD_INSPECTOR_CLASS(SiPixelFullLorentzAngleMapByROC);
1092 PAYLOAD_INSPECTOR_CLASS(SiPixelLorentzAngleFullPixelMap);
1093 PAYLOAD_INSPECTOR_CLASS(SiPixelLorentzAngleFullMapCompareSingleTag);
1094 PAYLOAD_INSPECTOR_CLASS(SiPixelLorentzAngleFullMapCompareTwoTags);
1095 }