Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-07-26 17:49:08

0001 // -*- C++ -*-
0002 // Package:    SiPixelPhase1RawDataErrorComparator
0003 // Class:      SiPixelPhase1RawDataErrorComparator
0004 //
0005 /**\class SiPixelPhase1RawDataErrorComparator SiPixelPhase1RawDataErrorComparator.cc
0006 */
0007 //
0008 // Author: Marco Musich
0009 //
0010 #include "CondFormats/DataRecord/interface/SiPixelFedCablingMapRcd.h"
0011 #include "CondFormats/SiPixelObjects/interface/SiPixelFedCablingMap.h"
0012 #include "CondFormats/SiPixelObjects/interface/SiPixelFedCablingTree.h"
0013 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0014 #include "DQMServices/Core/interface/DQMStore.h"
0015 #include "DQMServices/Core/interface/MonitorElement.h"
0016 #include "DataFormats/Common/interface/DetSetVector.h"
0017 #include "DataFormats/Common/interface/Handle.h"
0018 #include "DataFormats/FEDRawData/interface/FEDNumbering.h"
0019 #include "DataFormats/SiPixelDetId/interface/PixelFEDChannel.h"
0020 #include "DataFormats/SiPixelDigi/interface/PixelDigi.h"
0021 #include "EventFilter/SiPixelRawToDigi/interface/PixelDataFormatter.h"
0022 #include "FWCore/Framework/interface/Event.h"
0023 #include "FWCore/Framework/interface/Frameworkfwd.h"
0024 #include "FWCore/Framework/interface/MakerMacros.h"
0025 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0026 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0027 #include "FWCore/Utilities/interface/InputTag.h"
0028 #include "Geometry/CommonTopologies/interface/SimplePixelTopology.h"
0029 // for string manipulations
0030 #include <fmt/printf.h>
0031 
0032 namespace {
0033   // same logic used for the MTV:
0034   // cf https://github.com/cms-sw/cmssw/blob/master/Validation/RecoTrack/src/MTVHistoProducerAlgoForTracker.cc
0035   typedef dqm::reco::DQMStore DQMStore;
0036 
0037   void setBinLog(TAxis* axis) {
0038     int bins = axis->GetNbins();
0039     float from = axis->GetXmin();
0040     float to = axis->GetXmax();
0041     float width = (to - from) / bins;
0042     std::vector<float> new_bins(bins + 1, 0);
0043     for (int i = 0; i <= bins; i++) {
0044       new_bins[i] = TMath::Power(10, from + i * width);
0045     }
0046     axis->Set(bins, new_bins.data());
0047   }
0048 
0049   void setBinLogX(TH1* h) {
0050     TAxis* axis = h->GetXaxis();
0051     setBinLog(axis);
0052   }
0053   void setBinLogY(TH1* h) {
0054     TAxis* axis = h->GetYaxis();
0055     setBinLog(axis);
0056   }
0057 
0058   template <typename... Args>
0059   dqm::reco::MonitorElement* make2DIfLog(DQMStore::IBooker& ibook, bool logx, bool logy, Args&&... args) {
0060     auto h = std::make_unique<TH2I>(std::forward<Args>(args)...);
0061     if (logx)
0062       setBinLogX(h.get());
0063     if (logy)
0064       setBinLogY(h.get());
0065     const auto& name = h->GetName();
0066     return ibook.book2I(name, h.release());
0067   }
0068 
0069   //errorType - a number (25-38) indicating the type of error recorded.
0070   enum SiPixelFEDErrorCodes {
0071     k_FED25 = 25,  // 25 indicates an invalid ROC of 25
0072     k_FED26 = 26,  // 26 indicates a gap word
0073     k_FED27 = 27,  // 27 indicates a dummy word
0074     k_FED28 = 28,  // 28 indicates a FIFO full error
0075     k_FED29 = 29,  // 29 indicates a timeout error
0076     k_FED30 = 30,  // 30 indicates a TBM error trailer
0077     k_FED31 = 31,  // 31 indicates an event number error (TBM and FED event number mismatch)
0078     k_FED32 = 32,  // 32 indicates an incorrectly formatted Slink Header
0079     k_FED33 = 33,  // 33 indicates an incorrectly formatted Slink Trailer
0080     k_FED34 = 34,  // 34 indicates the evt size encoded in Slink Trailer is different than size found at raw2digi
0081     k_FED35 = 35,  // 35 indicates an invalid FED channel number
0082     k_FED36 = 36,  // 36 indicates an invalid ROC value
0083     k_FED37 = 37,  // 37 indicates an invalid dcol or pixel value
0084     k_FED38 = 38   // 38 indicates the pixels on a ROC weren't read out from lowest to highest row and dcol value
0085   };
0086 
0087   using MapToCodes = std::map<SiPixelFEDErrorCodes, std::string>;
0088 
0089   const MapToCodes errorCodeToStringMap = {{k_FED25, "FED25 error"},
0090                                            {k_FED26, "FED26 error"},
0091                                            {k_FED27, "FED27 error"},
0092                                            {k_FED28, "FED28 error"},
0093                                            {k_FED29, "FED29 error"},
0094                                            {k_FED30, "FED30 error"},
0095                                            {k_FED31, "FED31 error"}};
0096 
0097   const MapToCodes errorCodeToTypeMap = {{k_FED25, "ROC of 25"},
0098                                          {k_FED26, "Gap word"},
0099                                          {k_FED27, "Dummy word"},
0100                                          {k_FED28, "FIFO full"},
0101                                          {k_FED29, "Timeout"},
0102                                          {k_FED30, "TBM error trailer"},
0103                                          {k_FED31, "Event number"},
0104                                          {k_FED32, "Slink header"},
0105                                          {k_FED33, "Slink trailer"},
0106                                          {k_FED34, "Event size"},
0107                                          {k_FED35, "Invalid channel#"},
0108                                          {k_FED36, "ROC value"},
0109                                          {k_FED37, "Dcol or pixel value"},
0110                                          {k_FED38, "Readout order"}};
0111 }  // namespace
0112 
0113 class SiPixelPhase1RawDataErrorComparator : public DQMEDAnalyzer {
0114 public:
0115   explicit SiPixelPhase1RawDataErrorComparator(const edm::ParameterSet&);
0116   ~SiPixelPhase1RawDataErrorComparator() override = default;
0117   void bookHistograms(DQMStore::IBooker& ibooker, edm::Run const& iRun, edm::EventSetup const& iSetup) override;
0118   void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override;
0119   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0120 
0121 private:
0122   const edm::InputTag pixelErrorSrcGPU_;
0123   const edm::InputTag pixelErrorSrcCPU_;
0124   const edm::EDGetTokenT<edm::DetSetVector<SiPixelRawDataError>> tokenErrorsGPU_;
0125   const edm::EDGetTokenT<edm::DetSetVector<SiPixelRawDataError>> tokenErrorsCPU_;
0126   const std::string topFolderName_;
0127 
0128   MonitorElement* h_totFEDErrors_;
0129   MonitorElement* h_FEDerrorVsFEDIdUnbalance_;
0130   std::unordered_map<SiPixelFEDErrorCodes, MonitorElement*> h_nFEDErrors_;
0131 
0132   // name of the plugins
0133   static constexpr const char* kName = "SiPixelPhase1RawDataErrorComparator";
0134 
0135   // Define the dimensions of the 2D array
0136   static constexpr int nFEDs = FEDNumbering::MAXSiPixeluTCAFEDID - FEDNumbering::MINSiPixeluTCAFEDID;
0137   static constexpr int nErrors = k_FED38 - k_FED25;
0138 };
0139 
0140 //
0141 // constructors
0142 //
0143 SiPixelPhase1RawDataErrorComparator::SiPixelPhase1RawDataErrorComparator(const edm::ParameterSet& iConfig)
0144     : pixelErrorSrcGPU_(iConfig.getParameter<edm::InputTag>("pixelErrorSrcGPU")),
0145       pixelErrorSrcCPU_(iConfig.getParameter<edm::InputTag>("pixelErrorSrcCPU")),
0146       tokenErrorsGPU_(consumes<edm::DetSetVector<SiPixelRawDataError>>(pixelErrorSrcGPU_)),
0147       tokenErrorsCPU_(consumes<edm::DetSetVector<SiPixelRawDataError>>(pixelErrorSrcCPU_)),
0148       topFolderName_(iConfig.getParameter<std::string>("topFolderName")) {}
0149 
0150 //
0151 // -- Analyze
0152 //
0153 void SiPixelPhase1RawDataErrorComparator::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0154   std::map<int, int> countsOnCPU;
0155   std::map<int, int> countsOnGPU;
0156 
0157   std::array<std::array<int, nErrors>, nFEDs> countsMatrixOnCPU;
0158   std::array<std::array<int, nErrors>, nFEDs> countsMatrixOnGPU;
0159 
0160   // initialize the counts for FED/error matrix
0161   for (int i = 0; i < nFEDs; i++) {
0162     for (int j = 0; j < nErrors; j++) {
0163       countsMatrixOnCPU[i][j] = 0;
0164       countsMatrixOnGPU[i][j] = 0;
0165     }
0166   }
0167 
0168   // initialize the counts for errors per type scatter plots
0169   for (unsigned int j = k_FED25; j <= k_FED31; j++) {
0170     countsOnCPU[j] = 0.;
0171     countsOnGPU[j] = 0.;
0172   }
0173 
0174   // check upfront if the error collection is present
0175   edm::Handle<edm::DetSetVector<SiPixelRawDataError>> inputFromCPU;
0176   iEvent.getByToken(tokenErrorsCPU_, inputFromCPU);
0177   if (!inputFromCPU.isValid()) {
0178     edm::LogError(kName) << "reference (cpu) SiPixelRawDataErrors collection (" << pixelErrorSrcCPU_.encode()
0179                          << ") not found; \n"
0180                          << "the comparison will not run.";
0181     return;
0182   }
0183 
0184   edm::Handle<edm::DetSetVector<SiPixelRawDataError>> inputFromGPU;
0185   iEvent.getByToken(tokenErrorsGPU_, inputFromGPU);
0186   if (!inputFromGPU.isValid()) {
0187     edm::LogError(kName) << "target (gpu) SiPixelRawDataErrors collection (" << pixelErrorSrcGPU_.encode()
0188                          << ") not found; \n"
0189                          << "the comparison will not run.";
0190     return;
0191   }
0192 
0193   // fill the counters on host
0194   uint errorsOnCPU{0};
0195   for (auto it = inputFromCPU->begin(); it != inputFromCPU->end(); ++it) {
0196     for (auto& siPixelRawDataError : *it) {
0197       int fed = siPixelRawDataError.getFedId();
0198       int type = siPixelRawDataError.getType();
0199       DetId id = it->detId();
0200 
0201       // fill the error matrices for CPU
0202       countsOnCPU[type] += 1;
0203       countsMatrixOnCPU[fed - FEDNumbering::MINSiPixeluTCAFEDID][type - k_FED25] += 1;
0204 
0205       LogDebug(kName) << " (on cpu) FED: " << fed << " detid: " << id.rawId() << " type:" << type;
0206       errorsOnCPU++;
0207     }
0208   }
0209 
0210   // fill the counters on device
0211   uint errorsOnGPU{0};
0212   for (auto it = inputFromGPU->begin(); it != inputFromGPU->end(); ++it) {
0213     for (auto& siPixelRawDataError : *it) {
0214       int fed = siPixelRawDataError.getFedId();
0215       int type = siPixelRawDataError.getType();
0216       DetId id = it->detId();
0217 
0218       // fill the error matrices for GPU
0219       countsOnGPU[type] += 1;
0220       countsMatrixOnGPU[fed - FEDNumbering::MINSiPixeluTCAFEDID][type - k_FED25] += 1;
0221 
0222       LogDebug(kName) << " (on gpu) FED: " << fed << " detid: " << id.rawId() << " type:" << type;
0223       errorsOnGPU++;
0224     }
0225   }
0226 
0227   LogDebug(kName) << " on gpu found: " << errorsOnGPU << " on cpu found: " << errorsOnCPU;
0228 
0229   h_totFEDErrors_->Fill(errorsOnCPU, errorsOnGPU);
0230 
0231   // fill the correlations per error type
0232   for (unsigned int j = k_FED25; j <= k_FED31; j++) {
0233     const SiPixelFEDErrorCodes code = static_cast<SiPixelFEDErrorCodes>(j);
0234     h_nFEDErrors_[code]->Fill(std::min(1000, countsOnCPU[j]), std::min(1000, countsOnGPU[j]));
0235   }
0236 
0237   // fill the error unbalance per FEDid per error type
0238   for (int i = 0; i < nFEDs; i++) {
0239     for (int j = 0; j < nErrors; j++) {
0240       if (countsMatrixOnGPU[i][j] != 0 || countsMatrixOnCPU[i][j] != 0) {
0241         edm::LogVerbatim(kName) << "FED: " << i + FEDNumbering::MINSiPixeluTCAFEDID << " error: " << j + k_FED25
0242                                 << " | GPU counts: " << countsMatrixOnGPU[i][j]
0243                                 << " CPU counts:" << countsMatrixOnCPU[i][j];
0244         h_FEDerrorVsFEDIdUnbalance_->Fill(
0245             j, i + FEDNumbering::MINSiPixeluTCAFEDID, countsMatrixOnGPU[i][j] - countsMatrixOnCPU[i][j]);
0246       }
0247     }
0248   }
0249 }
0250 
0251 //
0252 // -- Book Histograms
0253 //
0254 void SiPixelPhase1RawDataErrorComparator::bookHistograms(DQMStore::IBooker& iBook,
0255                                                          edm::Run const& iRun,
0256                                                          edm::EventSetup const& iSetup) {
0257   iBook.cd();
0258   iBook.setCurrentFolder(topFolderName_);
0259 
0260   h_FEDerrorVsFEDIdUnbalance_ =
0261       iBook.book2I("FEErrorVsFEDIdUnbalance",
0262                    "difference (GPU-CPU) of FED errors per FEDid per error type;;FED Id number;GPU counts - CPU counts",
0263                    nErrors,
0264                    -0.5,
0265                    nErrors - 0.5,
0266                    nFEDs,
0267                    static_cast<double>(FEDNumbering::MINSiPixeluTCAFEDID) - 0.5,
0268                    static_cast<double>(FEDNumbering::MAXSiPixeluTCAFEDID) - 0.5);
0269   for (int j = 0; j < nErrors; j++) {
0270     const auto& errorCode = static_cast<SiPixelFEDErrorCodes>(j + k_FED25);
0271     h_FEDerrorVsFEDIdUnbalance_->setBinLabel(j + 1, errorCodeToTypeMap.at(errorCode));
0272   }
0273 
0274   h_totFEDErrors_ = make2DIfLog(iBook,
0275                                 true,
0276                                 true,
0277                                 "nTotalFEDErrors",
0278                                 "n. of total Pixel FEDErrors per event; CPU; GPU",
0279                                 200,
0280                                 log10(0.5),
0281                                 log10(1000.),
0282                                 200,
0283                                 log10(0.5),
0284                                 log10(1000.));
0285 
0286   for (const auto& element : errorCodeToStringMap) {
0287     h_nFEDErrors_[element.first] = iBook.book2I(fmt::sprintf("nFED%i_Errors", static_cast<int>(element.first)),
0288                                                 fmt::sprintf("n. of %ss per event; CPU; GPU", element.second),
0289                                                 1000,
0290                                                 -0.5,
0291                                                 1000.5,
0292                                                 1000,
0293                                                 -0.5,
0294                                                 1000.5);
0295   }
0296 }
0297 
0298 void SiPixelPhase1RawDataErrorComparator::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0299   edm::ParameterSetDescription desc;
0300   desc.add<edm::InputTag>("pixelErrorSrcGPU", edm::InputTag("siPixelDigis@cuda"))
0301       ->setComment("input GPU SiPixel FED errors");
0302   desc.add<edm::InputTag>("pixelErrorSrcCPU", edm::InputTag("siPixelDigis@cpu"))
0303       ->setComment("input CPU SiPixel FED errors");
0304   desc.add<std::string>("topFolderName", "SiPixelHeterogeneous/PixelErrorCompareGPUvsCPU");
0305   descriptions.addWithDefaultLabel(desc);
0306 }
0307 
0308 DEFINE_FWK_MODULE(SiPixelPhase1RawDataErrorComparator);