CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DTDDUUnpacker.cc
Go to the documentation of this file.
1 
9 
12 
15 
18 
20 
21 #include <iostream>
22 
23 using namespace std;
24 using namespace edm;
25 
27 
28  // the ROS unpacker
30 
31  // parameters
32  localDAQ = dduPSet.getUntrackedParameter<bool>("localDAQ",false);
33  performDataIntegrityMonitor = dduPSet.getUntrackedParameter<bool>("performDataIntegrityMonitor",false);
34  debug = dduPSet.getUntrackedParameter<bool>("debug",false);
35 
36  // enable DQM if Service is available
37  if(performDataIntegrityMonitor) {
40  } else {
41  LogWarning("DTRawToDigi|DTDDUUnpacker") <<
42  "[DTDDUUnpacker] WARNING! Data Integrity Monitoring requested but no DTDataMonitorInterface Service available" << endl;
43  performDataIntegrityMonitor = false;
44  }
45  }
46 
47 }
48 
49 
51  delete ros25Unpacker;
52 }
53 
54 
55 void DTDDUUnpacker::interpretRawData(const unsigned int* index32, int datasize,
56  int dduID,
58  std::auto_ptr<DTDigiCollection>& detectorProduct,
59  std::auto_ptr<DTLocalTriggerCollection>& triggerProduct,
60  uint16_t rosList) {
61 
62  // Definitions
63  const int wordSize_32 = 4;
64  const int wordSize_64 = 8;
65 
66  int numberOf32Words = datasize/wordSize_32;
67 
68  const unsigned char* index8 = reinterpret_cast<const unsigned char*>(index32);
69 
70 
72  /* D D U d a t a */
74 
75  // DDU header
76  FEDHeader dduHeader(index8);
77  if (dduHeader.check()) {
78  if(debug) cout << "[DTDDUUnpacker] FED Header. BXID: "<<dduHeader.bxID()
79  << " L1ID: "<<dduHeader.lvl1ID() <<endl;
80  } else {
81  LogWarning("DTRawToDigi|DTDDUUnpacker") << "[DTDDUUnpacker] WARNING!, this is not a DDU Header, FED ID: "
82  << dduID << endl;
83  }
84 
85  // DDU trailer
86  // [BITS] stop before FED trailer := 8 bytes
87  FEDTrailer dduTrailer(index8 + datasize - 1*wordSize_64);
88 
89  if (dduTrailer.check()) {
90  if(debug) cout << "[DTDDUUnpacker] FED Trailer. Lenght of the DT event: "
91  << dduTrailer.lenght() << endl;
92  } else {
93  LogWarning("DTRawToDigi|DTDDUUnpacker") << "[DTDDUUnpacker] WARNING!, this is not a DDU Trailer, FED ID: "
94  << dduID << endl;
95  }
96 
97 
98  // Control DDU data
99  DTDDUData controlData(dduHeader,dduTrailer);
100  // check the CRC set in the FED trailer (FCRC errors)
101  controlData.checkCRCBit(index8 + datasize - 1*wordSize_64);
102 
103  // Check Status Words
104  vector<DTDDUFirstStatusWord> rosStatusWords;
105  // [BITS] 3 words of 8 bytes + "rosId" bytes
106  // In the case we are reading from DMA, the status word are swapped as the ROS data
107  if (localDAQ) {
108  // DDU channels from 1 to 4
109  for (int rosId = 0; rosId < 4; rosId++ ) {
110  int wordIndex8 = numberOf32Words*wordSize_32 - 3*wordSize_64 + wordSize_32 + rosId;
111  controlData.addROSStatusWord(DTDDUFirstStatusWord(index8[wordIndex8]));
112  }
113  // DDU channels from 5 to 8
114  for (int rosId = 0; rosId < 4; rosId++ ) {
115  int wordIndex8 = numberOf32Words*wordSize_32 - 3*wordSize_64 + rosId;
116  controlData.addROSStatusWord(DTDDUFirstStatusWord(index8[wordIndex8]));
117  }
118  // DDU channels from 9 to 12
119  for (int rosId = 0; rosId < 4; rosId++ ) {
120  int wordIndex8 = numberOf32Words*wordSize_32 - 2*wordSize_64 + wordSize_32 + rosId;
121  controlData.addROSStatusWord(DTDDUFirstStatusWord(index8[wordIndex8]));
122  }
123  }
124  else {
125  for (int rosId = 0; rosId < 12; rosId++ ) {
126  int wordIndex8 = numberOf32Words*wordSize_32 - 3*wordSize_64 + rosId;
127  controlData.addROSStatusWord(DTDDUFirstStatusWord(index8[wordIndex8]));
128  }
129  }
130 
131  int theROSList;
132  // [BITS] 2 words of 8 bytes + 4 bytes (half 64 bit word)
133  // In the case we are reading from DMA, the status word are swapped as the ROS data
134  if (localDAQ) {
135  DTDDUSecondStatusWord dduStatusWord(index32[numberOf32Words - 2*wordSize_64/wordSize_32]);
136  controlData.addDDUStatusWord(dduStatusWord);
137  theROSList = dduStatusWord.rosList();
138  }
139  else {
140  DTDDUSecondStatusWord dduStatusWord(index32[numberOf32Words - 2*wordSize_64/wordSize_32 + 1]);
141  controlData.addDDUStatusWord(dduStatusWord);
142  theROSList = dduStatusWord.rosList();
143  }
144 
145 
147  /* R O S d a t a */
149 
150  // Set the index to start looping on ROS data
151  // [BITS] one 8 bytes word
152  index32 += (wordSize_64)/wordSize_32;
153 
154  // Set the datasize to look only at ROS data
155  // [BITS] header, trailer, 2 status words
156  datasize -= 4*wordSize_64;
157 
158  // unpacking the ROS payload
159  ros25Unpacker->interpretRawData(index32, datasize, dduID, mapping, detectorProduct, triggerProduct, theROSList);
160 
161  // Perform DQM if requested
163  dataMonitor->processFED(controlData, ros25Unpacker->getROSsControlData(),dduID);
164 
165 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
void addDDUStatusWord(const DTDDUSecondStatusWord &word)
bool check()
Definition: FEDTrailer.cc:64
DTDataMonitorInterface * dataMonitor
Definition: DTDDUUnpacker.h:49
void addROSStatusWord(const DTDDUFirstStatusWord &word)
int rosList() const
Definition: DTDDUWords.h:901
virtual void interpretRawData(const unsigned int *index, int datasize, int dduID, edm::ESHandle< DTReadOutMapping > &mapping, std::auto_ptr< DTDigiCollection > &product, std::auto_ptr< DTLocalTriggerCollection > &product2, uint16_t rosList=0)
DTROS25Unpacker * ros25Unpacker
Definition: DTDDUUnpacker.h:47
bool isAvailable() const
Definition: Service.h:46
virtual void processFED(DTDDUData &dduData, const std::vector< DTROS25Data > &rosData, int ddu)=0
bool localDAQ
if data are read locally, status words are swapped
Definition: DTDDUUnpacker.h:40
virtual void interpretRawData(const unsigned int *index, int datasize, int dduID, edm::ESHandle< DTReadOutMapping > &mapping, std::auto_ptr< DTDigiCollection > &product, std::auto_ptr< DTLocalTriggerCollection > &product2, uint16_t rosList=0)
const edm::ParameterSet dduPSet
Definition: DTDDUUnpacker.h:37
void checkCRCBit(const unsigned char *trailer)
virtual ~DTDDUUnpacker()
Destructor.
bool check()
Check that the header is OK.
Definition: FEDHeader.cc:64
int lenght()
The length of the event fragment counted in 64-bit words including header and trailer.
Definition: FEDTrailer.cc:17
int bxID()
The bunch crossing number.
Definition: FEDHeader.cc:24
bool performDataIntegrityMonitor
perform DQM for DDU
Definition: DTDDUUnpacker.h:43
tuple cout
Definition: gather_cfg.py:121
int lvl1ID()
Level-1 event number generated by the TTC system.
Definition: FEDHeader.cc:20
DTDDUUnpacker(const edm::ParameterSet &ps)
Constructor.
const std::vector< DTROS25Data > & getROSsControlData() const