CMS 3D CMS Logo

SiPixelPhase1RawData.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: SiPixelPhase1RawData
4 // Class: SiPixelPhase1RawData
5 //
6 
7 // Original Author: Marcel Schneider
8 
13 
16 
17 namespace {
18 
19 class SiPixelPhase1RawData final : public SiPixelPhase1Base {
20  enum {
21  NERRORS,
22  FIFOFULL,
23  TBMMESSAGE,
24  TBMTYPE,
25  TYPE_NERRORS
26  };
27 
28  public:
29  explicit SiPixelPhase1RawData(const edm::ParameterSet& conf);
30  void analyze(const edm::Event&, const edm::EventSetup&) override;
31 
32  private:
34 };
35 
36 SiPixelPhase1RawData::SiPixelPhase1RawData(const edm::ParameterSet& iConfig) :
37  SiPixelPhase1Base(iConfig)
38 {
39  srcToken_ = consumes<edm::DetSetVector<SiPixelRawDataError>>(iConfig.getParameter<edm::InputTag>("src"));
40 }
41 
42 
44  if( !checktrigger(iEvent,iSetup,DCS) ) return;
45 
47  iEvent.getByToken(srcToken_, input);
48  if (!input.isValid()) return;
49 
50  for (auto it = input->begin(); it != input->end(); ++it) {
51  for (auto& siPixelRawDataError : *it) {
52  int fed = siPixelRawDataError.getFedId();
53  int type = siPixelRawDataError.getType();
54  DetId id = it->detId();
55 
56  // encoding of the channel number within the FED error word
57  const uint32_t LINK_bits = 6;
58  const uint32_t LINK_shift = 26;
59  const uint64_t LINK_mask = (1 << LINK_bits) - 1;
60 
61  uint64_t errorWord = 0;
62  // use 64bit word for some error types
63  // invalid header, invalid trailer, size mismatch
64  if (type == 32 || type == 33 || type == 34) {
65  errorWord = siPixelRawDataError.getWord64();
66  } else {
67  errorWord = siPixelRawDataError.getWord32();
68  }
69 
70  int32_t chanNmbr = (errorWord >> LINK_shift) & LINK_mask;
71  // timeout
72  if (type == 29) chanNmbr = -1; // TODO: different formula needed.
73 
74  uint32_t error_data = errorWord & 0xFF;
75 
76  if (type == 28) { // overflow.
77  for (uint32_t i = 0; i < 8; i++) {
78  if (error_data & (1 << i)) histo[FIFOFULL].fill(i, id, &iEvent, fed, chanNmbr);
79  }
80  }
81 
82  if (type == 30) { // TBM stuff.
83  uint32_t statemachine_state = errorWord >> 8 & 0xF; // next 4 bits after data
84  const uint32_t tbm_types[16] = {
85  0, 1, 2, 4,
86  2, 4, 2, 4,
87  3, 1, 4, 4,
88  4, 4, 4, 4
89  };
90 
91  histo[TBMTYPE].fill(tbm_types[statemachine_state], id, &iEvent, fed, chanNmbr);
92 
93  for (uint32_t i = 0; i < 8; i++) {
94  if (error_data & (1 << i)) histo[TBMMESSAGE].fill(i, id, &iEvent, fed, chanNmbr);
95  }
96  continue; // we don't really consider these as errors.
97  }
98 
99  // note that a DetId of 0xFFFFFFFF can mean 'no DetId'.
100  // We hijack column and row for FED and chan in this case,
101  // the GeometryInterface does understand that.
102 
103  histo[NERRORS].fill(id, &iEvent, fed, chanNmbr);
104  histo[TYPE_NERRORS].fill(type, id, &iEvent, fed, chanNmbr);
105  }
106  }
107 
108  histo[NERRORS].executePerEventHarvesting(&iEvent);
109 }
110 
111 } //namespace
112 
113 DEFINE_FWK_MODULE(SiPixelPhase1RawData);
114 
type
Definition: HCALResponse.h:21
T getParameter(std::string const &) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
example_stream void analyze(const edm::Event &, const edm::EventSetup &) override
static std::string const input
Definition: EdmProvDump.cc:48
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
bool isValid() const
Definition: HandleBase.h:74
void analyze(edm::Event const &e, edm::EventSetup const &eSetup) override=0
Definition: DetId.h:18
unsigned long long uint64_t
Definition: Time.h:15
std::vector< HistogramManager > histo