Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 10:46:20

0001 /*!
0002   \file SiPixelQuality_PayloadInspector
0003   \Payload Inspector Plugin for SiPixelQuality
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/SiPixelQuality.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 
0026 #include <memory>
0027 #include <sstream>
0028 #include <iostream>
0029 
0030 // include ROOT
0031 #include "TH2F.h"
0032 #include "TLegend.h"
0033 #include "TCanvas.h"
0034 #include "TLine.h"
0035 #include "TGraph.h"
0036 #include "TStyle.h"
0037 #include "TLatex.h"
0038 #include "TPave.h"
0039 #include "TPaveStats.h"
0040 
0041 namespace {
0042 
0043   using namespace cond::payloadInspector;
0044 
0045   /************************************************
0046     test class
0047   *************************************************/
0048 
0049   class SiPixelQualityTest : public Histogram1D<SiPixelQuality, SINGLE_IOV> {
0050   public:
0051     SiPixelQualityTest()
0052         : Histogram1D<SiPixelQuality, SINGLE_IOV>("SiPixelQuality test", "SiPixelQuality test", 10, 0.0, 10.0) {}
0053 
0054     bool fill() override {
0055       auto tag = PlotBase::getTag<0>();
0056       for (auto const& iov : tag.iovs) {
0057         std::shared_ptr<SiPixelQuality> payload = Base::fetchPayload(std::get<1>(iov));
0058         if (payload.get()) {
0059           fillWithValue(1.);
0060 
0061           auto theDisabledModules = payload->getBadComponentList();
0062           for (const auto& mod : theDisabledModules) {
0063             int BadRocCount(0);
0064             for (unsigned short n = 0; n < 16; n++) {
0065               unsigned short mask = 1 << n;  // 1 << n = 2^{n} using bitwise shift
0066               if (mod.BadRocs & mask)
0067                 BadRocCount++;
0068             }
0069             COUT << "detId:" << mod.DetID << " error type:" << mod.errorType << " BadRocs:" << BadRocCount << std::endl;
0070           }
0071         }  // payload
0072       }    // iovs
0073       return true;
0074     }  // fill
0075   };
0076 
0077   /************************************************
0078     Debugging class, to not be displayed
0079   *************************************************/
0080 
0081   class SiPixelQualityDebugger : public Histogram1D<SiPixelQuality, SINGLE_IOV> {
0082   public:
0083     SiPixelQualityDebugger()
0084         : Histogram1D<SiPixelQuality, SINGLE_IOV>("SiPixelQuality test", "SiPixelQuality test", 10, 0.0, 10.0) {}
0085 
0086     bool fill() override {
0087       auto tag = PlotBase::getTag<0>();
0088       for (auto const& iov : tag.iovs) {
0089         std::shared_ptr<SiPixelQuality> payload = Base::fetchPayload(std::get<1>(iov));
0090         if (payload.get()) {
0091           fillWithValue(1.);
0092 
0093           SiPixelPI::topolInfo t_info_fromXML;
0094           t_info_fromXML.init();
0095 
0096           auto theDisabledModules = payload->getBadComponentList();
0097           for (const auto& mod : theDisabledModules) {
0098             DetId detid(mod.DetID);
0099             auto PhInfo = SiPixelPI::PhaseInfo(SiPixelPI::phase1size);
0100             const char* path_toTopologyXML = PhInfo.pathToTopoXML();
0101             auto tTopo =
0102                 StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
0103             t_info_fromXML.fillGeometryInfo(detid, tTopo, PhInfo.phase());
0104             std::stringstream ss;
0105             t_info_fromXML.printAll(ss);
0106             std::bitset<16> bad_rocs(mod.BadRocs);
0107 
0108             if (t_info_fromXML.subDetId() == 1 && t_info_fromXML.layer() == 1) {
0109               std::cout << ss.str() << " s_module: " << SiPixelPI::signed_module(mod.DetID, tTopo, true)
0110                         << " s_ladder: " << SiPixelPI::signed_ladder(mod.DetID, tTopo, true)
0111                         << " error type:" << mod.errorType << " BadRocs: " << bad_rocs.to_string('O', 'X') << std::endl;
0112             }
0113           }
0114         }  // payload
0115       }    // iovs
0116       return true;
0117     }  // fill
0118   };
0119 
0120   /************************************************
0121     summary class
0122   *************************************************/
0123 
0124   class SiPixelQualityBadRocsSummary : public PlotImage<SiPixelQuality, MULTI_IOV, 1> {
0125   public:
0126     SiPixelQualityBadRocsSummary() : PlotImage<SiPixelQuality, MULTI_IOV, 1>("SiPixel Quality Summary") {}
0127 
0128     bool fill() override {
0129       auto tag = PlotBase::getTag<0>();
0130       for (const auto& iov : tag.iovs) {
0131         std::shared_ptr<SiPixelQuality> payload = fetchPayload(std::get<1>(iov));
0132         auto unpacked = SiPixelPI::unpack(std::get<0>(iov));
0133 
0134         COUT << "======================= " << unpacked.first << " : " << unpacked.second << std::endl;
0135         auto theDisabledModules = payload->getBadComponentList();
0136         for (const auto& mod : theDisabledModules) {
0137           COUT << "detId: " << mod.DetID << " |error type: " << mod.errorType << " |BadRocs: " << mod.BadRocs
0138                << std::endl;
0139         }
0140       }
0141 
0142       //=========================
0143       TCanvas canv("Partition summary", "partition summary", 1200, 1000);
0144       canv.SetBottomMargin(0.11);
0145       canv.SetLeftMargin(0.13);
0146       canv.SetRightMargin(0.05);
0147       canv.cd();
0148       SiPixelPI::displayNotSupported(canv, 0);
0149 
0150       std::string fileName(m_imageFileName);
0151       canv.SaveAs(fileName.c_str());
0152 
0153       return true;
0154     }
0155   };
0156 
0157   /************************************************
0158     time history class
0159   *************************************************/
0160 
0161   class SiPixelQualityBadRocsTimeHistory : public TimeHistoryPlot<SiPixelQuality, std::pair<double, double> > {
0162   public:
0163     SiPixelQualityBadRocsTimeHistory()
0164         : TimeHistoryPlot<SiPixelQuality, std::pair<double, double> >("bad ROCs count vs time", "bad ROCs count") {}
0165 
0166     std::pair<double, double> getFromPayload(SiPixelQuality& payload) override {
0167       return std::make_pair(extractBadRocCount(payload), 0.);
0168     }
0169 
0170     unsigned int extractBadRocCount(SiPixelQuality& payload) {
0171       unsigned int BadRocCount(0);
0172       auto theDisabledModules = payload.getBadComponentList();
0173       for (const auto& mod : theDisabledModules) {
0174         for (unsigned short n = 0; n < 16; n++) {
0175           unsigned short mask = 1 << n;  // 1 << n = 2^{n} using bitwise shift
0176           if (mod.BadRocs & mask)
0177             BadRocCount++;
0178         }
0179       }
0180       return BadRocCount;
0181     }
0182   };
0183 
0184   /************************************************
0185    occupancy style map whole Pixel
0186   *************************************************/
0187   template <SiPixelPI::DetType myType>
0188   class SiPixelQualityMap : public PlotImage<SiPixelQuality, SINGLE_IOV> {
0189   public:
0190     SiPixelQualityMap()
0191         : PlotImage<SiPixelQuality, SINGLE_IOV>("SiPixelQuality Pixel Map"),
0192           m_trackerTopo{StandaloneTrackerTopology::fromTrackerParametersXMLFile(
0193               edm::FileInPath("Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml").fullPath())} {}
0194 
0195     bool fill() override {
0196       auto tag = PlotBase::getTag<0>();
0197       auto iov = tag.iovs.front();
0198       auto tagname = tag.name;
0199       std::shared_ptr<SiPixelQuality> payload = fetchPayload(std::get<1>(iov));
0200 
0201       Phase1PixelROCMaps theMap("");
0202 
0203       auto theDisabledModules = payload->getBadComponentList();
0204       if (this->isPhase0(theDisabledModules)) {
0205         edm::LogError("SiPixelQuality_PayloadInspector")
0206             << "SiPixelQuality maps are not supported for non-Phase1 Pixel geometries !";
0207         TCanvas canvas("Canv", "Canv", 1200, 1000);
0208         SiPixelPI::displayNotSupported(canvas, 0);
0209         std::string fileName(m_imageFileName);
0210         canvas.SaveAs(fileName.c_str());
0211         return false;
0212       }
0213 
0214       for (const auto& mod : theDisabledModules) {
0215         int subid = DetId(mod.DetID).subdetId();
0216 
0217         if ((subid == PixelSubdetector::PixelBarrel && myType == SiPixelPI::t_barrel) ||
0218             (subid == PixelSubdetector::PixelEndcap && myType == SiPixelPI::t_forward) ||
0219             (myType == SiPixelPI::t_all)) {
0220           std::bitset<16> bad_rocs(mod.BadRocs);
0221           if (payload->IsModuleBad(mod.DetID)) {
0222             theMap.fillWholeModule(mod.DetID, 1.);
0223           } else {
0224             theMap.fillSelectedRocs(mod.DetID, bad_rocs, 1.);
0225           }
0226         }
0227       }
0228 
0229       gStyle->SetOptStat(0);
0230       //=========================
0231       TCanvas canvas("Summary", "Summary", 1200, k_height[myType]);
0232       canvas.cd();
0233 
0234       auto unpacked = SiPixelPI::unpack(std::get<0>(iov));
0235 
0236       std::string IOVstring = (unpacked.first == 0)
0237                                   ? std::to_string(unpacked.second)
0238                                   : (std::to_string(unpacked.first) + "," + std::to_string(unpacked.second));
0239 
0240       const auto headerText = fmt::sprintf("#color[4]{%s},  IOV: #color[4]{%s}", tagname, IOVstring);
0241 
0242       switch (myType) {
0243         case SiPixelPI::t_barrel:
0244           theMap.drawBarrelMaps(canvas, headerText);
0245           break;
0246         case SiPixelPI::t_forward:
0247           theMap.drawForwardMaps(canvas, headerText);
0248           break;
0249         case SiPixelPI::t_all:
0250           theMap.drawMaps(canvas, headerText);
0251           break;
0252         default:
0253           throw cms::Exception("SiPixelQualityMap") << "\nERROR: unrecognized Pixel Detector part " << std::endl;
0254       }
0255 
0256       std::string fileName(m_imageFileName);
0257       canvas.SaveAs(fileName.c_str());
0258 #ifdef MMDEBUG
0259       canvas.SaveAs("outAll.root");
0260 #endif
0261 
0262       return true;
0263     }
0264 
0265   private:
0266     static constexpr std::array<int, 3> k_height = {{1200, 600, 1600}};
0267     TrackerTopology m_trackerTopo;
0268 
0269     //_________________________________________________
0270     bool isPhase0(std::vector<SiPixelQuality::disabledModuleType> mods) {
0271       SiPixelDetInfoFileReader reader =
0272           SiPixelDetInfoFileReader(edm::FileInPath(SiPixelDetInfoFileReader::kPh0DefaultFile).fullPath());
0273       const auto& p0detIds = reader.getAllDetIds();
0274 
0275       std::vector<uint32_t> ownDetIds;
0276       std::transform(mods.begin(),
0277                      mods.end(),
0278                      std::back_inserter(ownDetIds),
0279                      [](SiPixelQuality::disabledModuleType d) -> uint32_t { return d.DetID; });
0280 
0281       for (const auto& det : ownDetIds) {
0282         // if found at least one phase-0 detId early return
0283         if (std::find(p0detIds.begin(), p0detIds.end(), det) != p0detIds.end()) {
0284           return true;
0285         }
0286       }
0287       return false;
0288     }
0289   };
0290 
0291   using SiPixelBPixQualityMap = SiPixelQualityMap<SiPixelPI::t_barrel>;
0292   using SiPixelFPixQualityMap = SiPixelQualityMap<SiPixelPI::t_forward>;
0293   using SiPixelFullQualityMap = SiPixelQualityMap<SiPixelPI::t_all>;
0294 
0295   /************************************************
0296    occupancy style map whole Pixel, difference of payloads
0297   *************************************************/
0298   template <SiPixelPI::DetType myType, IOVMultiplicity nIOVs, int ntags>
0299   class SiPixelQualityMapComparisonBase : public PlotImage<SiPixelQuality, nIOVs, ntags> {
0300   public:
0301     SiPixelQualityMapComparisonBase()
0302         : PlotImage<SiPixelQuality, nIOVs, ntags>(
0303               Form("SiPixelQuality %s Pixel Map", SiPixelPI::DetNames[myType].c_str())),
0304           m_trackerTopo{StandaloneTrackerTopology::fromTrackerParametersXMLFile(
0305               edm::FileInPath("Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml").fullPath())} {}
0306 
0307     bool fill() override {
0308       // trick to deal with the multi-ioved tag and two tag case at the same time
0309       auto theIOVs = PlotBase::getTag<0>().iovs;
0310       auto f_tagname = PlotBase::getTag<0>().name;
0311       std::string l_tagname = "";
0312       auto firstiov = theIOVs.front();
0313       std::tuple<cond::Time_t, cond::Hash> lastiov;
0314 
0315       // we don't support (yet) comparison with more than 2 tags
0316       assert(this->m_plotAnnotations.ntags < 3);
0317 
0318       if (this->m_plotAnnotations.ntags == 2) {
0319         auto tag2iovs = PlotBase::getTag<1>().iovs;
0320         l_tagname = PlotBase::getTag<1>().name;
0321         lastiov = tag2iovs.front();
0322       } else {
0323         lastiov = theIOVs.back();
0324       }
0325 
0326       std::shared_ptr<SiPixelQuality> last_payload = this->fetchPayload(std::get<1>(lastiov));
0327       std::shared_ptr<SiPixelQuality> first_payload = this->fetchPayload(std::get<1>(firstiov));
0328 
0329       if (this->isPhase0(first_payload) || this->isPhase0(last_payload)) {
0330         edm::LogError("SiPixelQuality_PayloadInspector")
0331             << "SiPixelQuality comparison maps are not supported for non-Phase1 Pixel geometries !";
0332         TCanvas canvas("Canv", "Canv", 1200, 1000);
0333         SiPixelPI::displayNotSupported(canvas, 0);
0334         std::string fileName(this->m_imageFileName);
0335         canvas.SaveAs(fileName.c_str());
0336         return false;
0337       }
0338 
0339       Phase1PixelROCMaps theMap("", "#Delta payload A - payload B");
0340 
0341       gStyle->SetOptStat(0);
0342       //=========================
0343       TCanvas canvas("Summary", "Summary", 1200, k_height[myType]);
0344       canvas.cd();
0345 
0346       auto f_unpacked = SiPixelPI::unpack(std::get<0>(firstiov));
0347       auto l_unpacked = SiPixelPI::unpack(std::get<0>(lastiov));
0348 
0349       std::string f_IOVstring = (f_unpacked.first == 0)
0350                                     ? std::to_string(f_unpacked.second)
0351                                     : (std::to_string(f_unpacked.first) + "," + std::to_string(f_unpacked.second));
0352 
0353       std::string l_IOVstring = (l_unpacked.first == 0)
0354                                     ? std::to_string(l_unpacked.second)
0355                                     : (std::to_string(l_unpacked.first) + "," + std::to_string(l_unpacked.second));
0356 
0357       std::string headerText;
0358 
0359       if (this->m_plotAnnotations.ntags == 2) {
0360         headerText = fmt::sprintf(
0361             "#Delta #color[2]{A: %s, %s} - #color[4]{B: %s, %s}", f_tagname, f_IOVstring, l_tagname, l_IOVstring);
0362       } else {
0363         headerText =
0364             fmt::sprintf("%s, #Delta IOV #color[2]{A: %s} - #color[4]{B: %s} ", f_tagname, f_IOVstring, l_IOVstring);
0365       }
0366 
0367       switch (myType) {
0368         case SiPixelPI::t_barrel:
0369           theMap.drawBarrelMaps(canvas, headerText);
0370           break;
0371         case SiPixelPI::t_forward:
0372           theMap.drawForwardMaps(canvas, headerText);
0373           break;
0374         case SiPixelPI::t_all:
0375           theMap.drawMaps(canvas, headerText);
0376           break;
0377         default:
0378           throw cms::Exception("SiPixelQualityMapComparison")
0379               << "\nERROR: unrecognized Pixel Detector part " << std::endl;
0380       }
0381 
0382       // first loop on the first payload (newest)
0383       fillTheMapFromPayload(theMap, first_payload, false);
0384 
0385       // then loop on the second payload (oldest)
0386       fillTheMapFromPayload(theMap, last_payload, true);  // true will subtract
0387 
0388       std::string fileName(this->m_imageFileName);
0389       canvas.SaveAs(fileName.c_str());
0390 #ifdef MMDEBUG
0391       canvas.SaveAs("outAll.root");
0392 #endif
0393 
0394       return true;
0395     }
0396 
0397   private:
0398     static constexpr std::array<int, 3> k_height = {{1200, 600, 1600}};
0399     TrackerTopology m_trackerTopo;
0400 
0401     //_________________________________________________
0402     bool isPhase0(const std::shared_ptr<SiPixelQuality>& payload) {
0403       const auto mods = payload->getBadComponentList();
0404       SiPixelDetInfoFileReader reader =
0405           SiPixelDetInfoFileReader(edm::FileInPath(SiPixelDetInfoFileReader::kPh0DefaultFile).fullPath());
0406       const auto& p0detIds = reader.getAllDetIds();
0407 
0408       std::vector<uint32_t> ownDetIds;
0409       std::transform(mods.begin(),
0410                      mods.end(),
0411                      std::back_inserter(ownDetIds),
0412                      [](SiPixelQuality::disabledModuleType d) -> uint32_t { return d.DetID; });
0413 
0414       for (const auto& det : ownDetIds) {
0415         // if found at least one phase-0 detId early return
0416         if (std::find(p0detIds.begin(), p0detIds.end(), det) != p0detIds.end()) {
0417           return true;
0418         }
0419       }
0420       return false;
0421     }
0422 
0423     //____________________________________________________________________________________________
0424     void fillTheMapFromPayload(Phase1PixelROCMaps& theMap,
0425                                const std::shared_ptr<SiPixelQuality>& payload,
0426                                bool subtract) {
0427       const auto theDisabledModules = payload->getBadComponentList();
0428       for (const auto& mod : theDisabledModules) {
0429         int subid = DetId(mod.DetID).subdetId();
0430         if ((subid == PixelSubdetector::PixelBarrel && myType == SiPixelPI::t_barrel) ||
0431             (subid == PixelSubdetector::PixelEndcap && myType == SiPixelPI::t_forward) ||
0432             (myType == SiPixelPI::t_all)) {
0433           std::bitset<16> bad_rocs(mod.BadRocs);
0434           if (payload->IsModuleBad(mod.DetID)) {
0435             theMap.fillWholeModule(mod.DetID, (subtract ? -1. : 1.));
0436           } else {
0437             theMap.fillSelectedRocs(mod.DetID, bad_rocs, (subtract ? -1. : 1.));
0438           }
0439         }
0440       }
0441     }
0442   };
0443 
0444   using SiPixelBPixQualityMapCompareSingleTag = SiPixelQualityMapComparisonBase<SiPixelPI::t_barrel, MULTI_IOV, 1>;
0445   using SiPixelFPixQualityMapCompareSingleTag = SiPixelQualityMapComparisonBase<SiPixelPI::t_forward, MULTI_IOV, 1>;
0446   using SiPixelFullQualityMapCompareSingleTag = SiPixelQualityMapComparisonBase<SiPixelPI::t_all, MULTI_IOV, 1>;
0447   using SiPixelBPixQualityMapCompareTwoTags = SiPixelQualityMapComparisonBase<SiPixelPI::t_barrel, SINGLE_IOV, 2>;
0448   using SiPixelFPixQualityMapCompareTwoTags = SiPixelQualityMapComparisonBase<SiPixelPI::t_forward, SINGLE_IOV, 2>;
0449   using SiPixelFullQualityMapCompareTwoTags = SiPixelQualityMapComparisonBase<SiPixelPI::t_all, SINGLE_IOV, 2>;
0450 
0451 }  // namespace
0452 
0453 // Register the classes as boost python plugin
0454 PAYLOAD_INSPECTOR_MODULE(SiPixelQuality) {
0455   PAYLOAD_INSPECTOR_CLASS(SiPixelQualityTest);
0456   PAYLOAD_INSPECTOR_CLASS(SiPixelQualityBadRocsSummary);
0457   PAYLOAD_INSPECTOR_CLASS(SiPixelQualityBadRocsTimeHistory);
0458   //PAYLOAD_INSPECTOR_CLASS(SiPixelQualityDebugger);
0459   PAYLOAD_INSPECTOR_CLASS(SiPixelBPixQualityMap);
0460   PAYLOAD_INSPECTOR_CLASS(SiPixelFPixQualityMap);
0461   PAYLOAD_INSPECTOR_CLASS(SiPixelFullQualityMap);
0462   PAYLOAD_INSPECTOR_CLASS(SiPixelBPixQualityMapCompareSingleTag);
0463   PAYLOAD_INSPECTOR_CLASS(SiPixelFPixQualityMapCompareSingleTag);
0464   PAYLOAD_INSPECTOR_CLASS(SiPixelFullQualityMapCompareSingleTag);
0465   PAYLOAD_INSPECTOR_CLASS(SiPixelBPixQualityMapCompareTwoTags);
0466   PAYLOAD_INSPECTOR_CLASS(SiPixelFPixQualityMapCompareTwoTags);
0467   PAYLOAD_INSPECTOR_CLASS(SiPixelFullQualityMapCompareTwoTags);
0468 }