Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2022-02-02 23:43:09

0001 /*!
0002   \file SiPixelLorentzAngle_PayloadInspector
0003   \Payload Inspector Plugin for SiPixel Lorentz angles
0004   \author M. Musich
0005   \version $Revision: 1.0 $
0006   \date $Date: 2019/06/20 10:59:56 $
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 // include ROOT
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     1d histogram of SiPixelLorentzAngle of 1 IOV 
0041   *************************************************/
0042 
0043   // inherit from one of the predefined plot class: Histogram1D
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         }  // payload
0061       }    // iovs
0062       return true;
0063     }  // fill
0064   };
0065 
0066   /************************************************
0067     1d histogram of SiPixelLorentzAngle of 1 IOV 
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");  // option "C" allows to center the header
0112       //legend.AddEntry(h1.get(), ("IOV: " + std::to_string(std::get<0>(iov))).c_str(), "PL");
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       //ltx.SetTextColor(kBlue);
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     1d histogram of SiPixelLorentzAngle of 1 IOV per region
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");  // option "C" allows to center the header
0184       //legend.AddEntry(h1.get(), ("IOV: " + std::to_string(std::get<0>(iov))).c_str(), "PL");
0185       legend.SetTextSize(0.025);
0186       legend.SetLineColor(10);
0187 
0188       unsigned int maxPads = canvas.GetListOfPrimitives()->GetSize();  //= isBarrel ? 4 : 12;
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     1d histogram of SiPixelLorentzAngle of 2 IOV per region
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       // trick to deal with the multi-ioved tag and two tag case at the same time
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       // we don't support (yet) comparison with more than 2 tags
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       // deal with last IOV
0288       const char *path_toTopologyXML = l_phaseInfo.pathToTopoXML();
0289 
0290       auto l_tTopo =
0291           StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
0292 
0293       auto l_myPlots = PixelRegions::PixelRegionContainers(&l_tTopo, l_phaseInfo.phase());
0294       l_myPlots.bookAll("SiPixel LA,last",
0295                         "SiPixel LorentzAngle #mu_{H}(tan#theta_{L}/B) [1/T]",
0296                         "#modules",
0297                         50,
0298                         min * 0.9,
0299                         max * 1.1);
0300 
0301       for (const auto &element : l_LAMap_) {
0302         l_myPlots.fill(element.first, element.second);
0303       }
0304 
0305       l_myPlots.beautify();
0306       l_myPlots.draw(canvas, isBarrel, "bar2", f_phaseInfo.isPhase1Comparison(l_phaseInfo));
0307 
0308       // deal with first IOV
0309       path_toTopologyXML = f_phaseInfo.pathToTopoXML();
0310 
0311       auto f_tTopo =
0312           StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
0313 
0314       auto f_myPlots = PixelRegions::PixelRegionContainers(&f_tTopo, f_phaseInfo.phase());
0315       f_myPlots.bookAll("SiPixel LA,first",
0316                         "SiPixel LorentzAngle #mu_{H}(tan#theta_{L}/B) [1/T]",
0317                         "#modules",
0318                         50,
0319                         min * 0.9,
0320                         max * 1.1);
0321 
0322       for (const auto &element : f_LAMap_) {
0323         f_myPlots.fill(element.first, element.second);
0324       }
0325 
0326       f_myPlots.beautify(kAzure, kBlue);
0327       f_myPlots.draw(canvas, isBarrel, "HISTsames", f_phaseInfo.isPhase1Comparison(l_phaseInfo));
0328 
0329       // rescale the y-axis ranges in order to fit the canvas
0330       l_myPlots.rescaleMax(f_myPlots);
0331 
0332       // done dealing with IOVs
0333 
0334       auto colorTag = PixelRegions::L1;  //: PixelRegions::Rm1l;
0335       std::unique_ptr<TLegend> legend;
0336       if (this->m_plotAnnotations.ntags == 2) {
0337         legend = std::make_unique<TLegend>(0.36, 0.86, 0.94, 0.92);
0338         legend->AddEntry(l_myPlots.getHistoFromMap(colorTag).get(), ("#color[2]{" + l_tagname + "}").c_str(), "F");
0339         legend->AddEntry(f_myPlots.getHistoFromMap(colorTag).get(), ("#color[4]{" + f_tagname + "}").c_str(), "F");
0340         legend->SetTextSize(0.024);
0341       } else {
0342         legend = std::make_unique<TLegend>(0.58, 0.80, 0.90, 0.92);
0343         legend->AddEntry(l_myPlots.getHistoFromMap(colorTag).get(), ("#color[2]{" + lastIOVsince + "}").c_str(), "F");
0344         legend->AddEntry(f_myPlots.getHistoFromMap(colorTag).get(), ("#color[4]{" + firstIOVsince + "}").c_str(), "F");
0345         legend->SetTextSize(0.040);
0346       }
0347       legend->SetLineColor(10);
0348 
0349       unsigned int maxPads = canvas.GetListOfPrimitives()->GetSize();  //= isBarrel ? 4 : 12;
0350       for (unsigned int c = 1; c <= maxPads; c++) {
0351         if (l_phaseInfo.phase() == SiPixelPI::phase::two && (c == 5 || c == 10))
0352           continue;
0353         canvas.cd(c);
0354         SiPixelPI::adjustCanvasMargins(canvas.cd(c), 0.06, 0.12, 0.12, 0.05);
0355         legend->Draw("same");
0356         canvas.cd(c)->Update();
0357       }
0358 
0359       f_myPlots.stats(0);
0360       l_myPlots.stats(1);
0361 
0362       auto ltx = TLatex();
0363       ltx.SetTextFont(62);
0364       ltx.SetTextSize(0.05);
0365       ltx.SetTextAlign(11);
0366 
0367       int index = 0;
0368       for (unsigned int c = 1; c <= maxPads; c++) {
0369         if (l_phaseInfo.phase() == SiPixelPI::phase::two && (c == 5 || c == 10))
0370           continue;
0371         canvas.cd(c);
0372 
0373         COUT << "c:" << c << " index:" << index << " : "
0374              << PixelRegions::getIDLabels(l_phaseInfo.phase(), isBarrel)[index] << "\n";
0375 
0376         ltx.DrawLatexNDC(
0377             gPad->GetLeftMargin(),
0378             1 - gPad->GetTopMargin() + 0.01,
0379             (PixelRegions::getIDLabels(l_phaseInfo.phase(), isBarrel)[index] + " : #color[4]{" +
0380              std::to_string(std::get<0>(firstiov)) + "} vs #color[2]{" + std::to_string(std::get<0>(lastiov)) + "}")
0381                 .c_str());
0382         index++;
0383       }
0384 
0385       std::string fileName(this->m_imageFileName);
0386       canvas.SaveAs(fileName.c_str());
0387 
0388 #ifdef MMDEBUG
0389       canvas.SaveAs("DEBUG.root");
0390 #endif
0391 
0392       return true;
0393     }
0394   };
0395 
0396   using SiPixelLorentzAngleValuesBarrelCompareSingleTag =
0397       SiPixelLorentzAngleValuesComparisonPerRegion<true, MULTI_IOV, 1>;
0398   using SiPixelLorentzAngleValuesEndcapCompareSingleTag =
0399       SiPixelLorentzAngleValuesComparisonPerRegion<false, MULTI_IOV, 1>;
0400 
0401   using SiPixelLorentzAngleValuesBarrelCompareTwoTags =
0402       SiPixelLorentzAngleValuesComparisonPerRegion<true, SINGLE_IOV, 2>;
0403   using SiPixelLorentzAngleValuesEndcapCompareTwoTags =
0404       SiPixelLorentzAngleValuesComparisonPerRegion<false, SINGLE_IOV, 2>;
0405 
0406   /************************************************
0407     1d histogram of SiPixelLorentzAngle comparisons
0408   *************************************************/
0409   template <IOVMultiplicity nIOVs, int ntags>
0410   class SiPixelLorentzAngleValueComparisonBase : public PlotImage<SiPixelLorentzAngle, nIOVs, ntags> {
0411   public:
0412     SiPixelLorentzAngleValueComparisonBase()
0413         : PlotImage<SiPixelLorentzAngle, nIOVs, ntags>(Form("SiPixelLorentzAngle Values Comparison %i tag(s)", ntags)) {
0414     }
0415 
0416     bool fill() override {
0417       TH1F::SetDefaultSumw2(true);
0418 
0419       // trick to deal with the multi-ioved tag and two tag case at the same time
0420       auto theIOVs = PlotBase::getTag<0>().iovs;
0421       auto f_tagname = PlotBase::getTag<0>().name;
0422       std::string l_tagname = "";
0423       auto firstiov = theIOVs.front();
0424       std::tuple<cond::Time_t, cond::Hash> lastiov;
0425 
0426       // we don't support (yet) comparison with more than 2 tags
0427       assert(this->m_plotAnnotations.ntags < 3);
0428 
0429       if (this->m_plotAnnotations.ntags == 2) {
0430         auto tag2iovs = PlotBase::getTag<1>().iovs;
0431         l_tagname = PlotBase::getTag<1>().name;
0432         lastiov = tag2iovs.front();
0433       } else {
0434         lastiov = theIOVs.back();
0435       }
0436 
0437       std::shared_ptr<SiPixelLorentzAngle> last_payload = this->fetchPayload(std::get<1>(lastiov));
0438       std::map<uint32_t, float> l_LAMap_ = last_payload->getLorentzAngles();
0439       auto l_extrema = SiPixelPI::findMinMaxInMap(l_LAMap_);
0440 
0441       std::shared_ptr<SiPixelLorentzAngle> first_payload = this->fetchPayload(std::get<1>(firstiov));
0442       std::map<uint32_t, float> f_LAMap_ = first_payload->getLorentzAngles();
0443       auto f_extrema = SiPixelPI::findMinMaxInMap(f_LAMap_);
0444 
0445       auto max = (l_extrema.second > f_extrema.second) ? l_extrema.second : f_extrema.second;
0446       auto min = (l_extrema.first < f_extrema.first) ? l_extrema.first : f_extrema.first;
0447 
0448       std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
0449       std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
0450 
0451       TCanvas canvas("Canv", "Canv", 1200, 1000);
0452       canvas.cd();
0453       auto hfirst =
0454           std::make_unique<TH1F>("value_first",
0455                                  "SiPixel LA value;SiPixel LorentzAngle #mu_{H}(tan#theta_{L}/B) [1/T];# modules",
0456                                  50,
0457                                  min * 0.9,
0458                                  max * 1.1);
0459       hfirst->SetStats(false);
0460 
0461       auto hlast =
0462           std::make_unique<TH1F>("value_last",
0463                                  "SiPixel LA value;SiPixel LorentzAngle #mu_{H}(tan#theta_{L}/B) [1/T];# modules",
0464                                  50,
0465                                  min * 0.9,
0466                                  max * 1.1);
0467       hlast->SetStats(false);
0468 
0469       SiPixelPI::adjustCanvasMargins(canvas.cd(), 0.06, 0.12, 0.12, 0.05);
0470       canvas.Modified();
0471 
0472       for (const auto &element : f_LAMap_) {
0473         hfirst->Fill(element.second);
0474       }
0475 
0476       for (const auto &element : l_LAMap_) {
0477         hlast->Fill(element.second);
0478       }
0479 
0480       auto extrema = SiPixelPI::getExtrema(hfirst.get(), hlast.get());
0481       hfirst->GetYaxis()->SetRangeUser(extrema.first, extrema.second * 1.10);
0482 
0483       hfirst->SetTitle("");
0484       hfirst->SetFillColor(kRed);
0485       hfirst->SetBarWidth(0.95);
0486       hfirst->Draw("histbar");
0487 
0488       hlast->SetTitle("");
0489       hlast->SetFillColorAlpha(kBlue, 0.20);
0490       hlast->SetBarWidth(0.95);
0491       hlast->Draw("histbarsame");
0492 
0493       SiPixelPI::makeNicePlotStyle(hfirst.get());
0494       SiPixelPI::makeNicePlotStyle(hlast.get());
0495 
0496       canvas.Update();
0497 
0498       TLegend legend = TLegend(0.30, 0.86, 0.95, 0.94);
0499       //legend.SetHeader("#font[22]{SiPixel Lorentz Angle Comparison}", "C");  // option "C" allows to center the header
0500       //legend.AddEntry(hfirst.get(), ("IOV: " + std::to_string(std::get<0>(firstiov))).c_str(), "FL");
0501       //legend.AddEntry(hlast.get(),  ("IOV: " + std::to_string(std::get<0>(lastiov))).c_str(), "FL");
0502       legend.AddEntry(hfirst.get(), ("payload: #color[2]{" + std::get<1>(firstiov) + "}").c_str(), "F");
0503       legend.AddEntry(hlast.get(), ("payload: #color[4]{" + std::get<1>(lastiov) + "}").c_str(), "F");
0504       legend.SetTextSize(0.025);
0505       legend.Draw("same");
0506 
0507       auto ltx = TLatex();
0508       ltx.SetTextFont(62);
0509       //ltx.SetTextColor(kBlue);
0510       ltx.SetTextSize(0.047);
0511       ltx.SetTextAlign(11);
0512       std::string ltxText;
0513       if (this->m_plotAnnotations.ntags == 2) {
0514         ltxText = fmt::sprintf("#color[2]{%s, %s} vs #color[4]{%s, %s}",
0515                                f_tagname,
0516                                std::to_string(std::get<0>(firstiov)),
0517                                l_tagname,
0518                                std::to_string(std::get<0>(lastiov)));
0519       } else {
0520         ltxText = fmt::sprintf("%s IOV: #color[2]{%s} vs IOV: #color[4]{%s}",
0521                                f_tagname,
0522                                std::to_string(std::get<0>(firstiov)),
0523                                std::to_string(std::get<0>(lastiov)));
0524       }
0525       ltx.DrawLatexNDC(gPad->GetLeftMargin(), 1 - gPad->GetTopMargin() + 0.01, ltxText.c_str());
0526 
0527       std::string fileName(this->m_imageFileName);
0528       canvas.SaveAs(fileName.c_str());
0529 
0530       return true;
0531     }
0532   };
0533 
0534   using SiPixelLorentzAngleValueComparisonSingleTag = SiPixelLorentzAngleValueComparisonBase<MULTI_IOV, 1>;
0535   using SiPixelLorentzAngleValueComparisonTwoTags = SiPixelLorentzAngleValueComparisonBase<SINGLE_IOV, 2>;
0536 
0537   /************************************************
0538    Summary Comparison per region of SiPixelLorentzAngle between 2 IOVs
0539   *************************************************/
0540   template <IOVMultiplicity nIOVs, int ntags>
0541   class SiPixelLorentzAngleByRegionComparisonBase : public PlotImage<SiPixelLorentzAngle, nIOVs, ntags> {
0542   public:
0543     SiPixelLorentzAngleByRegionComparisonBase()
0544         : PlotImage<SiPixelLorentzAngle, nIOVs, ntags>(
0545               Form("SiPixelLorentzAngle Comparison by Region %i tag(s)", ntags)) {}
0546 
0547     bool fill() override {
0548       gStyle->SetPaintTextFormat(".3f");
0549 
0550       // trick to deal with the multi-ioved tag and two tag case at the same time
0551       auto theIOVs = PlotBase::getTag<0>().iovs;
0552       auto f_tagname = PlotBase::getTag<0>().name;
0553       std::string l_tagname = "";
0554       auto firstiov = theIOVs.front();
0555       std::tuple<cond::Time_t, cond::Hash> lastiov;
0556 
0557       // we don't support (yet) comparison with more than 2 tags
0558       assert(this->m_plotAnnotations.ntags < 3);
0559 
0560       if (this->m_plotAnnotations.ntags == 2) {
0561         auto tag2iovs = PlotBase::getTag<1>().iovs;
0562         l_tagname = PlotBase::getTag<1>().name;
0563         lastiov = tag2iovs.front();
0564       } else {
0565         lastiov = theIOVs.back();
0566       }
0567 
0568       std::shared_ptr<SiPixelLorentzAngle> last_payload = this->fetchPayload(std::get<1>(lastiov));
0569       std::map<uint32_t, float> l_LAMap_ = last_payload->getLorentzAngles();
0570       std::shared_ptr<SiPixelLorentzAngle> first_payload = this->fetchPayload(std::get<1>(firstiov));
0571       std::map<uint32_t, float> f_LAMap_ = first_payload->getLorentzAngles();
0572 
0573       std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
0574       std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
0575 
0576       SiPixelPI::PhaseInfo l_phaseInfo(l_LAMap_.size());
0577       SiPixelPI::PhaseInfo f_phaseInfo(f_LAMap_.size());
0578 
0579       TCanvas canvas("Comparison", "Comparison", 1600, 800);
0580       std::map<SiPixelPI::regions, std::shared_ptr<TH1F>> FirstLA_spectraByRegion;
0581       std::map<SiPixelPI::regions, std::shared_ptr<TH1F>> LastLA_spectraByRegion;
0582       std::shared_ptr<TH1F> summaryFirst;
0583       std::shared_ptr<TH1F> summaryLast;
0584 
0585       // book the intermediate histograms
0586       for (int r = SiPixelPI::BPixL1o; r != SiPixelPI::NUM_OF_REGIONS; r++) {
0587         SiPixelPI::regions part = static_cast<SiPixelPI::regions>(r);
0588         std::string s_part = SiPixelPI::getStringFromRegionEnum(part);
0589 
0590         FirstLA_spectraByRegion[part] = std::make_shared<TH1F>(Form("hfirstLA_%s", s_part.c_str()),
0591                                                                Form(";%s #mu_{H} [1/T];n. of modules", s_part.c_str()),
0592                                                                1000,
0593                                                                0.,
0594                                                                1000.);
0595         LastLA_spectraByRegion[part] = std::make_shared<TH1F>(Form("hlastLA_%s", s_part.c_str()),
0596                                                               Form(";%s #mu_{H} [1/T];n. of modules", s_part.c_str()),
0597                                                               1000,
0598                                                               0.,
0599                                                               1000.);
0600       }
0601 
0602       summaryFirst = std::make_shared<TH1F>("first Summary",
0603                                             "Summary for #LT tan#theta_{L}/B #GT;;average LA #LT #mu_{H} #GT [1/T]",
0604                                             FirstLA_spectraByRegion.size(),
0605                                             0,
0606                                             FirstLA_spectraByRegion.size());
0607       summaryLast = std::make_shared<TH1F>("last Summary",
0608                                            "Summary for #LT tan#theta_{L}/B #GT;;average LA #LT #mu_{H}  #GT [1/T]",
0609                                            LastLA_spectraByRegion.size(),
0610                                            0,
0611                                            LastLA_spectraByRegion.size());
0612 
0613       // deal with first IOV
0614       const char *path_toTopologyXML = f_phaseInfo.pathToTopoXML();
0615 
0616       auto f_tTopo =
0617           StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
0618 
0619       // -------------------------------------------------------------------
0620       // loop on the first LA Map
0621       // -------------------------------------------------------------------
0622       for (const auto &it : f_LAMap_) {
0623         if (DetId(it.first).det() != DetId::Tracker) {
0624           edm::LogWarning("SiPixelLorentzAngle_PayloadInspector")
0625               << "Encountered invalid Tracker DetId:" << it.first << " - terminating ";
0626           return false;
0627         }
0628 
0629         SiPixelPI::topolInfo t_info_fromXML;
0630         t_info_fromXML.init();
0631         DetId detid(it.first);
0632         t_info_fromXML.fillGeometryInfo(detid, f_tTopo, f_phaseInfo.phase());
0633 
0634         SiPixelPI::regions thePart = t_info_fromXML.filterThePartition();
0635         if (thePart != SiPixelPI::NUM_OF_REGIONS) {
0636           FirstLA_spectraByRegion[thePart]->Fill(it.second);
0637         }
0638       }  // ends loop on the vector of error transforms
0639 
0640       // deal with last IOV
0641       path_toTopologyXML = l_phaseInfo.pathToTopoXML();
0642 
0643       auto l_tTopo =
0644           StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
0645 
0646       // -------------------------------------------------------------------
0647       // loop on the second LA Map
0648       // -------------------------------------------------------------------
0649       for (const auto &it : l_LAMap_) {
0650         if (DetId(it.first).det() != DetId::Tracker) {
0651           edm::LogWarning("SiPixelLorentzAngle_PayloadInspector")
0652               << "Encountered invalid Tracker DetId:" << it.first << " - terminating ";
0653           return false;
0654         }
0655 
0656         SiPixelPI::topolInfo t_info_fromXML;
0657         t_info_fromXML.init();
0658         DetId detid(it.first);
0659         t_info_fromXML.fillGeometryInfo(detid, l_tTopo, l_phaseInfo.phase());
0660 
0661         SiPixelPI::regions thePart = t_info_fromXML.filterThePartition();
0662         if (thePart != SiPixelPI::NUM_OF_REGIONS) {
0663           LastLA_spectraByRegion[thePart]->Fill(it.second);
0664         }
0665       }  // ends loop on the vector of error transforms
0666 
0667       // fill the summary plots
0668       int bin = 1;
0669       for (int r = SiPixelPI::BPixL1o; r != SiPixelPI::NUM_OF_REGIONS; r++) {
0670         SiPixelPI::regions part = static_cast<SiPixelPI::regions>(r);
0671 
0672         summaryFirst->GetXaxis()->SetBinLabel(bin, SiPixelPI::getStringFromRegionEnum(part).c_str());
0673         // avoid filling the histogram with numerical noise
0674         float f_mean =
0675             FirstLA_spectraByRegion[part]->GetMean() > 10.e-6 ? FirstLA_spectraByRegion[part]->GetMean() : 0.;
0676         summaryFirst->SetBinContent(bin, f_mean);
0677         //summaryFirst->SetBinError(bin,LA_spectraByRegion[hash]->GetRMS());
0678 
0679         summaryLast->GetXaxis()->SetBinLabel(bin, SiPixelPI::getStringFromRegionEnum(part).c_str());
0680         // avoid filling the histogram with numerical noise
0681         float l_mean = LastLA_spectraByRegion[part]->GetMean() > 10.e-6 ? LastLA_spectraByRegion[part]->GetMean() : 0.;
0682         summaryLast->SetBinContent(bin, l_mean);
0683         //summaryLast->SetBinError(bin,LA_spectraByRegion[hash]->GetRMS());
0684         bin++;
0685       }
0686 
0687       SiPixelPI::makeNicePlotStyle(summaryFirst.get());  //, kBlue);
0688       summaryFirst->SetMarkerColor(kRed);
0689       summaryFirst->GetXaxis()->LabelsOption("v");
0690       summaryFirst->GetXaxis()->SetLabelSize(0.05);
0691       summaryFirst->GetYaxis()->SetTitleOffset(0.9);
0692 
0693       SiPixelPI::makeNicePlotStyle(summaryLast.get());  //, kRed);
0694       summaryLast->SetMarkerColor(kBlue);
0695       summaryLast->GetYaxis()->SetTitleOffset(0.9);
0696       summaryLast->GetXaxis()->LabelsOption("v");
0697       summaryLast->GetXaxis()->SetLabelSize(0.05);
0698 
0699       canvas.cd()->SetGridy();
0700 
0701       canvas.SetBottomMargin(0.18);
0702       canvas.SetLeftMargin(0.11);
0703       canvas.SetRightMargin(0.02);
0704       canvas.Modified();
0705 
0706       summaryFirst->SetFillColor(kRed);
0707       summaryLast->SetFillColor(kBlue);
0708 
0709       summaryFirst->SetBarWidth(0.45);
0710       summaryFirst->SetBarOffset(0.1);
0711 
0712       summaryLast->SetBarWidth(0.4);
0713       summaryLast->SetBarOffset(0.55);
0714 
0715       summaryLast->SetMarkerSize(1.5);
0716       summaryFirst->SetMarkerSize(1.5);
0717 
0718       float max = (summaryFirst->GetMaximum() > summaryLast->GetMaximum()) ? summaryFirst->GetMaximum()
0719                                                                            : summaryLast->GetMaximum();
0720 
0721       summaryFirst->GetYaxis()->SetRangeUser(0., std::max(0., max * 1.40));
0722 
0723       summaryFirst->Draw("b text0");
0724       summaryLast->Draw("b text0 same");
0725 
0726       TLegend legend = TLegend(0.52, 0.80, 0.98, 0.9);
0727       legend.SetHeader("#mu_{H} value comparison", "C");  // option "C" allows to center the header
0728       std::string l_tagOrHash, f_tagOrHash;
0729       if (this->m_plotAnnotations.ntags == 2) {
0730         l_tagOrHash = l_tagname;
0731         f_tagOrHash = f_tagname;
0732       } else {
0733         l_tagOrHash = std::get<1>(lastiov);
0734         f_tagOrHash = std::get<1>(firstiov);
0735       }
0736 
0737       legend.AddEntry(
0738           summaryLast.get(),
0739           ("IOV: #scale[1.2]{" + std::to_string(std::get<0>(lastiov)) + "} | #color[4]{" + l_tagOrHash + "}").c_str(),
0740           "F");
0741       legend.AddEntry(
0742           summaryFirst.get(),
0743           ("IOV: #scale[1.2]{" + std::to_string(std::get<0>(firstiov)) + "} | #color[2]{" + f_tagOrHash + "}").c_str(),
0744           "F");
0745 
0746       legend.SetTextSize(0.025);
0747       legend.Draw("same");
0748 
0749       std::string fileName(this->m_imageFileName);
0750       canvas.SaveAs(fileName.c_str());
0751       return true;
0752     }
0753   };
0754 
0755   using SiPixelLorentzAngleByRegionComparisonSingleTag = SiPixelLorentzAngleByRegionComparisonBase<MULTI_IOV, 1>;
0756   using SiPixelLorentzAngleByRegionComparisonTwoTags = SiPixelLorentzAngleByRegionComparisonBase<SINGLE_IOV, 2>;
0757 
0758   /************************************************
0759    occupancy style map Pixel LA
0760   *************************************************/
0761 
0762   template <SiPixelPI::DetType myType>
0763   class SiPixelLorentzAngleMap : public PlotImage<SiPixelLorentzAngle, SINGLE_IOV> {
0764   public:
0765     SiPixelLorentzAngleMap()
0766         : PlotImage<SiPixelLorentzAngle, SINGLE_IOV>("SiPixelLorentzAngle Pixel Map"),
0767           m_trackerTopo{StandaloneTrackerTopology::fromTrackerParametersXMLFile(
0768               edm::FileInPath("Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml").fullPath())} {}
0769 
0770     bool fill() override {
0771       auto tag = PlotBase::getTag<0>();
0772       auto iov = tag.iovs.front();
0773       auto tagname = tag.name;
0774       std::shared_ptr<SiPixelLorentzAngle> payload = fetchPayload(std::get<1>(iov));
0775 
0776       Phase1PixelROCMaps thePixLAMap("");
0777 
0778       std::map<uint32_t, float> LAMap_ = payload->getLorentzAngles();
0779       if (LAMap_.size() != SiPixelPI::phase1size) {
0780         edm::LogError("SiPixelLorentzAngle_PayloadInspector")
0781             << "SiPixelLorentzAngle maps are not supported for non-Phase1 Pixel geometries !";
0782         TCanvas canvas("Canv", "Canv", 1200, 1000);
0783         SiPixelPI::displayNotSupported(canvas, LAMap_.size());
0784         std::string fileName(m_imageFileName);
0785         canvas.SaveAs(fileName.c_str());
0786         return false;
0787       }
0788 
0789       // hard-coded phase-I
0790       std::array<double, n_layers> b_minima = {{999., 999., 999., 999.}};
0791       std::array<double, n_rings> f_minima = {{999., 999.}};
0792 
0793       for (const auto &element : LAMap_) {
0794         int subid = DetId(element.first).subdetId();
0795         if (subid == PixelSubdetector::PixelBarrel) {
0796           auto layer = m_trackerTopo.pxbLayer(DetId(element.first));
0797           if (element.second < b_minima.at(layer - 1)) {
0798             b_minima.at(layer - 1) = element.second;
0799           }
0800         } else if (subid == PixelSubdetector::PixelEndcap) {
0801           auto ring = SiPixelPI::ring(DetId(element.first), m_trackerTopo, true);
0802           if (element.second < f_minima.at(ring - 1)) {
0803             f_minima.at(ring - 1) = element.second;
0804           }
0805         }
0806         thePixLAMap.fillWholeModule(element.first, element.second);
0807       }
0808 
0809       gStyle->SetOptStat(0);
0810       //=========================
0811       TCanvas canvas("Summary", "Summary", 1200, k_height[myType]);
0812       canvas.cd();
0813 
0814       auto unpacked = SiPixelPI::unpack(std::get<0>(iov));
0815 
0816       std::string IOVstring = (unpacked.first == 0)
0817                                   ? std::to_string(unpacked.second)
0818                                   : (std::to_string(unpacked.first) + "," + std::to_string(unpacked.second));
0819 
0820       const auto headerText = fmt::sprintf("#color[4]{%s},  IOV: #color[4]{%s}", tagname, IOVstring);
0821 
0822       switch (myType) {
0823         case SiPixelPI::t_barrel:
0824           thePixLAMap.drawBarrelMaps(canvas, headerText);
0825           break;
0826         case SiPixelPI::t_forward:
0827           thePixLAMap.drawForwardMaps(canvas, headerText);
0828           break;
0829         case SiPixelPI::t_all:
0830           thePixLAMap.drawMaps(canvas, headerText);
0831           break;
0832         default:
0833           throw cms::Exception("SiPixelLorentzAngleMap") << "\nERROR: unrecognized Pixel Detector part " << std::endl;
0834       }
0835 
0836       if (myType == SiPixelPI::t_barrel || myType == SiPixelPI::t_all) {
0837         // set the minima and maxima of the barrel
0838         for (unsigned int lay = 1; lay <= n_layers; lay++) {
0839           auto h_bpix_LA = thePixLAMap.getLayerMaps();
0840 
0841           COUT << " layer:" << lay << " max:" << h_bpix_LA[lay - 1]->GetMaximum() << " min: " << b_minima.at(lay - 1)
0842                << std::endl;
0843 
0844           h_bpix_LA[lay - 1]->GetZaxis()->SetRangeUser(b_minima.at(lay - 1) - 0.001,
0845                                                        h_bpix_LA[lay - 1]->GetMaximum() + 0.001);
0846         }
0847       }
0848 
0849       if (myType == SiPixelPI::t_forward || myType == SiPixelPI::t_all) {
0850         // set the minima and maxima of the endcaps
0851         for (unsigned int ring = 1; ring <= n_rings; ring++) {
0852           auto h_fpix_LA = thePixLAMap.getRingMaps();
0853 
0854           COUT << " ringer:" << ring << " max:" << h_fpix_LA[ring - 1]->GetMaximum()
0855                << " min: " << f_minima.at(ring - 1) << std::endl;
0856 
0857           h_fpix_LA[ring - 1]->GetZaxis()->SetRangeUser(f_minima.at(ring - 1) - 0.001,
0858                                                         h_fpix_LA[ring - 1]->GetMaximum() + 0.001);
0859         }
0860       }
0861 
0862       std::string fileName(m_imageFileName);
0863       canvas.SaveAs(fileName.c_str());
0864 #ifdef MMDEBUG
0865       canvas.SaveAs("outPixLA.root");
0866 #endif
0867 
0868       return true;
0869     }
0870 
0871   private:
0872     TrackerTopology m_trackerTopo;
0873     static constexpr std::array<int, 3> k_height = {{1200, 600, 1600}};
0874     static constexpr int n_layers = 4;
0875     static constexpr int n_rings = 2;
0876   };
0877 
0878   using SiPixelBPixLorentzAngleMap = SiPixelLorentzAngleMap<SiPixelPI::t_barrel>;
0879   using SiPixelFPixLorentzAngleMap = SiPixelLorentzAngleMap<SiPixelPI::t_forward>;
0880   using SiPixelFullLorentzAngleMapByROC = SiPixelLorentzAngleMap<SiPixelPI::t_all>;
0881 
0882   /************************************************
0883    Full Pixel Tracker Map class
0884   *************************************************/
0885   class SiPixelLorentzAngleFullPixelMap : public PlotImage<SiPixelLorentzAngle, SINGLE_IOV> {
0886   public:
0887     SiPixelLorentzAngleFullPixelMap() : PlotImage<SiPixelLorentzAngle, SINGLE_IOV>("SiPixelLorentzAngle Map") {
0888       label_ = "SiPixelLorentzAngleFullPixelMap";
0889       payloadString = "Lorentz Angle";
0890     }
0891 
0892     bool fill() override {
0893       gStyle->SetPalette(1);
0894       auto tag = PlotBase::getTag<0>();
0895       auto iov = tag.iovs.front();
0896       std::shared_ptr<SiPixelLorentzAngle> payload = this->fetchPayload(std::get<1>(iov));
0897 
0898       if (payload.get()) {
0899         Phase1PixelSummaryMap fullMap(
0900             "", fmt::sprintf("%s", payloadString), fmt::sprintf("%s #mu_{H} [1/T]", payloadString));
0901         fullMap.createTrackerBaseMap();
0902 
0903         std::map<uint32_t, float> LAMap_ = payload->getLorentzAngles();
0904 
0905         if (LAMap_.size() == SiPixelPI::phase0size || LAMap_.size() > SiPixelPI::phase1size) {
0906           edm::LogError(label_) << "There are " << LAMap_.size()
0907                                 << " DetIds in this payload. SiPixelLorentzAngleFullPixelMap maps are not supported "
0908                                    "for non-Phase1 Pixel geometries !";
0909           TCanvas canvas("Canv", "Canv", 1200, 1000);
0910           SiPixelPI::displayNotSupported(canvas, LAMap_.size());
0911           std::string fileName(this->m_imageFileName);
0912           canvas.SaveAs(fileName.c_str());
0913           return false;
0914         } else {
0915           if (LAMap_.size() < SiPixelPI::phase1size) {
0916             edm::LogWarning(label_) << "\n ********* WARNING! ********* \n There are " << LAMap_.size()
0917                                     << " DetIds in this payload !"
0918                                     << "\n **************************** \n";
0919           }
0920         }
0921 
0922         for (const auto &entry : LAMap_) {
0923           fullMap.fillTrackerMap(entry.first, entry.second);
0924         }
0925 
0926         TCanvas canvas("Canv", "Canv", 3000, 2000);
0927         fullMap.printTrackerMap(canvas);
0928 
0929         auto ltx = TLatex();
0930         ltx.SetTextFont(62);
0931         ltx.SetTextSize(0.025);
0932         ltx.SetTextAlign(11);
0933         ltx.DrawLatexNDC(
0934             gPad->GetLeftMargin() + 0.01,
0935             gPad->GetBottomMargin() + 0.01,
0936             ("#color[4]{" + tag.name + "}, IOV: #color[4]{" + std::to_string(std::get<0>(iov)) + "}").c_str());
0937 
0938         std::string fileName(this->m_imageFileName);
0939         canvas.SaveAs(fileName.c_str());
0940       }
0941       return true;
0942     }
0943 
0944   protected:
0945     std::string payloadString;
0946     std::string label_;
0947   };
0948 
0949   /************************************************
0950    Full Pixel Tracker Map for Comparison Base Class
0951   *************************************************/
0952   template <IOVMultiplicity nIOVs, int ntags>
0953   class SiPixelLorentzAngleFullMapCompareBase : public PlotImage<SiPixelLorentzAngle, nIOVs, ntags> {
0954   public:
0955     SiPixelLorentzAngleFullMapCompareBase() : PlotImage<SiPixelLorentzAngle, nIOVs, ntags>("SiPixelLorentzAngle Map") {
0956       label_ = "SiPixelLorentzAngleFullPixelMap";
0957       payloadString = "Lorentz Angle";
0958     }
0959 
0960     bool fill() override {
0961       gStyle->SetPalette(kBlueRedYellow);
0962 
0963       // trick to deal with the multi-ioved tag and two tag case at the same time
0964       auto theIOVs = PlotBase::getTag<0>().iovs;
0965       auto f_tagname = PlotBase::getTag<0>().name;
0966       std::string l_tagname = "";
0967       auto firstiov = theIOVs.front();
0968       std::tuple<cond::Time_t, cond::Hash> lastiov;
0969 
0970       // we don't support (yet) comparison with more than 2 tags
0971       assert(this->m_plotAnnotations.ntags < 3);
0972 
0973       if (this->m_plotAnnotations.ntags == 2) {
0974         auto tag2iovs = PlotBase::getTag<1>().iovs;
0975         l_tagname = PlotBase::getTag<1>().name;
0976         lastiov = tag2iovs.front();
0977       } else {
0978         lastiov = theIOVs.back();
0979       }
0980 
0981       std::shared_ptr<SiPixelLorentzAngle> last_payload = this->fetchPayload(std::get<1>(lastiov));
0982       std::shared_ptr<SiPixelLorentzAngle> first_payload = this->fetchPayload(std::get<1>(firstiov));
0983 
0984       if (first_payload.get() && last_payload.get()) {
0985         // creat the base map
0986         Phase1PixelSummaryMap fullMap("",
0987                                       fmt::sprintf("%s Diff", payloadString),
0988                                       fmt::sprintf("%s difference #Delta#mu_{H}/#mu_{H} [%%]", payloadString));
0989         fullMap.createTrackerBaseMap();
0990 
0991         std::map<uint32_t, float> l_LAMap_ = last_payload->getLorentzAngles();
0992         std::map<uint32_t, float> f_LAMap_ = first_payload->getLorentzAngles();
0993 
0994         std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
0995         std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
0996 
0997         if (l_LAMap_.size() != f_LAMap_.size()) {
0998           edm::LogError(label_) << "There are " << l_LAMap_.size() << "dets in payload" << std::get<1>(lastiov)
0999                                 << "and " << f_LAMap_.size() << " dets in payload" << std::get<1>(firstiov)
1000                                 << "display is not possible!";
1001 
1002           TCanvas canvas("Canv", "Canv", 1200, 1000);
1003           SiPixelPI::displayNotSupported(canvas, SiPixelPI::mismatched);
1004           std::string fileName(this->m_imageFileName);
1005           canvas.SaveAs(fileName.c_str());
1006           return false;
1007         }
1008 
1009         if (l_LAMap_.size() == SiPixelPI::phase0size || l_LAMap_.size() > SiPixelPI::phase1size) {
1010           edm::LogError(label_) << "There are " << l_LAMap_.size()
1011                                 << " DetIds in this payload. SiPixelLorentzAngleFullPixelMap maps are not supported "
1012                                    "for non-Phase1 Pixel geometries !";
1013           TCanvas canvas("Canv", "Canv", 1200, 1000);
1014           SiPixelPI::displayNotSupported(canvas, l_LAMap_.size());
1015           std::string fileName(this->m_imageFileName);
1016           canvas.SaveAs(fileName.c_str());
1017           return false;
1018         } else {
1019           if (l_LAMap_.size() < SiPixelPI::phase1size) {
1020             edm::LogWarning(label_) << "\n ********* WARNING! ********* \n There are " << l_LAMap_.size()
1021                                     << " DetIds in this payload !"
1022                                     << "\n **************************** \n";
1023           }
1024         }
1025 
1026         // fill the map
1027         for (const auto &[id, value] : l_LAMap_) {
1028           assert(value != 0.);  // do not divide by 0
1029           const auto &diff = (value - f_LAMap_[id]) * 100.f / value;
1030           fullMap.fillTrackerMap(id, diff);
1031         }
1032 
1033         // print the map
1034         TCanvas canvas("Canv", "Canv", 3000, 2000);
1035         fullMap.printTrackerMap(canvas, 0.03);  // is the top margin of the canvas
1036 
1037         // take care of the legend
1038         auto ltx = TLatex();
1039         ltx.SetTextFont(62);
1040         ltx.SetTextSize(0.025);
1041         ltx.SetTextAlign(11);
1042         std::string ltxText;
1043         if (this->m_plotAnnotations.ntags == 2) {
1044           ltxText = fmt::sprintf("#color[2]{%s, %s} vs #color[4]{%s, %s}",
1045                                  f_tagname,
1046                                  std::to_string(std::get<0>(firstiov)),
1047                                  l_tagname,
1048                                  std::to_string(std::get<0>(lastiov)));
1049         } else {
1050           ltxText = fmt::sprintf("%s IOV: #color[2]{%s} vs IOV: #color[4]{%s}",
1051                                  f_tagname,
1052                                  std::to_string(std::get<0>(firstiov)),
1053                                  std::to_string(std::get<0>(lastiov)));
1054         }
1055         ltx.DrawLatexNDC(gPad->GetLeftMargin(), 1 - gPad->GetTopMargin() + 0.01, ltxText.c_str());
1056 
1057         std::string fileName(this->m_imageFileName);
1058         canvas.SaveAs(fileName.c_str());
1059       }  // if has got the payloads
1060       return true;
1061     }
1062 
1063   protected:
1064     std::string payloadString;
1065     std::string label_;
1066   };
1067 
1068   using SiPixelLorentzAngleFullMapCompareSingleTag = SiPixelLorentzAngleFullMapCompareBase<MULTI_IOV, 1>;
1069   using SiPixelLorentzAngleFullMapCompareTwoTags = SiPixelLorentzAngleFullMapCompareBase<SINGLE_IOV, 2>;
1070 
1071 }  // namespace
1072 
1073 PAYLOAD_INSPECTOR_MODULE(SiPixelLorentzAngle) {
1074   PAYLOAD_INSPECTOR_CLASS(SiPixelLorentzAngleValue);
1075   PAYLOAD_INSPECTOR_CLASS(SiPixelLorentzAngleValues);
1076   PAYLOAD_INSPECTOR_CLASS(SiPixelLorentzAngleValuesBarrel);
1077   PAYLOAD_INSPECTOR_CLASS(SiPixelLorentzAngleValuesEndcap);
1078   PAYLOAD_INSPECTOR_CLASS(SiPixelLorentzAngleValuesBarrelCompareSingleTag);
1079   PAYLOAD_INSPECTOR_CLASS(SiPixelLorentzAngleValuesEndcapCompareSingleTag);
1080   PAYLOAD_INSPECTOR_CLASS(SiPixelLorentzAngleValuesBarrelCompareTwoTags);
1081   PAYLOAD_INSPECTOR_CLASS(SiPixelLorentzAngleValuesEndcapCompareTwoTags);
1082   PAYLOAD_INSPECTOR_CLASS(SiPixelLorentzAngleValueComparisonSingleTag);
1083   PAYLOAD_INSPECTOR_CLASS(SiPixelLorentzAngleValueComparisonTwoTags);
1084   PAYLOAD_INSPECTOR_CLASS(SiPixelLorentzAngleByRegionComparisonSingleTag);
1085   PAYLOAD_INSPECTOR_CLASS(SiPixelLorentzAngleByRegionComparisonTwoTags);
1086   PAYLOAD_INSPECTOR_CLASS(SiPixelBPixLorentzAngleMap);
1087   PAYLOAD_INSPECTOR_CLASS(SiPixelFPixLorentzAngleMap);
1088   PAYLOAD_INSPECTOR_CLASS(SiPixelFullLorentzAngleMapByROC);
1089   PAYLOAD_INSPECTOR_CLASS(SiPixelLorentzAngleFullPixelMap);
1090   PAYLOAD_INSPECTOR_CLASS(SiPixelLorentzAngleFullMapCompareSingleTag);
1091   PAYLOAD_INSPECTOR_CLASS(SiPixelLorentzAngleFullMapCompareTwoTags);
1092 }