Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:08:23

0001 // -*- C++ -*-
0002 //
0003 // Package:    SiPixelMonitorRawData
0004 // Class:      SiPixelRawDataErrorSource
0005 //
0006 /**\class
0007 
0008  Description:
0009  Produces histograms for error information generated at the raw2digi stage for
0010  the pixel tracker.
0011 
0012  Implementation:
0013  Takes a DetSetVector<SiPixelRawDataError> as input, and uses it to populate  a
0014  folder hierarchy (organized by detId) with histograms containing information
0015  about the errors.  Uses SiPixelRawDataErrorModule class to book and fill
0016  individual folders. Note that this source is different than other DQM sources
0017  in the creation of an unphysical detId folder (detId=0xffffffff) to hold
0018  information about errors for which there is no detId available (except the
0019  dummy detId given to it at raw2digi).
0020 
0021 */
0022 //
0023 // Original Author:  Andrew York
0024 //
0025 #include "DQM/SiPixelMonitorRawData/interface/SiPixelRawDataErrorSource.h"
0026 // Framework
0027 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0028 #include "FWCore/ServiceRegistry/interface/Service.h"
0029 #include "FWCore/Framework/interface/MakerMacros.h"
0030 // DQM Framework
0031 #include "DQM/SiPixelCommon/interface/SiPixelFolderOrganizer.h"
0032 #include "DQM/SiPixelCommon/interface/SiPixelHistogramId.h"
0033 #include "DQMServices/Core/interface/DQMStore.h"
0034 
0035 // Geometry
0036 #include "Geometry/CommonTopologies/interface/PixelTopology.h"
0037 #include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
0038 // DataFormats
0039 #include "DataFormats/DetId/interface/DetId.h"
0040 #include "DataFormats/FEDRawData/interface/FEDNumbering.h"
0041 #include "DataFormats/TrackerCommon/interface/PixelBarrelName.h"
0042 #include "DataFormats/SiPixelDetId/interface/PixelBarrelNameUpgrade.h"
0043 #include "DataFormats/TrackerCommon/interface/PixelEndcapName.h"
0044 #include "DataFormats/SiPixelDetId/interface/PixelEndcapNameUpgrade.h"
0045 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
0046 //
0047 #include <cstdlib>
0048 #include <string>
0049 
0050 using namespace std;
0051 using namespace edm;
0052 
0053 SiPixelRawDataErrorSource::SiPixelRawDataErrorSource(const edm::ParameterSet &iConfig)
0054     : conf_(iConfig),
0055       src_(consumes<DetSetVector<SiPixelRawDataError>>(conf_.getParameter<edm::InputTag>("src"))),
0056       saveFile(conf_.getUntrackedParameter<bool>("saveFile", false)),
0057       isPIB(conf_.getUntrackedParameter<bool>("isPIB", false)),
0058       slowDown(conf_.getUntrackedParameter<bool>("slowDown", false)),
0059       reducedSet(conf_.getUntrackedParameter<bool>("reducedSet", false)),
0060       modOn(conf_.getUntrackedParameter<bool>("modOn", true)),
0061       ladOn(conf_.getUntrackedParameter<bool>("ladOn", false)),
0062       bladeOn(conf_.getUntrackedParameter<bool>("bladeOn", false)),
0063       isUpgrade(conf_.getUntrackedParameter<bool>("isUpgrade", false)) {
0064   firstRun = true;
0065   LogInfo("PixelDQM") << "SiPixelRawDataErrorSource::SiPixelRawDataErrorSource:"
0066                          " Got DQM BackEnd interface"
0067                       << endl;
0068   topFolderName_ = conf_.getParameter<std::string>("TopFolderName");
0069   inputSourceToken_ = consumes<FEDRawDataCollection>(conf_.getUntrackedParameter<string>("inputSource", "source"));
0070   trackerTopoTokenBeginRun_ = esConsumes<TrackerTopology, TrackerTopologyRcd, edm::Transition::BeginRun>();
0071   trackerGeomTokenBeginRun_ = esConsumes<TrackerGeometry, TrackerDigiGeometryRecord, edm::Transition::BeginRun>();
0072 }
0073 
0074 SiPixelRawDataErrorSource::~SiPixelRawDataErrorSource() {
0075   // do anything here that needs to be done at desctruction time
0076   // (e.g. close files, deallocate resources etc.)
0077   LogInfo("PixelDQM") << "SiPixelRawDataErrorSource::~SiPixelRawDataErrorSource: Destructor" << endl;
0078 }
0079 
0080 void SiPixelRawDataErrorSource::dqmBeginRun(const edm::Run &r, const edm::EventSetup &iSetup) {
0081   LogInfo("PixelDQM") << " SiPixelRawDataErrorSource::beginRun - Initialisation ... " << std::endl;
0082   LogInfo("PixelDQM") << "Mod/Lad/Blade " << modOn << "/" << ladOn << "/" << bladeOn << std::endl;
0083 
0084   if (firstRun) {
0085     eventNo = 0;
0086 
0087     firstRun = false;
0088   }
0089 
0090   // Build map
0091   buildStructure(iSetup);
0092 }
0093 
0094 void SiPixelRawDataErrorSource::bookHistograms(DQMStore::IBooker &iBooker,
0095                                                edm::Run const &,
0096                                                edm::EventSetup const &iSetup) {
0097   // Book Monitoring Elements
0098   bookMEs(iBooker);
0099 }
0100 
0101 //------------------------------------------------------------------
0102 // Method called for every event
0103 //------------------------------------------------------------------
0104 void SiPixelRawDataErrorSource::analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) {
0105   eventNo++;
0106   // check feds in readout
0107   if (eventNo == 1) {
0108     // check if any Pixel FED is in readout:
0109     edm::Handle<FEDRawDataCollection> rawDataHandle;
0110     iEvent.getByToken(inputSourceToken_, rawDataHandle);
0111     if (!rawDataHandle.isValid()) {
0112       edm::LogInfo("SiPixelRawDataErrorSource") << "inputsource is empty";
0113     } else {
0114       const FEDRawDataCollection &rawDataCollection = *rawDataHandle;
0115       for (int i = 0; i != 40; i++) {
0116         if (rawDataCollection.FEDData(i).size() && rawDataCollection.FEDData(i).data())
0117           fedcounter->setBinContent(i + 1, 1);
0118       }
0119     }
0120   }
0121   // get input data
0122   edm::Handle<DetSetVector<SiPixelRawDataError>> input;
0123   iEvent.getByToken(src_, input);
0124   if (!input.isValid())
0125     return;
0126 
0127   int lumiSection = (int)iEvent.luminosityBlock();
0128 
0129   int nEventBPIXModuleErrors = 0;
0130   int nEventFPIXModuleErrors = 0;
0131   int nEventBPIXFEDErrors = 0;
0132   int nEventFPIXFEDErrors = 0;
0133   int nErrors = 0;
0134 
0135   std::map<uint32_t, SiPixelRawDataErrorModule *>::iterator struct_iter;
0136   std::map<uint32_t, SiPixelRawDataErrorModule *>::iterator struct_iter2;
0137   for (struct_iter = thePixelStructure.begin(); struct_iter != thePixelStructure.end(); struct_iter++) {
0138     int numberOfModuleErrors = (*struct_iter).second->fill(*input, &meMapFEDs_, modOn, ladOn, bladeOn);
0139     if (DetId((*struct_iter).first).subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel))
0140       nEventBPIXModuleErrors = nEventBPIXModuleErrors + numberOfModuleErrors;
0141     if (DetId((*struct_iter).first).subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap))
0142       nEventFPIXModuleErrors = nEventFPIXModuleErrors + numberOfModuleErrors;
0143     // cout<<"NErrors: "<<nEventBPIXModuleErrors<<"
0144     // "<<nEventFPIXModuleErrors<<endl;
0145     nErrors = nErrors + numberOfModuleErrors;
0146     // if(nErrors>0) cout<<"MODULES: nErrors: "<<nErrors<<endl;
0147   }
0148   for (struct_iter2 = theFEDStructure.begin(); struct_iter2 != theFEDStructure.end(); struct_iter2++) {
0149     int numberOfFEDErrors = (*struct_iter2).second->fillFED(*input, &meMapFEDs_);
0150     if ((*struct_iter2).first <= 31)
0151       nEventBPIXFEDErrors = nEventBPIXFEDErrors + numberOfFEDErrors;  // (*struct_iter2).first >= 0, since
0152                                                                       // (*struct_iter2).first is unsigned
0153     if ((*struct_iter2).first >= 32 && (*struct_iter2).first <= 39)
0154       nEventFPIXFEDErrors = nEventFPIXFEDErrors + numberOfFEDErrors;
0155     // cout<<"NFEDErrors: "<<nEventBPIXFEDErrors<<"
0156     // "<<nEventFPIXFEDErrors<<endl;
0157     nErrors = nErrors + numberOfFEDErrors;
0158     // if(nErrors>0) cout<<"FEDS: nErrors: "<<nErrors<<endl;
0159   }
0160   if (byLumiErrors) {
0161     byLumiErrors->setBinContent(0, eventNo);
0162     // cout<<"NErrors: "<<nEventBPIXModuleErrors<<" "<<nEventFPIXModuleErrors<<"
0163     // "<<nEventBPIXFEDErrors<<" "<<nEventFPIXFEDErrors<<endl;
0164     if (nEventBPIXModuleErrors + nEventBPIXFEDErrors > 0)
0165       byLumiErrors->Fill(0, 1.);
0166     if (nEventFPIXModuleErrors + nEventFPIXFEDErrors > 0)
0167       byLumiErrors->Fill(1, 1.);
0168     // cout<<"histo: "<<byLumiErrors->getBinContent(0)<<"
0169     // "<<byLumiErrors->getBinContent(1)<<"
0170     // "<<byLumiErrors->getBinContent(2)<<endl;
0171   }
0172 
0173   // Rate of errors per lumi section:
0174   if (errorRate)
0175     errorRate->Fill(lumiSection, nErrors);
0176 
0177   // slow down...
0178   if (slowDown)
0179     usleep(100000);
0180 }
0181 
0182 //------------------------------------------------------------------
0183 // Build data structure
0184 //------------------------------------------------------------------
0185 void SiPixelRawDataErrorSource::buildStructure(const edm::EventSetup &iSetup) {
0186   LogInfo("PixelDQM") << " SiPixelRawDataErrorSource::buildStructure";
0187 
0188   edm::ESHandle<TrackerGeometry> pDD = iSetup.getHandle(trackerGeomTokenBeginRun_);
0189 
0190   edm::ESHandle<TrackerTopology> tTopoHandle = iSetup.getHandle(trackerTopoTokenBeginRun_);
0191   const TrackerTopology *pTT = tTopoHandle.product();
0192 
0193   LogVerbatim("PixelDQM") << " *** Geometry node for TrackerGeom is  " << &(*pDD) << std::endl;
0194   LogVerbatim("PixelDQM") << " *** I have " << pDD->detsPXB().size() << " barrel pixel detectors" << std::endl;
0195   LogVerbatim("PixelDQM") << " *** I have " << pDD->detsPXF().size() << " endcap pixel detectors" << std::endl;
0196   // LogVerbatim ("PixelDQM") << " *** I have " << pDD->detTypes().size() <<"
0197   // types"<<std::endl;
0198 
0199   for (TrackerGeometry::DetContainer::const_iterator it = pDD->detsPXB().begin(); it != pDD->detsPXB().end(); it++) {
0200     const GeomDetUnit *geoUnit = dynamic_cast<const GeomDetUnit *>(*it);
0201     // check if it is a detUnit
0202     if (geoUnit == nullptr)
0203       LogError("PixelDQM") << "Pixel GeomDet is not a GeomDetUnit!" << std::endl;
0204     const PixelGeomDetUnit *pixDet = dynamic_cast<const PixelGeomDetUnit *>(geoUnit);
0205     int nrows = (pixDet->specificTopology()).nrows();
0206     int ncols = (pixDet->specificTopology()).ncolumns();
0207 
0208     if (isPIB)
0209       continue;
0210     DetId detId = (*it)->geographicalId();
0211     LogDebug("PixelDQM") << " ---> Adding Barrel Module " << detId.rawId() << endl;
0212     uint32_t id = detId();
0213     SiPixelRawDataErrorModule *theModule = new SiPixelRawDataErrorModule(id, ncols, nrows);
0214     thePixelStructure.insert(pair<uint32_t, SiPixelRawDataErrorModule *>(id, theModule));
0215   }
0216 
0217   for (TrackerGeometry::DetContainer::const_iterator it = pDD->detsPXF().begin(); it != pDD->detsPXF().end(); it++) {
0218     const GeomDetUnit *geoUnit = dynamic_cast<const GeomDetUnit *>(*it);
0219     // check if it is a detUnit
0220     if (geoUnit == nullptr)
0221       LogError("PixelDQM") << "Pixel GeomDet is not a GeomDetUnit!" << std::endl;
0222     const PixelGeomDetUnit *pixDet = dynamic_cast<const PixelGeomDetUnit *>(geoUnit);
0223     int nrows = (pixDet->specificTopology()).nrows();
0224     int ncols = (pixDet->specificTopology()).ncolumns();
0225 
0226     DetId detId = (*it)->geographicalId();
0227     LogDebug("PixelDQM") << " ---> Adding Endcap Module " << detId.rawId() << endl;
0228     uint32_t id = detId();
0229     SiPixelRawDataErrorModule *theModule = new SiPixelRawDataErrorModule(id, ncols, nrows);
0230 
0231     PixelEndcapName::HalfCylinder side = PixelEndcapName(DetId(id), pTT, isUpgrade).halfCylinder();
0232     int disk = PixelEndcapName(DetId(id), pTT, isUpgrade).diskName();
0233     int blade = PixelEndcapName(DetId(id), pTT, isUpgrade).bladeName();
0234     int panel = PixelEndcapName(DetId(id), pTT, isUpgrade).pannelName();
0235     int module = PixelEndcapName(DetId(id), pTT, isUpgrade).plaquetteName();
0236 
0237     char sside[80];
0238     sprintf(sside, "HalfCylinder_%i", side);
0239     char sdisk[80];
0240     sprintf(sdisk, "Disk_%i", disk);
0241     char sblade[80];
0242     sprintf(sblade, "Blade_%02i", blade);
0243     char spanel[80];
0244     sprintf(spanel, "Panel_%i", panel);
0245     char smodule[80];
0246     sprintf(smodule, "Module_%i", module);
0247     std::string side_str = sside;
0248     std::string disk_str = sdisk;
0249     bool mask = side_str.find("HalfCylinder_1") != string::npos || side_str.find("HalfCylinder_2") != string::npos ||
0250                 side_str.find("HalfCylinder_4") != string::npos || disk_str.find("Disk_2") != string::npos;
0251     // clutch to take all of FPIX, but no BPIX:
0252     mask = false;
0253     if (isPIB && mask)
0254       continue;
0255 
0256     thePixelStructure.insert(pair<uint32_t, SiPixelRawDataErrorModule *>(id, theModule));
0257   }
0258 
0259   LogDebug("PixelDQM") << " ---> Adding Module for Additional Errors " << endl;
0260   pair<int, int> fedIds(FEDNumbering::MINSiPixelFEDID, FEDNumbering::MAXSiPixelFEDID);
0261 
0262   fedIds.first = 0;
0263   fedIds.second = 39;
0264 
0265   for (int fedId = fedIds.first; fedId <= fedIds.second; fedId++) {
0266     // std::cout<<"Adding FED module: "<<fedId<<std::endl;
0267     uint32_t id = static_cast<uint32_t>(fedId);
0268     SiPixelRawDataErrorModule *theModule = new SiPixelRawDataErrorModule(id);
0269     theFEDStructure.insert(pair<uint32_t, SiPixelRawDataErrorModule *>(id, theModule));
0270   }
0271 
0272   LogInfo("PixelDQM") << " *** Pixel Structure Size " << thePixelStructure.size() << endl;
0273 }
0274 //------------------------------------------------------------------
0275 // Book MEs
0276 //------------------------------------------------------------------
0277 void SiPixelRawDataErrorSource::bookMEs(DQMStore::IBooker &iBooker) {
0278   iBooker.setCurrentFolder(topFolderName_ + "/EventInfo/DAQContents");
0279   char title0[80];
0280   sprintf(title0, "FED isPresent;FED ID;isPresent");
0281   fedcounter = iBooker.book1D("fedcounter", title0, 40, -0.5, 39.5);
0282   iBooker.setCurrentFolder(topFolderName_ + "/AdditionalPixelErrors");
0283   char title[80];
0284   sprintf(title, "By-LumiSection Error counters");
0285   {
0286     auto scope = DQMStore::IBooker::UseLumiScope(iBooker);
0287     byLumiErrors = iBooker.book1D("byLumiErrors", title, 2, 0., 2.);
0288   }
0289   char title1[80];
0290   sprintf(title1, "Errors per LumiSection;LumiSection;NErrors");
0291   errorRate = iBooker.book1D("errorRate", title1, 5000, 0., 5000.);
0292 
0293   std::map<uint32_t, SiPixelRawDataErrorModule *>::iterator struct_iter;
0294   std::map<uint32_t, SiPixelRawDataErrorModule *>::iterator struct_iter2;
0295 
0296   SiPixelFolderOrganizer theSiPixelFolder(false);
0297 
0298   for (struct_iter = thePixelStructure.begin(); struct_iter != thePixelStructure.end(); struct_iter++) {
0299     /// Create folder tree and book histograms
0300 
0301     if (modOn) {
0302       if (!theSiPixelFolder.setModuleFolder(iBooker, (*struct_iter).first, 0, isUpgrade)) {
0303         // std::cout<<"PIB! not booking histograms for non-PIB
0304         // modules!"<<std::endl;
0305         if (!isPIB)
0306           throw cms::Exception("LogicError") << "[SiPixelRawDataErrorSource::bookMEs] Creation of DQM folder "
0307                                                 "failed";
0308       }
0309     }
0310 
0311     if (ladOn) {
0312       if (!theSiPixelFolder.setModuleFolder(iBooker, (*struct_iter).first, 1, isUpgrade)) {
0313         LogDebug("PixelDQM") << "PROBLEM WITH LADDER-FOLDER\n";
0314       }
0315     }
0316 
0317     if (bladeOn) {
0318       if (!theSiPixelFolder.setModuleFolder(iBooker, (*struct_iter).first, 4, isUpgrade)) {
0319         LogDebug("PixelDQM") << "PROBLEM WITH BLADE-FOLDER\n";
0320       }
0321     }
0322 
0323   }  // for loop
0324 
0325   for (struct_iter2 = theFEDStructure.begin(); struct_iter2 != theFEDStructure.end(); struct_iter2++) {
0326     /// Create folder tree for errors without detId and book histograms
0327     if (!theSiPixelFolder.setFedFolder(iBooker, (*struct_iter2).first)) {
0328       throw cms::Exception("LogicError") << "[SiPixelRawDataErrorSource::bookMEs] Creation of DQM folder "
0329                                             "failed";
0330     }
0331   }
0332 
0333   // Booking FED histograms
0334   std::string hid;
0335   // Get collection name and instantiate Histo Id builder
0336   edm::InputTag src = conf_.getParameter<edm::InputTag>("src");
0337   SiPixelHistogramId *theHistogramId = new SiPixelHistogramId(src.label());
0338 
0339   for (uint32_t id = 0; id < 40; id++) {
0340     char temp[50];
0341     sprintf(temp, (topFolderName_ + "/AdditionalPixelErrors/FED_%d").c_str(), id);
0342     iBooker.cd(temp);
0343     // Types of errors
0344     hid = theHistogramId->setHistoId("errorType", id);
0345     meErrorType_[id] = iBooker.book1D(hid, "Type of errors", 15, 24.5, 39.5);
0346     meErrorType_[id]->setAxisTitle("Type of errors", 1);
0347     // Number of errors
0348     hid = theHistogramId->setHistoId("NErrors", id);
0349     meNErrors_[id] = iBooker.book1D(hid, "Number of errors", 36, 0., 36.);
0350     meNErrors_[id]->setAxisTitle("Number of errors", 1);
0351     // Type of FIFO full (errorType = 28).  FIFO 1 is 1-5 (where fullType =
0352     // channel of FIFO 1), fullType = 6 signifies FIFO 2 nearly full, 7
0353     // signifies trigger FIFO nearly full, 8 indicates an unexpected result
0354     hid = theHistogramId->setHistoId("fullType", id);
0355     meFullType_[id] = iBooker.book1D(hid, "Type of FIFO full", 7, 0.5, 7.5);
0356     meFullType_[id]->setAxisTitle("FIFO type", 1);
0357     // For error type 30, the type of problem encoded in the TBM trailer
0358     // 0 = stack full, 1 = Pre-cal issued, 2 = clear trigger counter, 3 = sync
0359     // trigger, 4 = sync trigger error, 5 = reset ROC, 6 = reset TBM, 7 = no
0360     // token bit pass
0361     hid = theHistogramId->setHistoId("TBMMessage", id);
0362     meTBMMessage_[id] = iBooker.book1D(hid, "TBM trailer message", 8, -0.5, 7.5);
0363     meTBMMessage_[id]->setAxisTitle("TBM message", 1);
0364     // For error type 30, the type of problem encoded in the TBM error trailer 0
0365     // = none 1 = data stream too long, 2 = FSM errors, 3 = invalid # of ROCs, 4
0366     // = multiple
0367     hid = theHistogramId->setHistoId("TBMType", id);
0368     meTBMType_[id] = iBooker.book1D(hid, "Type of TBM trailer", 5, -0.5, 4.5);
0369     meTBMType_[id]->setAxisTitle("TBM Type", 1);
0370     // For error type 31, the event number of the TBM header with the error
0371     hid = theHistogramId->setHistoId("EvtNbr", id);
0372     meEvtNbr_[id] = iBooker.book1D(hid, "Event number", 1, 0, 1);
0373     // For errorType = 34, datastream size according to error word
0374     hid = theHistogramId->setHistoId("evtSize", id);
0375     meEvtSize_[id] = iBooker.book1D(hid, "Event size", 1, 0, 1);
0376     //
0377     hid = theHistogramId->setHistoId("FedChNErr", id);
0378     meFedChNErr_[id] = iBooker.book1D(hid, "Number of errors per FED channel", 37, 0, 37);
0379     meFedChNErr_[id]->setAxisTitle("FED channel", 1);
0380     //
0381     hid = theHistogramId->setHistoId("FedChLErr", id);
0382     meFedChLErr_[id] = iBooker.book1D(hid, "Last error per FED channel", 37, 0, 37);
0383     meFedChLErr_[id]->setAxisTitle("FED channel", 1);
0384     //
0385     hid = theHistogramId->setHistoId("FedETypeNErr", id);
0386     meFedETypeNErr_[id] = iBooker.book1D(hid, "Number of errors per type", 21, 0, 21);
0387     meFedETypeNErr_[id]->setAxisTitle("Error type", 1);
0388   }
0389   // Add the booked histograms to the histogram map for booking
0390   meMapFEDs_["meErrorType_"] = meErrorType_;
0391   meMapFEDs_["meNErrors_"] = meNErrors_;
0392   meMapFEDs_["meFullType_"] = meFullType_;
0393   meMapFEDs_["meTBMMessage_"] = meTBMMessage_;
0394   meMapFEDs_["meTBMType_"] = meTBMType_;
0395   meMapFEDs_["meEvtNbr_"] = meEvtNbr_;
0396   meMapFEDs_["meEvtSize_"] = meEvtSize_;
0397   meMapFEDs_["meFedChNErr_"] = meFedChNErr_;
0398   meMapFEDs_["meFedChLErr_"] = meFedChLErr_;
0399   meMapFEDs_["meFedETypeNErr_"] = meFedETypeNErr_;
0400 
0401   // cout<<"...leaving SiPixelRawDataErrorSource::bookMEs now! "<<endl;
0402 }
0403 
0404 DEFINE_FWK_MODULE(SiPixelRawDataErrorSource);