Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-09-07 04:35:31

0001 #ifndef CONDCORE_SIPIXELPLUGINS_SIPIXELGAINCALIBHELPER_H
0002 #define CONDCORE_SIPIXELPLUGINS_SIPIXELGAINCALIBHELPER_H
0003 
0004 #include "CondCore/Utilities/interface/PayloadInspectorModule.h"
0005 #include "CondCore/Utilities/interface/PayloadInspector.h"
0006 #include "CondCore/CondDB/interface/Time.h"
0007 #include "FWCore/ParameterSet/interface/FileInPath.h"
0008 #include "FWCore/Utilities/interface/isFinite.h"
0009 #include "CalibTracker/StandaloneTrackerTopology/interface/StandaloneTrackerTopology.h"
0010 #include "CondCore/SiPixelPlugins/interface/SiPixelPayloadInspectorHelper.h"
0011 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
0012 #include "CondFormats/SiPixelObjects/interface/SiPixelGainCalibrationOffline.h"
0013 #include "CondFormats/SiPixelObjects/interface/SiPixelGainCalibrationForHLT.h"
0014 #include "CondCore/SiPixelPlugins/interface/PixelRegionContainers.h"
0015 #include "DQM/TrackerRemapper/interface/Phase1PixelROCMaps.h"
0016 
0017 #include <type_traits>
0018 #include <memory>
0019 #include <sstream>
0020 #include <fmt/printf.h>
0021 
0022 // include ROOT
0023 #include "TH2F.h"
0024 #include "TH1F.h"
0025 #include "TLegend.h"
0026 #include "TCanvas.h"
0027 #include "TLine.h"
0028 #include "TStyle.h"
0029 #include "TLatex.h"
0030 #include "TPave.h"
0031 #include "TPaveStats.h"
0032 #include "TGaxis.h"
0033 
0034 namespace gainCalibHelper {
0035 
0036   using AvgMap = std::map<uint32_t, float>;
0037 
0038   namespace gainCalibPI {
0039 
0040     enum type { t_gain = 0, t_pedestal = 1, t_correlation = 2 };
0041     static const std::array<std::string, 3> t_titles = {{"gain", "pedestal", "correlation"}};
0042 
0043     //===========================================================================
0044     // helper method to fill the ratio and diff distributions
0045     template <typename PayloadType>
0046     static std::array<std::shared_ptr<TH1F>, 2> fillDiffAndRatio(const std::shared_ptr<PayloadType>& payload_A,
0047                                                                  const std::shared_ptr<PayloadType>& payload_B,
0048                                                                  const gainCalibPI::type& theType) {
0049       std::array<std::shared_ptr<TH1F>, 2> arr;
0050 
0051       std::vector<uint32_t> detids_A;
0052       payload_A->getDetIds(detids_A);
0053       std::vector<uint32_t> detids_B;
0054       payload_B->getDetIds(detids_B);
0055 
0056       std::vector<float> v_ratios;
0057       std::vector<float> v_diffs;
0058 
0059       if (detids_A != detids_B) {
0060         edm::LogError("fillDiffAndRatio") << "the list of DetIds for the two payloads are not equal"
0061                                           << " cannot make any comparison!" << std::endl;
0062       }
0063 
0064       assert(detids_A == detids_B);
0065       for (const auto& d : detids_A) {
0066         auto range = payload_A->getRange(d);
0067         int numberOfRowsToAverageOver = payload_A->getNumberOfRowsToAverageOver();
0068         int ncols = payload_A->getNCols(d);
0069         int nRocsInRow = (range.second - range.first) / ncols / numberOfRowsToAverageOver;
0070         unsigned int nRowsForHLT = 1;
0071         int nrows = std::max((payload_A->getNumberOfRowsToAverageOver() * nRocsInRow),
0072                              nRowsForHLT);  // dirty trick to make it work for the HLT payload
0073 
0074         auto rAndCol_A = payload_A->getRangeAndNCols(d);
0075         auto rAndCol_B = payload_B->getRangeAndNCols(d);
0076         bool isDeadCol, isNoisyCol;
0077 
0078         float ratio(-.1), diff(-1.);
0079         for (int col = 0; col < ncols; col++) {
0080           for (int row = 0; row < nrows; row++) {
0081             switch (theType) {
0082               case gainCalibPI::t_gain: {
0083                 auto gainA = payload_A->getGain(col, row, rAndCol_A.first, rAndCol_A.second, isDeadCol, isNoisyCol);
0084                 auto gainB = payload_B->getGain(col, row, rAndCol_B.first, rAndCol_B.second, isDeadCol, isNoisyCol);
0085                 ratio = gainB != 0 ? (gainA / gainB) : -1.;  // if the ratio does not make sense the default is -1
0086                 diff = gainA - gainB;
0087                 break;
0088               }
0089               case gainCalibPI::t_pedestal: {
0090                 auto pedA = payload_A->getPed(col, row, rAndCol_A.first, rAndCol_A.second, isDeadCol, isNoisyCol);
0091                 auto pedB = payload_B->getPed(col, row, rAndCol_B.first, rAndCol_B.second, isDeadCol, isNoisyCol);
0092                 ratio = pedB != 0 ? (pedA / pedB) : -1.;  // if the ratio does not make sense the default is -1
0093                 diff = pedA - pedB;
0094                 break;
0095               }
0096               default:
0097                 edm::LogError("gainCalibPI::fillTheHisto") << "Unrecognized type " << theType << std::endl;
0098                 break;
0099             }
0100             // fill the containers
0101             v_diffs.push_back(diff);
0102             v_ratios.push_back(ratio);
0103           }  // loop on rows
0104         }  // loop on cols
0105       }  // loop on detids
0106 
0107       std::sort(v_diffs.begin(), v_diffs.end());
0108       std::sort(v_ratios.begin(), v_ratios.end());
0109 
0110       double minDiff = v_diffs.front();
0111       double maxDiff = v_diffs.back();
0112       double minRatio = v_ratios.front();
0113       double maxRatio = v_ratios.back();
0114 
0115       arr[0] = std::make_shared<TH1F>("h_Ratio", "", 50, minRatio - 1., maxRatio + 1.);
0116       arr[1] = std::make_shared<TH1F>("h_Diff", "", 50, minDiff - 1., maxDiff + 1.);
0117 
0118       for (const auto& r : v_ratios) {
0119         arr[0]->Fill(r);
0120       }
0121       for (const auto& d : v_diffs) {
0122         arr[1]->Fill(d);
0123       }
0124       return arr;
0125     }
0126 
0127     //============================================================================
0128     // helper method to fill the gain / pedestals distributions
0129     template <typename PayloadType>
0130     static void fillTheHisto(const std::shared_ptr<PayloadType>& payload,
0131                              std::shared_ptr<TH1F> h1,
0132                              gainCalibPI::type theType,
0133                              const std::vector<uint32_t>& wantedIds = {}) {
0134       std::vector<uint32_t> detids;
0135       if (wantedIds.empty()) {
0136         payload->getDetIds(detids);
0137       } else {
0138         detids.assign(wantedIds.begin(), wantedIds.end());
0139       }
0140 
0141       for (const auto& d : detids) {
0142         // skip the special case used to signal there are no attached dets
0143         if (d == 0xFFFFFFFF)
0144           continue;
0145 
0146         auto range = payload->getRange(d);
0147         int numberOfRowsToAverageOver = payload->getNumberOfRowsToAverageOver();
0148         int ncols = payload->getNCols(d);
0149         int nRocsInRow = (range.second - range.first) / ncols / numberOfRowsToAverageOver;
0150         unsigned int nRowsForHLT = 1;
0151         int nrows = std::max((payload->getNumberOfRowsToAverageOver() * nRocsInRow),
0152                              nRowsForHLT);  // dirty trick to make it work for the HLT payload
0153 
0154         auto rangeAndCol = payload->getRangeAndNCols(d);
0155         bool isDeadColumn;
0156         bool isNoisyColumn;
0157 
0158         COUT << "NCOLS: " << payload->getNCols(d) << " " << rangeAndCol.second << " NROWS:" << nrows
0159              << ", RANGES: " << rangeAndCol.first.second - rangeAndCol.first.first
0160              << ", Ratio: " << float(rangeAndCol.first.second - rangeAndCol.first.first) / rangeAndCol.second
0161              << std::endl;
0162 
0163         float quid(-99999.);
0164 
0165         for (int col = 0; col < ncols; col++) {
0166           for (int row = 0; row < nrows; row++) {
0167             switch (theType) {
0168               case gainCalibPI::t_gain:
0169                 quid = payload->getGain(col, row, rangeAndCol.first, rangeAndCol.second, isDeadColumn, isNoisyColumn);
0170                 break;
0171               case gainCalibPI::t_pedestal:
0172                 quid = payload->getPed(col, row, rangeAndCol.first, rangeAndCol.second, isDeadColumn, isNoisyColumn);
0173                 break;
0174               default:
0175                 edm::LogError("gainCalibPI::fillTheHisto") << "Unrecognized type " << theType << std::endl;
0176                 break;
0177             }
0178             h1->Fill(quid);
0179           }  // loop on rows
0180         }  // loop on cols
0181       }  // loop on detids
0182     }  // fillTheHisto
0183 
0184     //============================================================================
0185     // helper method to fill the gain / pedestal averages per module maps
0186     template <typename PayloadType>
0187     static void fillThePerModuleMap(const std::shared_ptr<PayloadType>& payload,
0188                                     AvgMap& map,
0189                                     gainCalibPI::type theType) {
0190       std::vector<uint32_t> detids;
0191       payload->getDetIds(detids);
0192 
0193       for (const auto& d : detids) {
0194         auto range = payload->getRange(d);
0195         int numberOfRowsToAverageOver = payload->getNumberOfRowsToAverageOver();
0196         int ncols = payload->getNCols(d);
0197         int nRocsInRow = (range.second - range.first) / ncols / numberOfRowsToAverageOver;
0198         unsigned int nRowsForHLT = 1;
0199         int nrows = std::max((payload->getNumberOfRowsToAverageOver() * nRocsInRow),
0200                              nRowsForHLT);  // dirty trick to make it work for the HLT payload
0201 
0202         auto rangeAndCol = payload->getRangeAndNCols(d);
0203         bool isDeadColumn;
0204         bool isNoisyColumn;
0205 
0206         float sumOfX(0.);
0207         int nPixels(0);
0208         for (int col = 0; col < ncols; col++) {
0209           for (int row = 0; row < nrows; row++) {
0210             nPixels++;
0211             switch (theType) {
0212               case gainCalibPI::t_gain:
0213                 sumOfX +=
0214                     payload->getGain(col, row, rangeAndCol.first, rangeAndCol.second, isDeadColumn, isNoisyColumn);
0215                 break;
0216               case gainCalibPI::t_pedestal:
0217                 sumOfX += payload->getPed(col, row, rangeAndCol.first, rangeAndCol.second, isDeadColumn, isNoisyColumn);
0218                 break;
0219               default:
0220                 edm::LogError("gainCalibPI::fillThePerModuleMap") << "Unrecognized type " << theType << std::endl;
0221                 break;
0222             }  // switch on the type
0223           }  // rows
0224         }  // columns
0225         // fill the return value map
0226         map[d] = sumOfX / nPixels;
0227       }  // loop on the detId
0228     }  // fillThePerModuleMap
0229 
0230     //============================================================================
0231     // helper method to fill the gain / pedestals distributions
0232     template <typename PayloadType>
0233     static void fillTheHistos(const std::shared_ptr<PayloadType>& payload,
0234                               std::shared_ptr<TH1> hBPix,
0235                               std::shared_ptr<TH1> hFPix,
0236                               gainCalibPI::type theType) {
0237       std::vector<uint32_t> detids;
0238       payload->getDetIds(detids);
0239 
0240       bool isCorrelation_ = hBPix.get()->InheritsFrom(TH2::Class()) && (theType == gainCalibPI::t_correlation);
0241 
0242       for (const auto& d : detids) {
0243         int subid = DetId(d).subdetId();
0244         auto range = payload->getRange(d);
0245         int numberOfRowsToAverageOver = payload->getNumberOfRowsToAverageOver();
0246         int ncols = payload->getNCols(d);
0247         int nRocsInRow = (range.second - range.first) / ncols / numberOfRowsToAverageOver;
0248         unsigned int nRowsForHLT = 1;
0249         int nrows = std::max((payload->getNumberOfRowsToAverageOver() * nRocsInRow),
0250                              nRowsForHLT);  // dirty trick to make it work for the HLT payload
0251 
0252         auto rangeAndCol = payload->getRangeAndNCols(d);
0253         bool isDeadColumn;
0254         bool isNoisyColumn;
0255 
0256         for (int col = 0; col < ncols; col++) {
0257           for (int row = 0; row < nrows; row++) {
0258             float gain = payload->getGain(col, row, rangeAndCol.first, rangeAndCol.second, isDeadColumn, isNoisyColumn);
0259             float ped = payload->getPed(col, row, rangeAndCol.first, rangeAndCol.second, isDeadColumn, isNoisyColumn);
0260 
0261             switch (subid) {
0262               case PixelSubdetector::PixelBarrel: {
0263                 if (isCorrelation_) {
0264                   hBPix->Fill(gain, ped);
0265                 } else {
0266                   hBPix->Fill((theType == gainCalibPI::t_gain ? gain : ped));
0267                 }
0268                 break;
0269               }
0270               case PixelSubdetector::PixelEndcap: {
0271                 if (isCorrelation_) {
0272                   hFPix->Fill(gain, ped);
0273                 } else {
0274                   hFPix->Fill((theType == gainCalibPI::t_gain ? gain : ped));
0275                 }
0276                 break;
0277               }
0278               default:
0279                 edm::LogError("gainCalibPI::fillTheHistos") << d << " is not a Pixel DetId" << std::endl;
0280                 break;
0281             }  // switch on subid
0282           }  // row loop
0283         }  // column loop
0284       }  // detid loop
0285     }  // filltheHistos
0286   }  // namespace gainCalibPI
0287 
0288   constexpr char const* TypeName[2] = {"Gains", "Pedestals"};
0289 
0290   /*******************************************************************
0291     1d histogram of SiPixelGainCalibration for Gains of 1 IOV 
0292   ********************************************************************/
0293   template <gainCalibPI::type myType, class PayloadType>
0294   class SiPixelGainCalibrationValues
0295       : public cond::payloadInspector::PlotImage<PayloadType, cond::payloadInspector::SINGLE_IOV> {
0296   public:
0297     SiPixelGainCalibrationValues()
0298         : cond::payloadInspector::PlotImage<PayloadType, cond::payloadInspector::SINGLE_IOV>(
0299               Form("SiPixelGainCalibration %s Values", TypeName[myType])) {
0300       if constexpr (std::is_same_v<PayloadType, SiPixelGainCalibrationOffline>) {
0301         isForHLT_ = false;
0302         label_ = "SiPixelGainCalibrationOffline_PayloadInspector";
0303       } else {
0304         isForHLT_ = true;
0305         label_ = "SiPixelGainCalibrationForHLT_PayloadInspector";
0306       }
0307     }
0308 
0309     bool fill() override {
0310       auto tag = cond::payloadInspector::PlotBase::getTag<0>();
0311       auto iov = tag.iovs.front();
0312 
0313       std::shared_ptr<PayloadType> payload = this->fetchPayload(std::get<1>(iov));
0314 
0315       gStyle->SetOptStat("emr");
0316 
0317       float minimum(9999.);
0318       float maximum(-9999.);
0319 
0320       switch (myType) {
0321         case gainCalibPI::t_gain:
0322           maximum = payload->getGainHigh();
0323           minimum = payload->getGainLow();
0324           break;
0325         case gainCalibPI::t_pedestal:
0326           maximum = payload->getPedHigh();
0327           minimum = payload->getPedLow();
0328           break;
0329         default:
0330           edm::LogError(label_) << "Unrecognized type " << myType << std::endl;
0331           break;
0332       }
0333 
0334       TCanvas canvas("Canv", "Canv", 1200, 1000);
0335       auto h1 = std::make_shared<TH1F>(Form("%s values", TypeName[myType]),
0336                                        Form("SiPixel Gain Calibration %s - %s;per %s %s;# %ss",
0337                                             (isForHLT_ ? "ForHLT" : "Offline"),
0338                                             TypeName[myType],
0339                                             (isForHLT_ ? "Column" : "Pixel"),
0340                                             TypeName[myType],
0341                                             (isForHLT_ ? "column" : "pixel")),
0342                                        200,
0343                                        minimum,
0344                                        maximum);
0345       canvas.cd();
0346       SiPixelPI::adjustCanvasMargins(canvas.cd(), 0.06, 0.12, 0.12, 0.05);
0347       canvas.Modified();
0348 
0349       // fill the histogram
0350       gainCalibPI::fillTheHisto(payload, h1, myType);
0351 
0352       canvas.cd()->SetLogy();
0353       h1->SetTitle("");
0354       h1->GetYaxis()->SetRangeUser(0.1, h1->GetMaximum() * 10.);
0355       h1->SetFillColor(kBlue);
0356       h1->SetMarkerStyle(20);
0357       h1->SetMarkerSize(1);
0358       h1->Draw("bar2");
0359 
0360       SiPixelPI::makeNicePlotStyle(h1.get());
0361       h1->SetStats(true);
0362 
0363       canvas.Update();
0364 
0365       TLegend legend = TLegend(0.40, 0.88, 0.94, 0.93);
0366       legend.SetHeader(("Payload hash: #bf{" + (std::get<1>(iov)) + "}").c_str(),
0367                        "C");  // option "C" allows to center the header
0368       //legend.AddEntry(h1.get(), ("IOV: " + std::to_string(std::get<0>(iov))).c_str(), "PL");
0369       legend.SetLineColor(10);
0370       legend.SetTextSize(0.025);
0371       legend.Draw("same");
0372 
0373       TPaveStats* st = (TPaveStats*)h1->FindObject("stats");
0374       st->SetTextSize(0.03);
0375       SiPixelPI::adjustStats(st, 0.15, 0.83, 0.39, 0.93);
0376 
0377       auto ltx = TLatex();
0378       ltx.SetTextFont(62);
0379       //ltx.SetTextColor(kBlue);
0380       ltx.SetTextSize(0.05);
0381       ltx.SetTextAlign(11);
0382       ltx.DrawLatexNDC(gPad->GetLeftMargin() + 0.1,
0383                        1 - gPad->GetTopMargin() + 0.01,
0384                        ("SiPixel Gain Calibration IOV:" + std::to_string(std::get<0>(iov))).c_str());
0385 
0386       std::string fileName(this->m_imageFileName);
0387       canvas.SaveAs(fileName.c_str());
0388 
0389       return true;
0390     }
0391 
0392   protected:
0393     bool isForHLT_;
0394     std::string label_;
0395   };
0396 
0397   /*******************************************************************
0398     1d histograms per region of SiPixelGainCalibration for Gains of 1 IOV
0399   ********************************************************************/
0400   template <bool isBarrel, gainCalibPI::type myType, class PayloadType>
0401   class SiPixelGainCalibrationValuesPerRegion
0402       : public cond::payloadInspector::PlotImage<PayloadType, cond::payloadInspector::SINGLE_IOV> {
0403   public:
0404     SiPixelGainCalibrationValuesPerRegion()
0405         : cond::payloadInspector::PlotImage<PayloadType, cond::payloadInspector::SINGLE_IOV>(
0406               Form("SiPixelGainCalibration %s Values Per Region", TypeName[myType])) {
0407       cond::payloadInspector::PlotBase::addInputParam("SetLog");
0408 
0409       if constexpr (std::is_same_v<PayloadType, SiPixelGainCalibrationOffline>) {
0410         isForHLT_ = false;
0411         label_ = "SiPixelGainCalibrationOffline_PayloadInspector";
0412       } else {
0413         isForHLT_ = true;
0414         label_ = "SiPixelGainCalibrationForHLT_PayloadInspector";
0415       }
0416     }
0417 
0418     bool fill() override {
0419       gStyle->SetOptStat("mr");
0420 
0421       auto tag = cond::payloadInspector::PlotBase::getTag<0>();
0422       auto iov = tag.iovs.front();
0423 
0424       // parse first if log
0425       bool setLog(true);
0426       auto paramValues = cond::payloadInspector::PlotBase::inputParamValues();
0427       auto ip = paramValues.find("SetLog");
0428       if (ip != paramValues.end()) {
0429         auto answer = ip->second;
0430         if (!SiPixelPI::checkAnswerOK(answer, setLog)) {
0431           throw cms::Exception(label_)
0432               << "\nERROR: " << answer
0433               << " is not a valid setting for this parameter, please use True,False,1,0,Yes,No \n\n";
0434         }
0435       }
0436 
0437       std::shared_ptr<PayloadType> payload = this->fetchPayload(std::get<1>(iov));
0438 
0439       std::vector<uint32_t> detids;
0440       payload->getDetIds(detids);
0441 
0442       float minimum(9999.);
0443       float maximum(-9999.);
0444 
0445       switch (myType) {
0446         case gainCalibPI::t_gain:
0447           maximum = payload->getGainHigh();
0448           minimum = payload->getGainLow();
0449           break;
0450         case gainCalibPI::t_pedestal:
0451           maximum = payload->getPedHigh();
0452           minimum = payload->getPedLow();
0453           break;
0454         default:
0455           edm::LogError(label_) << "Unrecognized type " << myType << std::endl;
0456           break;
0457       }
0458 
0459       TCanvas canvas("Canv", "Canv", isBarrel ? 1400 : 1800, 1200);
0460       if (detids.size() > SiPixelPI::phase1size) {
0461         SiPixelPI::displayNotSupported(canvas, detids.size());
0462         std::string fileName(this->m_imageFileName);
0463         canvas.SaveAs(fileName.c_str());
0464         return false;
0465       }
0466 
0467       canvas.cd();
0468 
0469       SiPixelPI::PhaseInfo phaseInfo(detids.size());
0470       const char* path_toTopologyXML = phaseInfo.pathToTopoXML();
0471 
0472       TrackerTopology tTopo =
0473           StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
0474 
0475       auto myPlots = PixelRegions::PixelRegionContainers(&tTopo, phaseInfo.phase());
0476       myPlots.bookAll(Form("SiPixel Gain Calibration %s - %s", (isForHLT_ ? "ForHLT" : "Offline"), TypeName[myType]),
0477                       Form("per %s %s", (isForHLT_ ? "Column" : "Pixel"), TypeName[myType]),
0478                       Form("# %ss", (isForHLT_ ? "column" : "pixel")),
0479                       200,
0480                       minimum,
0481                       maximum);
0482 
0483       canvas.Modified();
0484 
0485       // fill the histograms
0486       for (const auto& pixelId : PixelRegions::PixelIDs) {
0487         auto wantedDets = PixelRegions::attachedDets(pixelId, &tTopo, phaseInfo.phase());
0488         gainCalibPI::fillTheHisto(payload, myPlots.getHistoFromMap(pixelId), myType, wantedDets);
0489       }
0490 
0491       if (setLog) {
0492         myPlots.setLogScale();
0493       }
0494       myPlots.beautify(kBlue, -1);
0495       myPlots.draw(canvas, isBarrel, "HIST");
0496 
0497       TLegend legend = TLegend(0.45, 0.88, 0.91, 0.92);
0498       legend.SetHeader(("hash: #bf{" + (std::get<1>(iov)) + "}").c_str(),
0499                        "C");  // option "C" allows to center the header
0500       //legend.AddEntry(h1.get(), ("IOV: " + std::to_string(std::get<0>(iov))).c_str(), "PL");
0501       legend.SetLineColor(10);
0502       legend.SetTextSize(0.025);
0503       legend.Draw("same");
0504 
0505       unsigned int maxPads = isBarrel ? 4 : 12;
0506       for (unsigned int c = 1; c <= maxPads; c++) {
0507         canvas.cd(c);
0508         SiPixelPI::adjustCanvasMargins(canvas.cd(c), 0.06, 0.12, 0.12, 0.05);
0509         legend.Draw("same");
0510         canvas.cd(c)->Update();
0511       }
0512 
0513       myPlots.stats();
0514 
0515       auto ltx = TLatex();
0516       ltx.SetTextFont(62);
0517       ltx.SetTextSize(0.05);
0518       ltx.SetTextAlign(11);
0519 
0520       for (unsigned int c = 1; c <= maxPads; c++) {
0521         auto index = isBarrel ? c - 1 : c + 3;
0522         canvas.cd(c);
0523         auto leftX = setLog ? 0. : 0.1;
0524         ltx.DrawLatexNDC(gPad->GetLeftMargin() + leftX,
0525                          1 - gPad->GetTopMargin() + 0.01,
0526                          (PixelRegions::IDlabels.at(index) + ", IOV:" + std::to_string(std::get<0>(iov))).c_str());
0527       }
0528 
0529       std::string fileName(this->m_imageFileName);
0530       canvas.SaveAs(fileName.c_str());
0531 
0532       return true;
0533     }
0534 
0535   protected:
0536     bool isForHLT_;
0537     std::string label_;
0538   };
0539 
0540   /*******************************************************************
0541     1d histograms comparison per region of SiPixelGainCalibration for Gains of 2 IOV
0542   ********************************************************************/
0543   template <bool isBarrel,
0544             gainCalibPI::type myType,
0545             cond::payloadInspector::IOVMultiplicity nIOVs,
0546             int ntags,
0547             class PayloadType>
0548   class SiPixelGainCalibrationValuesComparisonPerRegion
0549       : public cond::payloadInspector::PlotImage<PayloadType, nIOVs, ntags> {
0550   public:
0551     SiPixelGainCalibrationValuesComparisonPerRegion()
0552         : cond::payloadInspector::PlotImage<PayloadType, nIOVs, ntags>(
0553               Form("SiPixelGainCalibration %s Values Per Region %i tag(s)", TypeName[myType], ntags)) {
0554       cond::payloadInspector::PlotBase::addInputParam("SetLog");
0555 
0556       if constexpr (std::is_same_v<PayloadType, SiPixelGainCalibrationOffline>) {
0557         isForHLT_ = false;
0558         label_ = "SiPixelGainCalibrationOffline_PayloadInspector";
0559       } else {
0560         isForHLT_ = true;
0561         label_ = "SiPixelGainCalibrationForHLT_PayloadInspector";
0562       }
0563     }
0564 
0565     bool fill() override {
0566       gStyle->SetOptStat("mr");
0567 
0568       COUT << "ntags: " << ntags << " this->m_plotAnnotations.ntags: " << this->m_plotAnnotations.ntags << std::endl;
0569 
0570       // trick to deal with the multi-ioved tag and two tag case at the same time
0571       auto theIOVs = cond::payloadInspector::PlotBase::getTag<0>().iovs;
0572       auto f_tagname = cond::payloadInspector::PlotBase::getTag<0>().name;
0573       std::string l_tagname = "";
0574       auto firstiov = theIOVs.front();
0575       std::tuple<cond::Time_t, cond::Hash> lastiov;
0576 
0577       // we don't support (yet) comparison with more than 2 tags
0578       assert(this->m_plotAnnotations.ntags < 3);
0579 
0580       if (this->m_plotAnnotations.ntags == 2) {
0581         auto tag2iovs = cond::payloadInspector::PlotBase::getTag<1>().iovs;
0582         l_tagname = cond::payloadInspector::PlotBase::getTag<1>().name;
0583         lastiov = tag2iovs.front();
0584       } else {
0585         lastiov = theIOVs.back();
0586       }
0587 
0588       // parse first if log
0589       bool setLog(true);
0590       auto paramValues = cond::payloadInspector::PlotBase::inputParamValues();
0591       auto ip = paramValues.find("SetLog");
0592       if (ip != paramValues.end()) {
0593         auto answer = ip->second;
0594         if (!SiPixelPI::checkAnswerOK(answer, setLog)) {
0595           throw cms::Exception(label_)
0596               << "\nERROR: " << answer
0597               << " is not a valid setting for this parameter, please use True,False,1,0,Yes,No \n\n";
0598         }
0599       }
0600 
0601       std::shared_ptr<PayloadType> last_payload = this->fetchPayload(std::get<1>(lastiov));
0602       std::shared_ptr<PayloadType> first_payload = this->fetchPayload(std::get<1>(firstiov));
0603 
0604       std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
0605       std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
0606 
0607       std::vector<uint32_t> f_detids, l_detids;
0608       last_payload->getDetIds(l_detids);
0609       first_payload->getDetIds(f_detids);
0610 
0611       float minimum(9999.);
0612       float maximum(-9999.);
0613 
0614       switch (myType) {
0615         case gainCalibPI::t_gain:
0616           maximum = std::max(last_payload->getGainHigh(), first_payload->getGainHigh());
0617           minimum = std::min(last_payload->getGainLow(), first_payload->getGainLow());
0618           break;
0619         case gainCalibPI::t_pedestal:
0620           maximum = std::max(last_payload->getPedHigh(), first_payload->getPedHigh());
0621           minimum = std::min(last_payload->getPedLow(), first_payload->getPedLow());
0622           break;
0623         default:
0624           edm::LogError(label_) << "Unrecognized type " << myType << std::endl;
0625           break;
0626       }
0627 
0628       TCanvas canvas("Canv", "Canv", isBarrel ? 1400 : 1800, 1200);
0629       if (std::max(l_detids.size(), f_detids.size()) > SiPixelPI::phase1size) {
0630         SiPixelPI::displayNotSupported(canvas, std::max(f_detids.size(), l_detids.size()));
0631         std::string fileName(this->m_imageFileName);
0632         canvas.SaveAs(fileName.c_str());
0633         return false;
0634       }
0635 
0636       canvas.cd();
0637 
0638       SiPixelPI::PhaseInfo l_phaseInfo(l_detids.size());
0639       SiPixelPI::PhaseInfo f_phaseInfo(f_detids.size());
0640       const char* path_toTopologyXML = l_phaseInfo.pathToTopoXML();
0641 
0642       auto l_tTopo =
0643           StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
0644 
0645       auto l_myPlots = PixelRegions::PixelRegionContainers(&l_tTopo, l_phaseInfo.phase());
0646       l_myPlots.bookAll(
0647           Form("Last SiPixel Gain Calibration %s - %s", (isForHLT_ ? "ForHLT" : "Offline"), TypeName[myType]),
0648           Form("per %s %s", (isForHLT_ ? "Column" : "Pixel"), TypeName[myType]),
0649           Form("# %ss", (isForHLT_ ? "column" : "pixel")),
0650           200,
0651           minimum,
0652           maximum);
0653 
0654       path_toTopologyXML = f_phaseInfo.pathToTopoXML();
0655 
0656       auto f_tTopo =
0657           StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
0658 
0659       auto f_myPlots = PixelRegions::PixelRegionContainers(&f_tTopo, f_phaseInfo.phase());
0660       f_myPlots.bookAll(
0661           Form("First SiPixel Gain Calibration %s - %s", (isForHLT_ ? "ForHLT" : "Offline"), TypeName[myType]),
0662           Form("per %s %s", (isForHLT_ ? "Column" : "Pixel"), TypeName[myType]),
0663           Form("# %ss", (isForHLT_ ? "column" : "pixel")),
0664           200,
0665           minimum,
0666           maximum);
0667 
0668       // fill the histograms
0669       for (const auto& pixelId : PixelRegions::PixelIDs) {
0670         auto f_wantedDets = PixelRegions::attachedDets(pixelId, &f_tTopo, f_phaseInfo.phase());
0671         auto l_wantedDets = PixelRegions::attachedDets(pixelId, &l_tTopo, l_phaseInfo.phase());
0672         gainCalibPI::fillTheHisto(first_payload, f_myPlots.getHistoFromMap(pixelId), myType, f_wantedDets);
0673         gainCalibPI::fillTheHisto(last_payload, l_myPlots.getHistoFromMap(pixelId), myType, l_wantedDets);
0674       }
0675 
0676       if (setLog) {
0677         f_myPlots.setLogScale();
0678         l_myPlots.setLogScale();
0679       }
0680 
0681       l_myPlots.beautify(kRed, -1);
0682       f_myPlots.beautify(kAzure, -1);
0683 
0684       l_myPlots.draw(canvas, isBarrel, "HIST", f_phaseInfo.isPhase1Comparison(l_phaseInfo));
0685       f_myPlots.draw(canvas, isBarrel, "HISTsames", f_phaseInfo.isPhase1Comparison(l_phaseInfo));
0686 
0687       // rescale the y-axis ranges in order to fit the canvas
0688       l_myPlots.rescaleMax(f_myPlots);
0689 
0690       // done dealing with IOVs
0691       auto colorTag = isBarrel ? PixelRegions::L1 : PixelRegions::Rm1l;
0692       std::unique_ptr<TLegend> legend;
0693       if (this->m_plotAnnotations.ntags == 2) {
0694         legend = std::make_unique<TLegend>(0.36, 0.86, 0.94, 0.92);
0695         legend->AddEntry(l_myPlots.getHistoFromMap(colorTag).get(), ("#color[2]{" + l_tagname + "}").c_str(), "F");
0696         legend->AddEntry(f_myPlots.getHistoFromMap(colorTag).get(), ("#color[4]{" + f_tagname + "}").c_str(), "F");
0697         legend->SetTextSize(0.024);
0698       } else {
0699         legend = std::make_unique<TLegend>(0.58, 0.80, 0.90, 0.92);
0700         legend->AddEntry(l_myPlots.getHistoFromMap(colorTag).get(), ("#color[2]{" + lastIOVsince + "}").c_str(), "F");
0701         legend->AddEntry(f_myPlots.getHistoFromMap(colorTag).get(), ("#color[4]{" + firstIOVsince + "}").c_str(), "F");
0702         legend->SetTextSize(0.040);
0703       }
0704       legend->SetLineColor(10);
0705 
0706       unsigned int maxPads = isBarrel ? 4 : 12;
0707       for (unsigned int c = 1; c <= maxPads; c++) {
0708         canvas.cd(c);
0709         SiPixelPI::adjustCanvasMargins(canvas.cd(c), 0.06, 0.12, 0.12, 0.05);
0710         legend->Draw("same");
0711         canvas.cd(c)->Update();
0712       }
0713 
0714       f_myPlots.stats(0);
0715       l_myPlots.stats(1);
0716 
0717       auto ltx = TLatex();
0718       ltx.SetTextFont(62);
0719       ltx.SetTextSize(0.05);
0720       ltx.SetTextAlign(11);
0721 
0722       for (unsigned int c = 1; c <= maxPads; c++) {
0723         auto index = isBarrel ? c - 1 : c + 3;
0724         canvas.cd(c);
0725         auto leftX = setLog ? 0. : 0.1;
0726         ltx.DrawLatexNDC(gPad->GetLeftMargin() + leftX,
0727                          1 - gPad->GetTopMargin() + 0.01,
0728                          (PixelRegions::IDlabels.at(index) + " : #color[4]{" + std::to_string(std::get<0>(firstiov)) +
0729                           "} vs #color[2]{" + std::to_string(std::get<0>(lastiov)) + "}")
0730                              .c_str());
0731       }
0732 
0733       std::string fileName(this->m_imageFileName);
0734       canvas.SaveAs(fileName.c_str());
0735 
0736       return true;
0737     }
0738 
0739   protected:
0740     bool isForHLT_;
0741     std::string label_;
0742   };
0743 
0744   /*******************************************************************
0745     1d histogram of SiPixelGainCalibration for Gain/Pedestals
0746     correlation of 1 IOV
0747   ********************************************************************/
0748   template <class PayloadType>
0749   class SiPixelGainCalibrationCorrelations
0750       : public cond::payloadInspector::PlotImage<PayloadType, cond::payloadInspector::SINGLE_IOV> {
0751   public:
0752     SiPixelGainCalibrationCorrelations()
0753         : cond::payloadInspector::PlotImage<PayloadType, cond::payloadInspector::SINGLE_IOV>(
0754               "SiPixelGainCalibration gain/pedestal correlations") {
0755       if constexpr (std::is_same_v<PayloadType, SiPixelGainCalibrationOffline>) {
0756         isForHLT_ = false;
0757         label_ = "SiPixelGainCalibrationOffline_PayloadInspector";
0758       } else {
0759         isForHLT_ = true;
0760         label_ = "SiPixelGainCalibrationForHLT_PayloadInspector";
0761       }
0762     }
0763 
0764     bool fill() override {
0765       auto tag = cond::payloadInspector::PlotBase::getTag<0>();
0766       auto iov = tag.iovs.front();
0767 
0768       gStyle->SetOptStat("emr");
0769       gStyle->SetPalette(1);
0770 
0771       std::shared_ptr<PayloadType> payload = this->fetchPayload(std::get<1>(iov));
0772 
0773       TCanvas canvas("Canv", "Canv", 1400, 800);
0774       canvas.Divide(2, 1);
0775       canvas.cd();
0776 
0777       auto hBPix = std::make_shared<TH2F>("Correlation BPIX",
0778                                           Form("SiPixel Gain Calibration %s BPIx;per %s gains;per %s pedestals",
0779                                                (isForHLT_ ? "ForHLT" : "Offline"),
0780                                                (isForHLT_ ? "column" : "pixel"),
0781                                                (isForHLT_ ? "column" : "pixel")),
0782                                           200,
0783                                           payload->getGainLow(),
0784                                           payload->getGainHigh(),
0785                                           200,
0786                                           payload->getPedLow(),
0787                                           payload->getPedHigh());
0788 
0789       auto hFPix = std::make_shared<TH2F>("Correlation FPIX",
0790                                           Form("SiPixel Gain Calibration %s FPix;per %s gains;per %s pedestals",
0791                                                (isForHLT_ ? "ForHLT" : "Offline"),
0792                                                (isForHLT_ ? "column" : "pixel"),
0793                                                (isForHLT_ ? "column" : "pixel")),
0794                                           200,
0795                                           payload->getGainLow(),
0796                                           payload->getGainHigh(),
0797                                           200,
0798                                           payload->getPedLow(),
0799                                           payload->getPedHigh());
0800 
0801       for (unsigned int i : {1, 2}) {
0802         SiPixelPI::adjustCanvasMargins(canvas.cd(i), 0.04, 0.12, 0.15, 0.13);
0803         canvas.cd(i)->Modified();
0804       }
0805 
0806       // actually fill the histograms
0807       fillTheHistos(payload, hBPix, hFPix, gainCalibPI::t_correlation);
0808 
0809       canvas.cd(1)->SetLogz();
0810       hBPix->SetTitle("");
0811       hBPix->Draw("colz");
0812 
0813       SiPixelPI::makeNicePlotStyle(hBPix.get());
0814       hBPix->GetYaxis()->SetTitleOffset(1.65);
0815 
0816       canvas.cd(2)->SetLogz();
0817       hFPix->SetTitle("");
0818       hFPix->Draw("colz");
0819 
0820       SiPixelPI::makeNicePlotStyle(hFPix.get());
0821       hFPix->GetYaxis()->SetTitleOffset(1.65);
0822       canvas.Update();
0823 
0824       TLegend legend = TLegend(0.3, 0.92, 0.70, 0.95);
0825       legend.SetHeader(("Payload hash: #bf{" + (std::get<1>(iov)) + "}").c_str(),
0826                        "C");  // option "C" allows to center the header
0827       legend.SetLineColor(10);
0828       legend.SetTextSize(0.025);
0829       canvas.cd(1);
0830       legend.Draw("same");
0831       canvas.cd(2);
0832       legend.Draw("same");
0833 
0834       auto ltx = TLatex();
0835       ltx.SetTextFont(62);
0836       //ltx.SetTextColor(kBlue);
0837       ltx.SetTextSize(0.045);
0838       ltx.SetTextAlign(11);
0839       canvas.cd(1);
0840       ltx.DrawLatexNDC(gPad->GetLeftMargin() + 0.01,
0841                        1 - gPad->GetTopMargin() + 0.01,
0842                        ("SiPixel Gain Calibration IOV:" + std::to_string(std::get<0>(iov))).c_str());
0843 
0844       ltx.DrawLatexNDC(0.75, 0.15, "BPIX");
0845 
0846       canvas.cd(2);
0847       ltx.DrawLatexNDC(gPad->GetLeftMargin() + 0.01,
0848                        1 - gPad->GetTopMargin() + 0.01,
0849                        ("SiPixel Gain Calibration IOV:" + std::to_string(std::get<0>(iov))).c_str());
0850 
0851       ltx.DrawLatexNDC(0.75, 0.15, "FPIX");
0852 
0853       std::string fileName(this->m_imageFileName);
0854       canvas.SaveAs(fileName.c_str());
0855 #ifdef MMDEBUG
0856       canvas.SaveAs("out.root");
0857 #endif
0858       return true;
0859     }
0860 
0861   protected:
0862     bool isForHLT_;
0863     std::string label_;
0864   };
0865 
0866   /*******************************************************************
0867     1d histogram of SiPixelGainCalibration for Pedestals of 1 IOV
0868   ********************************************************************/
0869   template <gainCalibPI::type myType, class PayloadType>
0870   class SiPixelGainCalibrationValuesByPart
0871       : public cond::payloadInspector::PlotImage<PayloadType, cond::payloadInspector::SINGLE_IOV> {
0872   public:
0873     SiPixelGainCalibrationValuesByPart()
0874         : cond::payloadInspector::PlotImage<PayloadType, cond::payloadInspector::SINGLE_IOV>(
0875               Form("SiPixelGainCalibrationOffline %s Values By Partition", TypeName[myType])) {
0876       if constexpr (std::is_same_v<PayloadType, SiPixelGainCalibrationOffline>) {
0877         isForHLT_ = false;
0878         label_ = "SiPixelGainCalibrationOffline_PayloadInspector";
0879       } else {
0880         isForHLT_ = true;
0881         label_ = "SiPixelGainCalibrationForHLT_PayloadInspector";
0882       }
0883     }
0884 
0885     bool fill() override {
0886       auto tag = cond::payloadInspector::PlotBase::getTag<0>();
0887       auto iov = tag.iovs.front();
0888 
0889       gStyle->SetOptStat("emr");
0890 
0891       std::shared_ptr<PayloadType> payload = this->fetchPayload(std::get<1>(iov));
0892 
0893       TCanvas canvas("Canv", "Canv", 1400, 800);
0894       canvas.Divide(2, 1);
0895       canvas.cd();
0896 
0897       float minimum(9999.);
0898       float maximum(-9999.);
0899 
0900       switch (myType) {
0901         case gainCalibPI::t_gain:
0902           maximum = payload->getGainHigh();
0903           minimum = payload->getGainLow();
0904           break;
0905         case gainCalibPI::t_pedestal:
0906           maximum = payload->getPedHigh();
0907           minimum = payload->getPedLow();
0908           break;
0909         default:
0910           edm::LogError(label_) << "Unrecognized type " << myType << std::endl;
0911           break;
0912       }
0913 
0914       auto hBPix = std::make_shared<TH1F>(Form("%s BPIX", TypeName[myType]),
0915                                           Form("SiPixel Gain Calibration %s BPIx -%s;per %s %s (BPix);# %ss",
0916                                                (isForHLT_ ? "ForHLT" : "Offline"),
0917                                                TypeName[myType],
0918                                                (isForHLT_ ? "Column" : "Pixel"),
0919                                                TypeName[myType],
0920                                                (isForHLT_ ? "column" : "pixel")),
0921                                           200,
0922                                           minimum,
0923                                           maximum);
0924 
0925       auto hFPix = std::make_shared<TH1F>(Form("%s FPIX", TypeName[myType]),
0926                                           Form("SiPixel Gain Calibration %s FPix -%s;per %s %s (FPix);# %ss",
0927                                                (isForHLT_ ? "ForHLT" : "Offline"),
0928                                                TypeName[myType],
0929                                                (isForHLT_ ? "Column" : "Pixel"),
0930                                                TypeName[myType],
0931                                                (isForHLT_ ? "column" : "pixel")),
0932                                           200,
0933                                           minimum,
0934                                           maximum);
0935 
0936       for (unsigned int i : {1, 2}) {
0937         SiPixelPI::adjustCanvasMargins(canvas.cd(i), 0.04, 0.12, 0.12, 0.02);
0938         canvas.cd(i)->Modified();
0939       }
0940 
0941       // actually fill the histograms
0942       fillTheHistos(payload, hBPix, hFPix, myType);
0943 
0944       canvas.cd(1)->SetLogy();
0945       hBPix->SetTitle("");
0946       hBPix->GetYaxis()->SetRangeUser(0.1, hBPix->GetMaximum() * 10);
0947       hBPix->SetFillColor(kBlue);
0948       hBPix->SetMarkerStyle(20);
0949       hBPix->SetMarkerSize(1);
0950       hBPix->Draw("hist");
0951 
0952       SiPixelPI::makeNicePlotStyle(hBPix.get());
0953       hBPix->SetStats(true);
0954 
0955       canvas.cd(2)->SetLogy();
0956       hFPix->SetTitle("");
0957       hFPix->GetYaxis()->SetRangeUser(0.1, hFPix->GetMaximum() * 10);
0958       hFPix->SetFillColor(kBlue);
0959       hFPix->SetMarkerStyle(20);
0960       hFPix->SetMarkerSize(1);
0961       hFPix->Draw("hist");
0962 
0963       SiPixelPI::makeNicePlotStyle(hFPix.get());
0964       hFPix->SetStats(true);
0965 
0966       canvas.Update();
0967 
0968       TLegend legend = TLegend(0.32, 0.92, 0.97, 0.95);
0969       legend.SetHeader(("Payload hash: #bf{" + (std::get<1>(iov)) + "}").c_str(),
0970                        "C");  // option "C" allows to center the header
0971       //legend.AddEntry(h1.get(), ("IOV: " + std::to_string(std::get<0>(iov))).c_str(), "PL");
0972       legend.SetLineColor(10);
0973       legend.SetTextSize(0.025);
0974       canvas.cd(1);
0975       legend.Draw("same");
0976       canvas.cd(2);
0977       legend.Draw("same");
0978 
0979       canvas.cd(1);
0980       TPaveStats* st1 = (TPaveStats*)hBPix->FindObject("stats");
0981       st1->SetTextSize(0.03);
0982       SiPixelPI::adjustStats(st1, 0.13, 0.815, 0.44, 0.915);
0983 
0984       canvas.cd(2);
0985       TPaveStats* st2 = (TPaveStats*)hFPix->FindObject("stats");
0986       st2->SetTextSize(0.03);
0987       SiPixelPI::adjustStats(st2, 0.14, 0.815, 0.44, 0.915);
0988 
0989       auto ltx = TLatex();
0990       ltx.SetTextFont(62);
0991       //ltx.SetTextColor(kBlue);
0992       ltx.SetTextSize(0.045);
0993       ltx.SetTextAlign(11);
0994       canvas.cd(1);
0995       ltx.DrawLatexNDC(gPad->GetLeftMargin() + 0.01,
0996                        1 - gPad->GetTopMargin() + 0.01,
0997                        ("SiPixel Gain Calibration IOV:" + std::to_string(std::get<0>(iov))).c_str());
0998 
0999       canvas.cd(2);
1000       ltx.DrawLatexNDC(gPad->GetLeftMargin() + 0.01,
1001                        1 - gPad->GetTopMargin() + 0.01,
1002                        ("SiPixel Gain Calibration IOV:" + std::to_string(std::get<0>(iov))).c_str());
1003 
1004       std::string fileName(this->m_imageFileName);
1005       canvas.SaveAs(fileName.c_str());
1006 
1007       return true;
1008     }
1009 
1010   protected:
1011     bool isForHLT_;
1012     std::string label_;
1013   };
1014 
1015   /************************************************
1016     1d histogram comparison of SiPixelGainCalibration
1017   *************************************************/
1018   template <gainCalibPI::type myType, class PayloadType>
1019   class SiPixelGainCalibrationValueComparisonBase : public cond::payloadInspector::PlotImage<PayloadType> {
1020   public:
1021     SiPixelGainCalibrationValueComparisonBase()
1022         : cond::payloadInspector::PlotImage<PayloadType>(
1023               Form("SiPixelGainCalibration %s Values Comparison", TypeName[myType])) {
1024       if constexpr (std::is_same_v<PayloadType, SiPixelGainCalibrationOffline>) {
1025         isForHLT_ = false;
1026       } else {
1027         isForHLT_ = true;
1028       }
1029     }
1030     bool fill(const std::vector<std::tuple<cond::Time_t, cond::Hash>>& iovs) override {
1031       gStyle->SetOptStat("emr");
1032       TGaxis::SetExponentOffset(-0.1, 0.01, "y");  // Y offset
1033       TH1F::SetDefaultSumw2(true);
1034 
1035       std::vector<std::tuple<cond::Time_t, cond::Hash>> sorted_iovs = iovs;
1036       // make absolute sure the IOVs are sortd by since
1037       std::sort(begin(sorted_iovs), end(sorted_iovs), [](auto const& t1, auto const& t2) {
1038         return std::get<0>(t1) < std::get<0>(t2);
1039       });
1040       auto firstiov = sorted_iovs.front();
1041       auto lastiov = sorted_iovs.back();
1042 
1043       std::shared_ptr<PayloadType> last_payload = this->fetchPayload(std::get<1>(lastiov));
1044       std::shared_ptr<PayloadType> first_payload = this->fetchPayload(std::get<1>(firstiov));
1045 
1046       std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
1047       std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
1048 
1049       float minimum(9999.);
1050       float maximum(-9999.);
1051 
1052       switch (myType) {
1053         case gainCalibPI::t_gain:
1054           maximum = std::max(last_payload->getGainHigh(), first_payload->getGainHigh());
1055           minimum = std::min(last_payload->getGainLow(), first_payload->getGainLow());
1056           break;
1057         case gainCalibPI::t_pedestal:
1058           maximum = std::max(last_payload->getPedHigh(), first_payload->getPedHigh());
1059           minimum = std::min(last_payload->getPedLow(), first_payload->getPedLow());
1060           break;
1061         default:
1062           edm::LogError(label_) << "Unrecognized type " << myType << std::endl;
1063           break;
1064       }
1065 
1066       TCanvas canvas("Canv", "Canv", 1200, 1000);
1067       canvas.cd();
1068       auto hfirst = std::make_shared<TH1F>(Form("First, IOV %s", firstIOVsince.c_str()),
1069                                            Form("SiPixel Gain Calibration %s - %s;per %s %s;# %ss",
1070                                                 (isForHLT_ ? "ForHLT" : "Offline"),
1071                                                 TypeName[myType],
1072                                                 (isForHLT_ ? "Column" : "Pixel"),
1073                                                 TypeName[myType],
1074                                                 (isForHLT_ ? "column" : "pixel")),
1075                                            200,
1076                                            minimum,
1077                                            maximum);
1078 
1079       auto hlast = std::make_shared<TH1F>(Form("Last, IOV %s", lastIOVsince.c_str()),
1080                                           Form("SiPixel Gain Calibration %s - %s;per %s %s;# %ss",
1081                                                (isForHLT_ ? "ForHLT" : "Offline"),
1082                                                TypeName[myType],
1083                                                (isForHLT_ ? "Column" : "Pixel"),
1084                                                TypeName[myType],
1085                                                (isForHLT_ ? "column" : "pixel")),
1086                                           200,
1087                                           minimum,
1088                                           maximum);
1089 
1090       SiPixelPI::adjustCanvasMargins(canvas.cd(), 0.05, 0.12, 0.12, 0.03);
1091       canvas.Modified();
1092 
1093       gainCalibPI::fillTheHisto(first_payload, hfirst, myType);
1094       gainCalibPI::fillTheHisto(last_payload, hlast, myType);
1095 
1096       canvas.cd()->SetLogy();
1097       auto extrema = SiPixelPI::getExtrema(hfirst.get(), hlast.get());
1098       //hfirst->GetYaxis()->SetRangeUser(extrema.first, extrema.second * 1.10);
1099       hfirst->GetYaxis()->SetRangeUser(1., extrema.second * 10);
1100 
1101       hfirst->SetTitle("");
1102       hfirst->SetLineColor(kRed);
1103       hfirst->SetBarWidth(0.95);
1104       hfirst->Draw("hist");
1105 
1106       hlast->SetTitle("");
1107       hlast->SetFillColorAlpha(kBlue, 0.20);
1108       hlast->SetBarWidth(0.95);
1109       hlast->Draw("histsames");
1110 
1111       SiPixelPI::makeNicePlotStyle(hfirst.get());
1112       hfirst->SetStats(true);
1113       SiPixelPI::makeNicePlotStyle(hlast.get());
1114       hlast->SetStats(true);
1115 
1116       canvas.Update();
1117 
1118       TLegend legend = TLegend(0.45, 0.86, 0.74, 0.94);
1119       //legend.SetHeader("#font[22]{SiPixel Offline Gain Calibration Comparison}", "C");  // option "C" allows to center the header
1120       //legend.AddEntry(hfirst.get(), ("IOV: " + std::to_string(std::get<0>(firstiov))).c_str(), "FL");
1121       //legend.AddEntry(hlast.get(),  ("IOV: " + std::to_string(std::get<0>(lastiov))).c_str(), "FL");
1122       legend.AddEntry(hfirst.get(), ("payload: #color[2]{" + std::get<1>(firstiov) + "}").c_str(), "F");
1123       legend.AddEntry(hlast.get(), ("payload: #color[4]{" + std::get<1>(lastiov) + "}").c_str(), "F");
1124       legend.SetTextSize(0.022);
1125       legend.SetLineColor(10);
1126       legend.Draw("same");
1127 
1128       TPaveStats* st1 = (TPaveStats*)hfirst->FindObject("stats");
1129       st1->SetTextSize(0.021);
1130       st1->SetLineColor(kRed);
1131       st1->SetTextColor(kRed);
1132       SiPixelPI::adjustStats(st1, 0.13, 0.86, 0.31, 0.94);
1133 
1134       TPaveStats* st2 = (TPaveStats*)hlast->FindObject("stats");
1135       st2->SetTextSize(0.021);
1136       st2->SetLineColor(kBlue);
1137       st2->SetTextColor(kBlue);
1138       SiPixelPI::adjustStats(st2, 0.13, 0.77, 0.31, 0.85);
1139 
1140       auto ltx = TLatex();
1141       ltx.SetTextFont(62);
1142       //ltx.SetTextColor(kBlue);
1143       ltx.SetTextSize(0.047);
1144       ltx.SetTextAlign(11);
1145       ltx.DrawLatexNDC(gPad->GetLeftMargin(),
1146                        1 - gPad->GetTopMargin() + 0.01,
1147                        ("SiPixel Gain Calibration IOV:#color[2]{" + std::to_string(std::get<0>(firstiov)) +
1148                         "} vs IOV:#color[4]{" + std::to_string(std::get<0>(lastiov)) + "}")
1149                            .c_str());
1150 
1151       std::string fileName(this->m_imageFileName);
1152       canvas.SaveAs(fileName.c_str());
1153 #ifdef MMDEBUG
1154       canvas.SaveAs("out.root");
1155 #endif
1156 
1157       return true;
1158     }
1159 
1160   protected:
1161     bool isForHLT_;
1162     std::string label_;
1163   };
1164 
1165   template <gainCalibPI::type myType, class PayloadType>
1166   class SiPixelGainCalibrationValueComparisonSingleTag
1167       : public SiPixelGainCalibrationValueComparisonBase<myType, PayloadType> {
1168   public:
1169     SiPixelGainCalibrationValueComparisonSingleTag()
1170         : SiPixelGainCalibrationValueComparisonBase<myType, PayloadType>() {
1171       this->setSingleIov(false);
1172     }
1173   };
1174 
1175   template <gainCalibPI::type myType, class PayloadType>
1176   class SiPixelGainCalibrationValueComparisonTwoTags
1177       : public SiPixelGainCalibrationValueComparisonBase<myType, PayloadType> {
1178   public:
1179     SiPixelGainCalibrationValueComparisonTwoTags() : SiPixelGainCalibrationValueComparisonBase<myType, PayloadType>() {
1180       this->setTwoTags(true);
1181     }
1182   };
1183 
1184   /************************************************
1185     Diff and Ratio histograms of 2 IOVs
1186   *************************************************/
1187   template <gainCalibPI::type myType, cond::payloadInspector::IOVMultiplicity nIOVs, int ntags, class PayloadType>
1188   class SiPixelGainCalibDiffAndRatioBase : public cond::payloadInspector::PlotImage<PayloadType, nIOVs, ntags> {
1189   public:
1190     SiPixelGainCalibDiffAndRatioBase()
1191         : cond::payloadInspector::PlotImage<PayloadType, nIOVs, ntags>(
1192               Form("SiPixelGainCalibration %s Diff and Ratio %i tag(s)", TypeName[myType], ntags)) {
1193       if constexpr (std::is_same_v<PayloadType, SiPixelGainCalibrationOffline>) {
1194         isForHLT_ = false;
1195       } else {
1196         isForHLT_ = true;
1197       }
1198     }
1199 
1200     bool fill() override {
1201       gStyle->SetOptStat("emr");
1202       TGaxis::SetExponentOffset(-0.1, 0.01, "y");  // Y offset
1203       TH1F::SetDefaultSumw2(true);
1204 
1205       COUT << "ntags: " << ntags << " this->m_plotAnnotations.ntags: " << this->m_plotAnnotations.ntags << std::endl;
1206 
1207       // trick to deal with the multi-ioved tag and two tag case at the same time
1208       auto theIOVs = cond::payloadInspector::PlotBase::getTag<0>().iovs;
1209       auto f_tagname = cond::payloadInspector::PlotBase::getTag<0>().name;
1210       std::string l_tagname = "";
1211       auto firstiov = theIOVs.front();
1212       std::tuple<cond::Time_t, cond::Hash> lastiov;
1213 
1214       // we don't support (yet) comparison with more than 2 tags
1215       assert(this->m_plotAnnotations.ntags < 3);
1216 
1217       if (this->m_plotAnnotations.ntags == 2) {
1218         auto tag2iovs = cond::payloadInspector::PlotBase::getTag<1>().iovs;
1219         l_tagname = cond::payloadInspector::PlotBase::getTag<1>().name;
1220         lastiov = tag2iovs.front();
1221       } else {
1222         lastiov = theIOVs.back();
1223       }
1224 
1225       std::shared_ptr<PayloadType> last_payload = this->fetchPayload(std::get<1>(lastiov));
1226       std::shared_ptr<PayloadType> first_payload = this->fetchPayload(std::get<1>(firstiov));
1227 
1228       std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
1229       std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
1230 
1231       TCanvas canvas("Canv", "Canv", 1300, 800);
1232       canvas.Divide(2, 1);
1233       canvas.cd();
1234 
1235       SiPixelPI::adjustCanvasMargins(canvas.cd(1), 0.05, 0.12, 0.12, 0.04);
1236       SiPixelPI::adjustCanvasMargins(canvas.cd(2), 0.05, 0.12, 0.12, 0.04);
1237       canvas.Modified();
1238 
1239       auto array = gainCalibPI::fillDiffAndRatio(first_payload, last_payload, myType);
1240 
1241       array[0]->SetTitle(Form("SiPixel Gain Calibration %s - %s;per %s %s ratio;# %ss",
1242                               (isForHLT_ ? "ForHLT" : "Offline"),
1243                               TypeName[myType],
1244                               (isForHLT_ ? "Column" : "Pixel"),
1245                               TypeName[myType],
1246                               (isForHLT_ ? "column" : "pixel")));
1247 
1248       array[1]->SetTitle(Form("SiPixel Gain Calibration %s - %s;per %s %s difference;# %ss",
1249                               (isForHLT_ ? "ForHLT" : "Offline"),
1250                               TypeName[myType],
1251                               (isForHLT_ ? "Column" : "Pixel"),
1252                               TypeName[myType],
1253                               (isForHLT_ ? "column" : "pixel")));
1254 
1255       canvas.cd(1)->SetLogy();
1256       array[0]->SetTitle("");
1257       array[0]->SetLineColor(kBlack);
1258       array[0]->SetFillColor(kRed);
1259       array[0]->SetBarWidth(0.90);
1260       array[0]->SetMaximum(array[0]->GetMaximum() * 10);
1261       array[0]->Draw("bar");
1262       SiPixelPI::makeNicePlotStyle(array[0].get());
1263       array[0]->SetStats(true);
1264 
1265       canvas.cd(2)->SetLogy();
1266       array[1]->SetTitle("");
1267       array[1]->SetLineColor(kBlack);
1268       array[1]->SetFillColor(kBlue);
1269       array[1]->SetBarWidth(0.90);
1270       array[1]->SetMaximum(array[1]->GetMaximum() * 10);
1271       array[1]->Draw("bar");
1272       SiPixelPI::makeNicePlotStyle(array[1].get());
1273       array[1]->SetStats(true);
1274 
1275       canvas.Update();
1276 
1277       canvas.cd(1);
1278       TLatex latex;
1279       latex.SetTextSize(0.024);
1280       latex.SetTextAlign(13);  //align at top
1281       latex.DrawLatexNDC(
1282           .41,
1283           .94,
1284           fmt::sprintf("#scale[1.2]{SiPixelGainCalibration%s Ratio}", (isForHLT_ ? "ForHLT" : "Offline")).c_str());
1285       if (this->m_plotAnnotations.ntags == 2) {
1286         latex.DrawLatexNDC(
1287             .41, .91, ("#splitline{#font[12]{" + f_tagname + "}}{ / #font[12]{" + l_tagname + "}}").c_str());
1288       } else {
1289         latex.DrawLatexNDC(.41, .91, (firstIOVsince + " / " + lastIOVsince).c_str());
1290       }
1291 
1292       canvas.cd(2);
1293       TLatex latex2;
1294       latex2.SetTextSize(0.024);
1295       latex2.SetTextAlign(13);  //align at top
1296       latex2.DrawLatexNDC(
1297           .41,
1298           .94,
1299           fmt::sprintf("#scale[1.2]{SiPixelGainCalibration%s Diff}", (isForHLT_ ? "ForHLT" : "Offline")).c_str());
1300       if (this->m_plotAnnotations.ntags == 2) {
1301         latex2.DrawLatexNDC(
1302             .41, .91, ("#splitline{#font[12]{" + f_tagname + "}}{ - #font[12]{" + l_tagname + "}}").c_str());
1303       } else {
1304         latex2.DrawLatexNDC(.41, .91, (firstIOVsince + " - " + lastIOVsince).c_str());
1305       }
1306 
1307       TPaveStats* st1 = (TPaveStats*)array[0]->FindObject("stats");
1308       st1->SetTextSize(0.027);
1309       st1->SetLineColor(kRed);
1310       st1->SetTextColor(kRed);
1311       SiPixelPI::adjustStats(st1, 0.13, 0.84, 0.40, 0.94);
1312 
1313       TPaveStats* st2 = (TPaveStats*)array[1]->FindObject("stats");
1314       st2->SetTextSize(0.027);
1315       st2->SetLineColor(kBlue);
1316       st2->SetTextColor(kBlue);
1317       SiPixelPI::adjustStats(st2, 0.13, 0.84, 0.40, 0.94);
1318 
1319       auto ltx = TLatex();
1320       ltx.SetTextFont(62);
1321       //ltx.SetTextColor(kBlue);
1322       ltx.SetTextSize(0.040);
1323       ltx.SetTextAlign(11);
1324       canvas.cd(1);
1325       ltx.DrawLatexNDC(
1326           gPad->GetLeftMargin(),
1327           1 - gPad->GetTopMargin() + 0.01,
1328           fmt::sprintf("SiPixel %s Ratio, IOV %s / %s", TypeName[myType], firstIOVsince, lastIOVsince).c_str());
1329 
1330       canvas.cd(2);
1331       ltx.DrawLatexNDC(
1332           gPad->GetLeftMargin(),
1333           1 - gPad->GetTopMargin() + 0.01,
1334           fmt::sprintf("SiPixel %s Diff, IOV %s - %s", TypeName[myType], firstIOVsince, lastIOVsince).c_str());
1335 
1336       std::string fileName(this->m_imageFileName);
1337       canvas.SaveAs(fileName.c_str());
1338 #ifdef MMDEBUG
1339       canvas.SaveAs("out.root");
1340 #endif
1341 
1342       return true;
1343     }
1344 
1345   protected:
1346     bool isForHLT_;
1347     std::string label_;
1348   };
1349 
1350   // 2D MAPS
1351 
1352   /************************************************
1353    occupancy style map BPix
1354   *************************************************/
1355   template <gainCalibPI::type myType, class PayloadType, SiPixelPI::DetType myDetType>
1356   class SiPixelGainCalibrationMap
1357       : public cond::payloadInspector::PlotImage<PayloadType, cond::payloadInspector::SINGLE_IOV> {
1358   public:
1359     SiPixelGainCalibrationMap()
1360         : cond::payloadInspector::PlotImage<PayloadType, cond::payloadInspector::SINGLE_IOV>(
1361               Form("SiPixelGainCalibration %s Barrel Pixel Map", TypeName[myType])),
1362           m_trackerTopo{StandaloneTrackerTopology::fromTrackerParametersXMLFile(
1363               edm::FileInPath("Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml").fullPath())} {
1364       if constexpr (std::is_same_v<PayloadType, SiPixelGainCalibrationOffline>) {
1365         isForHLT_ = false;
1366         label_ = "SiPixelGainCalibrationOffline_PayloadInspector";
1367       } else {
1368         isForHLT_ = true;
1369         label_ = "SiPixelGainCalibrationForHLT_PayloadInspector";
1370       }
1371     };
1372 
1373     bool fill() override {
1374       TGaxis::SetMaxDigits(2);
1375 
1376       auto tag = cond::payloadInspector::PlotBase::getTag<0>();
1377       auto tagname = tag.name;
1378       auto iov = tag.iovs.front();
1379 
1380       std::shared_ptr<PayloadType> payload = this->fetchPayload(std::get<1>(iov));
1381 
1382       std::string ztitle =
1383           fmt::sprintf("average per %s %s", isForHLT_ ? "column" : "pixel", gainCalibPI::t_titles[myType]);
1384       Phase1PixelROCMaps theGainsMap("", ztitle);
1385 
1386       std::map<uint32_t, float> GainCalibMap_;
1387       gainCalibPI::fillThePerModuleMap(payload, GainCalibMap_, myType);
1388       if (GainCalibMap_.size() != SiPixelPI::phase1size) {
1389         edm::LogError(label_) << "SiPixelGainCalibration maps are not supported for non-Phase1 Pixel geometries !";
1390         TCanvas canvas("Canv", "Canv", 1200, 1000);
1391         SiPixelPI::displayNotSupported(canvas, GainCalibMap_.size());
1392         std::string fileName(this->m_imageFileName);
1393         canvas.SaveAs(fileName.c_str());
1394         return false;
1395       }
1396 
1397       // hard-coded phase-I
1398       std::array<double, n_layers> b_minima = {{999., 999., 999., 999.}};
1399       std::array<double, n_rings> f_minima = {{999., 999.}};
1400 
1401       for (const auto& element : GainCalibMap_) {
1402         int subid = DetId(element.first).subdetId();
1403         if (subid == PixelSubdetector::PixelBarrel) {
1404           auto layer = m_trackerTopo.pxbLayer(DetId(element.first));
1405           if (element.second < b_minima.at(layer - 1)) {
1406             b_minima.at(layer - 1) = element.second;
1407           }
1408           theGainsMap.fillWholeModule(element.first, element.second);
1409         } else if (subid == PixelSubdetector::PixelEndcap) {
1410           auto ring = SiPixelPI::ring(DetId(element.first), m_trackerTopo, true);
1411           if (element.second < f_minima.at(ring - 1)) {
1412             f_minima.at(ring - 1) = element.second;
1413           }
1414           theGainsMap.fillWholeModule(element.first, element.second);
1415         }
1416       }
1417 
1418       gStyle->SetOptStat(0);
1419       //=========================
1420       TCanvas canvas("Summary", "Summary", 1200, k_height[myDetType]);
1421       canvas.cd();
1422 
1423       auto unpacked = SiPixelPI::unpack(std::get<0>(iov));
1424 
1425       std::string IOVstring = (unpacked.first == 0)
1426                                   ? std::to_string(unpacked.second)
1427                                   : (std::to_string(unpacked.first) + "," + std::to_string(unpacked.second));
1428 
1429       const auto headerText = fmt::sprintf("#color[4]{%s},  IOV: #color[4]{%s}", tagname, IOVstring);
1430 
1431       switch (myDetType) {
1432         case SiPixelPI::t_barrel:
1433           theGainsMap.drawBarrelMaps(canvas, headerText);
1434           break;
1435         case SiPixelPI::t_forward:
1436           theGainsMap.drawForwardMaps(canvas, headerText);
1437           break;
1438         case SiPixelPI::t_all:
1439           theGainsMap.drawMaps(canvas, headerText);
1440           break;
1441         default:
1442           throw cms::Exception("SiPixelGainCalibrationMap")
1443               << "\nERROR: unrecognized Pixel Detector part " << std::endl;
1444       }
1445 
1446       if (myDetType == SiPixelPI::t_barrel || myDetType == SiPixelPI::t_all) {
1447         for (unsigned int lay = 1; lay <= n_layers; lay++) {
1448           //SiPixelPI::adjustCanvasMargins(canvas.cd(lay), -1, 0.08, 0.1, 0.13);
1449 
1450           auto h_bpix_Gains = theGainsMap.getLayerMaps();
1451 
1452           COUT << " layer:" << lay << " max:" << h_bpix_Gains[lay - 1]->GetMaximum() << " min: " << b_minima.at(lay - 1)
1453                << std::endl;
1454 
1455           h_bpix_Gains[lay - 1]->GetZaxis()->SetRangeUser(b_minima.at(lay - 1) - 0.001,
1456                                                           h_bpix_Gains[lay - 1]->GetMaximum() + 0.001);
1457         }
1458       }
1459 
1460       if (myDetType == SiPixelPI::t_forward || myDetType == SiPixelPI::t_all) {
1461         for (unsigned int ring = 1; ring <= n_rings; ring++) {
1462           //SiPixelPI::adjustCanvasMargins(canvas.cd(ring), -1, 0.08, 0.1, 0.13);
1463 
1464           auto h_fpix_Gains = theGainsMap.getRingMaps();
1465 
1466           COUT << " ring:" << ring << " max:" << h_fpix_Gains[ring - 1]->GetMaximum()
1467                << " min: " << f_minima.at(ring - 1) << std::endl;
1468 
1469           h_fpix_Gains[ring - 1]->GetZaxis()->SetRangeUser(f_minima.at(ring - 1) - 0.001,
1470                                                            h_fpix_Gains[ring - 1]->GetMaximum() + 0.001);
1471         }
1472       }
1473 
1474       std::string fileName(this->m_imageFileName);
1475       canvas.SaveAs(fileName.c_str());
1476 
1477       return true;
1478     }
1479 
1480   private:
1481     TrackerTopology m_trackerTopo;
1482     static constexpr std::array<int, 3> k_height = {{1200, 600, 1600}};
1483     static constexpr int n_layers = 4;
1484     static constexpr int n_rings = 2;
1485 
1486   protected:
1487     bool isForHLT_;
1488     std::string label_;
1489   };
1490 
1491   /************************************************
1492    Summary Comparison per region of SiPixelGainCalibration between 2 IOVs
1493   *************************************************/
1494   template <gainCalibPI::type myType, class PayloadType, cond::payloadInspector::IOVMultiplicity nIOVs, int ntags>
1495   class SiPixelGainCalibrationByRegionComparisonBase
1496       : public cond::payloadInspector::PlotImage<PayloadType, nIOVs, ntags> {
1497   public:
1498     SiPixelGainCalibrationByRegionComparisonBase()
1499         : cond::payloadInspector::PlotImage<PayloadType, nIOVs, ntags>(
1500               Form("SiPixelGainCalibration %s Comparison by Region", TypeName[myType])) {
1501       if constexpr (std::is_same_v<PayloadType, SiPixelGainCalibrationOffline>) {
1502         isForHLT_ = false;
1503         label_ = "SiPixelGainCalibrationOffline_PayloadInspector";
1504       } else {
1505         isForHLT_ = true;
1506         label_ = "SiPixelGainCalibrationForHLT_PayloadInspector";
1507       }
1508     }
1509 
1510     bool fill() override {
1511       gStyle->SetPaintTextFormat(".3f");
1512 
1513       // trick to deal with the multi-ioved tag and two tag case at the same time
1514       auto theIOVs = cond::payloadInspector::PlotBase::getTag<0>().iovs;
1515       auto f_tagname = cond::payloadInspector::PlotBase::getTag<0>().name;
1516       std::string l_tagname = "";
1517       auto firstiov = theIOVs.front();
1518       std::tuple<cond::Time_t, cond::Hash> lastiov;
1519 
1520       // we don't support (yet) comparison with more than 2 tags
1521       assert(this->m_plotAnnotations.ntags < 3);
1522 
1523       if (this->m_plotAnnotations.ntags == 2) {
1524         auto tag2iovs = cond::payloadInspector::PlotBase::getTag<1>().iovs;
1525         l_tagname = cond::payloadInspector::PlotBase::getTag<1>().name;
1526         lastiov = tag2iovs.front();
1527       } else {
1528         lastiov = theIOVs.back();
1529       }
1530 
1531       std::shared_ptr<PayloadType> last_payload = this->fetchPayload(std::get<1>(lastiov));
1532       std::shared_ptr<PayloadType> first_payload = this->fetchPayload(std::get<1>(firstiov));
1533 
1534       std::map<uint32_t, float> f_GainsMap_;
1535       gainCalibPI::fillThePerModuleMap(first_payload, f_GainsMap_, myType);
1536 
1537       std::map<uint32_t, float> l_GainsMap_;
1538       gainCalibPI::fillThePerModuleMap(last_payload, l_GainsMap_, myType);
1539 
1540       std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
1541       std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
1542 
1543       TCanvas canvas("Comparison", "Comparison", 1600, 800);
1544 
1545       SiPixelPI::PhaseInfo f_phaseInfo(f_GainsMap_.size());
1546       SiPixelPI::PhaseInfo l_phaseInfo(l_GainsMap_.size());
1547 
1548       std::map<SiPixelPI::regions, std::shared_ptr<TH1F>> FirstGains_spectraByRegion;
1549       std::map<SiPixelPI::regions, std::shared_ptr<TH1F>> LastGains_spectraByRegion;
1550       std::shared_ptr<TH1F> summaryFirst;
1551       std::shared_ptr<TH1F> summaryLast;
1552 
1553       float minimum(9999.);
1554       float maximum(-9999.);
1555 
1556       switch (myType) {
1557         case gainCalibPI::t_gain:
1558           maximum = std::max(last_payload->getGainHigh(), first_payload->getGainHigh());
1559           minimum = std::min(last_payload->getGainLow(), first_payload->getGainLow());
1560           break;
1561         case gainCalibPI::t_pedestal:
1562           maximum = std::max(last_payload->getPedHigh(), first_payload->getPedHigh());
1563           minimum = std::min(last_payload->getPedLow(), first_payload->getPedLow());
1564           break;
1565         default:
1566           edm::LogError(label_) << "Unrecognized type " << myType << std::endl;
1567           break;
1568       }
1569 
1570       // book the intermediate histograms
1571       for (int r = SiPixelPI::BPixL1o; r != SiPixelPI::NUM_OF_REGIONS; r++) {
1572         SiPixelPI::regions part = static_cast<SiPixelPI::regions>(r);
1573         std::string s_part = SiPixelPI::getStringFromRegionEnum(part);
1574 
1575         FirstGains_spectraByRegion[part] =
1576             std::make_shared<TH1F>(Form("hfirstGains_%s", s_part.c_str()),
1577                                    Form(";%s #LT %s #GT;n. of modules", s_part.c_str(), TypeName[myType]),
1578                                    200,
1579                                    minimum,
1580                                    maximum);
1581 
1582         LastGains_spectraByRegion[part] =
1583             std::make_shared<TH1F>(Form("hlastGains_%s", s_part.c_str()),
1584                                    Form(";%s #LT %s #GT;n. of modules", s_part.c_str(), TypeName[myType]),
1585                                    200,
1586                                    minimum,
1587                                    maximum);
1588       }
1589 
1590       summaryFirst = std::make_shared<TH1F>("first Summary",
1591                                             Form("Summary of #LT per %s %s #GT;;average %s",
1592                                                  (isForHLT_ ? "Column" : "Pixel"),
1593                                                  TypeName[myType],
1594                                                  TypeName[myType]),
1595                                             FirstGains_spectraByRegion.size(),
1596                                             0,
1597                                             FirstGains_spectraByRegion.size());
1598 
1599       summaryLast = std::make_shared<TH1F>("last Summary",
1600                                            Form("Summary of #LT per %s %s #GT;;average %s",
1601                                                 (isForHLT_ ? "Column" : "Pixel"),
1602                                                 TypeName[myType],
1603                                                 TypeName[myType]),
1604                                            LastGains_spectraByRegion.size(),
1605                                            0,
1606                                            LastGains_spectraByRegion.size());
1607 
1608       // deal with first IOV
1609       const char* path_toTopologyXML = f_phaseInfo.pathToTopoXML();
1610 
1611       auto f_tTopo =
1612           StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
1613 
1614       // -------------------------------------------------------------------
1615       // loop on the first Gains Map
1616       // -------------------------------------------------------------------
1617       for (const auto& it : f_GainsMap_) {
1618         if (DetId(it.first).det() != DetId::Tracker) {
1619           edm::LogWarning(label_) << "Encountered invalid Tracker DetId:" << it.first << " - terminating ";
1620           return false;
1621         }
1622 
1623         SiPixelPI::topolInfo t_info_fromXML;
1624         t_info_fromXML.init();
1625         DetId detid(it.first);
1626         t_info_fromXML.fillGeometryInfo(detid, f_tTopo, f_phaseInfo.phase());
1627 
1628         SiPixelPI::regions thePart = t_info_fromXML.filterThePartition();
1629         FirstGains_spectraByRegion[thePart]->Fill(it.second);
1630       }  // ends loop on the vector of error transforms
1631 
1632       // deal with last IOV
1633       path_toTopologyXML = l_phaseInfo.pathToTopoXML();
1634 
1635       auto l_tTopo =
1636           StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
1637 
1638       // -------------------------------------------------------------------
1639       // loop on the second Gains Map
1640       // -------------------------------------------------------------------
1641       for (const auto& it : l_GainsMap_) {
1642         if (DetId(it.first).det() != DetId::Tracker) {
1643           edm::LogWarning(label_) << "Encountered invalid Tracker DetId:" << it.first << " - terminating ";
1644           return false;
1645         }
1646 
1647         SiPixelPI::topolInfo t_info_fromXML;
1648         t_info_fromXML.init();
1649         DetId detid(it.first);
1650         t_info_fromXML.fillGeometryInfo(detid, l_tTopo, l_phaseInfo.phase());
1651 
1652         SiPixelPI::regions thePart = t_info_fromXML.filterThePartition();
1653         LastGains_spectraByRegion[thePart]->Fill(it.second);
1654       }  // ends loop on the vector of error transforms
1655 
1656       // fill the summary plots
1657       int bin = 1;
1658       for (int r = SiPixelPI::BPixL1o; r != SiPixelPI::NUM_OF_REGIONS; r++) {
1659         SiPixelPI::regions part = static_cast<SiPixelPI::regions>(r);
1660 
1661         summaryFirst->GetXaxis()->SetBinLabel(bin, SiPixelPI::getStringFromRegionEnum(part).c_str());
1662         // avoid filling the histogram with numerical noise
1663         float f_mean =
1664             FirstGains_spectraByRegion[part]->GetMean() > 10.e-6 ? FirstGains_spectraByRegion[part]->GetMean() : 0.;
1665         summaryFirst->SetBinContent(bin, f_mean);
1666         //summaryFirst->SetBinError(bin,Gains_spectraByRegion[hash]->GetRMS());
1667 
1668         summaryLast->GetXaxis()->SetBinLabel(bin, SiPixelPI::getStringFromRegionEnum(part).c_str());
1669         // avoid filling the histogram with numerical noise
1670         float l_mean =
1671             LastGains_spectraByRegion[part]->GetMean() > 10.e-6 ? LastGains_spectraByRegion[part]->GetMean() : 0.;
1672         summaryLast->SetBinContent(bin, l_mean);
1673         //summaryLast->SetBinError(bin,Gains_spectraByRegion[hash]->GetRMS());
1674         bin++;
1675       }
1676 
1677       SiPixelPI::makeNicePlotStyle(summaryFirst.get());  //, kBlue);
1678       summaryFirst->SetMarkerColor(kRed);
1679       summaryFirst->GetXaxis()->LabelsOption("v");
1680       summaryFirst->GetXaxis()->SetLabelSize(0.05);
1681       summaryFirst->GetYaxis()->SetTitleOffset(0.9);
1682 
1683       SiPixelPI::makeNicePlotStyle(summaryLast.get());  //, kRed);
1684       summaryLast->SetMarkerColor(kBlue);
1685       summaryLast->GetYaxis()->SetTitleOffset(0.9);
1686       summaryLast->GetXaxis()->LabelsOption("v");
1687       summaryLast->GetXaxis()->SetLabelSize(0.05);
1688 
1689       canvas.cd()->SetGridy();
1690 
1691       SiPixelPI::adjustCanvasMargins(canvas.cd(), -1, 0.18, 0.11, 0.02);
1692       canvas.Modified();
1693 
1694       summaryFirst->SetFillColor(kRed);
1695       summaryLast->SetFillColor(kBlue);
1696 
1697       summaryFirst->SetBarWidth(0.45);
1698       summaryFirst->SetBarOffset(0.1);
1699 
1700       summaryLast->SetBarWidth(0.4);
1701       summaryLast->SetBarOffset(0.55);
1702 
1703       summaryLast->SetMarkerSize(1.2);
1704       summaryFirst->SetMarkerSize(1.2);
1705 
1706       float max = (summaryFirst->GetMaximum() > summaryLast->GetMaximum()) ? summaryFirst->GetMaximum()
1707                                                                            : summaryLast->GetMaximum();
1708 
1709       summaryFirst->GetYaxis()->SetRangeUser(0., std::max(0., max * 1.40));
1710 
1711       summaryFirst->Draw("b text0");
1712       summaryLast->Draw("b text0 same");
1713 
1714       TLegend legend = TLegend(0.52, 0.80, 0.98, 0.9);
1715       legend.SetHeader(Form("#LT %s #GT value comparison", TypeName[myType]),
1716                        "C");  // option "C" allows to center the header
1717 
1718       legend.SetHeader("#mu_{H} value comparison", "C");  // option "C" allows to center the header
1719       std::string l_tagOrHash, f_tagOrHash;
1720       if (this->m_plotAnnotations.ntags == 2) {
1721         l_tagOrHash = l_tagname;
1722         f_tagOrHash = f_tagname;
1723       } else {
1724         l_tagOrHash = std::get<1>(lastiov);
1725         f_tagOrHash = std::get<1>(firstiov);
1726       }
1727 
1728       legend.AddEntry(
1729           summaryLast.get(),
1730           ("IOV: #scale[1.2]{" + std::to_string(std::get<0>(lastiov)) + "} | #color[4]{" + l_tagOrHash + "}").c_str(),
1731           "F");
1732       legend.AddEntry(
1733           summaryFirst.get(),
1734           ("IOV: #scale[1.2]{" + std::to_string(std::get<0>(firstiov)) + "} | #color[2]{" + f_tagOrHash + "}").c_str(),
1735           "F");
1736 
1737       legend.SetTextSize(0.025);
1738       legend.Draw("same");
1739 
1740       std::string fileName(this->m_imageFileName);
1741       canvas.SaveAs(fileName.c_str());
1742       return true;
1743     }
1744 
1745   protected:
1746     bool isForHLT_;
1747     std::string label_;
1748   };
1749 }  // namespace gainCalibHelper
1750 
1751 #endif