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:      SiPixelHLTSource
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 raw data and error data as input, and uses it to populate three
0014  histograms indexed by FED id.
0015 
0016 */
0017 //
0018 // Original Author:  Andrew York
0019 //
0020 // Framework
0021 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0022 #include "FWCore/ServiceRegistry/interface/Service.h"
0023 #include "FWCore/Framework/interface/MakerMacros.h"
0024 // DQM Framework
0025 #include "DQM/SiPixelCommon/interface/SiPixelFolderOrganizer.h"
0026 #include "DQM/SiPixelCommon/interface/SiPixelHistogramId.h"
0027 #include "DQM/SiPixelMonitorRawData/interface/SiPixelHLTSource.h"
0028 #include "DQMServices/Core/interface/DQMStore.h"
0029 // Geometry
0030 #include "Geometry/CommonDetUnit/interface/GeomDetEnumerators.h"
0031 #include "Geometry/CommonTopologies/interface/PixelTopology.h"
0032 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0033 #include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
0034 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0035 // DataFormats
0036 #include "DataFormats/DetId/interface/DetId.h"
0037 #include "DataFormats/FEDRawData/interface/FEDNumbering.h"
0038 #include "DataFormats/TrackerCommon/interface/PixelBarrelName.h"
0039 #include "DataFormats/SiPixelDetId/interface/PixelBarrelNameUpgrade.h"
0040 #include "DataFormats/TrackerCommon/interface/PixelEndcapName.h"
0041 #include "DataFormats/SiPixelDetId/interface/PixelEndcapNameUpgrade.h"
0042 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
0043 //
0044 #include <cstdlib>
0045 #include <string>
0046 
0047 using namespace std;
0048 using namespace edm;
0049 
0050 SiPixelHLTSource::SiPixelHLTSource(const edm::ParameterSet &iConfig)
0051     : conf_(iConfig),
0052       rawin_(consumes<FEDRawDataCollection>(conf_.getParameter<edm::InputTag>("RawInput"))),
0053       errin_(consumes<edm::DetSetVector<SiPixelRawDataError>>(conf_.getParameter<edm::InputTag>("ErrorInput"))),
0054       trackerGeomToken_(esConsumes<TrackerGeometry, TrackerDigiGeometryRecord>()),
0055       saveFile(conf_.getUntrackedParameter<bool>("saveFile", false)),
0056       slowDown(conf_.getUntrackedParameter<bool>("slowDown", false)),
0057       dirName_(conf_.getUntrackedParameter<string>("DirName", "Pixel/FEDIntegrity/")) {
0058   firstRun = true;
0059   LogInfo("PixelDQM") << "SiPixelHLTSource::SiPixelHLTSource: Got DQM BackEnd interface" << endl;
0060 }
0061 
0062 SiPixelHLTSource::~SiPixelHLTSource() {
0063   // do anything here that needs to be done at desctruction time
0064   // (e.g. close files, deallocate resources etc.)
0065   LogInfo("PixelDQM") << "SiPixelHLTSource::~SiPixelHLTSource: Destructor" << endl;
0066 }
0067 
0068 void SiPixelHLTSource::dqmBeginRun(const edm::Run &r, const edm::EventSetup &iSetup) {
0069   LogInfo("PixelDQM") << " SiPixelHLTSource::beginJob - Initialisation ... " << std::endl;
0070 
0071   if (firstRun) {
0072     eventNo = 0;
0073 
0074     firstRun = false;
0075   }
0076 }
0077 
0078 void SiPixelHLTSource::bookHistograms(DQMStore::IBooker &iBooker, edm::Run const &, edm::EventSetup const &) {
0079   // Book Monitoring Elements
0080   bookMEs(iBooker);
0081 }
0082 
0083 //------------------------------------------------------------------
0084 // Method called for every event
0085 //------------------------------------------------------------------
0086 void SiPixelHLTSource::analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) {
0087   eventNo++;
0088   edm::ESHandle<TrackerGeometry> pDD = iSetup.getHandle(trackerGeomToken_);
0089   // get raw input data
0090   edm::Handle<FEDRawDataCollection> rawinput;
0091   iEvent.getByToken(rawin_, rawinput);
0092   // get error input data
0093   edm::Handle<edm::DetSetVector<SiPixelRawDataError>> errorinput;
0094   iEvent.getByToken(errin_, errorinput);
0095   if (!errorinput.isValid())
0096     return;
0097 
0098   int fedId;
0099 
0100   for (fedId = 0; fedId <= 39; fedId++) {
0101     // get event data for this fed
0102     const FEDRawData &fedRawData = rawinput->FEDData(fedId);
0103     if (fedRawData.size() != 0)
0104       (meRawWords_)->Fill(fedId);
0105   }  // end for
0106 
0107   edm::DetSet<SiPixelRawDataError>::const_iterator di;
0108 
0109   for (TrackerGeometry::DetContainer::const_iterator it = pDD->dets().begin(); it != pDD->dets().end(); it++) {
0110     if (GeomDetEnumerators::isTrackerPixel((*it)->subDetector())) {
0111       uint32_t detId = (*it)->geographicalId();
0112       edm::DetSetVector<SiPixelRawDataError>::const_iterator isearch = errorinput->find(detId);
0113       if (isearch != errorinput->end()) {
0114         for (di = isearch->data.begin(); di != isearch->data.end(); di++) {
0115           fedId = di->getFedId();         // FED the error came from
0116           int errorType = di->getType();  // type of error
0117           switch (errorType) {
0118             case (35):
0119               (meNErrors_)->Fill(fedId);
0120               break;
0121             case (36):
0122               (meNErrors_)->Fill(fedId);
0123               break;
0124             case (37):
0125               (meNErrors_)->Fill(fedId);
0126               break;
0127             case (38):
0128               (meNErrors_)->Fill(fedId);
0129               break;
0130             default:
0131               break;
0132           };  // end switch
0133         }     // end for(di
0134       }       // end if( isearch
0135     }         // end if( ((*it)->subDetector()
0136   }           // for(TrackerGeometry
0137 
0138   edm::DetSetVector<SiPixelRawDataError>::const_iterator isearch = errorinput->find(0xffffffff);
0139 
0140   if (isearch != errorinput->end()) {  // Not at empty iterator
0141     for (di = isearch->data.begin(); di != isearch->data.end(); di++) {
0142       fedId = di->getFedId();         // FED the error came from
0143       int errorType = di->getType();  // type of error
0144       switch (errorType) {
0145         case (35):
0146           (meNErrors_)->Fill(fedId);
0147           break;
0148         case (36):
0149           (meNErrors_)->Fill(fedId);
0150           break;
0151         case (37):
0152           (meNErrors_)->Fill(fedId);
0153           break;
0154         case (38):
0155           (meNErrors_)->Fill(fedId);
0156           break;
0157         case (39):
0158           (meNCRCs_)->Fill(fedId);
0159           break;
0160         default:
0161           break;
0162       };  // end switch
0163     }     // end for(di
0164   }       // end if( isearch
0165   // slow down...
0166   if (slowDown)
0167     usleep(100000);
0168 }
0169 
0170 //------------------------------------------------------------------
0171 // Book MEs
0172 //------------------------------------------------------------------
0173 void SiPixelHLTSource::bookMEs(DQMStore::IBooker &iBooker) {
0174   iBooker.cd();
0175   iBooker.setCurrentFolder(dirName_);
0176 
0177   std::string rawhid;
0178   std::string errhid;
0179   // Get collection name and instantiate Histo Id builder
0180   edm::InputTag rawin = conf_.getParameter<edm::InputTag>("RawInput");
0181   SiPixelHistogramId *RawHistogramId = new SiPixelHistogramId(rawin.label());
0182   edm::InputTag errin = conf_.getParameter<edm::InputTag>("ErrorInput");
0183   SiPixelHistogramId *ErrorHistogramId = new SiPixelHistogramId(errin.label());
0184 
0185   // Is a FED sending raw data
0186   meRawWords_ = iBooker.book1D("FEDEntries", "Number of raw words", 40, -0.5, 39.5);
0187   meRawWords_->setAxisTitle("Number of raw words", 1);
0188 
0189   // Number of CRC errors
0190   meNCRCs_ = iBooker.book1D("FEDFatal", "Number of fatal errors", 40, -0.5, 39.5);
0191   meNCRCs_->setAxisTitle("Number of fatal errors", 1);
0192 
0193   // Number of translation error words
0194   meNErrors_ = iBooker.book1D("FEDNonFatal", "Number of non-fatal errors", 40, -0.5, 39.5);
0195   meNErrors_->setAxisTitle("Number of non-fatal errors", 1);
0196 
0197   delete RawHistogramId;
0198   delete ErrorHistogramId;
0199 }
0200 
0201 DEFINE_FWK_MODULE(SiPixelHLTSource);