Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 // -*- C++ -*-
0002 //
0003 // Package:     SiPixelPhase1RawData
0004 // Class:       SiPixelPhase1RawData
0005 //
0006 
0007 // Original Author: Marcel Schneider
0008 
0009 #include "DQM/SiPixelPhase1Common/interface/SiPixelPhase1Base.h"
0010 #include "DataFormats/Common/interface/DetSetVector.h"
0011 #include "DataFormats/SiPixelRawData/interface/SiPixelRawDataError.h"
0012 #include "FWCore/Framework/interface/MakerMacros.h"
0013 
0014 namespace {
0015 
0016   class SiPixelPhase1RawData final : public SiPixelPhase1Base {
0017     enum { NERRORS, FIFOFULL, TBMMESSAGE, TBMTYPE, TYPE_NERRORS };
0018 
0019   public:
0020     explicit SiPixelPhase1RawData(const edm::ParameterSet& conf);
0021     void analyze(const edm::Event&, const edm::EventSetup&) override;
0022 
0023   private:
0024     edm::EDGetTokenT<edm::DetSetVector<SiPixelRawDataError>> srcToken_;
0025   };
0026 
0027   SiPixelPhase1RawData::SiPixelPhase1RawData(const edm::ParameterSet& iConfig) : SiPixelPhase1Base(iConfig) {
0028     srcToken_ = consumes<edm::DetSetVector<SiPixelRawDataError>>(iConfig.getParameter<edm::InputTag>("src"));
0029   }
0030 
0031   void SiPixelPhase1RawData::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0032     if (!checktrigger(iEvent, iSetup, DCS))
0033       return;
0034 
0035     edm::Handle<edm::DetSetVector<SiPixelRawDataError>> input;
0036     iEvent.getByToken(srcToken_, input);
0037     if (!input.isValid())
0038       return;
0039 
0040     for (auto it = input->begin(); it != input->end(); ++it) {
0041       for (auto& siPixelRawDataError : *it) {
0042         int fed = siPixelRawDataError.getFedId();
0043         int type = siPixelRawDataError.getType();
0044         DetId id = it->detId();
0045 
0046         // encoding of the channel number within the FED error word
0047         const uint32_t LINK_bits = 6;
0048         const uint32_t LINK_shift = 26;
0049         const uint64_t LINK_mask = (1 << LINK_bits) - 1;
0050 
0051         uint64_t errorWord = 0;
0052         // use 64bit word for some error types
0053         // invalid header, invalid trailer, size mismatch
0054         if (type == 32 || type == 33 || type == 34) {
0055           errorWord = siPixelRawDataError.getWord64();
0056         } else {
0057           errorWord = siPixelRawDataError.getWord32();
0058         }
0059 
0060         int32_t chanNmbr = (errorWord >> LINK_shift) & LINK_mask;
0061         // timeout
0062         if (type == 29)
0063           chanNmbr = -1;  // TODO: different formula needed.
0064 
0065         uint32_t error_data = errorWord & 0xFF;
0066 
0067         if (type == 28) {  // overflow.
0068           for (uint32_t i = 0; i < 8; i++) {
0069             if (error_data & (1 << i))
0070               histo[FIFOFULL].fill(i, id, &iEvent, fed, chanNmbr);
0071           }
0072         }
0073 
0074         if (type == 30) {                                      // TBM stuff.
0075           uint32_t statemachine_state = errorWord >> 8 & 0xF;  // next 4 bits after data
0076           const uint32_t tbm_types[16] = {0, 1, 2, 4, 2, 4, 2, 4, 3, 1, 4, 4, 4, 4, 4, 4};
0077 
0078           histo[TBMTYPE].fill(tbm_types[statemachine_state], id, &iEvent, fed, chanNmbr);
0079 
0080           for (uint32_t i = 0; i < 8; i++) {
0081             if (error_data & (1 << i))
0082               histo[TBMMESSAGE].fill(i, id, &iEvent, fed, chanNmbr);
0083           }
0084           continue;  // we don't really consider these as errors.
0085         }
0086 
0087         // note that a DetId of 0xFFFFFFFF can mean 'no DetId'.
0088         // We hijack column and row for FED and chan in this case,
0089         // the GeometryInterface does understand that.
0090 
0091         histo[NERRORS].fill(id, &iEvent, fed, chanNmbr);
0092         histo[TYPE_NERRORS].fill(type, id, &iEvent, fed, chanNmbr);
0093       }
0094     }
0095 
0096     histo[NERRORS].executePerEventHarvesting(&iEvent);
0097   }
0098 
0099 }  //namespace
0100 
0101 DEFINE_FWK_MODULE(SiPixelPhase1RawData);