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