CbmRoot
Loading...
Searching...
No Matches
CbmMcbm2018UnpackerAlgoTof.cxx
Go to the documentation of this file.
1/* Copyright (C) 2019-2020 Facility for Antiproton and Ion Research in Europe, Darmstadt
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Pierre-Alain Loizeau [committer], Norbert Herrmann */
4
5// -----------------------------------------------------------------------------
6// ----- -----
7// ----- CbmMcbm2018UnpackerAlgoTof -----
8// ----- Created 10.02.2019 by P.-A. Loizeau -----
9// ----- -----
10// -----------------------------------------------------------------------------
11
13
16#include "CbmMcbm2018TofPar.h"
17#include "CbmTofAddress.h"
18#include "CbmTofDetectorId_v14a.h" // in cbmdata/tof
19
20#include <Logger.h>
21
22#include "TCanvas.h"
23#include "TH1.h"
24#include "TH2.h"
25#include "TList.h"
26#include "TProfile.h"
27#include "TROOT.h"
28#include "TString.h"
29
30#include <fstream>
31#include <iomanip>
32#include <iostream>
33#include <vector>
34
35#include <stdint.h>
36
37
38static Int_t NMappingWarnings = 100;
39/*
40static uint32_t pat_mess[8]={8*0};
41std::vector< std::vector <uint32_t> > Pat_Request;
42std::vector<Bool_t> bGdpbOK;
43static Bool_t bEnableOut=kFALSE;
44static Int_t fiMsGood=0;
45static Int_t fiMsBad=0;
46*/
47// -------------------------------------------------------------------------
50 ,
52 fbMonitorMode(kFALSE)
53 , fbDebugMonitorMode(kFALSE)
54 , fvbMaskedComponents()
55 , fUnpackPar(nullptr)
56 , fuNrOfGdpbs(0)
57 , fGdpbIdIndexMap()
58 , fuNrOfFeePerGdpb(0)
59 , fuNrOfGet4PerFee(0)
60 , fuNrOfChannelsPerGet4(0)
61 , fuNrOfChannelsPerFee(0)
62 , fuNrOfGet4(0)
63 , fuNrOfGet4PerGdpb(0)
64 , fuNrOfChannelsPerGdpb(0)
65 , fuNrOfGbtx(0)
66 , fuNrOfModules(0)
67 , fviNrOfRpc()
68 , fviRpcType()
69 , fviRpcSide()
70 , fviModuleId()
71 , fviRpcChUId()
72 , fdTimeOffsetNs(0.0)
73 , fuDiamondDpbIdx(99)
74 , fulCurrentTsIdx(0)
75 , fulCurrentMsIdx(0)
76 , fuCurrentMsSysId(0)
77 , fdTsStartTime(-1.0)
78 , fdTsStopTimeCore(-1.0)
79 , fdMsTime(-1.0)
80 , fuMsIndex(0)
81 , fuCurrentEquipmentId(0)
82 , fuCurrDpbId(0)
83 , fuCurrDpbIdx(0)
84 , fuGet4Id(0)
85 , fuGet4Nr(0)
86 , fvulCurrentEpoch()
87 , fvulCurrentEpochCycle()
88 , fvulCurrentEpochFull()
89 , fdStartTime(0.0)
90 , fdStartTimeMsSz(0.0)
91 , ftStartTimeUnix(std::chrono::steady_clock::now())
92 , fvvmEpSupprBuffer()
93 , fvulGdpbTsMsb()
94 , fvulGdpbTsLsb()
95 , fvulStarTsMsb()
96 , fvulStarTsMid()
97 , fvulGdpbTsFullLast()
98 , fvulStarTsFullLast()
99 , fvuStarTokenLast()
100 , fvuStarDaqCmdLast()
101 , fvuStarTrigCmdLast()
102 , fdRefTime(0.)
103 , fdLastDigiTime(0.)
104 , fdFirstDigiTimeDif(0.)
105 , fdEvTime0(0.)
106 , fhRawTDigEvBmon(nullptr)
107 , fhRawTDigRef0(nullptr)
108 , fhRawTDigRef(nullptr)
109 , fhRawTRefDig0(nullptr)
110 , fhRawTRefDig1(nullptr)
111 , fhRawDigiLastDigi(nullptr)
112 , fhRawTotCh()
113 , fhChCount()
114 , fhChCountRemap()
115 , fvbChanThere()
116 , fhChanCoinc()
117 , fhDetChanCoinc(nullptr)
118{
119}
121{
123 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
124 fvvmEpSupprBuffer[uGdpb].clear();
125 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
126 if (nullptr != fParCList) delete fParCList;
127 if (nullptr != fUnpackPar) delete fUnpackPar;
128}
129
130// -------------------------------------------------------------------------
132{
133 LOG(info) << "Initializing mCBM TOF 2019 unpacker algo";
134
135 return kTRUE;
136}
139{
140 /*
142 LOG(info)<<"<I> MS statistics - Good: " << fiMsGood <<", Bad: " << fiMsBad;
143*/
145}
146
147// -------------------------------------------------------------------------
149{
150 LOG(info) << "Init parameter containers for CbmMcbm2018UnpackerAlgoTof";
151
152 Bool_t initOK = ReInitContainers();
153
154 return initOK;
155}
157{
158 LOG(info) << "**********************************************";
159 LOG(info) << "ReInit parameter containers for CbmMcbm2018UnpackerAlgoTof";
160
161 fUnpackPar = (CbmMcbm2018TofPar*) fParCList->FindObject("CbmMcbm2018TofPar");
162 if (nullptr == fUnpackPar) {
163 LOG(error) << " CbmMcbm2018TofPar not found ";
164 return kFALSE;
165 }
166 Bool_t initOK = InitParameters();
167
168 return initOK;
169}
171{
172 if (nullptr == fParCList) fParCList = new TList();
173 fUnpackPar = new CbmMcbm2018TofPar("CbmMcbm2018TofPar");
174 fParCList->Add(fUnpackPar);
175
176 return fParCList;
177}
179{
180 LOG(info) << "InitParameters from " << fUnpackPar;
182
184 LOG(info) << "Nr. of Tof GDPBs: " << fuNrOfGdpbs;
185
187 LOG(info) << "Nr. of FEES per Tof GDPB: " << fuNrOfFeePerGdpb;
188
190 LOG(info) << "Nr. of GET4 per Tof FEE: " << fuNrOfGet4PerFee;
191
193 LOG(info) << "Nr. of channels per GET4: " << fuNrOfChannelsPerGet4;
194
196 LOG(info) << "Nr. of channels per FEE: " << fuNrOfChannelsPerFee;
197
199 LOG(info) << "Nr. of GET4s: " << fuNrOfGet4;
200
202 LOG(info) << "Nr. of GET4s per GDPB: " << fuNrOfGet4PerGdpb;
203
205 LOG(info) << "Nr. of channels per GDPB: " << fuNrOfChannelsPerGdpb;
206
207 fGdpbIdIndexMap.clear();
208 for (UInt_t i = 0; i < fuNrOfGdpbs; ++i) {
210 LOG(info) << "GDPB Id of TOF " << i << " : " << std::hex << fUnpackPar->GetGdpbId(i) << std::dec;
211 } // for( UInt_t i = 0; i < fuNrOfGdpbs; ++i )
212
214 LOG(info) << "Nr. of GBTx: " << fuNrOfGbtx;
215
216 fviRpcType.resize(fuNrOfGbtx);
217 fviModuleId.resize(fuNrOfGbtx);
218 fviNrOfRpc.resize(fuNrOfGbtx);
219 fviRpcSide.resize(fuNrOfGbtx);
220 for (UInt_t uGbtx = 0; uGbtx < fuNrOfGbtx; ++uGbtx) {
221 fviNrOfRpc[uGbtx] = fUnpackPar->GetNrOfRpc(uGbtx);
222 fviRpcType[uGbtx] = fUnpackPar->GetRpcType(uGbtx);
223 fviRpcSide[uGbtx] = fUnpackPar->GetRpcSide(uGbtx);
224 fviModuleId[uGbtx] = fUnpackPar->GetModuleId(uGbtx);
225 } // for( UInt_t uGbtx = 0; uGbtx < fuNrOfGbtx; ++uGbtx)
226
227 UInt_t uNrOfChannels = fuNrOfGet4 * fuNrOfChannelsPerGet4;
228 LOG(info) << "Nr. of possible Tof channels: " << uNrOfChannels;
229
230 // CbmTofDetectorId* fTofId = new CbmTofDetectorId_v14a();
231 fviRpcChUId.resize(uNrOfChannels);
232 UInt_t iCh = 0;
233 for (UInt_t iGbtx = 0; iGbtx < fuNrOfGbtx; ++iGbtx) {
234 Int_t iModuleIdMap = fviModuleId[iGbtx];
235 switch (fviRpcType[iGbtx]) {
236
237 case 0: // CBM modules
238 if (fviRpcSide[iGbtx] < 2) { // mTof modules
239 LOG(info) << " Map mTof box " << fviModuleId[iGbtx] << " at GBTX - iCh = " << iCh;
240 const Int_t RpcMap[5] = {4, 2, 0, 3, 1};
241 for (Int_t iRpc = 0; iRpc < fviNrOfRpc[iGbtx]; iRpc++) {
242 Int_t iStrMax = 32;
243 Int_t iChNext = 1;
244
245 for (Int_t iStr = 0; iStr < iStrMax; iStr++) {
246 Int_t iStrMap = iStr;
247 Int_t iRpcMap = RpcMap[iRpc];
248
249 if (fviRpcSide[iGbtx] == 0) iStrMap = 31 - iStr;
250 if (fviModuleId[iGbtx] > -1)
251 fviRpcChUId[iCh] = CbmTofAddress::GetUniqueAddress(fviModuleId[iGbtx], iRpcMap, iStrMap,
252 fviRpcSide[iGbtx], fviRpcType[iGbtx]);
253 else
254 fviRpcChUId[iCh] = 0;
255 // LOG(debug)<<Form("Map Ch %d to Address 0x%08x",iCh,fviRpcChUId[iCh]);
256 iCh += iChNext;
257 }
258 }
259 }
260 break;
261
262 case 1: // STAR eTOF modules
263 if (fviRpcSide[iGbtx] < 2) { // mTof modules
264 LOG(info) << "Start eTOF module side " << fviRpcSide[iGbtx] << " at " << iCh;
265 const Int_t RpcMap[3] = {0, 1, 2};
266 for (Int_t iRpc = 0; iRpc < fviNrOfRpc[iGbtx]; iRpc++) {
267 Int_t iStrMax = 32;
268 Int_t iChNext = 1;
269
270 for (Int_t iStr = 0; iStr < iStrMax; iStr++) {
271 Int_t iStrMap = iStr;
272 Int_t iRpcMap = RpcMap[iRpc];
273
274 if (fviRpcSide[iGbtx] == 0) iStrMap = 31 - iStr;
275 if (fviModuleId[iGbtx] > -1)
276 fviRpcChUId[iCh] = CbmTofAddress::GetUniqueAddress(fviModuleId[iGbtx], iRpcMap, iStrMap,
277 fviRpcSide[iGbtx], fviRpcType[iGbtx]);
278 else
279 fviRpcChUId[iCh] = 0;
280 // LOG(debug)<<Form("Map Ch %d to Address 0x%08x",iCh,fviRpcChUId[iCh]);
281 iCh += iChNext;
282 }
283 }
284 }
285 iCh += 64;
286 break;
287
289 case 5: {
290 LOG(info) << " Map diamond at GBTX - iCh = " << iCh;
291 for (UInt_t uFee = 0; uFee < fUnpackPar->GetNrOfFeePerGbtx(); ++uFee) {
292 for (UInt_t uCh = 0; uCh < fUnpackPar->GetNrOfChannelsPerFee(); ++uCh) {
293 /*
294 if( uFee < 4 && 0 == uCh ) {
295 fviRpcChUId[ iCh ] = CbmTofAddress::GetUniqueAddress(
296 fviModuleId[iGbtx],
297 0, uFee + 4 * fviRpcSide[iGbtx],
298 0, fviRpcType[iGbtx] );
299 LOG(info) << Form( "Map Bmon Ch %d to Address 0x%08x", iCh, fviRpcChUId[iCh] );
300 }
301 else fviRpcChUId[ iCh ] = 0;
302*/
303
305 if (0 == uFee && 1 == fviNrOfRpc[iGbtx]) {
306 switch (uCh % 8) {
307 case 0:
308 case 2:
309 case 4: {
314 UInt_t uChannelBmon = uCh / 8 + 4 * fviRpcSide[iGbtx];
315 fviRpcChUId[iCh] =
316 CbmTofAddress::GetUniqueAddress(fviModuleId[iGbtx], 0, uChannelBmon, 0, fviRpcType[iGbtx]);
317 LOG(info) << Form("Bmon channel: %u from GBTx %2u Fee %2u "
318 "Channel %2u, indx %d address %08x",
319 uChannelBmon, iGbtx, uFee, uCh, iCh, fviRpcChUId[iCh]);
320 break;
321 } // Valid Bmon channel
322 default: {
323 fviRpcChUId[iCh] = 0;
324 } // Invalid Bmon channel
325 } // switch( uCh % 4 )
326 } // if( 0 == uFee )
327
328 iCh++;
329 } // for( UInt_t uCh = 0; uCh < fUnpackPar->GetNrOfChannelsPerFee(); ++uCh )
330 } // for( UInt_t uFee = 0; uFee < fUnpackPar->GetNrOfFeePerGbtx(); ++uFee )
331 } // if( 5 == fviRpcType[iGbtx] )
332 break;
333
334 case 78: // cern-20-gap + ceramic module
335 {
336 LOG(info) << " Map CERN 20 gap at GBTX - iCh = " << iCh;
337 // clang-format off
338 const Int_t StrMap[32] = {0, 1, 2, 3, 4, 31, 5, 6, 7, 30, 8,
339 9, 10, 29, 11, 12, 13, 14, 28, 15, 16, 17,
340 18, 27, 26, 25, 24, 23, 22, 21, 20, 19};
341 // clang-format on
342 Int_t iModuleId = 0;
343 Int_t iModuleType = 7;
344 Int_t iRpcMap = 0;
345 for (Int_t iFeet = 0; iFeet < 2; iFeet++) {
346 for (Int_t iStr = 0; iStr < 32; iStr++) {
347 Int_t iStrMap = 31 - 12 - StrMap[iStr];
348 Int_t iSideMap = iFeet;
349 if (iStrMap < 20)
350 fviRpcChUId[iCh] = CbmTofAddress::GetUniqueAddress(iModuleId, iRpcMap, iStrMap, iSideMap, iModuleType);
351 else
352 fviRpcChUId[iCh] = 0;
353 iCh++;
354 }
355 }
356
357 LOG(info) << " Map end CERN 20 gap at GBTX - iCh = " << iCh;
358 }
359 [[fallthrough]]; // fall through is intended
360 case 8: // ceramics
361 {
362 Int_t iModuleId = 0;
363 Int_t iModuleType = 8;
364 for (Int_t iRpc = 0; iRpc < 8; iRpc++) {
365 fviRpcChUId[iCh] = CbmTofAddress::GetUniqueAddress(iModuleId, 7 - iRpc, 0, 0, iModuleType);
366 iCh++;
367 }
368 iCh += (24 + 2 * 32);
369 }
370
371 LOG(info) << " Map end ceramics box at GBTX - iCh = " << iCh;
372 break;
373
374 case 4: // intended fallthrough
375 [[fallthrough]];
376 case 9: // Star2 boxes
377 {
378 LOG(info) << " Map Star2 box " << fviModuleId[iGbtx] << " at GBTX - iCh = " << iCh;
379 const Int_t iRpc[5] = {1, -1, 1, 0, 0};
380 const Int_t iSide[5] = {1, -1, 0, 1, 0};
381 for (Int_t iFeet = 0; iFeet < 5; iFeet++) {
382 for (Int_t iStr = 0; iStr < 32; iStr++) {
383 Int_t iStrMap = iStr;
384 Int_t iRpcMap = iRpc[iFeet];
385 Int_t iSideMap = iSide[iFeet];
386 if (iSideMap == 0) iStrMap = 31 - iStr;
387 switch (fviRpcSide[iGbtx]) {
388 case 0:; break;
389 case 1:;
390 iRpcMap = 1 - iRpcMap; // swap counters
391 break;
392 case 2:
393 switch (iFeet) {
394 case 1:
395 iRpcMap = iRpc[4];
396 iSideMap = iSide[4];
397 iStrMap = 31 - iStrMap;
398 break;
399 case 4:
400 iRpcMap = iRpc[1];
401 iSideMap = iSide[1];
402 break;
403 default:;
404 }
405 break;
406 }
407 if (iSideMap > -1)
408 fviRpcChUId[iCh] =
409 CbmTofAddress::GetUniqueAddress(fviModuleId[iGbtx], iRpcMap, iStrMap, iSideMap, fviRpcType[iGbtx]);
410 else
411 fviRpcChUId[iCh] = 0;
412 iCh++;
413 }
414 }
415 } break;
416
417 case 6: // Buc box
418 {
419 LOG(info) << " Map Buc box at GBTX - iCh = " << iCh;
420 const Int_t iRpc[5] = {0, -1, 0, 1, 1};
421 const Int_t iSide[5] = {1, -1, 0, 1, 0};
422 for (Int_t iFeet = 0; iFeet < 5; iFeet++) {
423 for (Int_t iStr = 0; iStr < 32; iStr++) {
424 Int_t iStrMap = iStr;
425 Int_t iRpcMap = iRpc[iFeet];
426 Int_t iSideMap = iSide[iFeet];
427 switch (fviRpcSide[iGbtx]) {
428 case 0:; break;
429 case 1: // HD cosmic 2019, Buc2018, v18n
430 iStrMap = 31 - iStr;
431 iRpcMap = 1 - iRpcMap;
432 break;
433 case 2: // v18m_cosmicHD
434 // iStrMap=31-iStr;
435 iSideMap = 1 - iSideMap;
436 break;
437 case 3:
438 iStrMap = 31 - iStr;
439 iRpcMap = 1 - iRpcMap;
440 iSideMap = 1 - iSideMap;
441 break;
442 case 4: // HD cosmic 2019, Buc2018, v18o
443 iRpcMap = 1 - iRpcMap;
444 break;
445 case 5: // HD cosmic 2020, Buc2018, v20a
446 iStrMap = 31 - iStr;
447 break;
448 case 6: //BUC special
449 {
450 switch (fviModuleId[iGbtx]) {
451 case 0: iRpcMap = 0; break;
452 case 1: iRpcMap = 1; break;
453 }
454 if (iFeet % 2 == 1) iModuleIdMap = 1;
455 else
456 iModuleIdMap = 0;
457
458 switch (iFeet) {
459 case 0:
460 case 3: iSideMap = 0; break;
461 case 1:
462 case 2: iSideMap = 1; break;
463 }
464 } break;
465
466 case 7: {
467 // clang-format off
468 const Int_t iChMap[160]={
469 127, 126, 125, 124, 12, 13, 14, 15, 7, 6, 5, 4, 28, 29, 30, 31, 123, 122, 121, 120, 8, 9, 10, 11, 107, 106, 105, 104, 108, 109, 110, 111,
470 39, 38, 37, 36, 52, 53, 54, 55, 63, 62, 61, 60, 128, 129, 130, 131, 43, 42, 41, 40, 148, 149, 150, 151, 59, 58, 57, 56, 132, 133, 134, 135,
471 139, 138, 137, 136, 140, 141, 142, 143, 99, 98, 97, 96, 64, 65, 66, 67, 103, 102, 101, 100, 84, 85, 86, 87, 155, 154, 153, 152, 68, 69, 70, 71,
472 159, 158, 157, 156, 144, 145, 146, 147, 47, 46, 45, 44, 76, 77, 78, 79, 51, 50, 49, 48, 20, 21, 22, 23, 35, 34, 33, 32, 116, 117, 118, 119,
473 75, 74, 73, 72, 92, 93, 94, 95, 19, 18, 17, 16, 80, 81, 82, 83, 115, 114, 113, 112, 24, 25, 26, 27, 91, 90, 89, 88, 0, 1, 2, 3
474 };
475 // clang-format on
476 Int_t iInd = iFeet * 32 + iStr;
477 Int_t i = 0;
478 for (; i < 160; i++)
479 if (iInd == iChMap[i]) break;
480 iStrMap = i % 32;
481 Int_t iFeetInd = (i - iStrMap) / 32;
482 switch (iFeet) {
483 case 0:
484 iRpcMap = 0;
485 iSideMap = 1;
486 break;
487 case 1:
488 iRpcMap = 1;
489 iSideMap = 1;
490 break;
491 case 2:
492 iRpcMap = 0;
493 iSideMap = 0;
494 break;
495 case 3:
496 iRpcMap = 1;
497 iSideMap = 0;
498 break;
499 case 4: iSideMap = -1; break;
500 }
501 iModuleIdMap = fviModuleId[iGbtx];
502 LOG(info) << "Buc of GBTX " << iGbtx << " Ch " << iCh
503 << Form(", Feet %1d, Str %2d, Ind %3d, i %3d, FeetInd %1d, Rpc %1d, Side %1d, Str %2d ",
504 iFeet, iStr, iInd, i, iFeetInd, iRpcMap, iSideMap, iStrMap);
505 } break;
506 default:;
507 }
508 if (iSideMap > -1)
509 fviRpcChUId[iCh] =
510 CbmTofAddress::GetUniqueAddress(iModuleIdMap, iRpcMap, iStrMap, iSideMap, fviRpcType[iGbtx]);
511 else
512 fviRpcChUId[iCh] = 0;
513
514 iCh++;
515 }
516 }
517 } break;
518
519 case -1:
520 LOG(info) << " Found unused GBTX link at iCh = " << iCh;
521 iCh += 160;
522 break;
523
524 default: LOG(error) << "Invalid Tof Type specifier ";
525 }
526 }
527 TString sPrintout = "";
528 for (UInt_t uCh = 0; uCh < uNrOfChannels; ++uCh) {
529 if (0 == uCh % 8) sPrintout += "\n";
530 if (0 == uCh % fuNrOfChannelsPerGdpb) sPrintout += Form("\n Gdpb %u\n", uCh / fuNrOfChannelsPerGdpb);
531 sPrintout += Form(" 0x%08x", fviRpcChUId[uCh]);
532 } // for( UInt_t i = 0; i < uNrOfChannels; ++i)
533 LOG(info) << sPrintout;
534
535 // Request masks
536 /*
537 LOG(info) << " Load " << fUnpackPar->GetNrReqPattern() << " GET4 Request masks for " << fuNrOfGdpbs << " Gdpbs ";
538 if(fUnpackPar->GetNrReqPattern()>0){
539 bGdpbOK.resize(fuNrOfGdpbs);
540 Pat_Request.resize(fuNrOfGdpbs);
541 Int_t iInd=0;
542 for(Int_t iGdpb=0; iGdpb<fuNrOfGdpbs; iGdpb++) {
543 bGdpbOK[iGdpb]=kTRUE;
544 Pat_Request[iGdpb].resize(fUnpackPar->GetNrReqPattern());
545 for (Int_t iPat=0; iPat<fUnpackPar->GetNrReqPattern(); iPat++) {
546 UInt_t PatGet4=fUnpackPar->GetReqPattern(iInd++);
547 for( UInt_t uBit = 0; uBit < 32; ++uBit ) {
548 if( ( PatGet4 >> uBit ) & 0x1 ) {
549 UInt_t iGet4=iPat*32+uBit;
550 UInt_t uElink=fUnpackPar->Get4IdxToElinkIdx(iGet4);
551 UInt_t ubit=uElink%32;
552 UInt_t iEPat=(uElink-ubit)/32;
553 Pat_Request[iGdpb][iEPat] |= (0x1 << ubit);
554 }
555 }
556 }
557 }
558 }
559 */
561 fvulCurrentEpoch.resize(fuNrOfGdpbs, 0);
564
575 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
576 fvulGdpbTsMsb[uGdpb] = 0;
577 fvulGdpbTsLsb[uGdpb] = 0;
578 fvulStarTsMsb[uGdpb] = 0;
579 fvulStarTsMid[uGdpb] = 0;
580 fvulGdpbTsFullLast[uGdpb] = 0;
581 fvulStarTsFullLast[uGdpb] = 0;
582 fvuStarTokenLast[uGdpb] = 0;
583 fvuStarDaqCmdLast[uGdpb] = 0;
584 fvuStarTrigCmdLast[uGdpb] = 0;
585 } // for (Int_t iGdpb = 0; iGdpb < fuNrOfGdpbs; ++iGdpb)
586
589
590 return kTRUE;
591}
592// -------------------------------------------------------------------------
593
594void CbmMcbm2018UnpackerAlgoTof::AddMsComponentToList(size_t component, UShort_t usDetectorId)
595{
597 for (UInt_t uCompIdx = 0; uCompIdx < fvMsComponentsList.size(); ++uCompIdx)
598 if (component == fvMsComponentsList[uCompIdx]) return;
599
601 fvMsComponentsList.push_back(component);
602
603 LOG(info) << "CbmMcbm2018UnpackerAlgoTof::AddMsComponentToList => Component " << component << " with detector ID 0x"
604 << std::hex << usDetectorId << std::dec << " added to list";
605}
606// -------------------------------------------------------------------------
607
608Bool_t CbmMcbm2018UnpackerAlgoTof::ProcessTs(const fles::Timeslice& ts)
609{
610 fulCurrentTsIdx = ts.index();
611 fdTsStartTime = static_cast<Double_t>(ts.descriptor(0, 0).idx);
612 LOG(debug) << "ProcessTs " << fulCurrentTsIdx;
613
615 if (0 == fulCurrentTsIdx) { return kTRUE; } // if( 0 == fulCurrentTsIdx )
616
618 if (-1.0 == fdTsCoreSizeInNs) {
619 fuNbCoreMsPerTs = ts.num_core_microslices();
620 fuNbOverMsPerTs = ts.num_microslices(0) - ts.num_core_microslices();
623 LOG(info) << "Timeslice parameters: each TS has " << fuNbCoreMsPerTs << " Core MS and " << fuNbOverMsPerTs
624 << " Overlap MS, for a core duration of " << fdTsCoreSizeInNs << " ns and a full duration of "
625 << fdTsFullSizeInNs << " ns";
626
630 LOG(info) << "In each TS " << fuNbMsLoop << " MS will be looped over";
631 } // if( -1.0 == fdTsCoreSizeInNs )
632
635 // LOG(info) << Form( "TS %5d Start %12f Stop %12f", fulCurrentTsIdx, fdTsStartTime, fdTsStopTimeCore );
636
638 for (fuMsIndex = 0; fuMsIndex < fuNbMsLoop; fuMsIndex++) {
640 for (UInt_t uMsCompIdx = 0; uMsCompIdx < fvMsComponentsList.size(); ++uMsCompIdx) {
641 UInt_t uMsComp = fvMsComponentsList[uMsCompIdx];
642
643 if (kFALSE == ProcessMs(ts, uMsComp, fuMsIndex)) {
644 LOG(error) << "Failed to process ts " << fulCurrentTsIdx << " MS " << fuMsIndex << " for component " << uMsComp;
645 return kFALSE;
646 } // if( kFALSE == ProcessMs( ts, uMsCompIdx, fuMsIndex ) )
647 } // for( UInt_t uMsCompIdx = 0; uMsCompIdx < fvMsComponentsList.size(); ++uMsCompIdx )
648 } // for( fuMsIndex = 0; fuMsIndex < uNbMsLoop; fuMsIndex ++ )
649
651 if (fbMonitorMode) {
652 if (kFALSE == FillHistograms()) {
653 LOG(error) << "Failed to fill histos in ts " << fulCurrentTsIdx;
654 return kFALSE;
655 } // if( kFALSE == FillHistograms() )
656
658 fhVectorCapacity->Fill(fulCurrentTsIdx, fDigiVect.capacity());
659 } // if( fbMonitorMode )
660
661 if (fuTsMaxVectorSize < fDigiVect.size()) {
663 fDigiVect.shrink_to_fit();
665 } // if( fuTsMaxVectorSize < fDigiVect.size() )
666 /*
667 if(!bEnableOut) {
668 LOG(debug) << "Ts "<< fulCurrentTsIdx << " removed ";
669 fiMsBad++;
670 fDigiVect.clear();
671 }else {
672 LOG(debug) << "Ts "<< fulCurrentTsIdx << " accepted ";
673 fiMsGood++;
674 }
675*/
677 sort(fDigiVect.begin(), fDigiVect.end(),
678 [](const CbmTofDigi& a, const CbmTofDigi& b) -> bool { return a.GetTime() < b.GetTime(); });
679
680 return kTRUE;
681}
682
683Bool_t CbmMcbm2018UnpackerAlgoTof::ProcessMs(const fles::Timeslice& ts, size_t uMsCompIdx, size_t uMsIdx)
684{
685 auto msDescriptor = ts.descriptor(uMsCompIdx, uMsIdx);
686 fuCurrentEquipmentId = msDescriptor.eq_id;
687 const uint8_t* msContent = reinterpret_cast<const uint8_t*>(ts.content(uMsCompIdx, uMsIdx));
688
689 uint32_t uSize = msDescriptor.size;
690 fulCurrentMsIdx = msDescriptor.idx;
691 // Double_t dMsTime = (1e-9) * static_cast<double>(fulCurrentMsIdx);
692 LOG(debug) << "Microslice: " << fulCurrentMsIdx << " from EqId " << std::hex << fuCurrentEquipmentId << std::dec
693 << " has size: " << uSize;
694
695 if (0 == fvbMaskedComponents.size()) fvbMaskedComponents.resize(ts.num_components(), kFALSE);
696
697 fuCurrDpbId = static_cast<uint32_t>(fuCurrentEquipmentId & 0xFFFF);
698 // fuCurrDpbIdx = fDpbIdIndexMap[ fuCurrDpbId ];
699
701 auto it = fGdpbIdIndexMap.find(fuCurrDpbId);
702 if (it == fGdpbIdIndexMap.end()) {
703 if (kFALSE == fvbMaskedComponents[uMsCompIdx]) {
704 LOG(info) << "---------------------------------------------------------------";
705 LOG(info) << FormatMsHeaderPrintout(msDescriptor);
706 LOG(warning) << "Could not find the gDPB index for AFCK id 0x" << std::hex << fuCurrDpbId << std::dec
707 << " in timeslice " << fulCurrentTsIdx << " in microslice " << uMsIdx << " component " << uMsCompIdx
708 << "\n"
709 << "If valid this index has to be added in the TOF "
710 "parameter file in the DbpIdArray field";
711 fvbMaskedComponents[uMsCompIdx] = kTRUE;
712 } // if( kFALSE == fvbMaskedComponents[ uMsComp ] )
713 else
714 return kTRUE;
715
718
719 return kFALSE;
720 } // if( it == fGdpbIdIndexMap.end() )
721 else
723
724 fuCurrentMsSysId = static_cast<unsigned int>(msDescriptor.sys_id);
725
726 // If not integer number of message in input buffer, print warning/error
727 if (0 != (uSize % sizeof(gdpbv100::Message)))
728 LOG(error) << "The input microslice buffer does NOT "
729 << "contain only complete gDPB messages!";
730
731 // Compute the number of complete messages in the input microslice buffer
732 uint32_t uNbMessages = (uSize - (uSize % sizeof(gdpbv100::Message))) / sizeof(gdpbv100::Message);
733
734 // Prepare variables for the loop on contents
735 Int_t messageType = -111;
736 fbEpochFoundInThisMs = kFALSE;
737 const uint64_t* pInBuff = reinterpret_cast<const uint64_t*>(msContent); // for epoch cycle
738 const gdpbv100::Message* pMess = reinterpret_cast<const gdpbv100::Message*>(pInBuff);
739 for (uint32_t uIdx = 0; uIdx < uNbMessages; uIdx++) {
741 if (0 == uIdx) {
742 ProcessEpochCycle(pInBuff[uIdx]);
743 continue;
744 } // if( 0 == uIdx )
745
747 messageType = pMess[uIdx].getMessageType();
748
749 fuGet4Id = fUnpackPar->ElinkIdxToGet4Idx(pMess[uIdx].getGdpbGenChipId());
750 if (fuDiamondDpbIdx == fuCurrDpbIdx || 0x90 == fuCurrentMsSysId) fuGet4Id = pMess[uIdx].getGdpbGenChipId();
752
753 if (fuNrOfGet4PerGdpb <= fuGet4Id && !pMess[uIdx].isStarTrigger() && (gdpbv100::kuChipIdMergedEpoch != fuGet4Id))
754 LOG(warning) << "Message with Get4 ID too high: " << fuGet4Id << " VS " << fuNrOfGet4PerGdpb << " for GdpbIdx "
755 << fuCurrDpbIdx << " set in parameters.";
756
757 /*
758 if( 1 == uIdx && gdpbv100::MSG_EPOCH != messageType )
759 LOG(warning) << " in timeslice " << fulCurrentTsIdx
760 << " in microslice " << fuMsIndex
761 << " component " << uMsCompIdx
762 << " first message is not an epoch: type " << messageType;
763
764 if( uNbMessages - 1 == uIdx && gdpbv100::MSG_EPOCH != messageType )
765 LOG(warning) << " in timeslice " << fulCurrentTsIdx
766 << " in microslice " << fuMsIndex
767 << " component " << uMsCompIdx
768 << " last message is not an epoch: type " << messageType;
769*/
776 if (1 == uIdx && gdpbv100::MSG_EPOCH != messageType) {
777 LOG(debug) << " CbmMcbm2018UnpackerAlgoTof ==> In timeslice " << fulCurrentTsIdx << " in microslice " << fuMsIndex
778 << " component " << uMsCompIdx << " first message is not an epoch: type " << messageType
779 << " -> It will be ignored! ";
780 continue;
781 } // if( 1 == uIdx && gdpbv100::MSG_EPOCH != messageType )
782
783 switch (messageType) {
784 case gdpbv100::MSG_HIT: {
785 if (pMess[uIdx].getGdpbHitIs24b()) {
786 LOG(error) << "This event builder does not support 24b hit message!!!"
787 << " Message " << uIdx << "/" << uNbMessages << " 0x"
788 << FormatHexPrintout(pMess[uIdx].getData(), '0', 16);
789 continue;
790 } // if( getGdpbHitIs24b() )
791 else {
792 fvvmEpSupprBuffer[fuCurrDpbIdx].push_back(pMess[uIdx]);
793 } // else of if( getGdpbHitIs24b() )
794 break;
795 } // case gdpbv100::MSG_HIT:
796 case gdpbv100::MSG_EPOCH: {
798 fbEpochFoundInThisMs = kTRUE;
799 ProcessEpoch(pMess[uIdx], uIdx);
800 } // if this epoch message is a merged one valid for all chips
801 else {
803 LOG(debug2) << "This event builder does not support unmerged epoch "
804 "messages!!!.";
805 continue;
806 } // if single chip epoch message
807 break;
808 } // case gdpbv100::MSG_EPOCH:
809 case gdpbv100::MSG_SLOWC: {
810 fvvmEpSupprBuffer[fuCurrDpbIdx].push_back(pMess[uIdx]);
811 break;
812 } // case gdpbv100::MSG_SLOWC:
813 case gdpbv100::MSG_SYST: {
814 fvvmEpSupprBuffer[fuCurrDpbIdx].push_back(pMess[uIdx]);
815 break;
816 } // case gdpbv100::MSG_SYST:
821 ProcessStarTrigger(pMess[uIdx]);
822
825 /*
826 if( gdpbv100::MSG_STAR_TRI_A == messageType )
827 {
828 } // if( gdpbv100::MSG_STAR_TRI_A == messageType )
829*/
830 break;
831 } // case gdpbv100::MSG_STAR_TRI_A-D
832 default:
833 LOG(error) << "Message type " << std::hex << std::setw(2) << static_cast<uint16_t>(messageType)
834 << " not included in Get4 unpacker.";
835 } // switch( mess.getMessageType() )
836 } // for (uint32_t uIdx = 0; uIdx < uNbMessages; uIdx ++)
837
843 if (kTRUE == fbEpochFoundInThisMs) {
846 }
847 else {
849 }
850
851 return kTRUE;
852}
853
854// -------------------------------------------------------------------------
855void CbmMcbm2018UnpackerAlgoTof::ProcessEpochCycle(const uint64_t& ulCycleData)
856{
857 ULong64_t ulEpochCycleVal = ulCycleData & gdpbv100::kulEpochCycleFieldSz;
858
859 if (!(ulEpochCycleVal == fvulCurrentEpochCycle[fuCurrDpbIdx]
860 || ulEpochCycleVal == fvulCurrentEpochCycle[fuCurrDpbIdx] + 1)
861 && 0 < fulCurrentMsIdx) {
862 LOG(warning) << "CbmMcbm2018UnpackerAlgoTof::ProcessEpochCycle => "
863 << " Missmatch in epoch cycles detected for Gdpb " << fuCurrDpbIdx
864 << ", probably fake cycles due to epoch index corruption! "
865 << Form(" Current cycle 0x%09llX New cycle 0x%09llX", fvulCurrentEpochCycle[fuCurrDpbIdx],
866 ulEpochCycleVal);
867 } // if epoch cycle did not stay constant or increase by exactly 1, except if first MS of the TS
868 if (ulEpochCycleVal != fvulCurrentEpochCycle[fuCurrDpbIdx]) {
869 LOG(info) << "CbmMcbm2018UnpackerAlgoTof::ProcessEpochCycle => "
870 << " New epoch cycle for Gdpb " << fuCurrDpbIdx
871 << Form(": Current cycle 0x%09llX New cycle 0x%09llX", fvulCurrentEpochCycle[fuCurrDpbIdx],
872 ulEpochCycleVal);
873 } // if( ulEpochCycleVal != fvulCurrentEpochCycle[fuCurrDpbIdx] )
874 fvulCurrentEpochCycle[fuCurrDpbIdx] = ulEpochCycleVal;
875
876 return;
877}
878void CbmMcbm2018UnpackerAlgoTof::ProcessEpoch(const gdpbv100::Message& mess, uint32_t /*uMesgIdx*/)
879{
880 ULong64_t ulEpochNr = mess.getGdpbEpEpochNb();
881 /*
882 if( ( 6 == fulCurrentTsIdx && 2 == fuCurrDpbIdx ) ||
883 ( 57 == fulCurrentTsIdx && 0 == fuCurrDpbIdx ) )
884 LOG(info) << Form( "Gdpb %d TS %3llu MS %3u Msg %4u Ep %08llx",
885 fuCurrDpbIdx, fulCurrentTsIdx, fuMsIndex, uMesgIdx, ulEpochNr );
886
887 if( !( ulEpochNr == fvulCurrentEpoch[ fuCurrDpbIdx ] + 1 ||
888// ulEpochNr == fvulCurrentEpoch[ fuCurrDpbIdx ] || // For the fake "closing epoch"
889 ( 0 == ulEpochNr && gdpbv100::kuEpochCounterSz == fvulCurrentEpoch[ fuCurrDpbIdx ] )
890 )
891 )
892 {
893 Int_t iDiff = ulEpochNr;
894 iDiff -= fvulCurrentEpoch[ fuCurrDpbIdx ];
895 LOG(info) << "CbmMcbm2018UnpackerAlgoTof::ProcessEpoch => "
896 << " Non consecutive epochs for Gdpb " << fuCurrDpbIdx
897 << " in TS " << fulCurrentTsIdx
898 << " MS " << fuMsIndex
899 << " Msg " << uMesgIdx
900 << Form( " : old Ep %08llx new Ep %08llx Diff %d ", fvulCurrentEpoch[ fuCurrDpbIdx ], ulEpochNr, iDiff )
901 << std::endl;
902 } // if epoch neither +1 nor equal nor overflow
903*/
904 /*
907 if( 0 < fvulCurrentEpoch[ fuCurrDpbIdx ] && ulEpochNr < fvulCurrentEpoch[ fuCurrDpbIdx ] &&
908 1 < uMesgIdx )
909 {
910 LOG(info) << "CbmMcbm2018UnpackerAlgoTof::ProcessEpoch => "
911 << " New epoch cycle for Gdpb " << fuCurrDpbIdx
912 << Form( ": Current cycle 0x%09llX New cycle 0x%09llX", fvulCurrentEpochCycle[fuCurrDpbIdx], fvulCurrentEpochCycle[fuCurrDpbIdx] + 1 )
913 << Form( "(old Ep %08llx new Ep %08llx)", fvulCurrentEpoch[ fuCurrDpbIdx ], ulEpochNr )
914 << std::endl;
915 fvulCurrentEpochCycle[ fuCurrDpbIdx ]++;
916 } // if( 0 < fvulCurrentEpoch[ fuCurrDpbIdx ] && ulEpochNr < fvulCurrentEpoch[ fuCurrDpbIdx ] && 1 < uMesgIdx)
917*/
918 fvulCurrentEpoch[fuCurrDpbIdx] = ulEpochNr;
921
922 /*
925 if( 0 < ulEpochNr )
926 mess.setGdpbEpEpochNb( ulEpochNr - 1 );
927 else mess.setGdpbEpEpochNb( gdpbv100::kuEpochCounterSz );
928*/
931}
933{
935 ULong64_t ulEpochNr = (fvulCurrentEpoch[fuCurrDpbIdx] + 1) % gdpbv100::kuEpochCounterSz;
936
938 LOG(info) << "CbmMcbm2018UnpackerAlgoTof::ProcessEndOfMsEpoch => "
939 << " New epoch cycle for Gdpb " << fuCurrDpbIdx
940 << Form(": Current cycle 0x%09llX New cycle 0x%09llX", fvulCurrentEpochCycle[fuCurrDpbIdx],
942 << Form("(old Ep %08llx new Ep %08llx)", fvulCurrentEpoch[fuCurrDpbIdx], ulEpochNr);
944 } // if( 0 < fvulCurrentEpoch[ fuCurrDpbIdx ] && ulEpochNr < fvulCurrentEpoch[ fuCurrDpbIdx ] )
945 /*
946 fvulCurrentEpoch[ fuCurrDpbIdx ] = ulEpochNr;
947 fvulCurrentEpochFull[ fuCurrDpbIdx ] = ulEpochNr + ( gdpbv100::kuEpochCounterSz + 1 ) * fvulCurrentEpochCycle[ fuCurrDpbIdx ];
948*/
951
954}
956{
957 Int_t iMsgIndex = mess.getStarTrigMsgIndex();
958
959 switch (iMsgIndex) {
960 case 0: fvulGdpbTsMsb[fuCurrDpbIdx] = mess.getGdpbTsMsbStarA(); break;
961 case 1:
964 break;
965 case 2: fvulStarTsMid[fuCurrDpbIdx] = mess.getStarTsMidStarC(); break;
966 case 3: {
967 ULong64_t ulNewGdpbTsFull = (fvulGdpbTsMsb[fuCurrDpbIdx] << 24) + (fvulGdpbTsLsb[fuCurrDpbIdx]);
968 ULong64_t ulNewStarTsFull =
970 UInt_t uNewToken = mess.getStarTokenStarD();
971 UInt_t uNewDaqCmd = mess.getStarDaqCmdStarD();
972 UInt_t uNewTrigCmd = mess.getStarTrigCmdStarD();
973
974 if ((uNewToken == fvuStarTokenLast[fuCurrDpbIdx]) && (ulNewGdpbTsFull == fvulGdpbTsFullLast[fuCurrDpbIdx])
975 && (ulNewStarTsFull == fvulStarTsFullLast[fuCurrDpbIdx]) && (uNewDaqCmd == fvuStarDaqCmdLast[fuCurrDpbIdx])
976 && (uNewTrigCmd == fvuStarTrigCmdLast[fuCurrDpbIdx])) {
977 UInt_t uTrigWord = ((fvuStarTrigCmdLast[fuCurrDpbIdx] & 0x00F) << 16)
978 + ((fvuStarDaqCmdLast[fuCurrDpbIdx] & 0x00F) << 12)
979 + ((fvuStarTokenLast[fuCurrDpbIdx] & 0xFFF));
980 LOG(warning) << "Possible error: identical STAR tokens found twice in "
981 "a row => ignore 2nd! "
982 << " TS " << fulCurrentTsIdx << " gDBB #" << fuCurrDpbIdx << " "
983 << Form("token = %5u ", fvuStarTokenLast[fuCurrDpbIdx])
984 << Form("gDPB ts = %12llu ", fvulGdpbTsFullLast[fuCurrDpbIdx])
985 << Form("STAR ts = %12llu ", fvulStarTsFullLast[fuCurrDpbIdx])
986 << Form("DAQ cmd = %2u ", fvuStarDaqCmdLast[fuCurrDpbIdx])
987 << Form("TRG cmd = %2u ", fvuStarTrigCmdLast[fuCurrDpbIdx]) << Form("TRG Wrd = %5x ", uTrigWord);
988 return;
989 } // if exactly same message repeated
990
991 // GDPB TS counter reset detection
992 if (ulNewGdpbTsFull < fvulGdpbTsFullLast[fuCurrDpbIdx])
993 LOG(debug) << "Probable reset of the GDPB TS: old = " << Form("%16llu", fvulGdpbTsFullLast[fuCurrDpbIdx])
994 << " new = " << Form("%16llu", ulNewGdpbTsFull) << " Diff = -"
995 << Form("%8llu", fvulGdpbTsFullLast[fuCurrDpbIdx] - ulNewGdpbTsFull) << " GDPB #"
996 << Form("%2u", fuCurrDpbIdx);
997
998 // STAR TS counter reset detection
999 if (ulNewStarTsFull < fvulStarTsFullLast[fuCurrDpbIdx])
1000 LOG(debug) << "Probable reset of the STAR TS: old = " << Form("%16llu", fvulStarTsFullLast[fuCurrDpbIdx])
1001 << " new = " << Form("%16llu", ulNewStarTsFull) << " Diff = -"
1002 << Form("%8llu", fvulStarTsFullLast[fuCurrDpbIdx] - ulNewStarTsFull) << " GDPB #"
1003 << Form("%2u", fuCurrDpbIdx);
1004
1007 fvulGdpbTsFullLast[fuCurrDpbIdx] = ulNewGdpbTsFull;
1008 fvulStarTsFullLast[fuCurrDpbIdx] = ulNewStarTsFull;
1009 fvuStarTokenLast[fuCurrDpbIdx] = uNewToken;
1010 fvuStarDaqCmdLast[fuCurrDpbIdx] = uNewDaqCmd;
1011 fvuStarTrigCmdLast[fuCurrDpbIdx] = uNewTrigCmd;
1012 } // if( fuCurrentMs < fuNbCoreMsPerTs )
1013
1015 /*
1016 Double_t dTot = 1.;
1017 Double_t dTime = fulGdpbTsFullLast * 6.25;
1018 if( 0. == fdFirstDigiTimeDif && 0. != fdLastDigiTime )
1019 {
1020 fdFirstDigiTimeDif = dTime - fdLastDigiTime;
1021 LOG(info) << "Reference fake digi time shift initialized to " << fdFirstDigiTimeDif;
1022 } // if( 0. == fdFirstDigiTimeDif && 0. != fdLastDigiTime )
1023
1024 // dTime -= fdFirstDigiTimeDif;
1025
1026 LOG(debug) << "Insert fake digi with time " << dTime << ", Tot " << dTot;
1027 fhRawTRefDig0->Fill( dTime - fdLastDigiTime);
1028 fhRawTRefDig1->Fill( dTime - fdLastDigiTime);
1029
1030 fDigi = new CbmTofDigiExp(0x00005006, dTime, dTot); // fake start counter signal
1031 fBuffer->InsertData<CbmTofDigi>(fDigi);
1032*/
1033 break;
1034 } // case 3
1035 default: LOG(error) << "Unknown Star Trigger messageindex: " << iMsgIndex;
1036 } // switch( iMsgIndex )
1037}
1038// -------------------------------------------------------------------------
1040{
1041 Int_t iBufferSize = fvvmEpSupprBuffer[fuCurrDpbIdx].size();
1042
1043 if (0 == iBufferSize) return;
1044
1045 LOG(debug) << "Now processing stored messages for for gDPB " << fuCurrDpbIdx << " with epoch number "
1047
1050 // std::stable_sort( fvvmEpSupprBuffer[ fuCurrDpbIdx ].begin(), fvvmEpSupprBuffer[ fuCurrDpbIdx ].end() );
1051
1053 ULong64_t ulCurEpochGdpbGet4 = fvulCurrentEpochFull[fuCurrDpbIdx];
1054
1056 if (0 == ulCurEpochGdpbGet4) return;
1057
1059 ulCurEpochGdpbGet4--;
1060
1061 Int_t messageType = -111;
1062 for (Int_t iMsgIdx = 0; iMsgIdx < iBufferSize; iMsgIdx++) {
1063 messageType = fvvmEpSupprBuffer[fuCurrDpbIdx][iMsgIdx].getMessageType();
1064
1065 fuGet4Id = fUnpackPar->ElinkIdxToGet4Idx(fvvmEpSupprBuffer[fuCurrDpbIdx][iMsgIdx].getGdpbGenChipId());
1067 fuGet4Id = fvvmEpSupprBuffer[fuCurrDpbIdx][iMsgIdx].getGdpbGenChipId();
1069
1071 gdpbv100::FullMessage fullMess(fvvmEpSupprBuffer[fuCurrDpbIdx][iMsgIdx], ulCurEpochGdpbGet4);
1072
1074 switch (messageType) {
1075 case gdpbv100::MSG_HIT: {
1076 ProcessHit(fullMess);
1077 break;
1078 } // case gdpbv100::MSG_HIT:
1079 case gdpbv100::MSG_SLOWC: {
1080 ProcessSlCtrl(fullMess);
1081 break;
1082 } // case gdpbv100::MSG_SLOWC:
1083 case gdpbv100::MSG_SYST: {
1084 ProcessSysMess(fullMess);
1085 break;
1086 } // case gdpbv100::MSG_SYST:
1093 break;
1094 default:
1095 LOG(error) << "Message type " << std::hex << std::setw(2) << static_cast<uint16_t>(messageType)
1096 << " not included in Get4 unpacker.";
1097 } // switch( mess.getMessageType() )
1098 } // for( Int_t iMsgIdx = 0; iMsgIdx < iBufferSize; iMsgIdx++ )
1099
1101}
1102
1104// -------------------------------------------------------------------------
1106{
1107
1108 if (messlast == mess) return;
1109 messlast = mess;
1110
1111 UInt_t uChannel = mess.getGdpbHitChanId();
1112 UInt_t uTot = mess.getGdpbHit32Tot();
1113
1114 // In 32b mode the coarse counter is already computed back to 112 FTS bins
1115 // => need to hide its contribution from the Finetime
1116 // => FTS = Fullt TS modulo 112
1117 // UInt_t uFts = mess.getGdpbHitFullTs() % 112;
1118
1119 UInt_t uChannelNr = fuGet4Id * fuNrOfChannelsPerGet4 + uChannel;
1120 UInt_t uChannelNrInFee = (fuGet4Id % fuNrOfGet4PerFee) * fuNrOfChannelsPerGet4 + uChannel;
1121 UInt_t uFeeNr = (fuGet4Id / fuNrOfGet4PerFee);
1122 UInt_t uFeeNrInSys = fuCurrDpbIdx * fuNrOfFeePerGdpb + uFeeNr;
1123 UInt_t uRemappedChannelNr = uFeeNr * fuNrOfChannelsPerFee + fUnpackPar->Get4ChanToPadiChan(uChannelNrInFee);
1124 // UInt_t uGbtxNr = (uFeeNr / fUnpackPar->GetNrOfFeePerGbtx());
1125 // UInt_t uFeeInGbtx = (uFeeNr % fUnpackPar->GetNrOfFeePerGbtx());
1126 // UInt_t uGbtxNrInSys = fuCurrDpbIdx * fUnpackPar->GetNrOfGbtxPerGdpb() + uGbtxNr;
1127
1128 // UInt_t uChanInSyst = fuCurrDpbIdx * fuNrOfChannelsPerGdpb + uChannelNr;
1129 UInt_t uRemappedChannelNrInSys = fuCurrDpbIdx * fuNrOfChannelsPerGdpb + uFeeNr * fuNrOfChannelsPerFee
1130 + fUnpackPar->Get4ChanToPadiChan(uChannelNrInFee);
1132 if (fuDiamondDpbIdx == fuCurrDpbIdx || 0x90 == fuCurrentMsSysId) {
1133 uRemappedChannelNr = uChannelNr;
1134 uRemappedChannelNrInSys = fuCurrDpbIdx * fUnpackPar->GetNrOfChannelsPerGdpb() + uChannelNr;
1135 } // if( fuDiamondDpbIdx == fuCurrDpbIdx || 0x90 == fuCurrentMsSysId )
1136
1137 // ULong_t ulHitTime = mess.getMsgFullTime( mess.getExtendedEpoch() );
1138 Double_t dHitTime = mess.GetFullTimeNs();
1139 Double_t dHitTot = uTot; // in bins
1140
1141 // Histograms filling
1142 if (kTRUE == fbMonitorMode) {
1143 fhRawTotCh[fuCurrDpbIdx]->Fill(uRemappedChannelNr, dHitTot);
1144 fhChCount[fuCurrDpbIdx]->Fill(uChannelNr);
1145 fhChCountRemap[fuCurrDpbIdx]->Fill(uRemappedChannelNr);
1146 } // if( kTRUE == fbMonitorMode )
1147
1148 /*
1149 if( fUnpackPar->GetNumberOfChannels() < uChanInSyst )
1150 {
1151 LOG(fatal) << "Invalid mapping index " << uChanInSyst
1152 << " VS " << fUnpackPar->GetNumberOfChannels()
1153 <<", from " << fuCurrDpbIdx
1154 <<", " << fuGet4Id
1155 <<", " << uChannel;
1156 return;
1157 } // if( fUnpackPar->GetNumberOfChannels() < uChanUId )
1158*/
1159 if (fviRpcChUId.size() < uRemappedChannelNrInSys) {
1160 LOG(fatal) << "Invalid mapping index " << uRemappedChannelNrInSys << " VS " << fviRpcChUId.size()
1161 << ", from GdpbNr " << fuCurrDpbIdx << ", Get4 " << fuGet4Id << ", Ch " << uChannel << ", ChNr "
1162 << uChannelNr << ", ChNrIF " << uChannelNrInFee << ", FiS " << uFeeNrInSys;
1163 return;
1164 } // if( fviRpcChUId.size() < uRemappedChannelNrInSys )
1165
1166 UInt_t uChanUId = fviRpcChUId[uRemappedChannelNrInSys];
1167 /*
1168 if( 5 == fviRpcType[uGbtxNrInSys] )
1169 LOG(info) << "Bmon mapping index " << uRemappedChannelNrInSys
1170 << " UID " << std::hex << std::setw(8) << uChanUId << std::dec
1171 << ", from GdpbNr " << fuCurrDpbIdx
1172 << ", Get4 " << fuGet4Id
1173 << ", Ch " << uChannel
1174 << ", ChNr " << uChannelNr
1175 << ", ChNrIF " << uChannelNrInFee
1176 << ", FiS " << uFeeNrInSys
1177 << ", GBTx " << uGbtxNrInSys;
1178*/
1179
1180 if (0 == uChanUId) {
1181 if (0 < NMappingWarnings--)
1182 LOG(warning) << "Unused data item at " << uRemappedChannelNrInSys << ", from GdpbNr " << fuCurrDpbIdx << ", Get4 "
1183 << fuGet4Id << ", Ch " << uChannel << ", ChNr " << uChannelNr << ", ChNrIF " << uChannelNrInFee
1184 << ", FiS " << uFeeNrInSys;
1185 return; // Hit not mapped to digi
1186 }
1187
1189 if (0x90 != fuCurrentMsSysId) dHitTime -= fdTimeOffsetNs;
1190
1191 LOG(debug) << Form("Insert 0x%08x digi with time ", uChanUId) << dHitTime
1192 << Form(", Tot %4.0f", dHitTot)
1193 // << " into buffer with " << fBuffer->GetSize() << " data from "
1194 // << Form("%11.1f to %11.1f ", fBuffer->GetTimeFirst(), fBuffer->GetTimeLast())
1195
1196 << " at epoch " << mess.getExtendedEpoch();
1197
1198 // CbmTofDigi digi( uChanUId, dHitTime, dHitTot );
1199 // fDigiVect.emplace_back( digi );
1200 fDigiVect.emplace_back(uChanUId, dHitTime, dHitTot);
1201}
1202// -------------------------------------------------------------------------
1204// -------------------------------------------------------------------------
1206{
1207 switch (mess.getGdpbSysSubType()) {
1209 ProcessError(mess);
1210 break;
1211 } // case gdpbv100::SYSMSG_GET4_EVENT
1213 LOG(debug) << "Unknown GET4 message, data: " << std::hex << std::setw(8) << mess.getGdpbSysUnkwData() << std::dec
1214 << " Full message: " << std::hex << std::setw(16) << mess.getData() << std::dec;
1215 break;
1216 } // case gdpbv100::SYS_GDPB_UNKWN:
1218 if (mess.getGdpbSysFwErrResync())
1219 LOG(info) << Form("GET4 Resynchronization: Get4:0x%04x ", mess.getGdpbGenChipId()) << fuCurrDpbIdx;
1220 else
1221 LOG(info) << "GET4 synchronization pulse missing in gDPB " << fuCurrDpbIdx;
1222 break;
1223 } // case gdpbv100::SYS_GET4_SYNC_MISS:
1224 case gdpbv100::SYS_PATTERN: {
1225 ProcessPattern(mess);
1226 break;
1227 } // case gdpbv100::SYS_PATTERN:
1228 default: {
1229 LOG(info) << "Crazy system message, subtype " << mess.getGdpbSysSubType();
1230 break;
1231 } // default
1232 } // switch( mess.getGdpbSysSubType() )
1233}
1275
1277{
1278 uint16_t usType = mess.getGdpbSysPattType();
1279 uint16_t usIndex = mess.getGdpbSysPattIndex();
1280 uint32_t uPattern = mess.getGdpbSysPattPattern();
1281
1282 switch (usType) {
1284 LOG(debug) << Form("Missmatch pattern message => Type %u, Index %2d, Pattern 0x%08X", usType, usIndex, uPattern);
1285 /*
1286 if(usIndex==7) {
1287 TString Tok;
1288 if (bEnableOut) Tok="Ena";
1289 else Tok="Dis";
1290 LOG(debug) << Form( "Mismatch pat in TS %llu, MS %llu, Gdpb %u, T %u, Pattern 0x%08X %08X %08X %08X %08X %08X %08X %08X ",
1291 fulCurrentTsIdx, fulCurrentMsIdx, fuCurrDpbIdx, usType,
1292 pat_mess[0], pat_mess[1], pat_mess[2], pat_mess[3], pat_mess[4], pat_mess[5], pat_mess[6], pat_mess[7] )
1293 << Tok;
1294 }
1295 */
1296 break;
1297 } // case gdpbv100::PATT_MISSMATCH:
1298
1299 case gdpbv100::PATT_ENABLE: {
1300 LOG(debug2) << Form("Enable pattern message => Type %d, Index %2d, Pattern 0x%08X", usType, usIndex, uPattern);
1301 /*
1302 for( UInt_t uBit = 0; uBit < 32; ++uBit )
1303 if( ( uPattern >> uBit ) & 0x1 )
1304 {
1305 fhPatternEnable->Fill( 32 * usIndex + uBit, fuCurrDpbIdx );
1306 fvhGdpbPatternEnableEvo[ fuCurrDpbIdx ]->Fill( fulCurrentTsIndex, 32 * usIndex + uBit );
1307 } // if( ( uPattern >> uBit ) & 0x1 )
1308 */
1309 break;
1310 } // case gdpbv100::PATT_ENABLE:
1311
1312 case gdpbv100::PATT_RESYNC: {
1313 LOG(debug) << Form("RESYNC pattern message => Type %d, Index %2d, Pattern 0x%08X", usType, usIndex, uPattern);
1314
1315 break;
1316 } // case gdpbv100::PATT_RESYNC:
1317
1318 default: {
1319 LOG(debug) << "Crazy pattern message, subtype " << usType;
1320 break;
1321 } // default
1322 } // switch( usType )
1323 /*
1324 if(usIndex==7) {
1325 bEnableOut=kTRUE;
1326// for(Int_t iGdpb=0; iGdpb<fUnpackPar->GetNrOfGdpbs(); iGdpb++) {
1327// bEnableOut &= bGdpbOK[iGdpb];
1328// }
1329 }
1330*/
1331 return;
1332}
1333// -------------------------------------------------------------------------
1334
1336{
1337 std::string sFolder = "Tof_Raw_gDPB";
1338
1339 LOG(info) << "create Histos for " << fuNrOfGdpbs << " gDPBs ";
1340
1342 new TH1F(Form("Raw_TDig-EvBmon"), Form("Raw digi time difference to 1st digi ; time [ns]; cts"), 500, 0, 100.);
1343
1345 new TH1F(Form("Raw_TDig-Ref0"), Form("Raw digi time difference to Ref ; time [ns]; cts"), 6000, -10000, 50000);
1346
1347 fhRawTDigRef =
1348 new TH1F(Form("Raw_TDig-Ref"), Form("Raw digi time difference to Ref ; time [ns]; cts"), 6000, -1000, 5000);
1349
1350 fhRawTRefDig0 = new TH1F(Form("Raw_TRef-Dig0"), Form("Raw Ref time difference to last digi ; time [ns]; cts"), 9999,
1351 -50000, 50000);
1352
1354 new TH1F(Form("Raw_TRef-Dig1"), Form("Raw Ref time difference to last digi ; time [ns]; cts"), 9999, -5000, 5000);
1355
1356 fhRawDigiLastDigi = new TH1F(Form("Raw_Digi-LastDigi"),
1357 Form("Raw Digi time difference to last digi ; time [ns]; cts"), 9999, -5000, 5000);
1358
1366
1367 fhRawTotCh.resize(fuNrOfGdpbs);
1368 fhChCount.resize(fuNrOfGdpbs);
1370 fhChanCoinc.resize(fuNrOfGdpbs);
1371 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; uGdpb++) {
1372 fhRawTotCh[uGdpb] = new TH2F(Form("Raw_Tot_gDPB_%02u", uGdpb), Form("Raw TOT gDPB %02u; channel; TOT [bin]", uGdpb),
1373 fuNrOfGet4PerGdpb * 4, 0., fuNrOfGet4PerGdpb * 4, 256, 0., 256.);
1374
1375 fhChCount[uGdpb] =
1376 new TH1I(Form("ChCount_gDPB_%02u", uGdpb), Form("Channel counts gDPB %02u; channel; Hits", uGdpb),
1378
1379 fhChCountRemap[uGdpb] = new TH1I(Form("ChCountRemap_gDPB_%02u", uGdpb),
1380 Form("Remapped channel counts gDPB %02u; Remapped channel; Hits", uGdpb),
1382
1383 fhChanCoinc[uGdpb] =
1384 new TH2F(Form("fhChanCoinc_%02u", uGdpb), Form("Channels Coincidence %02u; Left; Right", uGdpb),
1386
1388 AddHistoToVector(fhRawTotCh[uGdpb], sFolder);
1389 AddHistoToVector(fhChCount[uGdpb], sFolder);
1390 AddHistoToVector(fhChCountRemap[uGdpb], sFolder);
1391 AddHistoToVector(fhChanCoinc[uGdpb], sFolder);
1392 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; uGdpb ++)
1393
1394 fhVectorSize = new TH1I("fhVectorSize", "Size of the vector VS TS index; TS index; Size [bytes]", 10000, 0., 10000.);
1396 new TH1I("fhVectorCapacity", "Size of the vector VS TS index; TS index; Size [bytes]", 10000, 0., 10000.);
1399
1400 /*
1402 Double_t w = 10;
1403 Double_t h = 10;
1404
1406 fcEventBuildStats = new TCanvas( "cEvtBuildStats", "Event building statistics", w, h);
1407 if( kTRUE == fbDebugMonitorMode )
1408 fcEventBuildStats->Divide( 2, 3 );
1409 else fcEventBuildStats->Divide( 2, 2 );
1410
1411 fcEventBuildStats->cd( 1 );
1412 gPad->SetGridx();
1413 gPad->SetGridy();
1414 gPad->SetLogy();
1415 fhEventNbPerTs->Draw();
1416
1417 fcEventBuildStats->cd( 2 );
1418 gPad->SetGridx();
1419 gPad->SetGridy();
1420 gPad->SetLogy();
1421 fhEventSizeDistribution->Draw();
1422
1423 fcEventBuildStats->cd( 3 );
1424 gPad->SetGridx();
1425 gPad->SetGridy();
1426 gPad->SetLogy();
1427 fhEventSizeEvolution->Draw();
1428
1429 fcEventBuildStats->cd( 4 );
1430 gPad->SetGridx();
1431 gPad->SetGridy();
1432 gPad->SetLogy();
1433 fhEventNbEvolution->Draw();
1434
1435 if( kTRUE == fbDebugMonitorMode )
1436 {
1437 fcEventBuildStats->cd( 5 );
1438 gPad->SetGridx();
1439 gPad->SetGridy();
1440 gPad->SetLogy();
1441 fhEventNbDistributionInTs->Draw();
1442
1443 fcEventBuildStats->cd( 6 );
1444 gPad->SetGridx();
1445 gPad->SetGridy();
1446 gPad->SetLogy();
1447 fhEventSizeDistributionInTs->Draw();
1448 } // if( kTRUE == fbDebugMonitorMode )
1449
1450 AddCanvasToVector( fcEventBuildStats, "canvases" );
1451*/
1452 return kTRUE;
1453}
1455{
1456 /*
1457 UInt_t uNbEvents = fvEventsBuffer.size();
1458 fhEventNbPerTs->Fill( uNbEvents );
1459
1460 for( UInt_t uEvent = 0; uEvent < uNbEvents; ++uEvent )
1461 {
1462 UInt_t uEventSize = fvEventsBuffer[ uEvent ].GetEventSize();
1463 Double_t dEventTimeSec = fvEventsBuffer[ uEvent ].GetEventTimeSec();
1464 Double_t dEventTimeMin = dEventTimeSec / 60.0;
1465
1466 fhEventSizeDistribution->Fill( uEventSize );
1467 fhEventSizeEvolution->Fill( dEventTimeMin, uEventSize );
1468 fhEventNbEvolution->Fill( dEventTimeMin );
1469
1470 if( kTRUE == fbDebugMonitorMode )
1471 {
1472 Double_t dEventTimeInTs = ( fvEventsBuffer[ uEvent ].GetTrigger().GetFullGdpbTs() * gdpbv100::kdClockCycleSizeNs
1473 - fdTsStartTime ) / 1000.0;
1474
1475 fhEventNbDistributionInTs->Fill( dEventTimeInTs );
1476 fhEventSizeDistributionInTs->Fill( dEventTimeInTs, uEventSize );
1477 } // if( kTRUE == fbDebugMonitorMode )
1478 } // for( UInt_t uEvent = 0; uEvent < uNbEvents; ++uEvent )
1479*/
1480 return kTRUE;
1481}
1483{
1484 /*
1485 for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1486 {
1487 fvhHitsTimeToTriggerRaw[ uGdpb ]->Reset();
1488 fvhHitsTimeToTriggerSel[ uGdpb ]->Reset();
1489
1490 if( kTRUE == fbDebugMonitorMode )
1491 {
1492 fvhHitsTimeToTriggerSelVsDaq[ uGdpb ]->Reset();
1493 fvhHitsTimeToTriggerSelVsTrig[ uGdpb ]->Reset();
1494 fvhTriggerDt[ uGdpb ]->Reset();
1495 fvhTriggerDistributionInTs[ uGdpb ]->Reset();
1496 fvhTriggerDistributionInMs[ uGdpb ]->Reset();
1497 fvhMessDistributionInMs[ uGdpb ]->Reset();
1498 } // if( kTRUE == fbDebugMonitorMode )
1499 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1500
1502 fhEventNbPerTs->Reset();
1503 fhEventSizeDistribution->Reset();
1504 fhEventSizeEvolution->Reset();
1505 fhEventNbEvolution->Reset();
1506
1507 if( kTRUE == fbDebugMonitorMode )
1508 {
1509 fhEventNbDistributionInTs->Reset();
1510 fhEventSizeDistributionInTs->Reset();
1511 fhRawTriggersStats->Reset();
1512 fhMissingTriggersEvolution->Reset();
1513 } // if( kTRUE == fbDebugMonitorMode )
1514*/
1515 return kTRUE;
1516}
1517// -------------------------------------------------------------------------
std::string FormatHexPrintout(uint64_t ulVal, char cFill, uint uWidth, bool bUppercase)
std::string FormatMsHeaderPrintout(const fles::MicrosliceDescriptor &msDescriptor)
static gdpbv100::FullMessage messlast
static Int_t NMappingWarnings
Int_t GetModuleId(Int_t i)
static constexpr UInt_t GetNrOfFeePerGbtx()
Int_t GetRpcSide(Int_t i)
Int_t GetNrOfGdpbs()
FIXME: replace with method returning the correspondign constants! see Star2019 parameter.
Int_t GetNrOfRpc(Int_t i)
Int_t GetGdpbId(Int_t i)
Int_t Get4ChanToPadiChan(UInt_t uChannelInFee)
Int_t GetRpcType(Int_t i)
static constexpr UInt_t GetNrOfChannelsPerGdpb()
Int_t ElinkIdxToGet4Idx(UInt_t uElink)
static constexpr UInt_t GetNrOfChannelsPerFee()
UInt_t fuMsIndex
Start Time in ns of current MS from its index field in header.
Double_t fdTimeOffsetNs
User settings: Data correction parameters.
void ProcessEpoch(const gdpbv100::Message &mess, uint32_t uMesgIdx)
UInt_t fuNrOfFeePerGdpb
gDPB ID to index map
std::vector< Bool_t > fvbMaskedComponents
Switch ON the filling of a additional set of histograms.
std::vector< ULong64_t > fvulStarTsFullLast
size_t fuCurrentMsSysId
Idx of the current MS in TS (0 to fuTotalMsNb)
Bool_t ProcessTs(const fles::Timeslice &ts)
std::vector< ULong64_t > fvulStarTsMsb
Double_t fdTsStopTimeCore
Time in ns of current TS from the index of the first MS first component.
void ProcessError(const gdpbv100::FullMessage &mess)
UInt_t fuNrOfChannelsPerFee
Number of channels in each GET4.
UInt_t fuGet4Nr
running number (0 to fuNrOfGet4PerGdpb) of the Get4 chip of a unique GDPB for current message
UInt_t fuNrOfGet4PerGdpb
Total number of Get4 chips in the system.
void ProcessSlCtrl(const gdpbv100::FullMessage &mess)
UInt_t fuNrOfGet4
Number of channels in each FEE.
std::vector< ULong64_t > fvulCurrentEpochCycle
Current epoch index, per DPB.
std::map< UInt_t, UInt_t > fGdpbIdIndexMap
Total number of GDPBs in the system.
void ProcessStarTrigger(const gdpbv100::Message &mess)
void ProcessPattern(const gdpbv100::FullMessage &mess)
std::vector< ULong64_t > fvulCurrentEpochFull
Epoch cycle from the Ms Start message and Epoch counter flip.
UInt_t fuGet4Id
Index of the DPB from which the MS currently unpacked is coming.
UInt_t fuNrOfGet4PerFee
Number of FEBs per GDPB.
UInt_t fuCurrDpbId
Current equipment ID, tells from which DPB the current MS is originating.
void AddMsComponentToList(size_t component, UShort_t usDetectorId)
void ProcessEpochCycle(const uint64_t &ulCycleData)
std::vector< ULong64_t > fvulGdpbTsLsb
void ProcessSysMess(const gdpbv100::FullMessage &mess)
UInt_t fuNrOfChannelsPerGet4
Number of GET4s per FEE.
CbmMcbm2018TofPar * fUnpackPar
Settings from parameter file.
UInt_t fuNrOfChannelsPerGdpb
Number of GET4s per GDPB.
std::vector< ULong64_t > fvulGdpbTsMsb
STAR TRIGGER detection.
Bool_t fbEpochFoundInThisMs
[DPB], FIXME: dimension to be removed as not needed with real MS
std::vector< ULong64_t > fvulCurrentEpoch
Data format control: Current time references for each GDPB: merged epoch marker, epoch cycle,...
UInt_t fuCurrDpbIdx
Temp holder until Current equipment ID is properly filled in MS.
std::vector< ULong64_t > fvulStarTsMid
void ProcessHit(const gdpbv100::FullMessage &mess)
Double_t fdTsStartTime
SysId of the current MS in TS (0 to fuTotalMsNb)
Bool_t ProcessMs(const fles::Timeslice &ts, size_t uMsCompIdx, size_t uMsIdx)
std::vector< std::vector< gdpbv100::Message > > fvvmEpSupprBuffer
Buffers.
ULong64_t fulCurrentMsIdx
Idx of the current TS.
std::vector< ULong64_t > fvulGdpbTsFullLast
void AddHistoToVector(TNamed *pointer, std::string sFolder="")
std::vector< size_t > fvMsComponentsList
std::vector< CbmTofDigi > fDigiVect
static uint32_t GetUniqueAddress(uint32_t Sm, uint32_t Rpc, uint32_t Channel, uint32_t Side=0, uint32_t SmType=0, uint32_t RpcType=0)
Data class for expanded digital TOF information.
Definition CbmTofDigi.h:47
double GetFullTimeNs() const
uint64_t getExtendedEpoch() const
uint16_t getStarTrigMsgIndex() const
uint32_t getGdpbEpEpochNb() const
uint16_t getGdpbHit32Tot() const
uint16_t getGdpbSysPattType() const
uint64_t getStarTsMidStarC() const
uint16_t getGdpbSysSubType() const
uint64_t getGdpbTsLsbStarB() const
uint16_t getGdpbSysPattIndex() const
uint32_t getGdpbSysFwErrResync() const
uint64_t getStarTsMsbStarB() const
uint32_t getStarTrigCmdStarD() const
uint32_t getGdpbSysPattPattern() const
uint32_t getStarDaqCmdStarD() const
uint64_t getGdpbTsMsbStarA() const
uint16_t getGdpbGenChipId() const
uint32_t getStarTokenStarD() const
uint64_t getStarTsLsbStarD() const
uint8_t getMessageType() const
Returns the message type. Valid for all message types. 4 bit.
uint32_t getGdpbSysUnkwData() const
uint64_t getData() const
uint16_t getGdpbSysErrData() const
uint16_t getGdpbHitChanId() const
@ GET4_V2X_ERR_LOST_EVT
@ GET4_V2X_ERR_TOT_RANGE
@ GET4_V2X_ERR_FIFO_WRITE
@ GET4_V2X_ERR_DLL_LOCK
@ GET4_V2X_ERR_EPOCH_OVERF
@ GET4_V2X_ERR_EP_CNT_SYNC
@ GET4_V2X_ERR_UNKNOWN
@ GET4_V2X_ERR_READ_INIT
@ GET4_V2X_ERR_DLL_RESET
@ GET4_V2X_ERR_ADD_RIS_EDG
@ GET4_V2X_ERR_TOK_RING_ST
@ GET4_V2X_ERR_TOKEN
@ GET4_V2X_ERR_TOT_OVERWRT
@ GET4_V2X_ERR_READOUT_ERR
@ GET4_V2X_ERR_EVT_DISCARD
@ GET4_V2X_ERR_UNPAIR_FALL
@ GET4_V2X_ERR_CHAN_STATE
@ GET4_V2X_ERR_SYNC
@ GET4_V2X_ERR_SEQUENCE_ER
const uint32_t kuChipIdMergedEpoch
const uint32_t kuEpochCounterSz
const uint64_t kulEpochCycleFieldSz
@ SYS_GET4_SYNC_MISS
Hash for CbmL1LinkKey.