Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2021-09-07 22:41:38

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 "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
0032 #include "DataFormats/DetId/interface/DetId.h"
0033 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0034 #include "DQM/TrackerRemapper/interface/Phase1PixelROCMaps.h"
0035 
0036 #include <memory>
0037 #include <sstream>
0038 #include <iostream>
0039 #include <fmt/printf.h>
0040 
0041 // include ROOT
0042 #include "TH2F.h"
0043 #include "TLegend.h"
0044 #include "TCanvas.h"
0045 #include "TLine.h"
0046 #include "TGraph.h"
0047 #include "TGaxis.h"
0048 #include "TStyle.h"
0049 #include "TLatex.h"
0050 #include "TPave.h"
0051 #include "TPaveStats.h"
0052 
0053 namespace {
0054 
0055   using namespace cond::payloadInspector;
0056 
0057   /************************************************
0058   1d histogram of SiPixelFEDChannelContainer of 1 IOV 
0059   *************************************************/
0060 
0061   template <SiPixelPI::DetType myType>
0062   class SiPixelFEDChannelContainerMap : public PlotImage<SiPixelFEDChannelContainer, SINGLE_IOV> {
0063   public:
0064     SiPixelFEDChannelContainerMap()
0065         : PlotImage<SiPixelFEDChannelContainer, SINGLE_IOV>("SiPixelFEDChannelContainer scenarios count"),
0066           m_trackerTopo{StandaloneTrackerTopology::fromTrackerParametersXMLFile(
0067               edm::FileInPath("Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml").fullPath())} {
0068       // for inputs
0069       PlotBase::addInputParam("Scenarios");
0070 
0071       // hardcoded connection to the MC cabling tag, though luck
0072       m_condDbCabling = "frontier://FrontierProd/CMS_CONDITIONS";
0073       m_CablingTagName = "SiPixelFedCablingMap_phase1_v7";
0074 
0075       m_connectionPool.setParameters(m_connectionPset);
0076       m_connectionPool.configure();
0077     }
0078 
0079     bool fill() override {
0080       std::vector<std::string> the_scenarios = {};
0081 
0082       auto paramValues = PlotBase::inputParamValues();
0083       auto ip = paramValues.find("Scenarios");
0084       if (ip != paramValues.end()) {
0085         auto input = ip->second;
0086         typedef boost::tokenizer<boost::char_separator<char>> tokenizer;
0087         boost::char_separator<char> sep{","};
0088         tokenizer tok{input, sep};
0089         for (const auto& t : tok) {
0090           the_scenarios.push_back(t);
0091         }
0092       } else {
0093         edm::LogWarning("SiPixelFEDChannelContainerTest")
0094             << "\n WARNING!!!! \n The needed parameter Scenarios has not been passed. Will use all the scenarios in "
0095                "the file!"
0096             << "\n Buckle your seatbelts... this might take a while... \n\n";
0097         the_scenarios.push_back("all");
0098       }
0099 
0100       Phase1PixelROCMaps theROCMap("");
0101 
0102       auto tag = PlotBase::getTag<0>();
0103       auto tagname = tag.name;
0104       auto iov = tag.iovs.front();
0105 
0106       // open db session for the cabling map
0107       edm::LogPrint("SiPixelFEDChannelContainerTest") << "[SiPixelFEDChannelContainerTest::" << __func__ << "] "
0108                                                       << "Query the condition database " << m_condDbCabling;
0109 
0110       cond::persistency::Session condDbSession = m_connectionPool.createSession(m_condDbCabling);
0111       condDbSession.transaction().start(true);
0112 
0113       // query the database
0114       edm::LogPrint("SiPixelFEDChannelContainerTest") << "[SiPixelFEDChannelContainerTest::" << __func__ << "] "
0115                                                       << "Reading IOVs from tag " << m_CablingTagName;
0116 
0117       const auto MIN_VAL = cond::timeTypeSpecs[cond::runnumber].beginValue;
0118       const auto MAX_VAL = cond::timeTypeSpecs[cond::runnumber].endValue;
0119 
0120       // get the list of payloads for the Cabling Map
0121       std::vector<std::tuple<cond::Time_t, cond::Hash>> m_cabling_iovs;
0122       condDbSession.readIov(m_CablingTagName).selectRange(MIN_VAL, MAX_VAL, m_cabling_iovs);
0123 
0124       std::vector<unsigned int> listOfCablingIOVs;
0125       std::transform(m_cabling_iovs.begin(),
0126                      m_cabling_iovs.end(),
0127                      std::back_inserter(listOfCablingIOVs),
0128                      [](std::tuple<cond::Time_t, cond::Hash> myIOV2) -> unsigned int { return std::get<0>(myIOV2); });
0129 
0130       edm::LogPrint("SiPixelFEDChannelContainerTest")
0131           << " Number of SiPixelFedCablngMap payloads: " << listOfCablingIOVs.size() << std::endl;
0132 
0133       auto it = std::find(
0134           listOfCablingIOVs.begin(), listOfCablingIOVs.end(), closest_from_below(listOfCablingIOVs, std::get<0>(iov)));
0135       int index = std::distance(listOfCablingIOVs.begin(), it);
0136 
0137       edm::LogPrint("SiPixelFEDChannelContainerTest")
0138           << " using the SiPixelFedCablingMap with hash: " << std::get<1>(m_cabling_iovs.at(index)) << std::endl;
0139 
0140       auto theCablingMap = condDbSession.fetchPayload<SiPixelFedCablingMap>(std::get<1>(m_cabling_iovs.at(index)));
0141       theCablingMap->initializeRocs();
0142       // auto theCablingTree = (*theCablingMap).cablingTree();
0143 
0144       //auto map = theCablingMap->det2fedMap();
0145       //for (const auto &element : map){
0146       //    std::cout << element.first << " " << element.second << std::endl;
0147       //}
0148 
0149       std::shared_ptr<SiPixelFEDChannelContainer> payload = fetchPayload(std::get<1>(iov));
0150       const auto& scenarioMap = payload->getScenarioMap();
0151 
0152       auto pIndexConverter = PixelIndices(numColumns, numRows);
0153 
0154       for (const auto& scenario : scenarioMap) {
0155         std::string scenName = scenario.first;
0156 
0157         if (std::find_if(the_scenarios.begin(), the_scenarios.end(), compareKeys(scenName)) != the_scenarios.end() ||
0158             the_scenarios[0] == "all") {
0159           edm::LogPrint("SiPixelFEDChannelContainerTest") << "\t Found Scenario: " << scenName << " ==> dumping it";
0160         } else {
0161           continue;
0162         }
0163 
0164         //if (strcmp(scenName.c_str(),"320824_103") != 0) continue;
0165 
0166         const auto& theDetSetBadPixelFedChannels = payload->getDetSetBadPixelFedChannels(scenName);
0167         for (const auto& disabledChannels : *theDetSetBadPixelFedChannels) {
0168           const auto t_detid = disabledChannels.detId();
0169           int subid = DetId(t_detid).subdetId();
0170           LogDebug("SiPixelFEDChannelContainerTest") << fmt::sprintf("DetId : %i \n", t_detid) << std::endl;
0171 
0172           std::bitset<16> badRocsFromFEDChannels;
0173 
0174           for (const auto& ch : disabledChannels) {
0175             std::string toOut_ = fmt::sprintf("fed : %i | link : %2i | roc_first : %2i | roc_last: %2i \n",
0176                                               ch.fed,
0177                                               ch.link,
0178                                               ch.roc_first,
0179                                               ch.roc_last);
0180 
0181             LogDebug("SiPixelFEDChannelContainerTest") << toOut_ << std::endl;
0182             const std::vector<sipixelobjects::CablingPathToDetUnit>& path =
0183                 theCablingMap->pathToDetUnit(disabledChannels.detId());
0184             for (unsigned int i_roc = ch.roc_first; i_roc <= ch.roc_last; ++i_roc) {
0185               for (const auto p : path) {
0186                 const sipixelobjects::PixelROC* myroc = theCablingMap->findItem(p);
0187                 if (myroc->idInDetUnit() == static_cast<unsigned int>(i_roc)) {
0188                   sipixelobjects::LocalPixel::RocRowCol local = {39, 25};  //corresponding to center of ROC row,col
0189                   sipixelobjects::GlobalPixel global = myroc->toGlobal(sipixelobjects::LocalPixel(local));
0190                   int chipIndex(0), colROC(0), rowROC(0);
0191 
0192                   pIndexConverter.transformToROC(global.col, global.row, chipIndex, colROC, rowROC);
0193 
0194                   LogDebug("SiPixelFEDChannelContainerTest")
0195                       << " => i_roc:" << i_roc << "  " << global.col << "-" << global.row << " | => " << chipIndex
0196                       << " : (" << colROC << "," << rowROC << ")" << std::endl;
0197 
0198                   badRocsFromFEDChannels[chipIndex] = true;
0199                 }
0200               }
0201             }
0202           }
0203 
0204           LogDebug("SiPixelFEDChannelContainerTest") << badRocsFromFEDChannels << std::endl;
0205 
0206           auto myDetId = DetId(t_detid);
0207 
0208           if (subid == PixelSubdetector::PixelBarrel) {
0209             theROCMap.fillSelectedRocs(myDetId, badRocsFromFEDChannels, 1.);
0210           }  // if it's barrel
0211           else if (subid == PixelSubdetector::PixelEndcap) {
0212             theROCMap.fillSelectedRocs(myDetId, badRocsFromFEDChannels, 1.);
0213           }  // if it's endcap
0214           else {
0215             throw cms::Exception("LogicError") << "Unknown Pixel SubDet ID " << std::endl;
0216           }  // else nonsense
0217         }    // loop on the channels
0218       }      // loop on the scenarios
0219 
0220       gStyle->SetOptStat(0);
0221       //=========================
0222       TCanvas canvas("Summary", "Summary", 1200, k_height[myType]);
0223       canvas.cd();
0224 
0225       auto unpacked = SiPixelPI::unpack(std::get<0>(iov));
0226 
0227       std::string IOVstring = (unpacked.first == 0)
0228                                   ? std::to_string(unpacked.second)
0229                                   : (std::to_string(unpacked.first) + "," + std::to_string(unpacked.second));
0230 
0231       const auto headerText = fmt::sprintf("#color[4]{%s},  IOV: #color[4]{%s}", tagname, IOVstring);
0232 
0233       switch (myType) {
0234         case SiPixelPI::t_barrel:
0235           theROCMap.drawBarrelMaps(canvas, headerText);
0236           break;
0237         case SiPixelPI::t_forward:
0238           theROCMap.drawForwardMaps(canvas, headerText);
0239           break;
0240         case SiPixelPI::t_all:
0241           theROCMap.drawMaps(canvas, headerText);
0242           break;
0243         default:
0244           throw cms::Exception("SiPixelQualityMap") << "\nERROR: unrecognized Pixel Detector part " << std::endl;
0245       }
0246 
0247       std::string fileName(m_imageFileName);
0248       canvas.SaveAs(fileName.c_str());
0249 
0250       // close the DB session
0251       condDbSession.transaction().commit();
0252 
0253       return true;
0254     }
0255 
0256   public:
0257     inline unsigned int closest_from_above(std::vector<unsigned int> const& vec, unsigned int value) {
0258       auto const it = std::lower_bound(vec.begin(), vec.end(), value);
0259       return vec.at(it - vec.begin() - 1);
0260     }
0261 
0262     inline unsigned int closest_from_below(std::vector<unsigned int> const& vec, unsigned int value) {
0263       auto const it = std::upper_bound(vec.begin(), vec.end(), value);
0264       return vec.at(it - vec.begin() - 1);
0265     }
0266 
0267     // auxilliary check
0268     struct compareKeys {
0269       std::string key;
0270       compareKeys(std::string const& i) : key(i) {}
0271 
0272       bool operator()(std::string const& i) { return (key == i); }
0273     };
0274 
0275   private:
0276     // tough luck, we can only do phase-1...
0277     static constexpr int numColumns = 416;
0278     static constexpr int numRows = 160;
0279     static constexpr int n_rings = 2;
0280     static constexpr int n_layers = 4;
0281 
0282     // graphics
0283     static constexpr std::array<int, 3> k_height = {{1200, 600, 1600}};
0284 
0285     TrackerTopology m_trackerTopo;
0286     edm::ParameterSet m_connectionPset;
0287     cond::persistency::ConnectionPool m_connectionPool;
0288     std::string m_CablingTagName;
0289     std::string m_condDbCabling;
0290   };
0291 
0292   using SiPixelBPixFEDChannelContainerMap = SiPixelFEDChannelContainerMap<SiPixelPI::t_barrel>;
0293   using SiPixelFPixFEDChannelContainerMap = SiPixelFEDChannelContainerMap<SiPixelPI::t_forward>;
0294   using SiPixelFullFEDChannelContainerMap = SiPixelFEDChannelContainerMap<SiPixelPI::t_all>;
0295 
0296   /************************************************
0297   1d histogram of SiPixelFEDChannelContainer of 1 IOV
0298   *************************************************/
0299 
0300   class SiPixelFEDChannelContainerScenarios : public PlotImage<SiPixelFEDChannelContainer, SINGLE_IOV> {
0301   public:
0302     SiPixelFEDChannelContainerScenarios()
0303         : PlotImage<SiPixelFEDChannelContainer, SINGLE_IOV>("SiPixelFEDChannelContainer scenarios count") {}
0304     bool fill() override {
0305       auto tag = PlotBase::getTag<0>();
0306       auto tagname = tag.name;
0307       auto iov = tag.iovs.front();
0308 
0309       TGaxis::SetMaxDigits(3);
0310 
0311       std::shared_ptr<SiPixelFEDChannelContainer> payload = fetchPayload(std::get<1>(iov));
0312       std::vector<std::string> scenarios = payload->getScenarioList();
0313       sort(scenarios.begin(), scenarios.end());
0314 
0315       TCanvas canvas("Canv", "Canv", 1200, 1000);
0316       canvas.cd();
0317       canvas.SetGrid();
0318       auto h1 = std::make_unique<TH1F>("Count",
0319                                        "SiPixelFEDChannelContainer Bad Roc count;Scenario index;n. of bad ROCs",
0320                                        scenarios.size(),
0321                                        1,
0322                                        scenarios.size());
0323       h1->SetStats(false);
0324 
0325       canvas.SetTopMargin(0.06);
0326       canvas.SetBottomMargin(0.12);
0327       canvas.SetLeftMargin(0.12);
0328       canvas.SetRightMargin(0.05);
0329       canvas.Modified();
0330 
0331       int scenarioIndex = 0;
0332       for (const auto& scenario : scenarios) {
0333         scenarioIndex++;
0334         int badRocCount = 0;
0335         LogDebug("SiPixelFEDChannelContainerScenarios") << scenario << std::endl;
0336         auto badChannelCollection = payload->getDetSetBadPixelFedChannels(scenario);
0337         for (const auto& disabledChannels : *badChannelCollection) {
0338           for (const auto& ch : disabledChannels) {
0339             int local_bad_rocs = ch.roc_last - ch.roc_first;
0340             badRocCount += local_bad_rocs;
0341           }  // loop on the channels
0342         }    // loop on the DetSetVector
0343 
0344         h1->SetBinContent(scenarioIndex, badRocCount);
0345       }  // loop on scenarios
0346 
0347       TGaxis::SetExponentOffset(-0.1, 0.01, "y");    // Y offset
0348       TGaxis::SetExponentOffset(-0.03, -0.10, "x");  // Y and Y offset for X axis
0349 
0350       h1->SetTitle("");
0351       h1->GetYaxis()->SetRangeUser(0., h1->GetMaximum() * 1.30);
0352       h1->SetFillColor(kRed);
0353       h1->SetMarkerStyle(20);
0354       h1->SetMarkerSize(1);
0355       h1->Draw("bar2");
0356 
0357       SiPixelPI::makeNicePlotStyle(h1.get());
0358 
0359       canvas.Update();
0360 
0361       TLegend legend = TLegend(0.30, 0.88, 0.95, 0.94);
0362       //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
0363 
0364       legend.SetHeader(fmt::sprintf("Payload hash: #bf{%s}", std::get<1>(iov)).c_str(), "C");
0365       legend.AddEntry(h1.get(), fmt::sprintf("total scenarios: #bf{%s}", std::to_string(scenarioIndex)).c_str(), "F");
0366       legend.SetTextSize(0.025);
0367       legend.Draw("same");
0368 
0369       auto ltx = TLatex();
0370       ltx.SetTextFont(62);
0371       //ltx.SetTextColor(kBlue);
0372       //ltx.SetTextAlign(11);
0373       ltx.SetTextSize(0.040);
0374       ltx.DrawLatexNDC(
0375           gPad->GetLeftMargin(),
0376           1 - gPad->GetTopMargin() + 0.01,
0377           fmt::sprintf("#color[4]{%s} IOV: #color[4]{%s}", tagname, std::to_string(std::get<0>(iov))).c_str());
0378 
0379       std::string fileName(m_imageFileName);
0380       canvas.SaveAs(fileName.c_str());
0381 
0382       return true;
0383 
0384     }  // fill
0385   };
0386 
0387 }  // namespace
0388 
0389 // Register the classes as boost python plugin
0390 PAYLOAD_INSPECTOR_MODULE(SiPixelFEDChannelContainer) {
0391   PAYLOAD_INSPECTOR_CLASS(SiPixelBPixFEDChannelContainerMap);
0392   PAYLOAD_INSPECTOR_CLASS(SiPixelFPixFEDChannelContainerMap);
0393   PAYLOAD_INSPECTOR_CLASS(SiPixelFullFEDChannelContainerMap);
0394   PAYLOAD_INSPECTOR_CLASS(SiPixelFEDChannelContainerScenarios);
0395 }