Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-05-05 02:47:07

0001 /*!
0002   \file SiPixelDynamicInefficiency_PayloadInspector
0003   \Payload Inspector Plugin for SiPixelDynamicInefficiency
0004   \author M. Musich
0005   \version $Revision: 1.0 $
0006   \date $Date: 2018/10/18 14:48:00 $
0007 */
0008 
0009 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0010 
0011 #include "CondCore/Utilities/interface/PayloadInspectorModule.h"
0012 #include "CondCore/Utilities/interface/PayloadInspector.h"
0013 #include "CondCore/CondDB/interface/Time.h"
0014 #include "CondCore/SiPixelPlugins/interface/SiPixelPayloadInspectorHelper.h"
0015 #include "FWCore/ParameterSet/interface/FileInPath.h"
0016 #include "CalibTracker/StandaloneTrackerTopology/interface/StandaloneTrackerTopology.h"
0017 #include "CalibTracker/SiPixelESProducers/interface/SiPixelDetInfoFileReader.h"
0018 
0019 // the data format of the condition to be inspected
0020 #include "CondFormats/SiPixelObjects/interface/SiPixelDynamicInefficiency.h"
0021 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
0022 #include "DataFormats/DetId/interface/DetId.h"
0023 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0024 #include "DQM/TrackerRemapper/interface/Phase1PixelROCMaps.h"
0025 #include "DQM/TrackerRemapper/interface/Phase1PixelSummaryMap.h"
0026 
0027 #include <cmath>
0028 #include <memory>
0029 #include <sstream>
0030 #include <iostream>
0031 
0032 // include ROOT
0033 #include "TCanvas.h"
0034 #include "TGraph.h"
0035 #include "TH2F.h"
0036 #include "TLatex.h"
0037 #include "TLegend.h"
0038 #include "TLine.h"
0039 #include "TPave.h"
0040 #include "TPaveStats.h"
0041 #include "TStyle.h"
0042 #include "TF1.h"
0043 #include "TMath.h"
0044 #include <Math/Polynomial.h>
0045 
0046 namespace {
0047 
0048   using namespace cond::payloadInspector;
0049   namespace SiPixDynIneff {
0050 
0051     // different types of geometrical inefficiency factors
0052     enum factor { geom = 0, colgeom = 1, chipgeom = 2, pu = 3, INVALID = 4 };
0053     const std::array<std::string, 5> factorString = {
0054         {"pixel geometry", "column geometry", "chip geometry", "PU", "invalid"}};
0055 
0056     using FactorMap = std::map<unsigned int, double>;
0057     using PUFactorMap = std::map<unsigned int, std::vector<double> >;
0058 
0059     // constants for ROC level simulation for Phase1
0060     enum shiftEnumerator { FPixRocIdShift = 3, BPixRocIdShift = 6 };
0061     const int rocIdMaskBits = 0x1F;
0062 
0063     struct packedBadRocFraction {
0064       std::vector<int> badRocNumber;
0065       std::vector<float> badRocFrac;
0066     };
0067 
0068     using BRFractions = std::unordered_map<uint32_t, packedBadRocFraction>;
0069 
0070     //_________________________________________________
0071     BRFractions pbrf(std::shared_ptr<SiPixelDynamicInefficiency> payload) {
0072       BRFractions f;
0073       const std::map<uint32_t, double>& PixelGeomFactorsDBIn = payload->getPixelGeomFactors();
0074 
0075       // first fill
0076       for (const auto db_factor : PixelGeomFactorsDBIn) {
0077         int subid = DetId(db_factor.first).subdetId();
0078         int shift = (subid == static_cast<int>(PixelSubdetector::PixelBarrel)) ? BPixRocIdShift : FPixRocIdShift;
0079         unsigned int rocMask = rocIdMaskBits << shift;
0080         unsigned int rocId = (((db_factor.first) & rocMask) >> shift);
0081         uint32_t rawid = db_factor.first & (~rocMask);
0082 
0083         if (f.find(rawid) == f.end()) {
0084           packedBadRocFraction p;
0085           f.insert(std::make_pair(rawid, p));
0086         }
0087 
0088         if (rocId != 0) {
0089           rocId--;
0090           double factor = db_factor.second;
0091           double badFraction = 1 - factor;
0092 
0093           f.at(rawid).badRocNumber.emplace_back(rocId);
0094           f.at(rawid).badRocFrac.emplace_back(badFraction);
0095         }
0096       }
0097       return f;
0098     }
0099 
0100     //_________________________________________________
0101     bool isPhase0(const BRFractions& fractions) {
0102       SiPixelDetInfoFileReader reader =
0103           SiPixelDetInfoFileReader(edm::FileInPath(SiPixelDetInfoFileReader::kPh0DefaultFile).fullPath());
0104       const auto& p0detIds = reader.getAllDetIds();
0105       std::vector<uint32_t> ownDetIds;
0106 
0107       std::transform(fractions.begin(),
0108                      fractions.end(),
0109                      std::back_inserter(ownDetIds),
0110                      [](std::pair<uint32_t, packedBadRocFraction> d) -> uint32_t { return d.first; });
0111 
0112       for (const auto& det : ownDetIds) {
0113         // if found at least one phase-0 detId early return
0114         if (std::find(p0detIds.begin(), p0detIds.end(), det) != p0detIds.end()) {
0115           return true;
0116         }
0117       }
0118       return false;
0119     }
0120 
0121     //_________________________________________________
0122     double getMatchingGeomFactor(const DetId& detid,
0123                                  const std::map<unsigned int, double>& map_geomfactor,
0124                                  const std::vector<uint32_t>& detIdmasks) {
0125       double geomfactor_db = 1;
0126       for (auto map_element : map_geomfactor) {
0127         const DetId mapid = DetId(map_element.first);
0128         if (mapid.subdetId() != detid.subdetId())
0129           continue;
0130         size_t __i = 0;
0131         for (; __i < detIdmasks.size(); __i++) {
0132           DetId maskid = DetId(detIdmasks.at(__i));
0133           if (maskid.subdetId() != mapid.subdetId())
0134             continue;
0135           if ((detid.rawId() & maskid.rawId()) != (mapid.rawId() & maskid.rawId()) &&
0136               (mapid.rawId() & maskid.rawId()) != DetId(mapid.det(), mapid.subdetId()).rawId())
0137             break;
0138         }
0139         if (__i != detIdmasks.size())
0140           continue;
0141         geomfactor_db *= map_element.second;
0142       }
0143       return geomfactor_db;
0144     }
0145 
0146     //_________________________________________________
0147     std::vector<double> getMatchingPUFactors(const DetId& detid,
0148                                              const std::map<unsigned int, std::vector<double> >& map_pufactory,
0149                                              const std::vector<uint32_t>& detIdmasks) {
0150       std::vector<double> pufactors_db;
0151       for (const auto& map_element : map_pufactory) {
0152         const DetId mapid = DetId(map_element.first);
0153         if (mapid.subdetId() != detid.subdetId())
0154           continue;
0155         size_t __i = 0;
0156         for (; __i < detIdmasks.size(); __i++) {
0157           DetId maskid = DetId(detIdmasks.at(__i));
0158           if (maskid.subdetId() != mapid.subdetId())
0159             continue;
0160           if ((detid.rawId() & maskid.rawId()) != (mapid.rawId() & maskid.rawId()) &&
0161               (mapid.rawId() & maskid.rawId()) != DetId(mapid.det(), mapid.subdetId()).rawId())
0162             break;
0163         }
0164         if (__i != detIdmasks.size())
0165           continue;
0166         pufactors_db = map_element.second;
0167       }
0168       return pufactors_db;
0169     }
0170 
0171     //(Not used for the moment)
0172     //_________________________________________________
0173     [[maybe_unused]] bool matches(const DetId& detid, const DetId& db_id, const std::vector<uint32_t>& DetIdmasks) {
0174       if (detid.subdetId() != db_id.subdetId())
0175         return false;
0176       for (size_t i = 0; i < DetIdmasks.size(); ++i) {
0177         DetId maskid = DetId(DetIdmasks.at(i));
0178         if (maskid.subdetId() != db_id.subdetId())
0179           continue;
0180         if ((detid.rawId() & maskid.rawId()) != (db_id.rawId() & maskid.rawId()) &&
0181             (db_id.rawId() & maskid.rawId()) != DetId(db_id.det(), db_id.subdetId()).rawId())
0182           return false;
0183       }
0184       return true;
0185     }
0186 
0187     //_________________________________________________
0188     bool checkPhase(const SiPixelPI::phase phase, const std::vector<uint32_t>& masks_db) {
0189       const char* inputFile;
0190       switch (phase) {
0191         case SiPixelPI::phase::zero:
0192           inputFile = "Geometry/TrackerCommonData/data/trackerParameters.xml";
0193           break;
0194         case SiPixelPI::phase::one:
0195           inputFile = "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
0196           break;
0197         case SiPixelPI::phase::two:
0198           inputFile = "Geometry/TrackerCommonData/data/PhaseII/trackerParameters.xml";
0199           break;
0200         default:
0201           throw cms::Exception("SiPixelDynamicInefficiency_PayloadInspector") << "checkPhase: unrecongnized phase!";
0202       }
0203 
0204       // create the standalone tracker topology
0205       const auto& tkTopo =
0206           StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(inputFile).fullPath());
0207 
0208       // Check what masks we would get using the current geometry
0209       // It has to match what is in the db content!!
0210 
0211       std::vector<uint32_t> masks_geom;
0212       uint32_t max = std::numeric_limits<uint32_t>::max();
0213 
0214       masks_geom.push_back(tkTopo.pxbDetId(max, 0, 0).rawId());
0215       masks_geom.push_back(tkTopo.pxbDetId(0, max, 0).rawId());
0216       masks_geom.push_back(tkTopo.pxbDetId(0, 0, max).rawId());
0217       masks_geom.push_back(tkTopo.pxfDetId(max, 0, 0, 0, 0).rawId());
0218       masks_geom.push_back(tkTopo.pxfDetId(0, max, 0, 0, 0).rawId());
0219       masks_geom.push_back(tkTopo.pxfDetId(0, 0, max, 0, 0).rawId());
0220       masks_geom.push_back(tkTopo.pxfDetId(0, 0, 0, max, 0).rawId());
0221       masks_geom.push_back(tkTopo.pxfDetId(0, 0, 0, 0, max).rawId());
0222 
0223       return (masks_geom.size() == masks_db.size() &&
0224               std::equal(masks_geom.begin(), masks_geom.end(), masks_db.begin()));
0225     }
0226 
0227     //_________________________________________________
0228     unsigned int maxDepthOfPUArray(const std::map<unsigned int, std::vector<double> >& map_pufactor) {
0229       unsigned int size{0};
0230       for (const auto& [id, vec] : map_pufactor) {
0231         if (vec.size() > size)
0232           size = vec.size();
0233       }
0234       return size;
0235     }
0236 
0237     //_________________________________________________
0238     std::pair<int, int> getClosestFactors(int input) {
0239       if ((input % 2 != 0) && input > 1) {
0240         input += 1;
0241       }
0242 
0243       int testNum = (int)sqrt(input);
0244       while (input % testNum != 0) {
0245         testNum--;
0246       }
0247       return std::make_pair(testNum, input / testNum);
0248     }
0249 
0250     /** Given an input list of std::string this
0251      *  finds the longest common substring
0252      */
0253     std::string findStem(const std::vector<std::string>& arr) {
0254       // Determine size of the array
0255       int n = arr.size();
0256 
0257       // Take first word from array as reference
0258       const std::string& s = arr[0];
0259       int len = s.length();
0260 
0261       std::string res = "";
0262 
0263       for (int i = 0; i < len; i++) {
0264         for (int j = i + 1; j <= len; j++) {
0265           // generating all possible substrings
0266           // of our reference string arr[0] i.e s
0267           std::string stem = s.substr(i, j);
0268           int k = 1;
0269           for (k = 1; k < n; k++) {
0270             // Check if the generated stem is
0271             // common to all words
0272             if (arr[k].find(stem) == std::string::npos)
0273               break;
0274           }
0275 
0276           // If current substring is present in
0277           // all strings and its length is greater
0278           // than current result
0279           if (k == n && res.length() < stem.length())
0280             res = stem;
0281         }
0282       }
0283       return res;
0284     }
0285 
0286     /** This function determines from the list of attached DetId which
0287      *  SiPixelPI::region represents them
0288      */
0289     std::string attachLocationLabel(const std::vector<uint32_t>& listOfDetIds, SiPixelPI::PhaseInfo& phInfo) {
0290       // collect all the regions in which it can be split
0291       std::vector<std::string> regions;
0292       for (const auto& rawId : listOfDetIds) {
0293         SiPixelPI::topolInfo t_info_fromXML;
0294         t_info_fromXML.init();
0295         DetId detid(rawId);
0296 
0297         const char* path_toTopologyXML = phInfo.pathToTopoXML();
0298         auto tTopo =
0299             StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
0300         t_info_fromXML.fillGeometryInfo(detid, tTopo, phInfo.phase());
0301         const auto& reg = SiPixelPI::getStringFromRegionEnum(t_info_fromXML.filterThePartition());
0302         if (!std::count(regions.begin(), regions.end(), reg)) {
0303           regions.push_back(reg);
0304         }
0305       }
0306 
0307       std::string retVal = "";
0308       // if perfect match (only one category)
0309       if (regions.size() == 1) {
0310         retVal = regions.front();
0311       } else {
0312         retVal = findStem(regions);
0313       }
0314 
0315       // if the last char is "/" strip it from the string
0316       if (retVal.back() == '/')
0317         retVal.pop_back();
0318 
0319       return retVal;
0320     }
0321 
0322     void fillParametrizations(std::vector<std::vector<double> >& listOfParametrizations,
0323                               std::vector<TF1*>& parametrizations,
0324                               std::vector<std::string>& formulas,
0325                               const std::vector<std::string>& namesOfParts) {
0326       static constexpr double xmin_ = 0.;   // x10e34 (min inst. lumi)
0327       static constexpr double xmax_ = 25.;  // x10e34 (max inst. lumi)
0328 
0329       // functional for polynomial of n-th degree
0330       auto func = [](double* x, double* p) {
0331         int n = p[0];
0332         double* params = p + 1;
0333         ROOT::Math::Polynomial pol(n);
0334         return pol(x, params);
0335       };
0336 
0337       int index{0};
0338       for (auto& params : listOfParametrizations) {
0339         index++;
0340         int n = params.size();
0341         int npar = n + 2;
0342         std::string str{namesOfParts[index - 1]};
0343         if (str.length() >= 2 && str.substr(str.length() - 2) == "/i") {
0344           str += "nner";
0345         } else if (str.length() >= 2 && str.substr(str.length() - 2) == "/o") {
0346           str += "uter";
0347         }
0348 
0349         TF1* f1 = new TF1((fmt::sprintf("region: #bf{%s}", str)).c_str(), func, xmin_, xmax_, npar);
0350 
0351         // push polynomial degree as first entry in the vector
0352         params.insert(params.begin(), n);
0353         // TF1::SetParameters needs a C-style array
0354         double* arr = params.data();
0355         f1->SetLineWidth(2);
0356 
0357         // fill in the parameter
0358         for (unsigned int j = 0; j < params.size(); j++) {
0359           f1->SetParameter(j, arr[j]);
0360         }
0361 
0362         parametrizations.push_back(f1);
0363 
0364         // build the formula to be displayed
0365         std::string formula;
0366         edm::LogVerbatim("fillParametrizations") << "index: " << index;
0367         for (unsigned int i = 1; i < params.size(); i++) {
0368           edm::LogVerbatim("fillParametrizations") << " " << params[i];
0369           formula += fmt::sprintf("%s%fx^{%i}", (i == 1 ? "" : (std::signbit(params[i]) ? "" : "+")), params[i], i - 1);
0370         }
0371         edm::LogVerbatim("fillParametrizations") << std::endl;
0372         formulas.push_back(formula);
0373       }
0374     }
0375   }  // namespace SiPixDynIneff
0376 
0377   /************************************************
0378     test class
0379   *************************************************/
0380 
0381   class SiPixelDynamicInefficiencyTest : public Histogram1D<SiPixelDynamicInefficiency, SINGLE_IOV> {
0382   public:
0383     SiPixelDynamicInefficiencyTest()
0384         : Histogram1D<SiPixelDynamicInefficiency, SINGLE_IOV>(
0385               "SiPixelDynamicInefficiency test", "SiPixelDynamicInefficiency test", 1, 0.0, 1.0) {}
0386 
0387     bool fill() override {
0388       auto tag = PlotBase::getTag<0>();
0389       for (auto const& iov : tag.iovs) {
0390         std::shared_ptr<SiPixelDynamicInefficiency> payload = Base::fetchPayload(std::get<1>(iov));
0391         if (payload.get()) {
0392           fillWithValue(1.);
0393 
0394           std::map<unsigned int, double> map_pixelgeomfactor = payload->getPixelGeomFactors();
0395           std::map<unsigned int, double> map_colgeomfactor = payload->getColGeomFactors();
0396           std::map<unsigned int, double> map_chipgeomfactor = payload->getChipGeomFactors();
0397           std::map<unsigned int, std::vector<double> > map_pufactor = payload->getPUFactors();
0398           std::vector<uint32_t> detIdmasks_db = payload->getDetIdmasks();
0399           double theInstLumiScaleFactor_db = payload->gettheInstLumiScaleFactor_();
0400 
0401           edm::LogPrint("SiPixelDynamicInefficiencyTest") << "-------------------------------------------------------";
0402           edm::LogPrint("SiPixelDynamicInefficiencyTest") << "Printing out DB content:\n";
0403 
0404           edm::LogPrint("SiPixelDynamicInefficiencyTest") << "  PixelGeomFactors:";
0405           for (auto pixel : map_pixelgeomfactor)
0406             edm::LogPrint("SiPixelDynamicInefficiencyTest")
0407                 << "    MapID = " << pixel.first << "\tFactor = " << pixel.second;
0408           edm::LogPrint("SiPixelDynamicInefficiencyTest");
0409 
0410           edm::LogPrint("SiPixelDynamicInefficiencyTest") << "  ColGeomFactors:";
0411           for (auto col : map_colgeomfactor)
0412             edm::LogPrint("SiPixelDynamicInefficiencyTest")
0413                 << "    MapID = " << col.first << "\tFactor = " << col.second;
0414           edm::LogPrint("SiPixelDynamicInefficiencyTest");
0415 
0416           edm::LogPrint("SiPixelDynamicInefficiencyTest") << "  ChipGeomFactors:";
0417           for (auto chip : map_chipgeomfactor)
0418             edm::LogPrint("SiPixelDynamicInefficiencyTest")
0419                 << "    MapID = " << chip.first << "\tFactor = " << chip.second;
0420           edm::LogPrint("SiPixelDynamicInefficiencyTest");
0421 
0422           edm::LogPrint("SiPixelDynamicInefficiencyTest") << "  PUFactors:";
0423           for (auto pu : map_pufactor) {
0424             edm::LogPrint("SiPixelDynamicInefficiencyTest")
0425                 << "    MapID = " << pu.first << "\t Factor" << (pu.second.size() > 1 ? "s" : "") << " = ";
0426             for (size_t i = 0, n = pu.second.size(); i < n; ++i)
0427               edm::LogPrint("SiPixelDynamicInefficiencyTest") << pu.second[i] << ((i == n - 1) ? "\n" : ", ");
0428           }
0429           edm::LogPrint("SiPixelDynamicInefficiencyTest");
0430 
0431           edm::LogPrint("SiPixelDynamicInefficiencyTest") << "  DetIdmasks:";
0432           for (auto mask : detIdmasks_db)
0433             edm::LogPrint("SiPixelDynamicInefficiencyTest") << "    MaskID = " << mask;
0434           edm::LogPrint("SiPixelDynamicInefficiencyTest");
0435 
0436           edm::LogPrint("SiPixelDynamicInefficiencyTest") << "  theInstLumiScaleFactor = " << theInstLumiScaleFactor_db;
0437 
0438         }  // payload
0439       }    // iovs
0440       return true;
0441     }  // fill
0442   };
0443 
0444   /************************************************
0445    occupancy style map whole Pixel of inefficient ROCs
0446   *************************************************/
0447   template <SiPixelPI::DetType myType>
0448   class SiPixelIneffROCfromDynIneffMap : public PlotImage<SiPixelDynamicInefficiency, SINGLE_IOV> {
0449   public:
0450     SiPixelIneffROCfromDynIneffMap()
0451         : PlotImage<SiPixelDynamicInefficiency, SINGLE_IOV>("SiPixel Inefficient ROC from Dyn Ineff Pixel Map"),
0452           m_trackerTopo{StandaloneTrackerTopology::fromTrackerParametersXMLFile(
0453               edm::FileInPath("Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml").fullPath())} {}
0454 
0455     bool fill() override {
0456       auto tag = PlotBase::getTag<0>();
0457       auto iov = tag.iovs.front();
0458       auto tagname = tag.name;
0459       std::shared_ptr<SiPixelDynamicInefficiency> payload = fetchPayload(std::get<1>(iov));
0460 
0461       const auto fr = SiPixDynIneff::pbrf(payload);
0462 
0463       if (SiPixDynIneff::isPhase0(fr)) {
0464         edm::LogError("SiPixelDynamicInefficiency_PayloadInspector")
0465             << "SiPixelIneffROCfromDynIneff maps are not supported for non-Phase1 Pixel geometries !";
0466         TCanvas canvas("Canv", "Canv", 1200, 1000);
0467         SiPixelPI::displayNotSupported(canvas, 0);
0468         std::string fileName(m_imageFileName);
0469         canvas.SaveAs(fileName.c_str());
0470         return false;
0471       }
0472 
0473       Phase1PixelROCMaps theMap("", "bad pixel fraction in ROC [%]");
0474 
0475       for (const auto& element : fr) {
0476         auto rawid = element.first;
0477         int subid = DetId(rawid).subdetId();
0478         auto packedinfo = element.second;
0479         auto badRocs = packedinfo.badRocNumber;
0480         auto badRocsF = packedinfo.badRocFrac;
0481 
0482         for (size_t i = 0; i < badRocs.size(); i++) {
0483           std::bitset<16> rocToMark;
0484           rocToMark.set(badRocs[i]);
0485           if ((subid == PixelSubdetector::PixelBarrel && myType == SiPixelPI::t_barrel) ||
0486               (subid == PixelSubdetector::PixelEndcap && myType == SiPixelPI::t_forward) ||
0487               (myType == SiPixelPI::t_all)) {
0488             theMap.fillSelectedRocs(rawid, rocToMark, badRocsF[i] * 100.f);
0489           }
0490         }
0491       }
0492 
0493       gStyle->SetOptStat(0);
0494       //=========================
0495       TCanvas canvas("Summary", "Summary", 1200, k_height[myType]);
0496       canvas.cd();
0497 
0498       auto unpacked = SiPixelPI::unpack(std::get<0>(iov));
0499 
0500       std::string IOVstring = (unpacked.first == 0)
0501                                   ? std::to_string(unpacked.second)
0502                                   : (std::to_string(unpacked.first) + "," + std::to_string(unpacked.second));
0503 
0504       const auto headerText = fmt::sprintf("#color[4]{%s},  IOV: #color[4]{%s}", tagname, IOVstring);
0505 
0506       switch (myType) {
0507         case SiPixelPI::t_barrel:
0508           theMap.drawBarrelMaps(canvas, headerText);
0509           break;
0510         case SiPixelPI::t_forward:
0511           theMap.drawForwardMaps(canvas, headerText);
0512           break;
0513         case SiPixelPI::t_all:
0514           theMap.drawMaps(canvas, headerText);
0515           break;
0516         default:
0517           throw cms::Exception("SiPixelIneffROCfromDynIneffMap") << "\nERROR: unrecognized Pixel Detector part ";
0518       }
0519 
0520       std::string fileName(m_imageFileName);
0521       canvas.SaveAs(fileName.c_str());
0522 #ifdef MMDEBUG
0523       canvas.SaveAs("outAll.root");
0524 #endif
0525 
0526       return true;
0527     }
0528 
0529   private:
0530     static constexpr std::array<int, 3> k_height = {{1200, 600, 1600}};
0531     TrackerTopology m_trackerTopo;
0532   };
0533 
0534   using SiPixelBPixIneffROCfromDynIneffMap = SiPixelIneffROCfromDynIneffMap<SiPixelPI::t_barrel>;
0535   using SiPixelFPixIneffROCfromDynIneffMap = SiPixelIneffROCfromDynIneffMap<SiPixelPI::t_forward>;
0536   using SiPixelFullIneffROCfromDynIneffMap = SiPixelIneffROCfromDynIneffMap<SiPixelPI::t_all>;
0537 
0538   /************************************************
0539    occupancy style map whole Pixel, difference of payloads
0540   *************************************************/
0541   template <SiPixelPI::DetType myType, IOVMultiplicity nIOVs, int ntags>
0542   class SiPixelIneffROCComparisonBase : public PlotImage<SiPixelDynamicInefficiency, nIOVs, ntags> {
0543   public:
0544     SiPixelIneffROCComparisonBase()
0545         : PlotImage<SiPixelDynamicInefficiency, nIOVs, ntags>(
0546               Form("SiPixelDynamicInefficiency %s Pixel Map", SiPixelPI::DetNames[myType].c_str())),
0547           m_trackerTopo{StandaloneTrackerTopology::fromTrackerParametersXMLFile(
0548               edm::FileInPath("Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml").fullPath())} {}
0549 
0550     bool fill() override {
0551       // trick to deal with the multi-ioved tag and two tag case at the same time
0552       auto theIOVs = PlotBase::getTag<0>().iovs;
0553       auto f_tagname = PlotBase::getTag<0>().name;
0554       std::string l_tagname = "";
0555       auto firstiov = theIOVs.front();
0556       std::tuple<cond::Time_t, cond::Hash> lastiov;
0557 
0558       // we don't support (yet) comparison with more than 2 tags
0559       assert(this->m_plotAnnotations.ntags < 3);
0560 
0561       if (this->m_plotAnnotations.ntags == 2) {
0562         auto tag2iovs = PlotBase::getTag<1>().iovs;
0563         l_tagname = PlotBase::getTag<1>().name;
0564         lastiov = tag2iovs.front();
0565       } else {
0566         lastiov = theIOVs.back();
0567       }
0568 
0569       std::shared_ptr<SiPixelDynamicInefficiency> last_payload = this->fetchPayload(std::get<1>(lastiov));
0570       std::shared_ptr<SiPixelDynamicInefficiency> first_payload = this->fetchPayload(std::get<1>(firstiov));
0571 
0572       const auto fp = SiPixDynIneff::pbrf(last_payload);
0573       const auto lp = SiPixDynIneff::pbrf(first_payload);
0574 
0575       if (SiPixDynIneff::isPhase0(fp) || SiPixDynIneff::isPhase0(lp)) {
0576         edm::LogError("SiPixelDynamicInefficiency_PayloadInspector")
0577             << "SiPixelDynamicInefficiency comparison maps are not supported for non-Phase1 Pixel geometries !";
0578         TCanvas canvas("Canv", "Canv", 1200, 1000);
0579         SiPixelPI::displayNotSupported(canvas, 0);
0580         std::string fileName(this->m_imageFileName);
0581         canvas.SaveAs(fileName.c_str());
0582         return false;
0583       }
0584 
0585       Phase1PixelROCMaps theMap("", "#Delta payload A - payload B");
0586 
0587       gStyle->SetOptStat(0);
0588       //=========================
0589       TCanvas canvas("Summary", "Summary", 1200, k_height[myType]);
0590       canvas.cd();
0591 
0592       auto f_unpacked = SiPixelPI::unpack(std::get<0>(firstiov));
0593       auto l_unpacked = SiPixelPI::unpack(std::get<0>(lastiov));
0594 
0595       std::string f_IOVstring = (f_unpacked.first == 0)
0596                                     ? std::to_string(f_unpacked.second)
0597                                     : (std::to_string(f_unpacked.first) + "," + std::to_string(f_unpacked.second));
0598 
0599       std::string l_IOVstring = (l_unpacked.first == 0)
0600                                     ? std::to_string(l_unpacked.second)
0601                                     : (std::to_string(l_unpacked.first) + "," + std::to_string(l_unpacked.second));
0602 
0603       std::string headerText;
0604 
0605       if (this->m_plotAnnotations.ntags == 2) {
0606         headerText =
0607             fmt::sprintf("#color[2]{A: %s, %s} - #color[4]{B: %s, %s}", f_tagname, f_IOVstring, l_tagname, l_IOVstring);
0608       } else {
0609         headerText = fmt::sprintf("%s,IOV #color[2]{A: %s} - #color[4]{B: %s} ", f_tagname, f_IOVstring, l_IOVstring);
0610       }
0611 
0612       switch (myType) {
0613         case SiPixelPI::t_barrel:
0614           theMap.drawBarrelMaps(canvas, headerText);
0615           break;
0616         case SiPixelPI::t_forward:
0617           theMap.drawForwardMaps(canvas, headerText);
0618           break;
0619         case SiPixelPI::t_all:
0620           theMap.drawMaps(canvas, headerText);
0621           break;
0622         default:
0623           throw cms::Exception("SiPixelDynamicInefficiencyMapComparison")
0624               << "\nERROR: unrecognized Pixel Detector part ";
0625       }
0626 
0627       // first loop on the first payload (newest)
0628       fillTheMapFromPayload(theMap, fp, false);
0629 
0630       // then loop on the second payload (oldest)
0631       fillTheMapFromPayload(theMap, lp, true);  // true will subtract
0632 
0633       std::string fileName(this->m_imageFileName);
0634       canvas.SaveAs(fileName.c_str());
0635 #ifdef MMDEBUG
0636       canvas.SaveAs("outAll.root");
0637 #endif
0638 
0639       return true;
0640     }
0641 
0642   private:
0643     static constexpr std::array<int, 3> k_height = {{1200, 600, 1600}};
0644     TrackerTopology m_trackerTopo;
0645 
0646     //____________________________________________________________________________________________
0647     void fillTheMapFromPayload(Phase1PixelROCMaps& theMap, const SiPixDynIneff::BRFractions& fr, bool subtract) {
0648       for (const auto& element : fr) {
0649         auto rawid = element.first;
0650         int subid = DetId(rawid).subdetId();
0651         auto packedinfo = element.second;
0652         auto badRocs = packedinfo.badRocNumber;
0653         auto badRocsF = packedinfo.badRocFrac;
0654 
0655         for (size_t i = 0; i < badRocs.size(); i++) {
0656           std::bitset<16> rocToMark;
0657           rocToMark.set(badRocs[i]);
0658           if ((subid == PixelSubdetector::PixelBarrel && myType == SiPixelPI::t_barrel) ||
0659               (subid == PixelSubdetector::PixelEndcap && myType == SiPixelPI::t_forward) ||
0660               (myType == SiPixelPI::t_all)) {
0661             theMap.fillSelectedRocs(rawid, rocToMark, badRocsF[i] * (subtract ? -1. : 1.));
0662           }
0663         }
0664       }
0665     }
0666   };
0667 
0668   /*
0669     These are not implemented for the time being, since the SiPixelDynamicInefficiency is a condition
0670     used only in simulation, hence there is no such thing as a multi-IoV Dynamic Inefficiency tag 
0671   */
0672 
0673   using SiPixelBPixIneffROCsMapCompareSingleTag = SiPixelIneffROCComparisonBase<SiPixelPI::t_barrel, MULTI_IOV, 1>;
0674   using SiPixelFPixIneffROCsMapCompareSingleTag = SiPixelIneffROCComparisonBase<SiPixelPI::t_forward, MULTI_IOV, 1>;
0675   using SiPixelFullIneffROCsMapCompareSingleTag = SiPixelIneffROCComparisonBase<SiPixelPI::t_all, MULTI_IOV, 1>;
0676 
0677   using SiPixelBPixIneffROCsMapCompareTwoTags = SiPixelIneffROCComparisonBase<SiPixelPI::t_barrel, SINGLE_IOV, 2>;
0678   using SiPixelFPixIneffROCsMapCompareTwoTags = SiPixelIneffROCComparisonBase<SiPixelPI::t_forward, SINGLE_IOV, 2>;
0679   using SiPixelFullIneffROCsMapCompareTwoTags = SiPixelIneffROCComparisonBase<SiPixelPI::t_all, SINGLE_IOV, 2>;
0680 
0681   /************************************************
0682    Full Pixel Tracker Map class (for geometrical factors)
0683   *************************************************/
0684   template <SiPixDynIneff::factor theFactor>
0685   class SiPixelDynamicInefficiencyFullPixelMap : public PlotImage<SiPixelDynamicInefficiency, SINGLE_IOV> {
0686   public:
0687     SiPixelDynamicInefficiencyFullPixelMap()
0688         : PlotImage<SiPixelDynamicInefficiency, SINGLE_IOV>("SiPixelDynamicInefficiency Map") {
0689       label_ = "SiPixelDynamicInefficiencyFullPixelMap";
0690       payloadString = fmt::sprintf("%s Dynamic Inefficiency", SiPixDynIneff::factorString[theFactor]);
0691     }
0692 
0693     bool fill() override {
0694       gStyle->SetPalette(1);
0695       auto tag = PlotBase::getTag<0>();
0696       auto iov = tag.iovs.front();
0697       std::shared_ptr<SiPixelDynamicInefficiency> payload = this->fetchPayload(std::get<1>(iov));
0698 
0699       if (payload.get()) {
0700         Phase1PixelSummaryMap fullMap("", fmt::sprintf("%s", payloadString), fmt::sprintf("%s", payloadString));
0701         fullMap.createTrackerBaseMap();
0702 
0703         SiPixDynIneff::FactorMap theMap{};
0704         switch (theFactor) {
0705           case SiPixDynIneff::geom:
0706             theMap = payload->getPixelGeomFactors();
0707             break;
0708           case SiPixDynIneff::colgeom:
0709             theMap = payload->getColGeomFactors();
0710             break;
0711           case SiPixDynIneff::chipgeom:
0712             theMap = payload->getChipGeomFactors();
0713             break;
0714           default:
0715             throw cms::Exception(label_) << "\nERROR: unrecognized type of geometry factor ";
0716         }
0717 
0718         std::vector<uint32_t> detIdmasks_db = payload->getDetIdmasks();
0719 
0720         if (!SiPixDynIneff::checkPhase(SiPixelPI::phase::one, detIdmasks_db)) {
0721           edm::LogError(label_) << label_ << " maps are not supported for non-Phase1 Pixel geometries !";
0722           TCanvas canvas("Canv", "Canv", 1200, 1000);
0723           SiPixelPI::displayNotSupported(canvas, 0);
0724           std::string fileName(m_imageFileName);
0725           canvas.SaveAs(fileName.c_str());
0726           return false;
0727         }
0728 
0729         SiPixelDetInfoFileReader reader =
0730             SiPixelDetInfoFileReader(edm::FileInPath(SiPixelDetInfoFileReader::kPh1DefaultFile).fullPath());
0731         const auto& p1detIds = reader.getAllDetIds();
0732 
0733         for (const auto& det : p1detIds) {
0734           const auto& value = SiPixDynIneff::getMatchingGeomFactor(det, theMap, detIdmasks_db);
0735           fullMap.fillTrackerMap(det, value);
0736         }
0737 
0738         const auto& range = fullMap.getZAxisRange();
0739         if (range.first == range.second) {
0740           // in case the map is completely filled with one value;
0741           // set the z-axis to be meaningful
0742           fullMap.setZAxisRange(range.first - 0.01, range.second + 0.01);
0743         }
0744 
0745         TCanvas canvas("Canv", "Canv", 3000, 2000);
0746         fullMap.printTrackerMap(canvas);
0747 
0748         auto ltx = TLatex();
0749         ltx.SetTextFont(62);
0750         ltx.SetTextSize(0.025);
0751         ltx.SetTextAlign(11);
0752         ltx.DrawLatexNDC(
0753             gPad->GetLeftMargin() + 0.01,
0754             gPad->GetBottomMargin() + 0.01,
0755             ("#color[4]{" + tag.name + "}, IOV: #color[4]{" + std::to_string(std::get<0>(iov)) + "}").c_str());
0756 
0757         std::string fileName(this->m_imageFileName);
0758         canvas.SaveAs(fileName.c_str());
0759       }
0760       return true;
0761     }
0762 
0763   protected:
0764     std::string payloadString;
0765     std::string label_;
0766   };
0767 
0768   using SiPixelDynamicInefficiencyGeomFactorMap = SiPixelDynamicInefficiencyFullPixelMap<SiPixDynIneff::geom>;
0769   using SiPixelDynamicInefficiencyColGeomFactorMap = SiPixelDynamicInefficiencyFullPixelMap<SiPixDynIneff::colgeom>;
0770   using SiPixelDynamicInefficiencyChipGeomFactorMap = SiPixelDynamicInefficiencyFullPixelMap<SiPixDynIneff::chipgeom>;
0771 
0772   /************************************************
0773    Full Pixel Tracker Map class (for PU factors)
0774   *************************************************/
0775   class SiPixelDynamicInefficiencyPUPixelMaps : public PlotImage<SiPixelDynamicInefficiency, SINGLE_IOV> {
0776   public:
0777     SiPixelDynamicInefficiencyPUPixelMaps()
0778         : PlotImage<SiPixelDynamicInefficiency, SINGLE_IOV>("SiPixelDynamicInefficiency Map") {
0779       label_ = "SiPixelDynamicInefficiencyFullPixelMap";
0780       payloadString = fmt::sprintf("%s Dynamic Inefficiency", SiPixDynIneff::factorString[SiPixDynIneff::pu]);
0781     }
0782 
0783     bool fill() override {
0784       gStyle->SetPalette(1);
0785       auto tag = PlotBase::getTag<0>();
0786       auto iov = tag.iovs.front();
0787       std::shared_ptr<SiPixelDynamicInefficiency> payload = this->fetchPayload(std::get<1>(iov));
0788 
0789       if (payload.get()) {
0790         std::vector<Phase1PixelSummaryMap> maps;
0791 
0792         SiPixDynIneff::PUFactorMap theMap = payload->getPUFactors();
0793         std::vector<uint32_t> detIdmasks_db = payload->getDetIdmasks();
0794 
0795         if (!SiPixDynIneff::checkPhase(SiPixelPI::phase::one, detIdmasks_db)) {
0796           edm::LogError(label_) << label_ << " maps are not supported for non-Phase1 Pixel geometries !";
0797           TCanvas canvas("Canv", "Canv", 1200, 1000);
0798           SiPixelPI::displayNotSupported(canvas, 0);
0799           std::string fileName(m_imageFileName);
0800           canvas.SaveAs(fileName.c_str());
0801           return false;
0802         }
0803 
0804         unsigned int depth = SiPixDynIneff::maxDepthOfPUArray(theMap);
0805 
0806         // create the maps
0807         for (unsigned int i = 0; i < depth; i++) {
0808           maps.emplace_back(
0809               "", fmt::sprintf("%s, factor %i", payloadString, i), fmt::sprintf("%s, factor %i", payloadString, i));
0810           maps[i].createTrackerBaseMap();
0811         }
0812 
0813         // retrieve the list of phase1 detids
0814         const auto& reader =
0815             SiPixelDetInfoFileReader(edm::FileInPath(SiPixelDetInfoFileReader::kPh1DefaultFile).fullPath());
0816         const auto& p1detIds = reader.getAllDetIds();
0817 
0818         // fill the maps
0819         for (const auto& det : p1detIds) {
0820           const auto& values = SiPixDynIneff::getMatchingPUFactors(det, theMap, detIdmasks_db);
0821           int index = 0;
0822           for (const auto& value : values) {
0823             maps[index].fillTrackerMap(det, value);
0824             index++;
0825           }
0826         }
0827 
0828         // in case the map is completely filled with one value;
0829         // set the z-axis to be meaningful
0830         for (unsigned int i = 0; i < depth; i++) {
0831           const auto& range = maps[i].getZAxisRange();
0832           if (range.first == range.second) {
0833             maps[i].setZAxisRange(range.first - 0.01, range.second + 0.01);
0834           }
0835         }
0836 
0837         // determine how the plot will be paginated
0838         auto sides = SiPixDynIneff::getClosestFactors(depth);
0839         TCanvas canvas("Canv", "Canv", sides.second * 900, sides.first * 600);
0840         canvas.Divide(sides.second, sides.first);
0841 
0842         // print the sub-canvases
0843         for (unsigned int i = 0; i < depth; i++) {
0844           maps[i].printTrackerMap(canvas, 0.035, i + 1);
0845           auto ltx = TLatex();
0846           ltx.SetTextFont(62);
0847           ltx.SetTextSize(0.025);
0848           ltx.SetTextAlign(11);
0849           ltx.DrawLatexNDC(
0850               gPad->GetLeftMargin() + 0.01,
0851               gPad->GetBottomMargin() + 0.01,
0852               ("#color[4]{" + tag.name + "}, IOV: #color[4]{" + std::to_string(std::get<0>(iov)) + "}").c_str());
0853         }
0854 
0855         std::string fileName(this->m_imageFileName);
0856         canvas.SaveAs(fileName.c_str());
0857       }
0858       return true;
0859     }
0860 
0861   protected:
0862     std::string payloadString;
0863     std::string label_;
0864   };
0865 
0866   /************************************************
0867    Per sector plot of the inefficiency parameterization vs inst lumi (for PU factors)
0868   *************************************************/
0869   class SiPixelDynamicInefficiencyPUParametrization : public PlotImage<SiPixelDynamicInefficiency, SINGLE_IOV> {
0870   public:
0871     SiPixelDynamicInefficiencyPUParametrization()
0872         : PlotImage<SiPixelDynamicInefficiency, SINGLE_IOV>("SiPixelDynamicInefficiency PU parametrization") {
0873       label_ = "SiPixelDynamicInefficiencyPUParameterization";
0874       payloadString =
0875           fmt::sprintf("%s Dynamic Inefficiency parametrization", SiPixDynIneff::factorString[SiPixDynIneff::pu]);
0876     }
0877 
0878     bool fill() override {
0879       gStyle->SetPalette(1);
0880       auto tag = PlotBase::getTag<0>();
0881       auto iov = tag.iovs.front();
0882       std::shared_ptr<SiPixelDynamicInefficiency> payload = this->fetchPayload(std::get<1>(iov));
0883 
0884       if (payload.get()) {
0885         SiPixDynIneff::PUFactorMap theMap = payload->getPUFactors();
0886         std::vector<uint32_t> detIdmasks_db = payload->getDetIdmasks();
0887 
0888         if (!SiPixDynIneff::checkPhase(SiPixelPI::phase::one, detIdmasks_db)) {
0889           edm::LogError(label_) << label_ << " maps are not supported for non-Phase1 Pixel geometries !";
0890           TCanvas canvas("Canv", "Canv", 1200, 1000);
0891           SiPixelPI::displayNotSupported(canvas, 0);
0892           std::string fileName(m_imageFileName);
0893           canvas.SaveAs(fileName.c_str());
0894           return false;
0895         }
0896 
0897         std::vector<std::vector<double> > listOfParametrizations;
0898 
0899         // retrieve the list of phase1 detids
0900         const auto& reader =
0901             SiPixelDetInfoFileReader(edm::FileInPath(SiPixelDetInfoFileReader::kPh1DefaultFile).fullPath());
0902         auto p1detIds = reader.getAllDetIds();
0903 
0904         // follows hack to get an inner ladder module first
0905         // skip the first 8 dets since they lay on the same ladder
0906         std::swap(p1detIds[0], p1detIds[8]);
0907 
0908         std::map<unsigned int, std::vector<uint32_t> > modules_per_region;
0909 
0910         // fill the maps
0911         for (const auto& det : p1detIds) {
0912           const auto& values = SiPixDynIneff::getMatchingPUFactors(det, theMap, detIdmasks_db);
0913           // find the index in the vector
0914           auto pIndex = std::find(listOfParametrizations.begin(), listOfParametrizations.end(), values);
0915           // if it's not there push it back in the list of parametrizations and then insert in the map
0916           if (pIndex == listOfParametrizations.end()) {
0917             listOfParametrizations.push_back(values);
0918             std::vector<unsigned int> toInsert = {det};
0919             modules_per_region.insert(std::make_pair(listOfParametrizations.size() - 1, toInsert));
0920           } else {
0921             modules_per_region.at(pIndex - listOfParametrizations.begin()).push_back(det);
0922           }
0923         }
0924 
0925         unsigned int depth = listOfParametrizations.size();
0926 
0927         std::vector<std::string> namesOfParts;
0928         namesOfParts.reserve(depth);
0929         // fill in the (ordered) information about regions
0930         for (const auto& [index, modules] : modules_per_region) {
0931           auto PhInfo = SiPixelPI::PhaseInfo(SiPixelPI::phase1size);
0932           const auto& regName = SiPixDynIneff::attachLocationLabel(modules, PhInfo);
0933           namesOfParts.push_back(regName);
0934 
0935           /*
0936        *  The following is needed for internal
0937        *  debug / cross-check
0938        */
0939 
0940           std::stringstream ss;
0941           ss << "region name: " << regName << " has the following modules attached: [";
0942           for (const auto& module : modules) {
0943             ss << module << ", ";
0944           }
0945           ss.seekp(-2, std::ios_base::end);  // remove last two chars
0946           ss << "] "
0947              << " has representation (";
0948           for (const auto& param : listOfParametrizations[index]) {
0949             ss << param << ", ";
0950           }
0951           ss.seekp(-2, std::ios_base::end);  // remove last two chars
0952           ss << ") ";
0953           edm::LogPrint(label_) << ss.str() << "\n\n";
0954           ss.str(std::string()); /* clear the stringstream */
0955         }
0956 
0957         std::vector<TF1*> parametrizations;
0958         std::vector<std::string> formulas;
0959         parametrizations.reserve(depth);
0960         formulas.reserve(depth);
0961 
0962         // fill the parametrization plots
0963         SiPixDynIneff::fillParametrizations(listOfParametrizations, parametrizations, formulas, namesOfParts);
0964 
0965         // determine how the plot will be paginated
0966         auto sides = SiPixDynIneff::getClosestFactors(depth);
0967         TCanvas canvas("Canv", "Canv", sides.second * 900, sides.first * 600);
0968         canvas.Divide(sides.second, sides.first);
0969 
0970         // print the sub-canvases
0971         for (unsigned int i = 0; i < depth; i++) {
0972           canvas.cd(i + 1);
0973           gPad->SetGrid();
0974           SiPixelPI::adjustCanvasMargins(gPad, 0.07, 0.14, 0.14, 0.03);
0975           parametrizations[i]->Draw();
0976 
0977           // set axis
0978           TH1* h = parametrizations[i]->GetHistogram();
0979           TAxis* ax = h->GetXaxis();
0980           ax->SetTitle("Inst. luminosity [10^{33} cm^{-2}s^{-1}]");
0981 
0982           TAxis* ay = h->GetYaxis();
0983           ay->SetRangeUser(0.80, 1.00);  // force the scale
0984           ay->SetTitle("Double Column Efficiency parametrization");
0985 
0986           // beautify
0987           SiPixelPI::makeNicePlotStyle(h);
0988 
0989           auto ltx = TLatex();
0990           ltx.SetTextFont(62);
0991           ltx.SetTextSize(0.045);
0992           ltx.SetTextAlign(11);
0993           ltx.DrawLatexNDC(
0994               gPad->GetLeftMargin() + 0.03,
0995               gPad->GetBottomMargin() + 0.08,
0996               ("#color[2]{" + tag.name + "},IOV: #color[2]{" + std::to_string(std::get<0>(iov)) + "}").c_str());
0997 
0998           auto ltxForm = TLatex();
0999           ltxForm.SetTextFont(62);
1000           ltxForm.SetTextColor(kRed);
1001           ltxForm.SetTextSize(0.035);
1002           ltxForm.SetTextAlign(11);
1003           ltxForm.DrawLatexNDC(gPad->GetLeftMargin() + 0.03, gPad->GetBottomMargin() + 0.04, formulas[i].c_str());
1004 
1005           edm::LogPrint(label_) << namesOfParts[i] << " => " << formulas[i] << std::endl;
1006         }
1007 
1008         std::string fileName(this->m_imageFileName);
1009         canvas.SaveAs(fileName.c_str());
1010 
1011       }  // if payload.get()
1012       return true;
1013     }  // fill()
1014 
1015   protected:
1016     std::string payloadString;
1017     std::string label_;
1018   };
1019 
1020   template <IOVMultiplicity nIOVs, int ntags>
1021   class SiPixelDynamicInefficiencyPUParamComparisonBase : public PlotImage<SiPixelDynamicInefficiency, nIOVs, ntags> {
1022   public:
1023     SiPixelDynamicInefficiencyPUParamComparisonBase()
1024         : PlotImage<SiPixelDynamicInefficiency, nIOVs, ntags>(
1025               Form("SiPixelDynamic Inefficiency parameterization comparison by Region %i tag(s)", ntags)) {
1026       label_ = "SiPixelDynamicInefficiencyPUParameterization";
1027     }
1028 
1029     bool fill() override {
1030       gStyle->SetPalette(1);
1031 
1032       // trick to deal with the multi-ioved tag and two tag case at the same time
1033       auto theIOVs = PlotBase::getTag<0>().iovs;
1034       auto f_tagname = PlotBase::getTag<0>().name;
1035       std::string l_tagname = "";
1036       auto firstiov = theIOVs.front();
1037       std::tuple<cond::Time_t, cond::Hash> lastiov;
1038 
1039       // we don't support (yet) comparison with more than 2 tags
1040       assert(this->m_plotAnnotations.ntags < 3);
1041 
1042       if (this->m_plotAnnotations.ntags == 2) {
1043         auto tag2iovs = PlotBase::getTag<1>().iovs;
1044         l_tagname = PlotBase::getTag<1>().name;
1045         lastiov = tag2iovs.front();
1046       } else {
1047         lastiov = theIOVs.back();
1048       }
1049 
1050       std::shared_ptr<SiPixelDynamicInefficiency> last_payload = this->fetchPayload(std::get<1>(lastiov));
1051       std::shared_ptr<SiPixelDynamicInefficiency> first_payload = this->fetchPayload(std::get<1>(firstiov));
1052 
1053       if (first_payload.get() && last_payload.get()) {
1054         SiPixDynIneff::PUFactorMap f_theMap = first_payload->getPUFactors();
1055         std::vector<uint32_t> f_detIdmasks_db = first_payload->getDetIdmasks();
1056 
1057         SiPixDynIneff::PUFactorMap l_theMap = last_payload->getPUFactors();
1058         std::vector<uint32_t> l_detIdmasks_db = last_payload->getDetIdmasks();
1059 
1060         if (!SiPixDynIneff::checkPhase(SiPixelPI::phase::one, f_detIdmasks_db) ||
1061             !SiPixDynIneff::checkPhase(SiPixelPI::phase::one, l_detIdmasks_db)) {
1062           edm::LogError(label_) << label_ << " maps are not supported for non-Phase1 Pixel geometries !";
1063           TCanvas canvas("Canv", "Canv", 1200, 1000);
1064           SiPixelPI::displayNotSupported(canvas, 0);
1065           std::string fileName(this->m_imageFileName);
1066           canvas.SaveAs(fileName.c_str());
1067           return false;
1068         }
1069 
1070         std::vector<std::vector<double> > f_listOfParametrizations;
1071         std::vector<std::vector<double> > l_listOfParametrizations;
1072 
1073         // retrieve the list of phase1 detids
1074         const auto& reader =
1075             SiPixelDetInfoFileReader(edm::FileInPath(SiPixelDetInfoFileReader::kPh1DefaultFile).fullPath());
1076         auto p1detIds = reader.getAllDetIds();
1077 
1078         // follows hack to get an inner ladder module first
1079         // skip the first 8 dets since they lay on the same ladder
1080         std::swap(p1detIds[0], p1detIds[8]);
1081 
1082         std::map<unsigned int, std::vector<uint32_t> > modules_per_region;
1083 
1084         // fill the maps
1085         for (const auto& det : p1detIds) {
1086           const auto& f_values = SiPixDynIneff::getMatchingPUFactors(det, f_theMap, f_detIdmasks_db);
1087           const auto& l_values = SiPixDynIneff::getMatchingPUFactors(det, l_theMap, l_detIdmasks_db);
1088 
1089           // find the index in the vector
1090           auto fIndex = std::find(f_listOfParametrizations.begin(), f_listOfParametrizations.end(), f_values);
1091           auto lIndex = std::find(l_listOfParametrizations.begin(), l_listOfParametrizations.end(), l_values);
1092 
1093           // if it's not there push it back in the list of parametrizations and then insert in the map
1094           if (fIndex == f_listOfParametrizations.end()) {
1095             f_listOfParametrizations.push_back(f_values);
1096             std::vector<unsigned int> toInsert = {det};
1097             modules_per_region.insert(std::make_pair(f_listOfParametrizations.size() - 1, toInsert));
1098           } else {
1099             modules_per_region.at(fIndex - f_listOfParametrizations.begin()).push_back(det);
1100           }
1101 
1102           if (lIndex == l_listOfParametrizations.end()) {
1103             l_listOfParametrizations.push_back(l_values);
1104           }
1105         }
1106 
1107         unsigned int f_depth = f_listOfParametrizations.size();
1108         unsigned int l_depth = l_listOfParametrizations.size();
1109 
1110         if (l_depth != f_depth) {
1111           edm::LogError(label_) << label_
1112                                 << " trying to compare dynamic inefficiencys payload with different detid masks!";
1113           TCanvas canvas("Canv", "Canv", 1200, 1000);
1114           SiPixelPI::displayNotSupported(canvas, 0);
1115           std::string fileName(this->m_imageFileName);
1116           canvas.SaveAs(fileName.c_str());
1117           return false;
1118         }
1119 
1120         assert(f_depth == l_depth);
1121 
1122         std::vector<std::string> namesOfParts;
1123         namesOfParts.reserve(f_depth);
1124         // fill in the (ordered) information about regions
1125         for (const auto& [index, modules] : modules_per_region) {
1126           auto PhInfo = SiPixelPI::PhaseInfo(SiPixelPI::phase1size);
1127           const auto& regName = SiPixDynIneff::attachLocationLabel(modules, PhInfo);
1128           namesOfParts.push_back(regName);
1129 
1130           /*
1131        *  The following is needed for internal
1132        *  debug / cross-check
1133        */
1134 
1135           std::stringstream ss;
1136           ss << "region name: " << regName << " has the following modules attached: [";
1137           for (const auto& module : modules) {
1138             ss << module << ", ";
1139           }
1140           ss.seekp(-2, std::ios_base::end);  // remove last two chars
1141           ss << "] "
1142              << " has representation (";
1143           for (const auto& param : f_listOfParametrizations[index]) {
1144             ss << param << ", ";
1145           }
1146           ss.seekp(-2, std::ios_base::end);  // remove last two chars
1147           ss << ") ";
1148           edm::LogPrint(label_) << ss.str() << "\n\n";
1149           ss.str(std::string()); /* clear the stringstream */
1150         }
1151 
1152         // now fill the paramtrizations
1153         std::vector<TF1*> f_parametrizations;
1154         std::vector<std::string> f_formulas;
1155         f_parametrizations.reserve(f_depth);
1156         f_formulas.reserve(f_depth);
1157 
1158         std::vector<TF1*> l_parametrizations;
1159         std::vector<std::string> l_formulas;
1160         l_parametrizations.reserve(l_depth);
1161         l_formulas.reserve(l_depth);
1162 
1163         // fill the parametrization plots
1164         SiPixDynIneff::fillParametrizations(f_listOfParametrizations, f_parametrizations, f_formulas, namesOfParts);
1165         SiPixDynIneff::fillParametrizations(l_listOfParametrizations, l_parametrizations, l_formulas, namesOfParts);
1166 
1167         // determine how the plot will be paginated
1168         auto sides = SiPixDynIneff::getClosestFactors(f_depth);
1169         TCanvas canvas("Canv", "Canv", sides.second * 900, sides.first * 600);
1170         canvas.Divide(sides.second, sides.first);
1171 
1172         // print the sub-canvases
1173         for (unsigned int i = 0; i < f_depth; i++) {
1174           canvas.cd(i + 1);
1175           gPad->SetGrid();
1176           SiPixelPI::adjustCanvasMargins(gPad, 0.07, 0.14, 0.14, 0.03);
1177           f_parametrizations[i]->Draw();
1178           l_parametrizations[i]->SetLineColor(kBlue);
1179           l_parametrizations[i]->SetLineStyle(kDashed);
1180           l_parametrizations[i]->Draw("same");
1181 
1182           // set axis
1183           TH1* h = f_parametrizations[i]->GetHistogram();
1184           TAxis* ax = h->GetXaxis();
1185           ax->SetTitle("Inst. luminosity [10^{33} cm^{-2}s^{-1}]");
1186 
1187           TAxis* ay = h->GetYaxis();
1188           ay->SetRangeUser(0.80, 1.00);  // force the scale
1189           ay->SetTitle("Double Column Efficiency parametrization");
1190 
1191           // beautify
1192           SiPixelPI::makeNicePlotStyle(h);
1193 
1194           auto ltx = TLatex();
1195           ltx.SetTextFont(62);
1196           ltx.SetTextSize(0.045);
1197           ltx.SetTextAlign(11);
1198           ltx.DrawLatexNDC(
1199               gPad->GetLeftMargin() + 0.03,
1200               gPad->GetBottomMargin() + 0.16,
1201               ("#color[2]{" + f_tagname + "},IOV: #color[2]{" + std::to_string(std::get<0>(firstiov)) + "}").c_str());
1202           ltx.DrawLatexNDC(
1203               gPad->GetLeftMargin() + 0.03,
1204               gPad->GetBottomMargin() + 0.08,
1205               ("#color[4]{" + l_tagname + "},IOV: #color[4]{" + std::to_string(std::get<0>(lastiov)) + "}").c_str());
1206 
1207           auto ltxForm = TLatex();
1208           ltxForm.SetTextFont(62);
1209           ltxForm.SetTextColor(kRed);
1210           ltxForm.SetTextSize(0.035);
1211           ltxForm.SetTextAlign(11);
1212           ltxForm.DrawLatexNDC(gPad->GetLeftMargin() + 0.03, gPad->GetBottomMargin() + 0.12, f_formulas[i].c_str());
1213 
1214           ltxForm.SetTextColor(kBlue);
1215           ltxForm.DrawLatexNDC(gPad->GetLeftMargin() + 0.03, gPad->GetBottomMargin() + 0.04, l_formulas[i].c_str());
1216 
1217           edm::LogPrint(label_) << namesOfParts[i] << " => " << f_formulas[i] << std::endl;
1218           edm::LogPrint(label_) << namesOfParts[i] << " => " << l_formulas[i] << std::endl;
1219         }
1220 
1221         std::string fileName(this->m_imageFileName);
1222         canvas.SaveAs(fileName.c_str());
1223         //canvas.SaveAs((fileName+".root").c_str());
1224 
1225       }  // if payload.get()
1226       return true;
1227     }  // fill()
1228 
1229   protected:
1230     std::string label_;
1231   };
1232 
1233   using SiPixelDynamicInefficiencyPUParamComparisonTwoTags =
1234       SiPixelDynamicInefficiencyPUParamComparisonBase<SINGLE_IOV, 2>;
1235 
1236 }  // namespace
1237 
1238 // Register the classes as boost python plugin
1239 PAYLOAD_INSPECTOR_MODULE(SiPixelDynamicInefficiency) {
1240   PAYLOAD_INSPECTOR_CLASS(SiPixelDynamicInefficiencyTest);
1241   PAYLOAD_INSPECTOR_CLASS(SiPixelBPixIneffROCfromDynIneffMap);
1242   PAYLOAD_INSPECTOR_CLASS(SiPixelFPixIneffROCfromDynIneffMap);
1243   PAYLOAD_INSPECTOR_CLASS(SiPixelFullIneffROCfromDynIneffMap);
1244   PAYLOAD_INSPECTOR_CLASS(SiPixelBPixIneffROCsMapCompareTwoTags);
1245   PAYLOAD_INSPECTOR_CLASS(SiPixelFPixIneffROCsMapCompareTwoTags);
1246   PAYLOAD_INSPECTOR_CLASS(SiPixelFullIneffROCsMapCompareTwoTags);
1247   PAYLOAD_INSPECTOR_CLASS(SiPixelDynamicInefficiencyGeomFactorMap);
1248   PAYLOAD_INSPECTOR_CLASS(SiPixelDynamicInefficiencyColGeomFactorMap);
1249   PAYLOAD_INSPECTOR_CLASS(SiPixelDynamicInefficiencyChipGeomFactorMap);
1250   PAYLOAD_INSPECTOR_CLASS(SiPixelDynamicInefficiencyPUPixelMaps);
1251   PAYLOAD_INSPECTOR_CLASS(SiPixelDynamicInefficiencyPUParametrization);
1252   PAYLOAD_INSPECTOR_CLASS(SiPixelDynamicInefficiencyPUParamComparisonTwoTags);
1253 }