Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2021-07-15 22:42:25

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