CMS 3D CMS Logo

DTDDUUnpacker.cc
Go to the documentation of this file.
1 
13 
14 #include <iostream>
15 
16 using namespace std;
17 using namespace edm;
18 
20  // the ROS unpacker
22 
23  // parameters
24  localDAQ = dduPSet.getUntrackedParameter<bool>("localDAQ", false);
25  debug = dduPSet.getUntrackedParameter<bool>("debug", false);
26 }
27 
29 
30 void DTDDUUnpacker::interpretRawData(const unsigned int* index32,
31  int datasize,
32  int dduID,
34  std::unique_ptr<DTDigiCollection>& detectorProduct,
35  std::unique_ptr<DTLocalTriggerCollection>& triggerProduct,
36  uint16_t rosList) {
37  // Definitions
38  const int wordSize_32 = 4;
39  const int wordSize_64 = 8;
40 
41  int numberOf32Words = datasize / wordSize_32;
42 
43  const unsigned char* index8 = reinterpret_cast<const unsigned char*>(index32);
44 
46  /* D D U d a t a */
48 
49  // DDU header
50  FEDHeader dduHeader(index8);
51  if (dduHeader.check()) {
52  if (debug)
53  cout << "[DTDDUUnpacker] FED Header. BXID: " << dduHeader.bxID() << " L1ID: " << dduHeader.lvl1ID() << endl;
54  } else {
55  LogWarning("DTRawToDigi|DTDDUUnpacker")
56  << "[DTDDUUnpacker] WARNING!, this is not a DDU Header, FED ID: " << dduID << endl;
57  }
58 
59  // DDU trailer
60  // [BITS] stop before FED trailer := 8 bytes
61  FEDTrailer dduTrailer(index8 + datasize - 1 * wordSize_64);
62 
63  if (dduTrailer.check()) {
64  if (debug)
65  cout << "[DTDDUUnpacker] FED Trailer. Length of the DT event: " << dduTrailer.fragmentLength() << endl;
66  } else {
67  LogWarning("DTRawToDigi|DTDDUUnpacker")
68  << "[DTDDUUnpacker] WARNING!, this is not a DDU Trailer, FED ID: " << dduID << endl;
69  }
70 
71  // Initialize control data
73  controlData.addDDUHeader(dduHeader);
74  controlData.addDDUTrailer(dduTrailer);
75  // check the CRC set in the FED trailer (FCRC errors)
76  controlData.checkCRCBit(index8 + datasize - 1 * wordSize_64);
77 
78  // Check Status Words
79  vector<DTDDUFirstStatusWord> rosStatusWords;
80  // [BITS] 3 words of 8 bytes + "rosId" bytes
81  // In the case we are reading from DMA, the status word are swapped as the ROS data
82  if (localDAQ) {
83  // DDU channels from 1 to 4
84  for (int rosId = 0; rosId < 4; rosId++) {
85  int wordIndex8 = numberOf32Words * wordSize_32 - 3 * wordSize_64 + wordSize_32 + rosId;
87  }
88  // DDU channels from 5 to 8
89  for (int rosId = 0; rosId < 4; rosId++) {
90  int wordIndex8 = numberOf32Words * wordSize_32 - 3 * wordSize_64 + rosId;
92  }
93  // DDU channels from 9 to 12
94  for (int rosId = 0; rosId < 4; rosId++) {
95  int wordIndex8 = numberOf32Words * wordSize_32 - 2 * wordSize_64 + wordSize_32 + rosId;
97  }
98  } else {
99  for (int rosId = 0; rosId < 12; rosId++) {
100  int wordIndex8 = numberOf32Words * wordSize_32 - 3 * wordSize_64 + rosId;
102  }
103  }
104 
105  int theROSList;
106  // [BITS] 2 words of 8 bytes + 4 bytes (half 64 bit word)
107  // In the case we are reading from DMA, the status word are swapped as the ROS data
108  if (localDAQ) {
109  DTDDUSecondStatusWord dduStatusWord(index32[numberOf32Words - 2 * wordSize_64 / wordSize_32]);
110  controlData.addDDUStatusWord(dduStatusWord);
111  theROSList = dduStatusWord.rosList();
112  } else {
113  DTDDUSecondStatusWord dduStatusWord(index32[numberOf32Words - 2 * wordSize_64 / wordSize_32 + 1]);
114  controlData.addDDUStatusWord(dduStatusWord);
115  theROSList = dduStatusWord.rosList();
116  }
117 
119  /* R O S d a t a */
121 
122  // Set the index to start looping on ROS data
123  // [BITS] one 8 bytes word
124  index32 += (wordSize_64) / wordSize_32;
125 
126  // Set the datasize to look only at ROS data
127  // [BITS] header, trailer, 2 status words
128  datasize -= 4 * wordSize_64;
129 
130  // unpacking the ROS payload
131  ros25Unpacker->interpretRawData(index32, datasize, dduID, mapping, detectorProduct, triggerProduct, theROSList);
132 }
void interpretRawData(const unsigned int *index, int datasize, int dduID, edm::ESHandle< DTReadOutMapping > &mapping, std::unique_ptr< DTDigiCollection > &product, std::unique_ptr< DTLocalTriggerCollection > &product2, uint16_t rosList=0) override
void addDDUStatusWord(const DTDDUSecondStatusWord &word)
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
void addROSStatusWord(const DTDDUFirstStatusWord &word)
int rosList() const
Definition: DTDDUWords.h:690
bool check() const
Check that the header is OK.
Definition: FEDHeader.cc:44
uint32_t fragmentLength() const
The length of the event fragment counted in 64-bit words including header and trailer.
Definition: FEDTrailer.cc:13
void interpretRawData(const unsigned int *index, int datasize, int dduID, edm::ESHandle< DTReadOutMapping > &mapping, std::unique_ptr< DTDigiCollection > &product, std::unique_ptr< DTLocalTriggerCollection > &product2, uint16_t rosList=0) override
T getUntrackedParameter(std::string const &, T const &) const
DTROS25Unpacker * ros25Unpacker
Definition: DTDDUUnpacker.h:49
void addDDUHeader(const FEDHeader &word)
Setters.
void addDDUTrailer(const FEDTrailer &word)
bool localDAQ
if data are read locally, status words are swapped
Definition: DTDDUUnpacker.h:42
const edm::ParameterSet dduPSet
Definition: DTDDUUnpacker.h:39
void checkCRCBit(const unsigned char *trailer)
uint16_t bxID() const
The bunch crossing number.
Definition: FEDHeader.cc:17
void clean()
HLT enums.
~DTDDUUnpacker() override
Destructor.
uint32_t lvl1ID() const
Level-1 event number generated by the TTC system.
Definition: FEDHeader.cc:15
Log< level::Warning, false > LogWarning
DTDDUUnpacker(const edm::ParameterSet &ps)
Constructor.
DTDDUData controlData
Definition: DTDDUUnpacker.h:51
bool check() const
Check that the trailer is OK.
Definition: FEDTrailer.cc:45