Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-09-07 04:35:07

0001 // C++ standard
0002 #include <string>
0003 // ROOT
0004 #include "TMath.h"
0005 // CMSSW FW
0006 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0007 #include "FWCore/Framework/interface/MakerMacros.h"
0008 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0009 // CMSSW DataFormats
0010 #include "DataFormats/Common/interface/ConditionsInEdm.h"
0011 #include "DataFormats/Common/interface/DetSetVector.h"
0012 #include "DataFormats/FEDRawData/interface/FEDNumbering.h"
0013 #include "DataFormats/SiPixelCluster/interface/SiPixelCluster.h"
0014 #include "DataFormats/TrackerCommon/interface/PixelBarrelName.h"
0015 #include "DataFormats/TrackerCommon/interface/PixelEndcapName.h"
0016 #include "DataFormats/SiPixelDigi/interface/PixelDigi.h"
0017 // "FED error 25"
0018 #include "DataFormats/SiPixelDetId/interface/PixelFEDChannel.h"
0019 // CMSSW CondFormats
0020 #include "CondFormats/RunInfo/interface/RunSummary.h"
0021 #include "CondFormats/RunInfo/interface/RunInfo.h"
0022 #include "CondFormats/DataRecord/interface/RunSummaryRcd.h"
0023 #include "CondFormats/SiPixelObjects/interface/SiPixelFrameConverter.h"
0024 #include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
0025 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0026 
0027 // Header file
0028 #include "CalibTracker/SiPixelQuality/plugins/SiPixelStatusProducer.h"
0029 
0030 SiPixelStatusProducer::SiPixelStatusProducer(const edm::ParameterSet& iConfig, SiPixelStatusCache const* iCache) {
0031   //NOTE: Token for all stream replicas are identical and constructors for the replicas are called
0032   // sequentially so there is no race condition.
0033   iCache->trackerGeometryToken_ = esConsumes<edm::Transition::BeginRun>();
0034   iCache->trackerTopologyToken_ = esConsumes<edm::Transition::BeginRun>();
0035   iCache->siPixelFedCablingMapToken_ = esConsumes<edm::Transition::BeginRun>();
0036 
0037   /* badPixelFEDChannelCollections */
0038   std::vector<edm::InputTag> badPixelFEDChannelCollectionLabels =
0039       iConfig.getParameter<edm::ParameterSet>("SiPixelStatusProducerParameters")
0040           .getParameter<std::vector<edm::InputTag>>("badPixelFEDChannelCollections");
0041   for (auto& t : badPixelFEDChannelCollectionLabels)
0042     theBadPixelFEDChannelsTokens_.push_back(consumes<PixelFEDChannelCollection>(t));
0043 
0044   /* pixel clusters */
0045   fPixelClusterLabel_ = iConfig.getParameter<edm::ParameterSet>("SiPixelStatusProducerParameters")
0046                             .getUntrackedParameter<edm::InputTag>("pixelClusterLabel");
0047   fSiPixelClusterToken_ = consumes<edmNew::DetSetVector<SiPixelCluster>>(fPixelClusterLabel_);
0048 
0049   //debug_ = iConfig.getUntrackedParameter<bool>("debug");
0050 
0051   /* register products */
0052   produces<SiPixelDetectorStatus, edm::Transition::EndLuminosityBlock>("siPixelStatus");
0053 }
0054 
0055 SiPixelStatusProducer::~SiPixelStatusProducer() {}
0056 
0057 //--------------------------------------------------------------------------------------------------
0058 
0059 std::shared_ptr<SiPixelTopoFinder> SiPixelStatusProducer::globalBeginRun(edm::Run const& iRun,
0060                                                                          edm::EventSetup const& iSetup,
0061                                                                          GlobalCache const* iCache) {
0062   const TrackerGeometry* trackerGeometry = &iSetup.getData(iCache->trackerGeometryToken_);
0063   const TrackerTopology* trackerTopology = &iSetup.getData(iCache->trackerTopologyToken_);
0064   const SiPixelFedCablingMap* cablingMap = &iSetup.getData(iCache->siPixelFedCablingMapToken_);
0065 
0066   auto returnValue = std::make_shared<SiPixelTopoFinder>();
0067 
0068   returnValue->init(trackerGeometry, trackerTopology, cablingMap);
0069   return returnValue;
0070 }
0071 
0072 void SiPixelStatusProducer::beginRun(edm::Run const&, edm::EventSetup const&) {
0073   /*Is it good to pass the objects stored in runCache to set class private members values  *
0074      or just call runCahche every time in the calss function?*/
0075 
0076   edm::LogInfo("SiPixelStatusProducer") << "beginRun: update the std::map for pixel geo/topo " << std::endl;
0077   /* update the std::map for pixel geo/topo */
0078   /* vector of all <int> detIds */
0079   fDetIds_ = runCache()->getDetIds();  //getDetIds();
0080   /* ROC size (number of row, number of columns for each det id) */
0081   fSensors_ = runCache()->getSensors();
0082   /* the roc layout on a module */
0083   fSensorLayout_ = runCache()->getSensorLayout();
0084   /* fedId as a function of detId */
0085   fFedIds_ = runCache()->getFedIds();
0086   /* map the index ROC to rocId */
0087   fRocIds_ = runCache()->getRocIds();
0088 }
0089 
0090 void SiPixelStatusProducer::beginLuminosityBlock(edm::LuminosityBlock const&, edm::EventSetup const&) {
0091   edm::LogInfo("SiPixelStatusProducer") << "beginlumi instance" << std::endl;
0092 
0093   /* initialize fDet_ with a set of modules(detIds) and clean the fFEDerror25_ */
0094   fDet_ = SiPixelDetectorStatus();
0095   for (unsigned int itDetId = 0; itDetId < fDetIds_.size(); ++itDetId) {
0096     int detid = fDetIds_[itDetId];
0097     int nrocs = fSensorLayout_[detid].first * fSensorLayout_[detid].second;
0098 
0099     fDet_.addModule(detid, nrocs);
0100   }
0101 
0102   fFEDerror25_.clear();
0103   ftotalevents_ = 0;
0104 }
0105 
0106 void SiPixelStatusProducer::accumulate(edm::Event const& iEvent, edm::EventSetup const&) {
0107   edm::LogInfo("SiPixelStatusProducer") << "start cluster analyzer " << std::endl;
0108 
0109   /* count number of events for the current module instance in the luminosityBlock */
0110   ftotalevents_++;
0111 
0112   /* ----------------------------------------------------------------------
0113      -- Pixel cluster analysis
0114      ----------------------------------------------------------------------*/
0115 
0116   edm::Handle<edmNew::DetSetVector<SiPixelCluster>> hClusterColl;
0117   if (!iEvent.getByToken(fSiPixelClusterToken_, hClusterColl)) {
0118     edm::LogWarning("SiPixelStatusProducer")
0119         << " edmNew::DetSetVector<SiPixelCluster> " << fPixelClusterLabel_ << " does not exist!" << std::endl;
0120     return;
0121   }
0122 
0123   iEvent.getByToken(fSiPixelClusterToken_, hClusterColl);
0124 
0125   if (hClusterColl.isValid()) {
0126     for (const auto& clusters : *hClusterColl) { /*loop over different clusters in a clusters vector (module)*/
0127       for (const auto& clu : clusters) {         /*loop over cluster in a given detId (module)*/
0128         int detid = clusters.detId();
0129         int rowsperroc = fSensors_[detid].first;
0130         int colsperroc = fSensors_[detid].second;
0131 
0132         //int nROCrows    = fSensorLayout_[detid].first;
0133         int nROCcolumns = fSensorLayout_[detid].second;
0134 
0135         int roc(-1);
0136         std::map<int, int> rocIds_detid;
0137         if (fRocIds_.find(detid) != fRocIds_.end()) {
0138           rocIds_detid = fRocIds_[detid];
0139         }
0140 
0141         /* A module is made with a few ROCs
0142            Need to convert global row/column (on a module) to local row/column (on a ROC) */
0143         const std::vector<SiPixelCluster::Pixel>& pixvector = clu.pixels();
0144         for (unsigned int i = 0; i < pixvector.size(); ++i) {
0145           int mr0 = pixvector[i].x; /* constant column direction is along x-axis */
0146           int mc0 = pixvector[i].y; /* constant row direction is along y-axis */
0147 
0148           int irow = mr0 / rowsperroc;
0149           int icol = mc0 / colsperroc;
0150 
0151           int key = indexROC(irow, icol, nROCcolumns);
0152           if (rocIds_detid.find(key) != rocIds_detid.end()) {
0153             roc = rocIds_detid[key];
0154           }
0155 
0156           fDet_.fillDIGI(detid, roc);
0157 
0158         } /* loop over pixels in a cluster */
0159 
0160       } /* loop over cluster in a detId (module) */
0161 
0162     } /* loop over detId-grouped clusters in cluster detId-grouped clusters-vector* */
0163 
0164   } /* hClusterColl.isValid() */
0165   else {
0166     edm::LogWarning("SiPixelStatusProducer")
0167         << " edmNew::DetSetVector<SiPixelCluster> " << fPixelClusterLabel_ << " is NOT Valid!" << std::endl;
0168   }
0169 
0170   /*|||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||*/
0171 
0172   /* FEDerror25 per-ROC per-event */
0173   edm::Handle<PixelFEDChannelCollection> pixelFEDChannelCollectionHandle;
0174 
0175   /* look over different resouces of tokens */
0176   for (const edm::EDGetTokenT<PixelFEDChannelCollection>& tk : theBadPixelFEDChannelsTokens_) {
0177     if (!iEvent.getByToken(tk, pixelFEDChannelCollectionHandle)) {
0178       edm::LogWarning("SiPixelStatusProducer")
0179           << " PixelFEDChannelCollection with index " << tk.index() << " does NOT exist!" << std::endl;
0180       continue;
0181     }
0182 
0183     iEvent.getByToken(tk, pixelFEDChannelCollectionHandle);
0184     if (!pixelFEDChannelCollectionHandle.isValid()) {
0185       edm::LogWarning("SiPixelStatusProducer")
0186           << " PixelFEDChannelCollection with index " << tk.index() << " is NOT valid!" << std::endl;
0187       continue;
0188     }
0189 
0190     /* FEDerror channels for the current events */
0191     std::map<int, std::vector<PixelFEDChannel>> tmpFEDerror25;
0192     for (const auto& disabledChannels : *pixelFEDChannelCollectionHandle) {
0193       /* loop over different PixelFED in a PixelFED vector (module) */
0194       for (const auto& ch : disabledChannels) {
0195         DetId detId = disabledChannels.detId();
0196         int detid = detId.rawId();
0197 
0198         if (ftotalevents_ == 1) {
0199           /* FEDerror25 channels for the "first" event in the lumi section (first for the current instance of the module) */
0200           fFEDerror25_[detid].push_back(ch);
0201         } else
0202           tmpFEDerror25[detid].push_back(ch);
0203 
0204       } /* loop over different PixelFED in a PixelFED vector (different channel for a module) */
0205 
0206     } /* loop over different (different DetId) PixelFED vectors in PixelFEDChannelCollection */
0207 
0208     /* Compare the current FEDerror list with the first event's FEDerror list
0209  *        and save the common channels */
0210     if (!tmpFEDerror25.empty() && !fFEDerror25_.empty()) {
0211       std::map<int, std::vector<PixelFEDChannel>>::iterator itFEDerror25;
0212       for (itFEDerror25 = fFEDerror25_.begin(); itFEDerror25 != fFEDerror25_.end(); itFEDerror25++) {
0213         int detid = itFEDerror25->first;
0214         if (tmpFEDerror25.find(detid) != tmpFEDerror25.end()) {
0215           std::vector<PixelFEDChannel> chs = itFEDerror25->second;
0216           std::vector<PixelFEDChannel> chs_tmp = tmpFEDerror25[detid];
0217 
0218           std::vector<PixelFEDChannel> chs_common;
0219           for (unsigned int ich = 0; ich < chs.size(); ich++) {
0220             PixelFEDChannel ch = chs[ich];
0221             /* loop over the current FEDerror25 channels, save the common FED channels */
0222             for (unsigned int ich_tmp = 0; ich_tmp < chs_tmp.size(); ich_tmp++) {
0223               PixelFEDChannel ch_tmp = chs_tmp[ich_tmp];
0224               if ((ch.fed == ch_tmp.fed) && (ch.link == ch_tmp.link)) { /* the same FED channel */
0225                 chs_common.push_back(ch);
0226                 break;
0227               }
0228             }
0229           }
0230           /* remove the full module from FEDerror25 list if no common channels are left */
0231           if (chs_common.empty())
0232             fFEDerror25_.erase(itFEDerror25);
0233           /* otherwise replace with the common channels */
0234           else {
0235             fFEDerror25_[detid].clear();
0236             fFEDerror25_[detid] = chs_common;
0237           }
0238         } else { /* remove the full module from FEDerror25 list if the module doesn't appear in the current event's FEDerror25 list */
0239           fFEDerror25_.erase(itFEDerror25);
0240         }
0241 
0242       } /* loop over modules that have FEDerror25 in the first event in the lumi section */
0243 
0244     } /* non-empty FEDerror lists */
0245 
0246   } /* look over different resouces of takens */
0247 
0248   /* Caveat
0249      no per-event collection put into iEvent
0250      If use produce() function and no collection is put into iEvent, produce() will not run in unScheduled mode
0251      Now since CMSSW_10_1_X, the accumulate() function will run whatsoever in the unScheduled mode
0252      Accumulate() is NOT available for releases BEFORE CMSSW_10_1_X */
0253 }
0254 
0255 void SiPixelStatusProducer::endLuminosityBlock(edm::LuminosityBlock const&, edm::EventSetup const&) {
0256   /* set total number of events through ftotalevents_ */
0257   fDet_.setNevents(ftotalevents_);
0258 
0259   if (ftotalevents_ > 0) {
0260     /* Add FEDerror25 information into SiPixelDetectorStatus fDet_ for FED channels stored in fFEDerror25_ */
0261     if (!fFEDerror25_.empty()) {  // non-empty FEDerror25
0262       std::map<int, std::vector<PixelFEDChannel>>::iterator itFEDerror25;
0263       for (itFEDerror25 = fFEDerror25_.begin(); itFEDerror25 != fFEDerror25_.end();
0264            itFEDerror25++) {  // loop over detIds
0265         int detid = itFEDerror25->first;
0266         std::vector<PixelFEDChannel> chs = itFEDerror25->second;
0267         for (unsigned int ich = 0; ich < chs.size(); ich++) {
0268           fDet_.fillFEDerror25(detid, chs[ich]);
0269         }
0270       }  // loop over detIds
0271     }  // if non-empty FEDerror25
0272 
0273   }  // only for non-zero events
0274 }
0275 
0276 void SiPixelStatusProducer::endLuminosityBlockSummary(
0277     edm::LuminosityBlock const& iLumi,
0278     edm::EventSetup const&,
0279     std::vector<SiPixelDetectorStatus>* siPixelDetectorStatusVtr) const {
0280   /*add the Stream's partial information to the full information*/
0281 
0282   /* only save for the lumi sections with NON-ZERO events */
0283   if (ftotalevents_ > 0)
0284     siPixelDetectorStatusVtr->push_back(fDet_);
0285 }
0286 
0287 /* helper function */
0288 int SiPixelStatusProducer::indexROC(int irow, int icol, int nROCcolumns) {
0289   return int(icol + irow * nROCcolumns);
0290 
0291   /* generate the folling roc index that is going to map with ROC id as
0292      8  9  10 11 12 13 14 15
0293      0  1  2  3  4  5  6  7 */
0294 }
0295 
0296 DEFINE_FWK_MODULE(SiPixelStatusProducer);