Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:01:45

0001 /*!
0002   \file SiPixelFEDChannelContainer_PayloadInspector
0003   \Payload Inspector Plugin for SiPixelFEDChannelContainer
0004   \author M. Musich
0005   \version $Revision: 1.0 $
0006   \date $Date: 2020/02/22 10:00:00 $
0007 */
0008 
0009 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0010 #include "CondCore/CondDB/interface/ConnectionPool.h"
0011 
0012 #include "CondCore/Utilities/interface/PayloadInspectorModule.h"
0013 #include "CondCore/Utilities/interface/PayloadInspector.h"
0014 #include "CondCore/CondDB/interface/Time.h"
0015 #include "CondCore/SiPixelPlugins/interface/SiPixelPayloadInspectorHelper.h"
0016 #include "FWCore/ParameterSet/interface/FileInPath.h"
0017 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0018 #include "CalibTracker/StandaloneTrackerTopology/interface/StandaloneTrackerTopology.h"
0019 
0020 #include "CondFormats/SiPixelObjects/interface/CablingPathToDetUnit.h"
0021 #include "CondFormats/SiPixelObjects/interface/LocalPixel.h"
0022 #include "CondFormats/SiPixelObjects/interface/PixelIndices.h"
0023 #include "CondFormats/SiPixelObjects/interface/PixelROC.h"
0024 #include "CondFormats/SiPixelObjects/interface/SiPixelFedCabling.h"
0025 #include "CondFormats/SiPixelObjects/interface/SiPixelFedCablingMap.h"
0026 #include "CondFormats/SiPixelObjects/interface/SiPixelFedCablingTree.h"
0027 #include "CondFormats/SiPixelObjects/interface/SiPixelFrameConverter.h"
0028 
0029 // the data format of the condition to be inspected
0030 #include "CondFormats/SiPixelObjects/interface/SiPixelFEDChannelContainer.h"
0031 #include "CondFormats/SiPixelObjects/interface/SiPixelQualityProbabilities.h"  // to display aggregate probability
0032 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
0033 #include "DataFormats/DetId/interface/DetId.h"
0034 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0035 #include "DQM/TrackerRemapper/interface/Phase1PixelROCMaps.h"
0036 
0037 #include <memory>
0038 #include <sstream>
0039 #include <iostream>
0040 #include <fmt/printf.h>
0041 
0042 // include ROOT
0043 #include "TH2F.h"
0044 #include "TLegend.h"
0045 #include "TCanvas.h"
0046 #include "TLine.h"
0047 #include "TGraph.h"
0048 #include "TGaxis.h"
0049 #include "TStyle.h"
0050 #include "TLatex.h"
0051 #include "TPave.h"
0052 #include "TPaveStats.h"
0053 
0054 namespace {
0055 
0056   using namespace cond::payloadInspector;
0057 
0058   /************************************************
0059   1d histogram of SiPixelFEDChannelContainer of 1 IOV 
0060   *************************************************/
0061 
0062   template <SiPixelPI::DetType myType>
0063   class SiPixelFEDChannelContainerMap : public PlotImage<SiPixelFEDChannelContainer, SINGLE_IOV> {
0064   public:
0065     SiPixelFEDChannelContainerMap()
0066         : PlotImage<SiPixelFEDChannelContainer, SINGLE_IOV>(
0067               "SiPixelFEDChannelContainer Pixel Track Map of one (or more scenarios"),
0068           m_trackerTopo{StandaloneTrackerTopology::fromTrackerParametersXMLFile(
0069               edm::FileInPath("Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml").fullPath())} {
0070       // for inputs
0071       PlotBase::addInputParam("Scenarios");
0072 
0073       // hardcoded connection to the MC cabling tag, though luck
0074       m_condDbCabling = "frontier://FrontierProd/CMS_CONDITIONS";
0075       m_CablingTagName = "SiPixelFedCablingMap_phase1_v7";
0076 
0077       m_connectionPool.setParameters(m_connectionPset);
0078       m_connectionPool.configure();
0079     }
0080 
0081     bool fill() override {
0082       std::vector<std::string> the_scenarios = {};
0083 
0084       auto paramValues = PlotBase::inputParamValues();
0085       auto ip = paramValues.find("Scenarios");
0086       if (ip != paramValues.end()) {
0087         auto input = ip->second;
0088         typedef boost::tokenizer<boost::char_separator<char>> tokenizer;
0089         boost::char_separator<char> sep{","};
0090         tokenizer tok{input, sep};
0091         for (const auto& t : tok) {
0092           the_scenarios.push_back(t);
0093         }
0094       } else {
0095         edm::LogWarning(k_ClassName)
0096             << "\n WARNING!!!! \n The needed parameter Scenarios has not been passed. Will use all the scenarios in "
0097                "the file!"
0098             << "\n Buckle your seatbelts... this might take a while... \n\n";
0099         the_scenarios.push_back("all");
0100       }
0101 
0102       Phase1PixelROCMaps theROCMap("");
0103 
0104       auto tag = PlotBase::getTag<0>();
0105       auto tagname = tag.name;
0106       auto iov = tag.iovs.front();
0107 
0108       // open db session for the cabling map
0109       edm::LogPrint(k_ClassName) << "[SiPixelFEDChannelContainerTest::" << __func__ << "] "
0110                                  << "Query the condition database " << m_condDbCabling;
0111 
0112       cond::persistency::Session condDbSession = m_connectionPool.createSession(m_condDbCabling);
0113       condDbSession.transaction().start(true);
0114 
0115       // query the database
0116       edm::LogPrint(k_ClassName) << "[SiPixelFEDChannelContainerTest::" << __func__ << "] "
0117                                  << "Reading IOVs from tag " << m_CablingTagName;
0118 
0119       const auto MIN_VAL = cond::timeTypeSpecs[cond::runnumber].beginValue;
0120       const auto MAX_VAL = cond::timeTypeSpecs[cond::runnumber].endValue;
0121 
0122       // get the list of payloads for the Cabling Map
0123       std::vector<std::tuple<cond::Time_t, cond::Hash>> m_cabling_iovs;
0124       condDbSession.readIov(m_CablingTagName).selectRange(MIN_VAL, MAX_VAL, m_cabling_iovs);
0125 
0126       std::vector<unsigned int> listOfCablingIOVs;
0127       std::transform(m_cabling_iovs.begin(),
0128                      m_cabling_iovs.end(),
0129                      std::back_inserter(listOfCablingIOVs),
0130                      [](std::tuple<cond::Time_t, cond::Hash> myIOV2) -> unsigned int { return std::get<0>(myIOV2); });
0131 
0132       edm::LogPrint(k_ClassName) << " Number of SiPixelFedCablngMap payloads: " << listOfCablingIOVs.size()
0133                                  << std::endl;
0134 
0135       auto it = std::find(
0136           listOfCablingIOVs.begin(), listOfCablingIOVs.end(), closest_from_below(listOfCablingIOVs, std::get<0>(iov)));
0137       int index = std::distance(listOfCablingIOVs.begin(), it);
0138 
0139       edm::LogPrint(k_ClassName) << " using the SiPixelFedCablingMap with hash: "
0140                                  << std::get<1>(m_cabling_iovs.at(index)) << std::endl;
0141 
0142       auto theCablingMap = condDbSession.fetchPayload<SiPixelFedCablingMap>(std::get<1>(m_cabling_iovs.at(index)));
0143       theCablingMap->initializeRocs();
0144       // auto theCablingTree = (*theCablingMap).cablingTree();
0145 
0146       //auto map = theCablingMap->det2fedMap();
0147       //for (const auto &element : map){
0148       //    std::cout << element.first << " " << element.second << std::endl;
0149       //}
0150 
0151       std::shared_ptr<SiPixelFEDChannelContainer> payload = fetchPayload(std::get<1>(iov));
0152       const auto& scenarioMap = payload->getScenarioMap();
0153 
0154       auto pIndexConverter = PixelIndices(numColumns, numRows);
0155 
0156       for (const auto& scenario : scenarioMap) {
0157         std::string scenName = scenario.first;
0158 
0159         if (std::find_if(the_scenarios.begin(), the_scenarios.end(), compareKeys(scenName)) != the_scenarios.end() ||
0160             the_scenarios[0] == "all") {
0161           edm::LogPrint(k_ClassName) << "\t Found Scenario: " << scenName << " ==> dumping it";
0162         } else {
0163           continue;
0164         }
0165 
0166         //if (strcmp(scenName.c_str(),"320824_103") != 0) continue;
0167 
0168         const auto& theDetSetBadPixelFedChannels = payload->getDetSetBadPixelFedChannels(scenName);
0169         for (const auto& disabledChannels : *theDetSetBadPixelFedChannels) {
0170           const auto t_detid = disabledChannels.detId();
0171           int subid = DetId(t_detid).subdetId();
0172           LogDebug(k_ClassName) << fmt::sprintf("DetId : %i \n", t_detid) << std::endl;
0173 
0174           std::bitset<16> badRocsFromFEDChannels;
0175 
0176           for (const auto& ch : disabledChannels) {
0177             std::string toOut_ = fmt::sprintf("fed : %i | link : %2i | roc_first : %2i | roc_last: %2i \n",
0178                                               ch.fed,
0179                                               ch.link,
0180                                               ch.roc_first,
0181                                               ch.roc_last);
0182 
0183             LogDebug(k_ClassName) << toOut_ << std::endl;
0184             const std::vector<sipixelobjects::CablingPathToDetUnit>& path =
0185                 theCablingMap->pathToDetUnit(disabledChannels.detId());
0186             for (unsigned int i_roc = ch.roc_first; i_roc <= ch.roc_last; ++i_roc) {
0187               for (const auto p : path) {
0188                 const sipixelobjects::PixelROC* myroc = theCablingMap->findItem(p);
0189                 if (myroc->idInDetUnit() == static_cast<unsigned int>(i_roc)) {
0190                   sipixelobjects::LocalPixel::RocRowCol local = {39, 25};  //corresponding to center of ROC row,col
0191                   sipixelobjects::GlobalPixel global = myroc->toGlobal(sipixelobjects::LocalPixel(local));
0192                   int chipIndex(0), colROC(0), rowROC(0);
0193 
0194                   pIndexConverter.transformToROC(global.col, global.row, chipIndex, colROC, rowROC);
0195 
0196                   LogDebug(k_ClassName) << " => i_roc:" << i_roc << "  " << global.col << "-" << global.row << " | => "
0197                                         << chipIndex << " : (" << colROC << "," << rowROC << ")" << std::endl;
0198 
0199                   badRocsFromFEDChannels[chipIndex] = true;
0200                 }
0201               }
0202             }
0203           }
0204 
0205           LogDebug(k_ClassName) << badRocsFromFEDChannels << std::endl;
0206 
0207           auto myDetId = DetId(t_detid);
0208 
0209           if (subid == PixelSubdetector::PixelBarrel) {
0210             theROCMap.fillSelectedRocs(myDetId, badRocsFromFEDChannels, 1.);
0211           }  // if it's barrel
0212           else if (subid == PixelSubdetector::PixelEndcap) {
0213             theROCMap.fillSelectedRocs(myDetId, badRocsFromFEDChannels, 1.);
0214           }  // if it's endcap
0215           else {
0216             throw cms::Exception("LogicError") << "Unknown Pixel SubDet ID " << std::endl;
0217           }  // else nonsense
0218         }    // loop on the channels
0219       }      // loop on the scenarios
0220 
0221       gStyle->SetOptStat(0);
0222       //=========================
0223       TCanvas canvas("Summary", "Summary", 1200, k_height[myType]);
0224       canvas.cd();
0225 
0226       auto unpacked = SiPixelPI::unpack(std::get<0>(iov));
0227 
0228       std::string IOVstring = (unpacked.first == 0)
0229                                   ? std::to_string(unpacked.second)
0230                                   : (std::to_string(unpacked.first) + "," + std::to_string(unpacked.second));
0231 
0232       const auto headerText = fmt::sprintf("#color[4]{%s},  IOV: #color[4]{%s}", tagname, IOVstring);
0233 
0234       switch (myType) {
0235         case SiPixelPI::t_barrel:
0236           theROCMap.drawBarrelMaps(canvas, headerText);
0237           break;
0238         case SiPixelPI::t_forward:
0239           theROCMap.drawForwardMaps(canvas, headerText);
0240           break;
0241         case SiPixelPI::t_all:
0242           theROCMap.drawMaps(canvas, headerText);
0243           break;
0244         default:
0245           throw cms::Exception("LogicError") << "\nERROR: unrecognized Pixel Detector part " << std::endl;
0246       }
0247 
0248       std::string fileName(m_imageFileName);
0249       canvas.SaveAs(fileName.c_str());
0250 
0251       // close the DB session
0252       condDbSession.transaction().commit();
0253 
0254       return true;
0255     }
0256 
0257   public:
0258     inline unsigned int closest_from_above(std::vector<unsigned int> const& vec, unsigned int value) {
0259       auto const it = std::lower_bound(vec.begin(), vec.end(), value);
0260       return vec.at(it - vec.begin() - 1);
0261     }
0262 
0263     inline unsigned int closest_from_below(std::vector<unsigned int> const& vec, unsigned int value) {
0264       auto const it = std::upper_bound(vec.begin(), vec.end(), value);
0265       return vec.at(it - vec.begin() - 1);
0266     }
0267 
0268     // auxilliary check
0269     struct compareKeys {
0270       std::string key;
0271       compareKeys(std::string const& i) : key(i) {}
0272 
0273       bool operator()(std::string const& i) { return (key == i); }
0274     };
0275 
0276   private:
0277     // tough luck, we can only do phase-1...
0278     static constexpr int numColumns = 416;
0279     static constexpr int numRows = 160;
0280     static constexpr int n_rings = 2;
0281     static constexpr int n_layers = 4;
0282 
0283     // graphics
0284     static constexpr std::array<int, 3> k_height = {{1200, 600, 1600}};
0285     static constexpr const char* k_ClassName = "SiPixelFEDChannelContainerMap";
0286 
0287     TrackerTopology m_trackerTopo;
0288     edm::ParameterSet m_connectionPset;
0289     cond::persistency::ConnectionPool m_connectionPool;
0290     std::string m_CablingTagName;
0291     std::string m_condDbCabling;
0292   };
0293 
0294   /************************************************
0295   1d histogram of SiPixelFEDChannelContainer of 1 IOV
0296   *************************************************/
0297 
0298   template <SiPixelPI::DetType myType>
0299   class SiPixelFEDChannelContainerMapSimple : public PlotImage<SiPixelFEDChannelContainer, SINGLE_IOV> {
0300   public:
0301     SiPixelFEDChannelContainerMapSimple()
0302         : PlotImage<SiPixelFEDChannelContainer, SINGLE_IOV>(
0303               "SiPixelFEDChannelContainer Pixel Track Map of one (or more scenarios)") {
0304       // for inputs
0305       PlotBase::addInputParam("Scenarios");
0306     }
0307 
0308     bool fill() override {
0309       std::vector<std::string> the_scenarios = {};
0310 
0311       auto paramValues = PlotBase::inputParamValues();
0312       auto ip = paramValues.find("Scenarios");
0313       if (ip != paramValues.end()) {
0314         auto input = ip->second;
0315         typedef boost::tokenizer<boost::char_separator<char>> tokenizer;
0316         boost::char_separator<char> sep{","};
0317         tokenizer tok{input, sep};
0318         for (const auto& t : tok) {
0319           the_scenarios.push_back(t);
0320         }
0321       } else {
0322         edm::LogWarning(k_ClassName)
0323             << "\n WARNING!!!! \n The needed parameter Scenarios has not been passed. Will use all the scenarios in "
0324                "the file!"
0325             << "\n Buckle your seatbelts... this might take a while... \n\n";
0326         the_scenarios.push_back("all");
0327       }
0328 
0329       Phase1PixelROCMaps theROCMap("");
0330 
0331       auto tag = PlotBase::getTag<0>();
0332       auto tagname = tag.name;
0333       auto iov = tag.iovs.front();
0334 
0335       std::shared_ptr<SiPixelFEDChannelContainer> payload = fetchPayload(std::get<1>(iov));
0336       const auto& scenarioMap = payload->getScenarioMap();
0337 
0338       for (const auto& scenario : scenarioMap) {
0339         std::string scenName = scenario.first;
0340 
0341         if (std::find_if(the_scenarios.begin(), the_scenarios.end(), compareKeys(scenName)) != the_scenarios.end() ||
0342             the_scenarios[0] == "all") {
0343           edm::LogPrint(k_ClassName) << "\t Found Scenario: " << scenName << " ==> dumping it";
0344         } else {
0345           continue;
0346         }
0347 
0348         const auto& theDetSetBadPixelFedChannels = payload->getDetSetBadPixelFedChannels(scenName);
0349         for (const auto& disabledChannels : *theDetSetBadPixelFedChannels) {
0350           const auto t_detid = disabledChannels.detId();
0351           int subid = DetId(t_detid).subdetId();
0352           LogDebug(k_ClassName) << fmt::sprintf("DetId : %i \n", t_detid) << std::endl;
0353 
0354           std::bitset<16> badRocsFromFEDChannels;
0355 
0356           for (const auto& ch : disabledChannels) {
0357             std::string toOut_ = fmt::sprintf("fed : %i | link : %2i | roc_first : %2i | roc_last: %2i \n",
0358                                               ch.fed,
0359                                               ch.link,
0360                                               ch.roc_first,
0361                                               ch.roc_last);
0362 
0363             LogDebug(k_ClassName) << toOut_ << std::endl;
0364             for (unsigned int i_roc = ch.roc_first; i_roc <= ch.roc_last; ++i_roc) {
0365               badRocsFromFEDChannels.set(i_roc);
0366             }
0367           }
0368 
0369           LogDebug(k_ClassName) << badRocsFromFEDChannels << std::endl;
0370 
0371           const auto& myDetId = DetId(t_detid);
0372 
0373           if (subid == PixelSubdetector::PixelBarrel) {
0374             theROCMap.fillSelectedRocs(myDetId, badRocsFromFEDChannels, 1.);
0375           }  // if it's barrel
0376           else if (subid == PixelSubdetector::PixelEndcap) {
0377             theROCMap.fillSelectedRocs(myDetId, badRocsFromFEDChannels, 1.);
0378           }  // if it's endcap
0379           else {
0380             throw cms::Exception("LogicError") << "Unknown Pixel SubDet ID " << std::endl;
0381           }  // else nonsense
0382         }    // loop on the channels
0383       }      // loop on the scenarios
0384 
0385       gStyle->SetOptStat(0);
0386       //=========================
0387       TCanvas canvas("Summary", "Summary", 1200, k_height[myType]);
0388       canvas.cd();
0389 
0390       auto unpacked = SiPixelPI::unpack(std::get<0>(iov));
0391 
0392       std::string IOVstring = (unpacked.first == 0)
0393                                   ? std::to_string(unpacked.second)
0394                                   : (std::to_string(unpacked.first) + "," + std::to_string(unpacked.second));
0395 
0396       const auto headerText = fmt::sprintf("#color[4]{%s},  IOV: #color[4]{%s}", tagname, IOVstring);
0397 
0398       switch (myType) {
0399         case SiPixelPI::t_barrel:
0400           theROCMap.drawBarrelMaps(canvas, headerText);
0401           break;
0402         case SiPixelPI::t_forward:
0403           theROCMap.drawForwardMaps(canvas, headerText);
0404           break;
0405         case SiPixelPI::t_all:
0406           theROCMap.drawMaps(canvas, headerText);
0407           break;
0408         default:
0409           throw cms::Exception("LogicError") << "\nERROR: unrecognized Pixel Detector part " << std::endl;
0410       }
0411 
0412       // add list of scenarios watermark
0413       canvas.cd();
0414       auto ltx = TLatex();
0415       ltx.SetTextFont(62);
0416       //ltx.SetTextColor(kMagenta);
0417       ltx.SetTextSize(0.023);
0418       ltx.DrawLatexNDC(
0419           gPad->GetLeftMargin() - 0.09,
0420           gPad->GetBottomMargin() - 0.09,
0421           ("scenarios: #color[4]{" +
0422            std::accumulate(the_scenarios.begin(),
0423                            the_scenarios.end(),
0424                            std::string(),
0425                            [](const std::string& acc, const std::string& str) { return acc + " " + str; }) +
0426            "}")
0427               .c_str());
0428 
0429       std::string fileName(m_imageFileName);
0430       canvas.SaveAs(fileName.c_str());
0431       return true;
0432     }
0433 
0434   public:
0435     // auxilliary check
0436     struct compareKeys {
0437       std::string key;
0438       compareKeys(std::string const& i) : key(i) {}
0439 
0440       bool operator()(std::string const& i) { return (key == i); }
0441     };
0442 
0443   private:
0444     // graphics
0445     static constexpr std::array<int, 3> k_height = {{1200, 600, 1600}};
0446     static constexpr const char* k_ClassName = "SiPixelFEDChannelContainerMapSimple";
0447   };
0448 
0449   using SiPixelBPixFEDChannelContainerMap = SiPixelFEDChannelContainerMapSimple<SiPixelPI::t_barrel>;
0450   using SiPixelFPixFEDChannelContainerMap = SiPixelFEDChannelContainerMapSimple<SiPixelPI::t_forward>;
0451   using SiPixelFullFEDChannelContainerMap = SiPixelFEDChannelContainerMapSimple<SiPixelPI::t_all>;
0452 
0453   /*
0454     Produces an aggregate map of the masked components for all scenarios,
0455     weighted on the probability per PU unit from SiPixelQualityProbabilities
0456     assuming a flat PU profile in the range encoded in  SiPixelQualityProbabilities
0457     The SiPixelQualityProbabilities tag comes from user input
0458   */
0459   template <SiPixelPI::DetType myType>
0460   class SiPixelFEDChannelContainerMapWeigthed : public PlotImage<SiPixelFEDChannelContainer, SINGLE_IOV> {
0461   public:
0462     SiPixelFEDChannelContainerMapWeigthed()
0463         : PlotImage<SiPixelFEDChannelContainer, SINGLE_IOV>(
0464               "SiPixelFEDChannelContainer Pixel Track Map of one (or more scenarios)") {
0465       // for inputs
0466       PlotBase::addInputParam("SiPixelQualityProbabilitiesTag");
0467 
0468       // hardcoded connection to the SiPixelQualityProbability tag, though luck
0469       m_condSiPixelProb = "frontier://FrontierProd/CMS_CONDITIONS";
0470       m_connectionPool.setParameters(m_connectionPset);
0471       m_connectionPool.configure();
0472     }
0473 
0474     bool fill() override {
0475       auto paramValues = PlotBase::inputParamValues();
0476       auto ip = paramValues.find("SiPixelQualityProbabilitiesTag");
0477       if (ip != paramValues.end()) {
0478         m_SiPixelProbTagName = ip->second;
0479       } else {
0480         edm::LogWarning(k_ClassName) << "\n WARNING!!!! \n The needed parameter SiPixelQualityProbabilitiesTag was not "
0481                                         "inputed from the user \n Display will be aborted \n\n";
0482         return false;
0483       }
0484 
0485       Phase1PixelROCMaps theROCMap("", "Masking Probability [%]");
0486 
0487       auto tag = PlotBase::getTag<0>();
0488       auto tagname = tag.name;
0489       auto iov = tag.iovs.front();
0490 
0491       // open db session for the cabling map
0492       edm::LogPrint(k_ClassName) << "[SiPixelFEDChannelContainerTest::" << __func__ << "] "
0493                                  << "Query the condition database " << m_condSiPixelProb;
0494 
0495       cond::persistency::Session condDbSession = m_connectionPool.createSession(m_condSiPixelProb);
0496       condDbSession.transaction().start(true);
0497 
0498       // query the database
0499       edm::LogPrint(k_ClassName) << "[SiPixelFEDChannelContainerTest::" << __func__ << "] "
0500                                  << "Reading IOVs from tag " << m_SiPixelProbTagName;
0501 
0502       const auto MIN_VAL = cond::timeTypeSpecs[cond::runnumber].beginValue;
0503       const auto MAX_VAL = cond::timeTypeSpecs[cond::runnumber].endValue;
0504 
0505       // get the list of payloads for the Cabling Map
0506       std::vector<std::tuple<cond::Time_t, cond::Hash>> m_pixelProb_iovs;
0507       condDbSession.readIov(m_SiPixelProbTagName).selectRange(MIN_VAL, MAX_VAL, m_pixelProb_iovs);
0508 
0509       // in MC there should be only 1 IOV, oh well...
0510       edm::LogPrint(k_ClassName) << " using the SiPixelQualityProbabilities with hash: "
0511                                  << std::get<1>(m_pixelProb_iovs.front()) << std::endl;
0512 
0513       auto probabilitiesPayload =
0514           condDbSession.fetchPayload<SiPixelQualityProbabilities>(std::get<1>(m_pixelProb_iovs.front()));
0515 
0516       const auto& PUbins = probabilitiesPayload->getPileUpBins();
0517 
0518       SiPixelQualityProbabilities::probabilityMap m_probabilities = probabilitiesPayload->getProbability_Map();
0519 
0520       // find the PU-averaged (assuming flat PU in the range) probabilities for each scenario
0521       std::map<std::string, float> puAvgedProbabilities;
0522       for (const auto& [PUbin, ProbMap] : m_probabilities) {
0523         float totProbInPUbin = 0.f;
0524         for (const auto& [scenName, prob] : ProbMap) {
0525           totProbInPUbin += prob;
0526           if (puAvgedProbabilities.find(scenName) == puAvgedProbabilities.end()) {
0527             puAvgedProbabilities[scenName] += prob;
0528           } else {
0529             puAvgedProbabilities.insert({scenName, prob});
0530           }
0531         }
0532         LogDebug(k_ClassName) << "PU bin: " << PUbin << " tot probability " << totProbInPUbin << std::endl;
0533       }
0534 
0535       std::shared_ptr<SiPixelFEDChannelContainer> payload = fetchPayload(std::get<1>(iov));
0536       const auto& scenarioMap = payload->getScenarioMap();
0537 
0538       float totProb{0.f};
0539       for (const auto& [scenName, prob] : puAvgedProbabilities) {
0540         // only sum up the scenarios that are in the SiPixelFEDChannelContainer payload!
0541         if (scenarioMap.find(scenName) != scenarioMap.end()) {
0542           LogDebug(k_ClassName) << scenName << " : " << prob << std::endl;
0543           totProb += prob;
0544         }
0545       }
0546 
0547       LogDebug(k_ClassName) << "Total probability to normalize to: " << totProb << std::endl;
0548 
0549       //normalize the probabilities per scenario to the toal probability
0550       for (auto& pair : puAvgedProbabilities) {
0551         pair.second /= totProb;
0552       }
0553 
0554       for (const auto& scenario : scenarioMap) {
0555         std::string scenName = scenario.first;
0556         LogDebug(k_ClassName) << "\t Found Scenario: " << scenName << " ==> dumping it";
0557 
0558         // calculate the weight
0559         float w_frac = 0.f;
0560         if (puAvgedProbabilities.find(scenName) != puAvgedProbabilities.end()) {
0561           w_frac = puAvgedProbabilities[scenName];
0562         }
0563 
0564         // if scenario is not in the probability payload, continue
0565         if (w_frac == 0.f)
0566           continue;
0567 
0568         LogDebug(k_ClassName) << "scen: " << scenName << " weight: " << w_frac << " log(weight):" << log10(w_frac)
0569                               << std::endl;
0570 
0571         const auto& theDetSetBadPixelFedChannels = payload->getDetSetBadPixelFedChannels(scenName);
0572         for (const auto& disabledChannels : *theDetSetBadPixelFedChannels) {
0573           const auto t_detid = disabledChannels.detId();
0574           int subid = DetId(t_detid).subdetId();
0575           LogDebug(k_ClassName) << fmt::sprintf("DetId : %i \n", t_detid) << std::endl;
0576 
0577           std::bitset<16> badRocsFromFEDChannels;
0578 
0579           for (const auto& ch : disabledChannels) {
0580             std::string toOut_ = fmt::sprintf("fed : %i | link : %2i | roc_first : %2i | roc_last: %2i \n",
0581                                               ch.fed,
0582                                               ch.link,
0583                                               ch.roc_first,
0584                                               ch.roc_last);
0585 
0586             LogDebug(k_ClassName) << toOut_ << std::endl;
0587             for (unsigned int i_roc = ch.roc_first; i_roc <= ch.roc_last; ++i_roc) {
0588               badRocsFromFEDChannels.set(i_roc);
0589             }
0590           }
0591 
0592           LogDebug(k_ClassName) << badRocsFromFEDChannels << std::endl;
0593 
0594           const auto& myDetId = DetId(t_detid);
0595 
0596           if (subid == PixelSubdetector::PixelBarrel) {
0597             theROCMap.fillSelectedRocs(myDetId, badRocsFromFEDChannels, w_frac * 100);
0598           }  // if it's barrel
0599           else if (subid == PixelSubdetector::PixelEndcap) {
0600             theROCMap.fillSelectedRocs(myDetId, badRocsFromFEDChannels, w_frac * 100);
0601           }  // if it's endcap
0602           else {
0603             throw cms::Exception("LogicError") << "Unknown Pixel SubDet ID " << std::endl;
0604           }  // else nonsense
0605         }    // loop on the channels
0606       }      // loop on the scenarios
0607 
0608       gStyle->SetOptStat(0);
0609       //=========================
0610       TCanvas canvas("Summary", "Summary", 1200, k_height[myType]);
0611       canvas.cd();
0612 
0613       auto unpacked = SiPixelPI::unpack(std::get<0>(iov));
0614 
0615       std::string IOVstring = (unpacked.first == 0)
0616                                   ? std::to_string(unpacked.second)
0617                                   : (std::to_string(unpacked.first) + "," + std::to_string(unpacked.second));
0618 
0619       const auto headerText =
0620           fmt::sprintf("#bf{#scale[0.6]{#color[2]{%s}, #color[4]{%s}}}", tagname, m_SiPixelProbTagName);
0621 
0622       switch (myType) {
0623         case SiPixelPI::t_barrel:
0624           theROCMap.drawBarrelMaps(canvas, headerText);
0625           break;
0626         case SiPixelPI::t_forward:
0627           theROCMap.drawForwardMaps(canvas, headerText);
0628           break;
0629         case SiPixelPI::t_all:
0630           theROCMap.drawMaps(canvas, headerText);
0631           break;
0632         default:
0633           throw cms::Exception("LogicError") << "\nERROR: unrecognized Pixel Detector part " << std::endl;
0634       }
0635 
0636       // add list of scenarios watermark
0637       canvas.cd();
0638       auto ltx = TLatex();
0639       ltx.SetTextFont(62);
0640       //ltx.SetTextColor(kMagenta);
0641       ltx.SetTextSize(0.023);
0642       ltx.DrawLatexNDC(gPad->GetLeftMargin() - 0.09, gPad->GetBottomMargin() - 0.09, "");
0643 
0644       std::string fileName(m_imageFileName);
0645       canvas.SaveAs(fileName.c_str());
0646       return true;
0647     }
0648 
0649   private:
0650     // graphics
0651     static constexpr std::array<int, 3> k_height = {{1200, 600, 1600}};
0652     static constexpr const char* k_ClassName = "SiPixelFEDChannelContainerMapWeigthed";
0653 
0654     // parameters for auxilliary DB connection
0655     edm::ParameterSet m_connectionPset;
0656     cond::persistency::ConnectionPool m_connectionPool;
0657     std::string m_SiPixelProbTagName;
0658     std::string m_condSiPixelProb;
0659   };
0660 
0661   using SiPixelBPixFEDChannelContainerWeightedMap = SiPixelFEDChannelContainerMapWeigthed<SiPixelPI::t_barrel>;
0662   using SiPixelFPixFEDChannelContainerWeightedMap = SiPixelFEDChannelContainerMapWeigthed<SiPixelPI::t_forward>;
0663   using SiPixelFullFEDChannelContainerWeightedMap = SiPixelFEDChannelContainerMapWeigthed<SiPixelPI::t_all>;
0664 
0665   /************************************************
0666   1d histogram of number of SiPixelFEDChannelContainer scenarios
0667   *************************************************/
0668 
0669   class SiPixelFEDChannelContainerScenarios : public PlotImage<SiPixelFEDChannelContainer, SINGLE_IOV> {
0670   public:
0671     SiPixelFEDChannelContainerScenarios()
0672         : PlotImage<SiPixelFEDChannelContainer, SINGLE_IOV>("SiPixelFEDChannelContainer scenarios count") {}
0673     bool fill() override {
0674       auto tag = PlotBase::getTag<0>();
0675       auto tagname = tag.name;
0676       auto iov = tag.iovs.front();
0677 
0678       TGaxis::SetMaxDigits(3);
0679 
0680       std::shared_ptr<SiPixelFEDChannelContainer> payload = fetchPayload(std::get<1>(iov));
0681       std::vector<std::string> scenarios = payload->getScenarioList();
0682       sort(scenarios.begin(), scenarios.end());
0683 
0684       TCanvas canvas("Canv", "Canv", 1200, 1000);
0685       canvas.cd();
0686       canvas.SetGrid();
0687       auto h1 = std::make_unique<TH1F>("Count",
0688                                        "SiPixelFEDChannelContainer Bad Roc count;Scenario index;n. of bad ROCs",
0689                                        scenarios.size(),
0690                                        1,
0691                                        scenarios.size());
0692       h1->SetStats(false);
0693 
0694       canvas.SetTopMargin(0.06);
0695       canvas.SetBottomMargin(0.12);
0696       canvas.SetLeftMargin(0.12);
0697       canvas.SetRightMargin(0.05);
0698       canvas.Modified();
0699 
0700       int scenarioIndex = 0;
0701       for (const auto& scenario : scenarios) {
0702         scenarioIndex++;
0703         int badRocCount = 0;
0704         LogDebug("SiPixelFEDChannelContainerScenarios") << scenario << std::endl;
0705         auto badChannelCollection = payload->getDetSetBadPixelFedChannels(scenario);
0706         for (const auto& disabledChannels : *badChannelCollection) {
0707           for (const auto& ch : disabledChannels) {
0708             int local_bad_rocs = ch.roc_last - ch.roc_first;
0709             badRocCount += local_bad_rocs;
0710           }  // loop on the channels
0711         }    // loop on the DetSetVector
0712 
0713         h1->SetBinContent(scenarioIndex, badRocCount);
0714       }  // loop on scenarios
0715 
0716       TGaxis::SetExponentOffset(-0.1, 0.01, "y");    // Y offset
0717       TGaxis::SetExponentOffset(-0.03, -0.10, "x");  // Y and Y offset for X axis
0718 
0719       h1->SetTitle("");
0720       h1->GetYaxis()->SetRangeUser(0., h1->GetMaximum() * 1.30);
0721       h1->SetFillColor(kRed);
0722       h1->SetMarkerStyle(20);
0723       h1->SetMarkerSize(1);
0724       h1->Draw("bar2");
0725 
0726       SiPixelPI::makeNicePlotStyle(h1.get());
0727 
0728       canvas.Update();
0729 
0730       TLegend legend = TLegend(0.30, 0.88, 0.95, 0.94);
0731       //legend.SetHeader(("#splitline{Payload hash: #bf{" + (std::get<1>(iov)) + "}}{Total Scenarios:"+std::to_string(scenarioIndex)+"}").c_str(),"C");  // option "C" allows to center the header
0732 
0733       legend.SetHeader(fmt::sprintf("Payload hash: #bf{%s}", std::get<1>(iov)).c_str(), "C");
0734       legend.AddEntry(h1.get(), fmt::sprintf("total scenarios: #bf{%s}", std::to_string(scenarioIndex)).c_str(), "F");
0735       legend.SetTextSize(0.025);
0736       legend.Draw("same");
0737 
0738       auto ltx = TLatex();
0739       ltx.SetTextFont(62);
0740       //ltx.SetTextColor(kBlue);
0741       //ltx.SetTextAlign(11);
0742       ltx.SetTextSize(0.040);
0743       ltx.DrawLatexNDC(
0744           gPad->GetLeftMargin(),
0745           1 - gPad->GetTopMargin() + 0.01,
0746           fmt::sprintf("#color[4]{%s} IOV: #color[4]{%s}", tagname, std::to_string(std::get<0>(iov))).c_str());
0747 
0748       std::string fileName(m_imageFileName);
0749       canvas.SaveAs(fileName.c_str());
0750 
0751       return true;
0752 
0753     }  // fill
0754   };
0755 
0756 }  // namespace
0757 
0758 // Register the classes as boost python plugin
0759 PAYLOAD_INSPECTOR_MODULE(SiPixelFEDChannelContainer) {
0760   PAYLOAD_INSPECTOR_CLASS(SiPixelBPixFEDChannelContainerMap);
0761   PAYLOAD_INSPECTOR_CLASS(SiPixelFPixFEDChannelContainerMap);
0762   PAYLOAD_INSPECTOR_CLASS(SiPixelFullFEDChannelContainerMap);
0763   PAYLOAD_INSPECTOR_CLASS(SiPixelBPixFEDChannelContainerWeightedMap);
0764   PAYLOAD_INSPECTOR_CLASS(SiPixelFPixFEDChannelContainerWeightedMap);
0765   PAYLOAD_INSPECTOR_CLASS(SiPixelFullFEDChannelContainerWeightedMap);
0766   PAYLOAD_INSPECTOR_CLASS(SiPixelFEDChannelContainerScenarios);
0767 }