CMS 3D CMS Logo

RawDataUnpacker.cc
Go to the documentation of this file.
1 /****************************************************************************
2 *
3 * This is a part of TOTEM offline software.
4 * Authors:
5 * Jan Kašpar (jan.kaspar@gmail.com)
6 * Nicola Minafra
7 *
8 ****************************************************************************/
9 
13 
14 using namespace std;
15 using namespace edm;
16 using namespace pps;
17 
18 RawDataUnpacker::RawDataUnpacker(const edm::ParameterSet &iConfig)
19  : verbosity(iConfig.getUntrackedParameter<unsigned int>("verbosity", 0)) {}
20 
22  const FEDRawData &data,
23  vector<TotemFEDInfo> &fedInfoColl,
24  SimpleVFATFrameCollection &coll) const {
25  unsigned int size_in_words = data.size() / 8; // bytes -> words
26  if (size_in_words < 2) {
27  if (verbosity)
28  LogWarning("Totem") << "Error in RawDataUnpacker::run > "
29  << "Data in FED " << fedId << " too short (size = " << size_in_words << " words).";
30  return 1;
31  }
32 
33  fedInfoColl.emplace_back(fedId);
34 
35  return processOptoRxFrame((const word *)data.data(), size_in_words, fedInfoColl.back(), &coll);
36 }
37 
39  unsigned int frameSize,
40  TotemFEDInfo &fedInfo,
41  SimpleVFATFrameCollection *fc) const {
42  // get OptoRx metadata
43  unsigned long long head = buf[0];
44  unsigned long long foot = buf[frameSize - 1];
45 
46  fedInfo.setHeader(head);
47  fedInfo.setFooter(foot);
48 
49  unsigned int boe = (head >> 60) & 0xF;
50  unsigned int h0 = (head >> 0) & 0xF;
51 
52  unsigned long lv1 = (head >> 32) & 0xFFFFFF;
53  unsigned long bx = (head >> 20) & 0xFFF;
54  unsigned int optoRxId = (head >> 8) & 0xFFF;
55  unsigned int fov = (head >> 4) & 0xF;
56 
57  unsigned int eoe = (foot >> 60) & 0xF;
58  unsigned int f0 = (foot >> 0) & 0xF;
59  unsigned int fSize = (foot >> 32) & 0x3FF;
60 
61  // check header and footer structure
62  if (boe != 5 || h0 != 0 || eoe != 10 || f0 != 0 || fSize != frameSize) {
63  if (verbosity)
64  LogWarning("Totem") << "Error in RawDataUnpacker::processOptoRxFrame > "
65  << "Wrong structure of OptoRx header/footer: "
66  << "BOE=" << boe << ", H0=" << h0 << ", EOE=" << eoe << ", F0=" << f0
67  << ", size (OptoRx)=" << fSize << ", size (DATE)=" << frameSize << ". OptoRxID=" << optoRxId
68  << ". Skipping frame." << endl;
69  return 0;
70  }
71 
72  LogDebug("Totem") << "RawDataUnpacker::processOptoRxFrame: "
73  << "OptoRxId = " << optoRxId << ", BX = " << bx << ", LV1 = " << lv1
74  << ", frameSize = " << frameSize;
75 
78  processOptoRxFrameSampic(buf, frameSize, fedInfo, fc);
79  return 0;
80  }
81 
82  // parallel or serial transmission?
83  switch (fov) {
84  case 1:
85  return processOptoRxFrameSerial(buf, frameSize, fc);
86  case 2:
87  case 3:
88  return processOptoRxFrameParallel(buf, frameSize, fedInfo, fc);
89  default:
90  break;
91  }
92 
93  if (verbosity)
94  LogWarning("Totem") << "Error in RawDataUnpacker::processOptoRxFrame > "
95  << "Unknown FOV = " << fov << endl;
96 
97  return 0;
98 }
99 
101  unsigned int frameSize,
102  SimpleVFATFrameCollection *fc) const {
103  // get OptoRx metadata
104  unsigned int optoRxId = (buf[0] >> 8) & 0xFFF;
105 
106  // get number of subframes
107  unsigned int subFrames = (frameSize - 2) / 194;
108 
109  // process all sub-frames
110  unsigned int errorCounter = 0;
111  for (unsigned int r = 0; r < subFrames; ++r) {
112  for (unsigned int c = 0; c < 4; ++c) {
113  unsigned int head = (buf[1 + 194 * r] >> (16 * c)) & 0xFFFF;
114  unsigned int foot = (buf[194 + 194 * r] >> (16 * c)) & 0xFFFF;
115 
116  LogDebug("Totem") << "r = " << r << ", c = " << c << ": "
117  << "S = " << (head & 0x1) << ", BOF = " << (head >> 12) << ", EOF = " << (foot >> 12)
118  << ", ID = " << ((head >> 8) & 0xF) << ", ID' = " << ((foot >> 8) & 0xF);
119 
120  // stop if this GOH is NOT active
121  if ((head & 0x1) == 0)
122  continue;
123 
124  LogDebug("Totem") << "Header active (" << head << " -> " << (head & 0x1) << ").";
125 
126  // check structure
127  if (head >> 12 != 0x4 || foot >> 12 != 0xB || ((head >> 8) & 0xF) != ((foot >> 8) & 0xF)) {
128  std::ostringstream oss;
129  if (head >> 12 != 0x4)
130  oss << "\n\tHeader is not 0x4 as expected (0x" << std::hex << head << ").";
131  if (foot >> 12 != 0xB)
132  oss << "\n\tFooter is not 0xB as expected (0x" << std::hex << foot << ").";
133  if (((head >> 8) & 0xF) != ((foot >> 8) & 0xF))
134  oss << "\n\tIncompatible GOH IDs in header (0x" << std::hex << ((head >> 8) & 0xF) << ") and footer (0x"
135  << std::hex << ((foot >> 8) & 0xF) << ").";
136 
137  if (verbosity)
138  LogWarning("Totem") << "Error in RawDataUnpacker::processOptoRxFrame > "
139  << "Wrong payload structure (in GOH block row " << r << " and column " << c
140  << ") in OptoRx frame ID " << optoRxId << ". GOH block omitted." << oss.str() << endl;
141 
142  errorCounter++;
143  continue;
144  }
145 
146  // allocate memory for VFAT frames
147  unsigned int goh = (head >> 8) & 0xF;
148  vector<VFATFrame::word *> dataPtrs;
149  for (unsigned int fi = 0; fi < 16; fi++) {
150  TotemFramePosition fp(0, 0, optoRxId, goh, fi);
151  dataPtrs.push_back(fc->InsertEmptyFrame(fp)->getData());
152  }
153 
154  LogDebug("Totem").log([&](auto &l) {
155  l << "transposing GOH block at prefix: " << (optoRxId * 192 + goh * 16) << ", dataPtrs = ";
156  for (auto p : dataPtrs) {
157  l << p << " ";
158  }
159  });
160  // deserialization
161  for (int i = 0; i < 192; i++) {
162  int iword = 11 - i / 16; // number of current word (11...0)
163  int ibit = 15 - i % 16; // number of current bit (15...0)
164  unsigned int w = (buf[i + 2 + 194 * r] >> (16 * c)) & 0xFFFF;
165 
166  // Fill the current bit of the current word of all VFAT frames
167  for (int idx = 0; idx < 16; idx++) {
168  if (w & (1 << idx))
169  dataPtrs[idx][iword] |= (1 << ibit);
170  }
171  }
172  }
173  }
174 
175  return errorCounter;
176 }
177 
179  unsigned int frameSize,
180  TotemFEDInfo &fedInfo,
181  SimpleVFATFrameCollection *fc) const {
182  // get OptoRx metadata
183  unsigned long long head = buf[0];
184  unsigned int optoRxId = (head >> 8) & 0xFFF;
185 
186  // recast data as buffer or 16bit words, skip header
187  const uint16_t *payload = (const uint16_t *)(buf + 1);
188 
189  // read in OrbitCounter block
190  const uint32_t *ocPtr = (const uint32_t *)payload;
191  fedInfo.setOrbitCounter(*ocPtr);
192  payload += 2;
193 
194  // size in 16bit words, without header, footer and orbit counter block
195  unsigned int nWords = (frameSize - 2) * 4 - 2;
196 
197  // process all VFAT data
198  for (unsigned int offset = 0; offset < nWords;) {
199  unsigned int wordsProcessed = processVFATDataParallel(payload + offset, nWords, optoRxId, fc);
200  offset += wordsProcessed;
201  }
202 
203  return 0;
204 }
205 
207  unsigned int maxWords,
208  unsigned int optoRxId,
209  SimpleVFATFrameCollection *fc) const {
210  // start counting processed words
211  unsigned int wordsProcessed = 1;
212 
213  // padding word? skip it
214  if (buf[0] == 0xFFFF)
215  return wordsProcessed;
216 
217  // check header flag
218  unsigned int hFlag = (buf[0] >> 8) & 0xFF;
219  if (hFlag != vmCluster && hFlag != vmRaw && hFlag != vmDiamondCompact) {
220  if (verbosity)
221  LogWarning("Totem") << "Error in RawDataUnpacker::processVFATDataParallel > "
222  << "Unknown header flag " << hFlag << ". Skipping this word." << endl;
223  return wordsProcessed;
224  }
225 
226  // compile frame position
227  // NOTE: DAQ group uses terms GOH and fiber in the other way
228  unsigned int gohIdx = (buf[0] >> 4) & 0xF;
229  unsigned int fiberIdx = (buf[0] >> 0) & 0xF;
230  TotemFramePosition fp(0, 0, optoRxId, gohIdx, fiberIdx);
231 
232  // prepare temporary VFAT frame
233  VFATFrame f;
234  VFATFrame::word *fd = f.getData();
235 
236  // copy footprint, BC, EC, Flags, ID, if they exist
237  uint8_t presenceFlags = 0;
238 
239  if (((buf[wordsProcessed] >> 12) & 0xF) == 0xA) // BC
240  {
241  presenceFlags |= 0x1;
242  fd[11] = buf[wordsProcessed];
243  wordsProcessed++;
244  }
245 
246  if (((buf[wordsProcessed] >> 12) & 0xF) == 0xC) // EC, flags
247  {
248  presenceFlags |= 0x2;
249  fd[10] = buf[wordsProcessed];
250  wordsProcessed++;
251  }
252 
253  if (((buf[wordsProcessed] >> 12) & 0xF) == 0xE) // ID
254  {
255  presenceFlags |= 0x4;
256  fd[9] = buf[wordsProcessed];
257  wordsProcessed++;
258  }
259 
260  // save offset where channel data start
261  unsigned int dataOffset = wordsProcessed;
262 
263  // find trailer
264  switch (hFlag) {
265  case vmCluster: {
266  unsigned int nCl = 0;
267  while ((buf[wordsProcessed + nCl] >> 12) != 0xF && (wordsProcessed + nCl < maxWords))
268  nCl++;
269  wordsProcessed += nCl;
270  } break;
271  case vmRaw:
272  wordsProcessed += 9;
273  break;
274  case vmDiamondCompact: {
275  wordsProcessed--;
276  while ((buf[wordsProcessed] & 0xFFF0) != 0xF000 && (wordsProcessed < maxWords))
277  wordsProcessed++;
278  } break;
279  }
280 
281  // process trailer
282  unsigned int tSig = buf[wordsProcessed] >> 12;
283  unsigned int tErrFlags = (buf[wordsProcessed] >> 8) & 0xF;
284  unsigned int tSize = buf[wordsProcessed] & 0xFF;
285 
286  f.setDAQErrorFlags(tErrFlags);
287 
288  // consistency checks
289  bool skipFrame = false;
290  stringstream ess;
291 
292  if (tSig != 0xF) {
293  if (verbosity)
294  ess << " Wrong trailer signature (" << tSig << ")." << endl;
295  skipFrame = true;
296  }
297 
298  if (tErrFlags != 0) {
299  if (verbosity)
300  ess << " Error flags not zero (" << tErrFlags << ")." << endl;
301  skipFrame = true;
302  }
303 
304  wordsProcessed++;
305 
306  if (tSize != wordsProcessed) {
307  if (verbosity)
308  ess << " Trailer size (" << tSize << ") does not match with words processed (" << wordsProcessed << ")."
309  << endl;
310  skipFrame = true;
311  }
312 
313  if (skipFrame) {
314  if (verbosity)
315  LogWarning("Totem") << "Error in RawDataUnpacker::processVFATDataParallel > Frame at " << fp
316  << " has the following problems and will be skipped.\n"
317  << endl
318  << ess.rdbuf();
319 
320  return wordsProcessed;
321  }
322 
323  // get channel data - cluster mode
324  if (hFlag == vmCluster) {
325  for (unsigned int nCl = 0; (buf[dataOffset + nCl] >> 12) != 0xF && (dataOffset + nCl < maxWords); ++nCl) {
326  const uint16_t &w = buf[dataOffset + nCl];
327  unsigned int upperBlock = w >> 8;
328  unsigned int clSize = upperBlock & 0x7F;
329  unsigned int clPos = (w >> 0) & 0xFF;
330 
331  // special case: upperBlock=0xD0 => numberOfClusters
332  if (upperBlock == 0xD0) {
333  presenceFlags |= 0x10;
334  f.setNumberOfClusters(clPos);
335  continue;
336  }
337 
338  // special case: size=0 means chip full
339  if (clSize == 0)
340  clSize = 128;
341 
342  // activate channels
343  // convention - range <pos, pos-size+1>
344  signed int chMax = clPos;
345  signed int chMin = clPos - clSize + 1;
346  if (chMax < 0 || chMax > 127 || chMin < 0 || chMin > 127 || chMin > chMax) {
347  if (verbosity)
348  LogWarning("Totem") << "Error in RawDataUnpacker::processVFATDataParallel > "
349  << "Invalid cluster (pos=" << clPos << ", size=" << clSize << ", min=" << chMin
350  << ", max=" << chMax << ") at " << fp << ". Skipping this cluster." << endl;
351  continue;
352  }
353 
354  for (signed int ch = chMin; ch <= chMax; ch++) {
355  unsigned int wi = ch / 16;
356  unsigned int bi = ch % 16;
357  fd[wi + 1] |= (1 << bi);
358  }
359  }
360  }
361 
362  // get channel data and CRC - raw mode
363  if (hFlag == vmRaw) {
364  for (unsigned int i = 0; i < 8; i++)
365  fd[8 - i] = buf[dataOffset + i];
366 
367  // copy CRC
368  presenceFlags |= 0x8;
369  fd[0] = buf[dataOffset + 8];
370  }
371 
372  // get channel data for diamond compact mode
373  if (hFlag == vmDiamondCompact) {
374  for (unsigned int i = 1; (buf[i + 1] & 0xFFF0) != 0xF000 && (i + 1 < maxWords); i++) {
375  if ((buf[i] & 0xF000) == VFAT_HEADER_OF_EC) {
376  // Event Counter word is found
377  fd[10] = buf[i];
378  continue;
379  }
380  switch (buf[i] & 0xF800) {
382  // word 2 of the diamond VFAT frame is found
383  fd[1] = buf[i + 1];
384  fd[2] = buf[i];
385  break;
387  // word 3 of the diamond VFAT frame is found
388  fd[3] = buf[i];
389  fd[4] = buf[i - 1];
390  break;
392  // word 5 of the diamond VFAT frame is found
393  fd[5] = buf[i];
394  fd[6] = buf[i - 1];
395  break;
397  // word 7 of the diamond VFAT frame is found
398  fd[7] = buf[i];
399  fd[8] = buf[i - 1];
400  break;
401  default:
402  break;
403  }
404  presenceFlags |= 0x8;
405  }
406  }
407 
408  // save frame to output
409  f.setPresenceFlags(presenceFlags);
410  fc->Insert(fp, f);
411 
412  return wordsProcessed;
413 }
414 
416  unsigned int frameSize,
417  TotemFEDInfo &fedInfo,
418  SimpleVFATFrameCollection *fc) const {
419  unsigned int optoRxId = (buf[0] >> 8) & 0xFFF;
420 
421  LogDebug("RawDataUnpacker::processOptoRxFrameSampic")
422  << "Processing sampic frame: OptoRx " << optoRxId << " framesize: " << frameSize;
423 
424  unsigned int orbitCounterVFATFrameWords = 6;
425  unsigned int sizeofVFATPayload = 12;
426 
427  const VFATFrame::word *VFATFrameWordPtr = (const VFATFrame::word *)buf;
428  VFATFrameWordPtr += orbitCounterVFATFrameWords - 1;
429 
430  LogDebug("RawDataUnpacker::processOptoRxFrameSampic")
431  << "Framesize: " << frameSize << "\tframes: " << frameSize / (sizeofVFATPayload + 2);
432 
433  unsigned int nWords = (frameSize - 2) * 4 - 2;
434 
435  for (unsigned int i = 1; i * (sizeofVFATPayload + 2) < nWords; ++i) {
436  // compile frame position
437  // NOTE: DAQ group uses terms GOH and fiber in the other way
438  unsigned int fiberIdx = (*(++VFATFrameWordPtr)) & 0xF;
439  unsigned int gohIdx = (*VFATFrameWordPtr >> 4) & 0xF;
440  TotemFramePosition fp(0, 0, optoRxId, gohIdx, fiberIdx);
441 
442  LogDebug("RawDataUnpacker::processOptoRxFrameSampic")
443  << "OptoRx: " << optoRxId << " Goh: " << gohIdx << " Idx: " << fiberIdx;
444 
445  // prepare temporary VFAT frame
446  VFATFrame frame(++VFATFrameWordPtr);
447  VFATFrameWordPtr += sizeofVFATPayload;
448 
449  if (*(VFATFrameWordPtr) != 0xf00e) {
450  edm::LogError("RawDataUnpacker::processOptoRxFrameSampic") << "Wrong trailer " << *VFATFrameWordPtr;
451  continue;
452  }
453  // save frame to output
454  frame.setPresenceFlags(1);
455  fc->Insert(fp, frame);
456 
457  LogDebug("RawDataUnpacker::processOptoRxFrameSampic") << "Trailer: " << std::hex << *VFATFrameWordPtr;
458  }
459 
460  return 0;
461 }
OptoRx headers and footers.
Definition: TotemFEDInfo.h:17
static constexpr unsigned int VFAT_DIAMOND_HEADER_OF_WORD_7
static constexpr unsigned int VFAT_HEADER_OF_EC
VFATFrame::word * getData()
Definition: VFATFrame.h:38
int processOptoRxFrame(const word *buf, unsigned int frameSize, TotemFEDInfo &fedInfo, SimpleVFATFrameCollection *fc) const
Process one Opto-Rx (or LoneG) frame.
T w() const
uint16_t word
Definition: VFATFrame.h:21
int processOptoRxFrameSampic(const word *buffer, unsigned int frameSize, TotemFEDInfo &fedInfo, SimpleVFATFrameCollection *fc) const
Process one Opto-Rx frame that contains SAMPIC frames.
static constexpr unsigned int VFAT_DIAMOND_HEADER_OF_WORD_5
int processVFATDataParallel(const uint16_t *buf, unsigned int maxWords, unsigned int OptoRxId, SimpleVFATFrameCollection *fc) const
Process data from one VFAT in parallel (new) format.
Log< level::Error, false > LogError
static constexpr unsigned int VFAT_DIAMOND_HEADER_OF_WORD_3
void setFooter(uint64_t _f)
Definition: TotemFEDInfo.h:35
double f[11][100]
static constexpr unsigned int VFAT_DIAMOND_HEADER_OF_WORD_2
VFATFrame * InsertEmptyFrame(TotemFramePosition index)
inserts an empty (default) frame to the given position and returns pointer to the frame ...
unsigned char verbosity
int run(int fedId, const FEDRawData &data, std::vector< TotemFEDInfo > &fedInfoColl, SimpleVFATFrameCollection &coll) const
Unpack data from FED with fedId into ‘coll’ collection.
void setOrbitCounter(uint32_t oc)
Definition: TotemFEDInfo.h:32
HLT enums.
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:79
int processOptoRxFrameSerial(const word *buffer, unsigned int frameSize, SimpleVFATFrameCollection *fc) const
Process one Opto-Rx frame in serial (old) format.
Log< level::Warning, false > LogWarning
fd
Definition: ztee.py:136
int processOptoRxFrameParallel(const word *buffer, unsigned int frameSize, TotemFEDInfo &fedInfo, SimpleVFATFrameCollection *fc) const
Process one Opto-Rx frame in parallel (new) format.
void setHeader(uint64_t _h)
Definition: TotemFEDInfo.h:24
void Insert(const TotemFramePosition &index, const VFATFrame &frame)
#define LogDebug(id)