Back to home page

Project CMSSW displayed by LXR

 
 

    


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

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