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