CbmRoot
Loading...
Searching...
No Matches
CbmStar2019EventBuilderEtofAlgo.cxx
Go to the documentation of this file.
1/* Copyright (C) 2018-2021 Facility for Antiproton and Ion Research in Europe, Darmstadt
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Pierre-Alain Loizeau [committer] */
4
5// -----------------------------------------------------------------------------
6// ----- -----
7// ----- CbmStar2019EventBuilderEtofAlgo -----
8// ----- Created 03.11.2018 by P.-A. Loizeau -----
9// ----- -----
10// -----------------------------------------------------------------------------
11
13
14#include "CbmHistManager.h"
15#include "CbmStar2019TofPar.h"
16
17#include "FairRootManager.h"
18#include "FairRun.h"
19#include "FairRunOnline.h"
20#include "FairRuntimeDb.h"
21#include <Logger.h>
22
23#include "TCanvas.h"
24#include "TH1.h"
25#include "TH2.h"
26#include "TList.h"
27#include "TProfile.h"
28#include "TROOT.h"
29#include "TString.h"
30
31#include <fstream>
32#include <iomanip>
33#include <iostream>
34
35#include <stdint.h>
36
37// -------------------------------------------------------------------------
40 ,
42 fbMonitorMode(kFALSE)
43 , fbDebugMonitorMode(kFALSE)
44 , fbStoreLostEventMsg(kFALSE)
45 , fbAddStatusToEvent(kTRUE)
46 , fUnpackPar(nullptr)
47 , fuNrOfGdpbs(0)
48 , fGdpbIdIndexMap()
49 , fuNrOfFeePerGdpb(0)
50 , fuNrOfGet4PerFee(0)
51 , fuNrOfChannelsPerGet4(0)
52 , fuNrOfChannelsPerFee(0)
53 , fuNrOfGet4(0)
54 , fuNrOfGet4PerGdpb(0)
55 , fuNrOfChannelsPerGdpb(0)
56 , fuNrOfGbtx(0)
57 , fuNrOfModules(0)
58 , fviNrOfRpc()
59 , fviRpcType()
60 , fviRpcSide()
61 , fviModuleId()
62 , fdAllowedTriggersSpread(-1.0)
63 , fdStarTriggerDeadtime()
64 , fdStarTriggerDelay()
65 , fdStarTriggerWinSize()
66 , fuMinTotPulser(90)
67 , fuMaxTotPulser(110)
68 , fulCurrentTsIndex(0)
69 , fdTsStartTime(-1.0)
70 , fdTsStopTimeCore(-1.0)
71 , fuCurrentMs(0)
72 , fdMsTime(-1.0)
73 , fuMsIndex(0)
74 , fuGdpbId(0)
75 , fuGdpbNr(0)
76 , fuGet4Id(0)
77 , fuGet4Nr(0)
78 , fiEquipmentId(0)
79 , fvulCurrentEpoch()
80 , fvulCurrentEpochCycle()
81 , fvulCurrentEpochFull()
82 , fvvmEpSupprBuffer()
83 , fvvBufferMajorAsicErrors()
84 , fvvBufferMessages()
85 , fvvBufferTriggers()
86 , fbTriggerFoundA(kFALSE)
87 , fbTriggerFoundB(kFALSE)
88 , fbTriggerFoundC(kFALSE)
89 , fvulGdpbTsMsb()
90 , fvulGdpbTsLsb()
91 , fvulStarTsMsb()
92 , fvulStarTsMid()
93 , fvulGdpbTsFullLast()
94 , fvulStarTsFullLast()
95 , fvuStarTokenLast()
96 , fvuStarDaqCmdLast()
97 , fvuStarTrigCmdLast()
98 , fvulGdpbTsFullLastCore()
99 , fvulStarTsFullLastCore()
100 , fvuStarTokenLastCore()
101 , fvuStarDaqCmdLastCore()
102 , fvuStarTrigCmdLastCore()
103 , fvdMessCandidateTimeStart()
104 , fvdMessCandidateTimeStop()
105 , fvdTrigCandidateTimeStart()
106 , fvdTrigCandidateTimeStop()
107 , fbEpochAfterCandWinFound(kFALSE)
108 , fbTriggerAfterCandWinFound(kFALSE)
109 , fvbGdpbLastMissmatchPattern()
110 , fvbGdpbLastEnablePattern()
111 , fvbGdpbLastResyncPattern()
112 , fvSectorStatusPattern()
113 , fvhHitsTimeToTriggerRaw()
114 , fvhHitsTimeToTriggerRawPulser()
115 , fvhHitsTimeToTriggerSel()
116 , fvhHitsTimeToTriggerSelVsDaq()
117 , fvhHitsTimeToTriggerSelVsTrig()
118 , fvhTriggerDt()
119 , fvhTriggerDistributionInTs()
120 , fvhTriggerDistributionInMs()
121 , fvhMessDistributionInMs()
122 , fhEventSizeDistribution(nullptr)
123 , fhEventSizeEvolution(nullptr)
124 , fhEventNbEvolution(nullptr)
125 , fhEventNbDistributionInTs(nullptr)
126 , fhEventSizeDistributionInTs(nullptr)
127 , fhRawTriggersStats(nullptr)
128 , fhRawTriggersStatsCore(nullptr)
129 , fhRawTriggersStatsOver(nullptr)
130 , fhRawTriggersStatsSel(nullptr)
131 , fhMissingTriggersEvolution(nullptr)
132 , fcTimeToTrigRaw(nullptr)
133 , fcTimeToTrigSel(nullptr)
134 , fcTrigDistMs(nullptr)
135 , fcMessDistMs(nullptr)
136 , fcEventBuildStats(nullptr)
137 , fcTriggerStats(nullptr)
138{
139}
141{
143 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
144 fvvmEpSupprBuffer[uGdpb].clear();
145 fvvBufferMajorAsicErrors[uGdpb].clear();
146 fvvBufferMessages[uGdpb].clear();
147 fvvBufferTriggers[uGdpb].clear();
148 } // for (Int_t iGdpb = 0; iGdpb < fuNrOfGdpbs; ++iGdpb)
149}
150
151// -------------------------------------------------------------------------
153{
154 LOG(info) << "Initializing STAR eTOF 2019 event builder algo";
155
156 return kTRUE;
157}
165
166// -------------------------------------------------------------------------
168{
169 LOG(info) << "Init parameter containers for CbmStar2019EventBuilderEtofAlgo";
170 Bool_t initOK = ReInitContainers();
171
172 return initOK;
173}
175{
176 LOG(info) << "**********************************************";
177 LOG(info) << "ReInit parameter containers for CbmStar2019EventBuilderEtofAlgo";
178
179 fUnpackPar = (CbmStar2019TofPar*) fParCList->FindObject("CbmStar2019TofPar");
180 if (nullptr == fUnpackPar) return kFALSE;
181
182 Bool_t initOK = InitParameters();
183
184 return initOK;
185}
187{
188 if (nullptr == fParCList) fParCList = new TList();
189 fUnpackPar = new CbmStar2019TofPar("CbmStar2019TofPar");
190 fParCList->Add(fUnpackPar);
191
192 return fParCList;
193}
195{
197 // fbMonitorMode = fUnpackPar->GetMonitorMode();
198 LOG(info) << "Monitor mode: " << (fbMonitorMode ? "ON" : "OFF");
199
200 // fbDebugMonitorMode = fUnpackPar->GetDebugMonitorMode();
201 LOG(info) << "Debug Monitor mode: " << (fbDebugMonitorMode ? "ON" : "OFF");
202
203 LOG(info) << "Store Lost Event Msg in event: " << (fbStoreLostEventMsg ? "ON" : "OFF");
204
205 LOG(info) << "Add status pattern to event: " << (fbAddStatusToEvent ? "ON" : "OFF");
206
209 LOG(info) << "Nr. of Tof GDPBs: " << fuNrOfGdpbs;
210
212 LOG(info) << "Nr. of FEEs per Tof GDPB: " << fuNrOfFeePerGdpb;
213
215 LOG(info) << "Nr. of GET4 per Tof FEE: " << fuNrOfGet4PerFee;
216
218 LOG(info) << "Nr. of channels per GET4: " << fuNrOfChannelsPerGet4;
219
221 LOG(info) << "Nr. of channels per FEE: " << fuNrOfChannelsPerFee;
222
224 LOG(info) << "Nr. of GET4s: " << fuNrOfGet4;
225
227 LOG(info) << "Nr. of GET4s per GDPB: " << fuNrOfGet4PerGdpb;
228
230 LOG(info) << "Nr. of channels per GDPB: " << fuNrOfChannelsPerGdpb;
231
232 fGdpbIdIndexMap.clear();
233 for (UInt_t i = 0; i < fuNrOfGdpbs; ++i) {
235 LOG(info) << "GDPB Id of TOF " << Form("%02d", i) << " : " << std::hex << Form("0x%04x", fUnpackPar->GetGdpbId(i))
236 << std::dec;
237 } // for( UInt_t i = 0; i < fuNrOfGdpbs; ++i )
238
240 LOG(info) << "Nr. of GBTx: " << fuNrOfGbtx;
241
243 LOG(info) << "Nr. of GBTx: " << fuNrOfModules;
244
246 fviRpcType.resize(fuNrOfGbtx);
247 fviModuleId.resize(fuNrOfGbtx);
248 fviNrOfRpc.resize(fuNrOfGbtx);
249 fviRpcSide.resize(fuNrOfGbtx);
250 for (UInt_t uGbtx = 0; uGbtx < fuNrOfGbtx; ++uGbtx) {
251 fviNrOfRpc[uGbtx] = fUnpackPar->GetNrOfRpc(uGbtx);
252 fviRpcType[uGbtx] = fUnpackPar->GetRpcType(uGbtx);
253 fviRpcSide[uGbtx] = fUnpackPar->GetRpcSide(uGbtx);
254 fviModuleId[uGbtx] = fUnpackPar->GetModuleId(uGbtx);
255 } // for( UInt_t uGbtx = 0; uGbtx < uNrOfGbtx; ++uGbtx)
256
257 LOG(info) << "Nr. of RPCs per GBTx: ";
258 std::stringstream ss;
259 for (UInt_t uGbtx = 0; uGbtx < fuNrOfGbtx; ++uGbtx)
260 ss << Form(" %2d", fviNrOfRpc[uGbtx]);
261 LOG(info) << ss.str();
262
263 LOG(info) << "RPC type per GBTx: ";
264 ss.clear();
265 for (UInt_t uGbtx = 0; uGbtx < fuNrOfGbtx; ++uGbtx)
266 ss << Form(" %2d", fviRpcType[uGbtx]);
267 LOG(info) << ss.str();
268
269 LOG(info) << "RPC side per GBTx: ";
270 ss.clear();
271 for (UInt_t uGbtx = 0; uGbtx < fuNrOfGbtx; ++uGbtx)
272 ss << Form(" %2d", fviRpcSide[uGbtx]);
273 LOG(info) << ss.str();
274
275 LOG(info) << "Module ID per GBTx: ";
276 ss.clear();
277 for (UInt_t uGbtx = 0; uGbtx < fuNrOfGbtx; ++uGbtx)
278 ss << Form(" %2d", fviModuleId[uGbtx]);
279 LOG(info) << ss.str();
280
283 LOG(info) << "Timeslice parameters: each MS is " << fdMsSizeInNs << " ns";
284
288 fdStarTriggerDelay.resize(fuNrOfGdpbs, 0.0);
298 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
302 LOG(info) << Form("Trigger window parameters for gDPB %2u are: ", uGdpb) << fdStarTriggerDeadtime[uGdpb]
303 << " ns deadtime, " << fdStarTriggerDelay[uGdpb] << " ns delay, " << fdStarTriggerWinSize[uGdpb]
304 << " ns window width, " << fdAllowedTriggersSpread << " ns allowed spread";
305 } // for (Int_t iGdpb = 0; iGdpb < fuNrOfGdpbs; ++iGdpb)
306
308 fvulCurrentEpoch.resize(fuNrOfGdpbs, 0);
311
317
319 fvulGdpbTsMsb.resize(fuNrOfGdpbs, 0);
320 fvulGdpbTsLsb.resize(fuNrOfGdpbs, 0);
321 fvulStarTsMsb.resize(fuNrOfGdpbs, 0);
322 fvulStarTsMid.resize(fuNrOfGdpbs, 0);
325 fvuStarTokenLast.resize(fuNrOfGdpbs, 0);
333
334 return kTRUE;
335}
336// -------------------------------------------------------------------------
337
338void CbmStar2019EventBuilderEtofAlgo::AddMsComponentToList(size_t component, UShort_t usDetectorId)
339{
341 for (UInt_t uCompIdx = 0; uCompIdx < fvMsComponentsList.size(); ++uCompIdx)
342 if (component == fvMsComponentsList[uCompIdx]) return;
343
345 fvMsComponentsList.push_back(component);
346
347 LOG(info) << "CbmStar2019EventBuilderEtofAlgo::AddMsComponentToList => Component " << component
348 << " with detector ID 0x" << std::hex << usDetectorId << std::dec << " added to list";
349}
350// -------------------------------------------------------------------------
351
352Bool_t CbmStar2019EventBuilderEtofAlgo::ProcessTs(const fles::Timeslice& ts)
353{
354 fulCurrentTsIndex = ts.index();
355 fdTsStartTime = static_cast<Double_t>(ts.descriptor(0, 0).idx);
356
358 if (-1.0 == fdTsCoreSizeInNs) {
359 fuNbCoreMsPerTs = ts.num_core_microslices();
360 fuNbOverMsPerTs = ts.num_microslices(0) - ts.num_core_microslices();
363 LOG(info) << "Timeslice parameters: each TS has " << fuNbCoreMsPerTs << " Core MS and " << fuNbOverMsPerTs
364 << " Overlap MS, for a core duration of " << fdTsCoreSizeInNs << " ns and a full duration of "
365 << fdTsFullSizeInNs << " ns";
366
370 LOG(info) << "In each TS " << fuNbMsLoop << " MS will be looped over";
371 } // if( -1.0 == fdTsCoreSizeInNs )
372
375 // LOG(info) << Form( "TS %5d Start %12f Stop %12f", fulCurrentTsIndex, fdTsStartTime, fdTsStopTimeCore );
376
379 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
380 if (fdStarTriggerDelay[uGdpb] < 0.0) {
382 if (fdStarTriggerDelay[uGdpb] + fdStarTriggerWinSize[uGdpb] < 0.0) {
385 // << Accept more than needed as this should be safer and small amounts >>
387 + 2.0 * fdStarTriggerWinSize[uGdpb]; // + fdStarTriggerWinSize[ uGdpb ];
390
391 } // if( fdStarTriggerDelay[ uGdpb ] + fdStarTriggerWinSize[ uGdpb ] < 0.0 )
392 else {
395 // << Accept more than needed as this should be safer and small amounts >>
397 + 2.0 * fdStarTriggerWinSize[uGdpb]; // + fdStarTriggerDelay[ uGdpb ];
400 } // else of if( fdStarTriggerDelay[ uGdpb ] + fdStarTriggerWinSize[ uGdpb ] < 0.0 )
401 } // if( fdStarTriggerDeadtime[ uGdpb ] < 0.0 )
402 else {
404 // << Accept more than needed as this should be safer and small amounts >>
405 fvdMessCandidateTimeStart[uGdpb] = fdTsStartTime; // - fdAllowedTriggersSpread + fdStarTriggerDelay[ uGdpb ];
410 } // else of if( fdStarTriggerDelay[ uGdpb ] < 0.0 )
411 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
412
414 for (UInt_t uMsCompIdx = 0; uMsCompIdx < fvMsComponentsList.size(); ++uMsCompIdx) {
415 UInt_t uMsComp = fvMsComponentsList[uMsCompIdx];
416
418 fbTriggerFoundA = kFALSE;
419 fbTriggerFoundB = kFALSE;
420 fbTriggerFoundC = kFALSE;
421
425
427 for (fuMsIndex = 0; fuMsIndex < fuNbMsLoop; fuMsIndex++) {
428 if (kFALSE == ProcessMs(ts, uMsComp, fuMsIndex)) {
429 LOG(error) << "Failed to process ts " << fulCurrentTsIndex << " MS " << fuMsIndex << " for component "
430 << uMsComp;
431 return kFALSE;
432 } // if( kFALSE == ProcessMs( ts, uMsCompIdx, fuMsIndex ) )
433 } // for( fuMsIndex = 0; fuMsIndex < uNbMsLoop; fuMsIndex ++ )
434 } // for( UInt_t uMsCompIdx = 0; uMsCompIdx < fvMsComponentsList.size(); ++uMsCompIdx )
435
438 fvEventsBuffer.clear();
439
441 if (kFALSE == BuildEvents()) {
442 LOG(error) << "Failed to build events in ts " << fulCurrentTsIndex;
443 return kFALSE;
444 } // if( kFALSE == BuildEvents() )
445
447 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
448 fvvmEpSupprBuffer[uGdpb].clear();
449 fvvBufferMajorAsicErrors[uGdpb].clear();
450 fvvBufferMessages[uGdpb].clear();
451 fvvBufferTriggers[uGdpb].clear();
452
454 fvbGdpbLastEnablePattern[uGdpb] = 0;
455 fvbGdpbLastResyncPattern[uGdpb] = 0;
456 fvSectorStatusPattern[uGdpb].clear();
457 } // for (Int_t iGdpb = 0; iGdpb < fuNrOfGdpbs; ++iGdpb)
458
460 if (fbMonitorMode) {
461 if (kFALSE == FillHistograms()) {
462 LOG(error) << "Failed to fill histos in ts " << fulCurrentTsIndex;
463 return kFALSE;
464 } // if( kFALSE == FillHistograms() )
465 } // if( fbMonitorMode )
466
467 return kTRUE;
468}
469
470static Int_t iWarn = 0;
471static Int_t iWarnMess = 0;
472
473Bool_t CbmStar2019EventBuilderEtofAlgo::ProcessMs(const fles::Timeslice& ts, size_t uMsCompIdx, size_t uMsIdx)
474{
475 UInt_t uMsComp = fvMsComponentsList[uMsCompIdx];
476 auto msDescriptor = ts.descriptor(uMsComp, uMsIdx);
477
479 fuCurrentMs = uMsIdx;
480 fiEquipmentId = msDescriptor.eq_id;
481 fdMsTime = static_cast<double>(msDescriptor.idx);
482 uint32_t size = msDescriptor.size;
483 // fulLastMsIdx = msDescriptor.idx;
484 if (size > 0) LOG(debug) << "Microslice: " << msDescriptor.idx << " has size: " << size;
485
487 const uint8_t* msContent = reinterpret_cast<const uint8_t*>(ts.content(uMsComp, uMsIdx));
488
489 /*
490 if( fdStartTimeMsSz < 0 )
491 fdStartTimeMsSz = (1e-9) * fdMsTime;
492
493 fvhMsSzPerLink[ uMsComp ]->Fill(size);
494 if( 2 * fuHistoryHistoSize < (1e-9) * fdMsTime - fdStartTimeMsSz )
495 {
496 // Reset the evolution Histogram and the start time when we reach the end of the range
497 fvhMsSzTimePerLink[ uMsComp ]->Reset();
498 fdStartTimeMsSz = (1e-9) * fdMsTime;
499 } // if( 2 * fuHistoryHistoSize < (1e-9) * fdMsTime - fdStartTimeMsSz )
500 fvhMsSzTimePerLink[ uMsComp ]->Fill((1e-9) * fdMsTime - fdStartTimeMsSz, size);
501*/
502
504 if (0 != (size % fUnpackPar->GetNbByteMessage()))
505 LOG(error) << "The input microslice buffer does NOT "
506 << "contain only complete nDPB messages!";
507
509 uint32_t uNbMessages = (size - (size % fUnpackPar->GetNbByteMessage())) / fUnpackPar->GetNbByteMessage();
510
513
515 auto it = fGdpbIdIndexMap.find(fuGdpbId);
516 if (it == fGdpbIdIndexMap.end()) {
517 iWarn++;
518 LOG(warning) << "Could not find the gDPB index for AFCK id 0x" << std::hex << fuGdpbId << std::dec
519 << " in microslice " << fdMsTime << "\n"
520 << "If valid this index has to be added in the TOF parameter "
521 "file in the RocIdArray field";
522 if (iWarn == 100) LOG(fatal) << "Got max number of Warnings!";
523 return kFALSE;
524 } // if( it == fGdpbIdIndexMap.end() )
525 else
527
535 } // if( 0 < fuNbOverMsPerTs && fuNbCoreMsPerTs == fuCurrentMs )
536
538 if (0 < fuNbOverMsPerTs && 0 == fuCurrentMs) {
544 } // if( 0 < fuNbOverMsPerTs && 0 == fuCurrentMs )
545
546 // Prepare variables for the loop on contents
547 Int_t messageType = -111;
548 const uint64_t* pInBuff = reinterpret_cast<const uint64_t*>(msContent);
549 for (uint32_t uIdx = 0; uIdx < uNbMessages; uIdx++) {
552 break;
553
554 // Fill message
555 uint64_t ulData = static_cast<uint64_t>(pInBuff[uIdx]);
556
558 if (0 == uIdx) {
559 // std::cout << Form( "gDPB %2d", fuGdpbNr) << " Epoch cycle " << Form( "0x%012lx", ulData ) << std::endl;
560 ProcessEpochCycle(ulData);
561 continue;
562 } // if( 0 == uIdx )
563
564 gdpbv100::Message mess(ulData);
566 messageType = mess.getMessageType();
567 if (fUnpackPar->GetNrOfGet4PerGdpb() <= mess.getGdpbGenChipId() && 255 > mess.getGdpbGenChipId()) { //FIXME
568 if (iWarnMess < 100)
569 LOG(warn) << iWarnMess << ": Invalid ChipID " << mess.getGdpbGenChipId() << " in messageType " << messageType;
570 iWarnMess++;
571 continue;
572 }
573
576 /*
577 if( 5916 == fulCurrentTsIndex && gdpbv100::MSG_STAR_TRI_A <= messageType)
578 mess.printDataCout();
579*/
580 /*
581 if( gdpbv100::MSG_STAR_TRI_A <= messageType )
582 mess.printDataCout();
583 continue;
584*/
585
587 LOG(warning) << "Message with Get4 ID too high: " << fuGet4Id << " VS " << fuNrOfGet4PerGdpb
588 << " set in parameters.";
589
590 switch (messageType) {
591 case gdpbv100::MSG_HIT: {
592 if (mess.getGdpbHitIs24b()) {
593 LOG(error) << "This event builder does not support 24b hit message!!!.";
594 continue;
595 } // if( getGdpbHitIs24b() )
596 else {
597 fvvmEpSupprBuffer[fuGdpbNr].push_back(mess);
598 } // else of if( getGdpbHitIs24b() )
599 break;
600 } // case gdpbv100::MSG_HIT:
601 case gdpbv100::MSG_EPOCH: {
603 ProcessEpoch(mess);
604 /*
605 if( kTRUE == fbPrintAllEpochsEnable )
606 {
607 LOG(info) << "Epoch: " << Form("0x%08x ", fuGdpbId)
608 << ", Merg"
609 << ", Link " << std::setw(1) << mess.getGdpbEpLinkId()
610 << ", epoch " << std::setw(8) << mess.getGdpbEpEpochNb()
611 << ", Sync " << std::setw(1) << mess.getGdpbEpSync()
612 << ", Data loss " << std::setw(1) << mess.getGdpbEpDataLoss()
613 << ", Epoch loss " << std::setw(1) << mess.getGdpbEpEpochLoss()
614 << ", Epoch miss " << std::setw(1) << mess.getGdpbEpMissmatch();
615 } // if( kTRUE == fbPrintAllEpochsEnable )
616*/
617 } // if this epoch message is a merged one valid for all chips
618 else {
620 LOG(debug2) << "This event builder does not support unmerged epoch "
621 "messages!!!.";
622 continue;
623 } // if single chip epoch message
624 break;
625 } // case gdpbv100::MSG_EPOCH:
626 case gdpbv100::MSG_SLOWC: {
627 fvvmEpSupprBuffer[fuGdpbNr].push_back(mess);
628 break;
629 } // case gdpbv100::MSG_SLOWC:
630 case gdpbv100::MSG_SYST: {
631 fvvmEpSupprBuffer[fuGdpbNr].push_back(mess);
632 break;
633 } // case gdpbv100::MSG_SYST:
638 ProcessStarTrigger(mess);
639
642 /*
643 if( gdpbv100::MSG_STAR_TRI_A == messageType )
644 {
645 } // if( gdpbv100::MSG_STAR_TRI_A == messageType )
646*/
647 break;
648 } // case gdpbv100::MSG_STAR_TRI_A-D
649 default:
650 LOG(error) << "Message type " << std::hex << std::setw(2) << static_cast<uint16_t>(messageType)
651 << " not included in Get4 unpacker.";
652 } // switch( mess.getMessageType() )
653 } // for (uint32_t uIdx = 0; uIdx < uNbMessages; uIdx ++)
654
655 return kTRUE;
656}
657
658// -------------------------------------------------------------------------
660{
661 ULong64_t ulEpochCycleVal = ulCycleData & gdpbv100::kulEpochCycleFieldSz;
662 /*
663 if( fuRawDataPrintMsgIdx < fuRawDataPrintMsgNb || fair::Logger::Logging(fair::Severity::debug2) )
664 {
665 LOG(info) << "CbmMcbm2018MonitorTof::ProcessEpochCyle => "
666 << Form( " TS %5d MS %3d In data 0x%016X Cycle 0x%016X",
667 fulCurrentTsIndex, fuCurrentMs, ulCycleData, ulEpochCycleVal );
668 fuRawDataPrintMsgIdx ++;
669 }
670*/
671 if (!(ulEpochCycleVal == fvulCurrentEpochCycle[fuGdpbNr] || ulEpochCycleVal == fvulCurrentEpochCycle[fuGdpbNr] + 1)) {
672 LOG(warning) << "CbmStar2019EventBuilderEtofAlgo::ProcessEpochCycle => "
673 << " Missmatch in epoch cycles detected for Gdpb " << fuGdpbNr
674 << ", probably fake cycles due to epoch index corruption! "
675 << Form(" Current cycle 0x%09llX New cycle 0x%09llX", fvulCurrentEpochCycle[fuGdpbNr],
676 ulEpochCycleVal);
677 }
678 if (ulEpochCycleVal != fvulCurrentEpochCycle[fuGdpbNr]) {
679 LOG(info) << "CbmStar2019EventBuilderEtofAlgo::ProcessEpochCycle => "
680 << " New epoch cycle for Gdpb " << fuGdpbNr
681 << Form(": Current cycle 0x%09llX New cycle 0x%09llX", fvulCurrentEpochCycle[fuGdpbNr], ulEpochCycleVal);
682 } // if( ulEpochCycleVal != fvulCurrentEpochCycle[fuGdpbNr] )
683 fvulCurrentEpochCycle[fuGdpbNr] = ulEpochCycleVal;
684
685 return;
686}
687
689{
690 ULong64_t ulEpochNr = mess.getGdpbEpEpochNb();
691 /*
692 * /// FIXME: Need proper handling of overlap MS
693 if( 0 < fvulCurrentEpoch[ fuGdpbNr ] && ulEpochNr < fvulCurrentEpoch[ fuGdpbNr ] )
694 {
695 std::cout << Form( "gDPB %2d", fuGdpbNr) << " New Epoch cycle "
696 << Form( "0x%012llx old Ep %08llx new Ep %08llx", fvulCurrentEpochCycle[ fuGdpbNr ], fvulCurrentEpoch[ fuGdpbNr ], ulEpochNr )
697 << std::endl;
698 fvulCurrentEpochCycle[ fuGdpbNr ]++;
699 } // if( 0 < fvulCurrentEpoch[ fuGdpbNr ] && ulEpochNr < fvulCurrentEpoch[ fuGdpbNr ] )
700*/
701 fvulCurrentEpoch[fuGdpbNr] = ulEpochNr;
703
704 // fulCurrentEpochTime = mess.getMsgFullTime(ulEpochNr);
705
708 if (0 < ulEpochNr) mess.setGdpbEpEpochNb(ulEpochNr - 1);
709 else
711
714
719 } // if( fuNbCoreMsPerTs <= fuCurrentMs && fvdMessCandidateTimeStop[ fuGdpbNr ] < fvulCurrentEpochFull[ fuGdpbNr ] * gdpbv100::kdEpochInNs )
720}
721
723{
724 Int_t iMsgIndex = mess.getStarTrigMsgIndex();
725
726 switch (iMsgIndex) {
727 case 0: {
729
735 else
737
738 if (kTRUE == fbTriggerFoundA)
739 LOG(info) << Form("Problem: Overwritting a A trigger in TS %12llu MS %3u Core? %u", fulCurrentTsIndex,
741 << Form(" (%u %u %u)", fbTriggerFoundA, fbTriggerFoundB, fbTriggerFoundC);
742 } // if( fbMonitorMode && fbDebugMonitorMode )
743 fbTriggerFoundA = kTRUE;
744
745 break;
746 } // 1st trigger message
747 case 1: {
750
756 else
758
759 if (kFALSE == fbTriggerFoundA)
760 LOG(info) << Form("Problem: Found B trigger before A in TS %12llu MS %3u", fulCurrentTsIndex, fuMsIndex);
761
762 if (kTRUE == fbTriggerFoundB)
763 LOG(info) << Form("Problem: Overwritting a B trigger in TS %12llu MS %3u Core? %u", fulCurrentTsIndex,
765 << Form(" (%u %u %u)", fbTriggerFoundA, fbTriggerFoundB, fbTriggerFoundC);
766 } // if( fbMonitorMode && fbDebugMonitorMode )
767 fbTriggerFoundB = kTRUE;
768
769 break;
770 } // 2nd trigger message
771 case 2: {
773
779 else
781
782 if (kFALSE == fbTriggerFoundB || kFALSE == fbTriggerFoundA)
783 LOG(info) << Form("Problem: Found C trigger before A or B in TS %12llu MS %3u", fulCurrentTsIndex, fuMsIndex);
784
785 if (kTRUE == fbTriggerFoundC)
786 LOG(info) << Form("Problem: Overwritting a C trigger in TS %12llu MS %3u Core? %u", fulCurrentTsIndex,
788 << Form(" (%u %u %u)", fbTriggerFoundA, fbTriggerFoundB, fbTriggerFoundC);
789 } // if( fbMonitorMode && fbDebugMonitorMode )
790 fbTriggerFoundC = kTRUE;
791
792 break;
793 } // 3rd trigger message
794 case 3: {
795 ULong64_t ulNewGdpbTsFull = (fvulGdpbTsMsb[fuGdpbNr] << 24) + (fvulGdpbTsLsb[fuGdpbNr]);
796 ULong64_t ulNewStarTsFull =
798 UInt_t uNewToken = mess.getStarTokenStarD();
799 UInt_t uNewDaqCmd = mess.getStarDaqCmdStarD();
800 UInt_t uNewTrigCmd = mess.getStarTrigCmdStarD();
801
805 if (kFALSE == fbTriggerFoundC || kFALSE == fbTriggerFoundB || kFALSE == fbTriggerFoundA) {
807 LOG(info) << Form("Problem: Found D trigger before A or B or C in TS "
808 "%12llu MS %3u Core? %u",
810 << Form(" (%u %u %u) token %5u", fbTriggerFoundA, fbTriggerFoundB, fbTriggerFoundC, uNewToken);
811
813 fbTriggerFoundA = kFALSE;
814 fbTriggerFoundB = kFALSE;
815 fbTriggerFoundC = kFALSE;
816
818 if (0 == fuCurrentMs) return;
819 } // if( kFALSE == fbTriggerFoundC || kFALSE == fbTriggerFoundB || kFALSE == fbTriggerFoundA )
820
822 fbTriggerFoundA = kFALSE;
823 fbTriggerFoundB = kFALSE;
824 fbTriggerFoundC = kFALSE;
825
830 else
832 } // if( fbMonitorMode && fbDebugMonitorMode )
833
834 /*
835 UInt_t uNewTrigWord = ( (uNewTrigCmd & 0x00F) << 16 )
836 + ( (uNewDaqCmd & 0x00F) << 12 )
837 + ( (uNewToken & 0xFFF) );
838 LOG(info) << "New STAR trigger "
839 << " TS " << fulCurrentTsIndex
840 << " gDBB #" << fuGdpbNr << " "
841 << Form("token = %5u ", uNewToken )
842 << Form("gDPB ts = %12llu ", ulNewGdpbTsFull )
843 << Form("STAR ts = %12llu ", ulNewStarTsFull )
844 << Form("DAQ cmd = %2u ", uNewDaqCmd )
845 << Form("TRG cmd = %2u ", uNewTrigCmd )
846 << Form("TRG Wrd = %5x ", uNewTrigWord );
847*/
848
849 if ((uNewToken == fvuStarTokenLast[fuGdpbNr]) && (ulNewGdpbTsFull == fvulGdpbTsFullLast[fuGdpbNr])
850 && (ulNewStarTsFull == fvulStarTsFullLast[fuGdpbNr]) && (uNewDaqCmd == fvuStarDaqCmdLast[fuGdpbNr])
851 && (uNewTrigCmd == fvuStarTrigCmdLast[fuGdpbNr])) {
852 UInt_t uTrigWord = ((fvuStarTrigCmdLast[fuGdpbNr] & 0x00F) << 16)
853 + ((fvuStarDaqCmdLast[fuGdpbNr] & 0x00F) << 12) + ((fvuStarTokenLast[fuGdpbNr] & 0xFFF));
854 LOG(warning) << "Possible error: identical STAR tokens found twice in "
855 "a row => ignore 2nd! "
856 << " TS " << fulCurrentTsIndex << " gDBB #" << fuGdpbNr << " "
857 << Form("token = %5u ", fvuStarTokenLast[fuGdpbNr])
858 << Form("gDPB ts = %12llu ", fvulGdpbTsFullLast[fuGdpbNr])
859 << Form("STAR ts = %12llu ", fvulStarTsFullLast[fuGdpbNr])
860 << Form("DAQ cmd = %2u ", fvuStarDaqCmdLast[fuGdpbNr])
861 << Form("TRG cmd = %2u ", fvuStarTrigCmdLast[fuGdpbNr]) << Form("TRG Wrd = %5x ", uTrigWord);
862 return;
863 } // if exactly same message repeated
864
865 // GDPB TS counter reset detection
866 if (ulNewGdpbTsFull < fvulGdpbTsFullLast[fuGdpbNr])
867 LOG(debug) << "Probable reset of the GDPB TS: old = " << Form("%16llu", fvulGdpbTsFullLast[fuGdpbNr])
868 << " new = " << Form("%16llu", ulNewGdpbTsFull) << " Diff = -"
869 << Form("%8llu", fvulGdpbTsFullLast[fuGdpbNr] - ulNewGdpbTsFull) << " GDPB #"
870 << Form("%2u", fuGdpbNr);
871
872 // STAR TS counter reset detection
873 if (ulNewStarTsFull < fvulStarTsFullLast[fuGdpbNr])
874 LOG(debug) << "Probable reset of the STAR TS: old = " << Form("%16llu", fvulStarTsFullLast[fuGdpbNr])
875 << " new = " << Form("%16llu", ulNewStarTsFull) << " Diff = -"
876 << Form("%8llu", fvulStarTsFullLast[fuGdpbNr] - ulNewStarTsFull) << " GDPB #"
877 << Form("%2u", fuGdpbNr);
878
879 /*
880 LOG(info) << "Updating trigger token for " << fuGdpbNr
881 << " " << fvuStarTokenLast[fuGdpbNr] << " " << uNewToken;
882*/
883
885 CbmTofStarTrigger2019 newTrig(ulNewGdpbTsFull, ulNewStarTsFull, uNewToken, uNewDaqCmd, uNewTrigCmd, fuGdpbId);
886 Double_t dTriggerTime = newTrig.GetFullGdpbTs() * gdpbv100::kdClockCycleSizeNs;
887 if (fvdTrigCandidateTimeStart[fuGdpbNr] < dTriggerTime && dTriggerTime < fvdTrigCandidateTimeStop[fuGdpbNr]) {
888 fvvBufferTriggers[fuGdpbNr].push_back(newTrig);
889
894 else
896 } // if( fbMonitorMode && fbDebugMonitorMode )
897 } // if( fvdTrigCandidateTimeStart[ fuGdpbNr ] < dTriggerTime && dTriggerTime < fvdTrigCandidateTimeStop[ fuGdpbNr ] )
898 else if (fuNbCoreMsPerTs <= fuCurrentMs && fvdTrigCandidateTimeStop[fuGdpbNr] < dTriggerTime) {
901 } // else if( fuNbCoreMsPerTs <= fuCurrentMs && fvdTrigCandidateTimeStop[ fuGdpbNr ] < dTriggerTime )
902
903 // if( fuCurrentMs < fuNbCoreMsPerTs )
904 {
905 // ULong64_t ulGdpbTsDiff = ulNewGdpbTsFull - fvulGdpbTsFullLast[fuGdpbNr];
906 fvulGdpbTsFullLast[fuGdpbNr] = ulNewGdpbTsFull;
907 fvulStarTsFullLast[fuGdpbNr] = ulNewStarTsFull;
908 fvuStarTokenLast[fuGdpbNr] = uNewToken;
909 fvuStarDaqCmdLast[fuGdpbNr] = uNewDaqCmd;
910 fvuStarTrigCmdLast[fuGdpbNr] = uNewTrigCmd;
911 }
918 else
920 fvhTriggerDistributionInTs[fuGdpbNr]->Fill((dTriggerTime - fdTsStartTime) / 1000.0);
921 fvhTriggerDistributionInMs[fuGdpbNr]->Fill((dTriggerTime - fdMsTime) / 1000.0);
922
923 if ((dTriggerTime - fdMsTime) / 1000.0 < -1000) {
924 LOG(info) << Form("Trigger in wrong MS TS %6llu", fulCurrentTsIndex) << Form(" MS %3u ", fuMsIndex)
925 << Form(" Sector %3u ", fuGdpbNr + 13) << Form(" Ttrig %15.2f", dTriggerTime)
926 << Form(" Tms %15.2f", fdMsTime) << Form(" dT %15.5f", ((dTriggerTime - fdMsTime) / 1000.0));
927 LOG(info) << Form("Full token, gDPB TS LSB bits: 0x%16llx, STAR TS "
928 "MSB bits: 0x%16llx, token is %4u",
929 ulNewGdpbTsFull, ulNewStarTsFull, uNewToken);
930 } // if( (dTriggerTime - fdMsTime) / 1000.0 < -1000 )
931 } // if( fbMonitorMode && fbDebugMonitorMode )
932
933 // LOG(info) << "First full trigger in TS " << fulCurrentTsIndex;
934
935 break;
936 } // case 3
937 default: LOG(error) << "Unknown Star Trigger messageindex: " << iMsgIndex;
938 } // switch( iMsgIndex )
939}
940
941// -------------------------------------------------------------------------
943{
944 Int_t iBufferSize = fvvmEpSupprBuffer[fuGdpbNr].size();
945
946 if (0 == iBufferSize) return;
947
948 LOG(debug) << "Now processing stored messages for for gDPB " << fuGdpbNr << " with epoch number "
949 << (fvulCurrentEpoch[fuGdpbNr] - 1);
950
953 std::stable_sort(fvvmEpSupprBuffer[fuGdpbNr].begin(), fvvmEpSupprBuffer[fuGdpbNr].begin());
954
956 ULong64_t ulCurEpochGdpbGet4 = fvulCurrentEpochFull[fuGdpbNr];
957
959 if (0 == ulCurEpochGdpbGet4) return;
960
962 ulCurEpochGdpbGet4--;
963
964 Int_t messageType = -111;
965 for (Int_t iMsgIdx = 0; iMsgIdx < iBufferSize; iMsgIdx++) {
966
967 messageType = fvvmEpSupprBuffer[fuGdpbNr][iMsgIdx].getMessageType();
968
969 fuGet4Id = fUnpackPar->ElinkIdxToGet4Idx(fvvmEpSupprBuffer[fuGdpbNr][iMsgIdx].getGdpbGenChipId());
971
973 gdpbv100::FullMessage fullMess(fvvmEpSupprBuffer[fuGdpbNr][iMsgIdx], ulCurEpochGdpbGet4);
974
976 switch (messageType) {
977 case gdpbv100::MSG_HIT: {
978 // ProcessHit( fullMess, ulCurEpochGdpbGet4 );
980 break;
981 } // case gdpbv100::MSG_HIT:
982 case gdpbv100::MSG_SLOWC: {
983 // ProcessSlCtrl( fullMess, ulCurEpochGdpbGet4 );
985 break;
986 } // case gdpbv100::MSG_SLOWC:
987 case gdpbv100::MSG_SYST: {
988 if (gdpbv100::SYS_PATTERN == fullMess.getGdpbSysSubType()) {
989 if (fbAddStatusToEvent) ProcessPattern(fullMess, ulCurEpochGdpbGet4);
990 } // if( gdpbv100::SYS_PATTERN == fvvmEpSupprBuffer[ fuGdpbNr ][ iMsgIdx ].getGdpbSysSubType() )
991 /*
992 else ProcessSysMess( fullMess, ulCurEpochGdpbGet4 );
993*/
994 else
996
997 break;
998 } // case gdpbv100::MSG_SYST:
1003 case gdpbv100::MSG_STAR_TRI_D: break;
1004 default:
1005 LOG(error) << "Message type " << std::hex << std::setw(2) << static_cast<uint16_t>(messageType)
1006 << " not included in Get4 unpacker.";
1007 } // switch( mess.getMessageType() )
1008 } // for( Int_t iMsgIdx = 0; iMsgIdx < iBufferSize; iMsgIdx++ )
1009
1010 fvvmEpSupprBuffer[fuGdpbNr].clear();
1011}
1013{
1015 uint16_t usG4ErrorType = fullMess.getGdpbSysErrData();
1017 && ((usG4ErrorType < gdpbv100::GET4_V2X_ERR_TOT_OVERWRT && gdpbv100::GET4_V2X_ERR_LOST_EVT != usG4ErrorType)
1018 || usG4ErrorType > gdpbv100::GET4_V2X_ERR_SEQUENCE_ER)) {
1019 fvvBufferMajorAsicErrors[fuGdpbNr].push_back(fullMess);
1020 return;
1021 } // if GET4 error out of TOT/hit building error range
1022 else if (fbStoreLostEventMsg && gdpbv100::GET4_V2X_ERR_LOST_EVT == usG4ErrorType) {
1023 fvvBufferMajorAsicErrors[fuGdpbNr].push_back(fullMess);
1024 return;
1025 } // else if( fbStoreLostEventMsg && gdpbv100::GET4_V2X_ERR_LOST_EVT == usG4ErrorType )
1026
1028 fvhMessDistributionInMs[fuGdpbNr]->Fill((fullMess.GetFullTimeNs() - fdMsTime) / 1000.0);
1029 /*
1030 LOG(info) << Form( "Message Full Time ns: %f MS time ns: %f diff: %f Start %f Stop %f",
1031 fullMess.GetFullTimeNs(), fdMsTime, fullMess.GetFullTimeNs() - fdMsTime,
1032 fvdMessCandidateTimeStart[ fuGdpbNr ], fvdMessCandidateTimeStop[ fuGdpbNr ] );
1033 LOG(info) << Form( "Current epoch %llu Current cycle %llu Current Full epoch %llu",
1034 fvulCurrentEpoch[ fuGdpbNr ], fvulCurrentEpochCycle[ fuGdpbNr ], fvulCurrentEpochFull[ fuGdpbNr ] );
1035*/
1040 fvvBufferMessages[fuGdpbNr].push_back(fullMess);
1041 } // if in limits
1042}
1043
1044// -------------------------------------------------------------------------
1045void CbmStar2019EventBuilderEtofAlgo::ProcessHit(gdpbv100::Message mess, uint64_t /*ulCurEpochGdpbGet4*/)
1046{
1047 UInt_t uChannel = mess.getGdpbHitChanId();
1048 // UInt_t uTot = mess.getGdpbHit32Tot();
1049
1050 // In 32b mode the coarse counter is already computed back to 112 FTS bins
1051 // => need to hide its contribution from the Finetime
1052 // => FTS = Fullt TS modulo 112
1053 // UInt_t uFts = mess.getGdpbHitFullTs() % 112;
1054
1055 // UInt_t uChannelNr = fuGet4Id * fuNrOfChannelsPerGet4 + uChannel;
1056 // UInt_t uChannelNrInFee = (fuGet4Id % fuNrOfGet4PerFee) * fuNrOfChannelsPerGet4 + uChannel;
1057 // UInt_t uFeeNr = (fuGet4Id / fuNrOfGet4PerFee);
1058 // UInt_t uFeeNrInSys = fuGdpbNr * fuNrOfFeePerGdpb + uFeeNr;
1059 // UInt_t uRemappedChannelNr = uFeeNr * fuNrOfChannelsPerFee + fUnpackPar->Get4ChanToPadiChan( uChannelNrInFee );
1060 // UInt_t uGbtxNr = (uFeeNr / fUnpackPar->GetNrOfFeePerGbtx());
1061 // UInt_t uFeeInGbtx = (uFeeNr % fUnpackPar->GetNrOfFeePerGbtx());
1062 // UInt_t uGbtxNrInSys = fuGdpbNr * fUnpackPar->GetNrOfGbtxPerGdpb() + uGbtxNr;
1063
1064 // ULong_t ulhitTime = mess.getMsgFullTime( ulCurEpochGdpbGet4 );
1065 // Double_t dHitTime = mess.getMsgFullTimeD( ulCurEpochGdpbGet4 );
1066 // Double_t dHitTot = uTot; // in bins
1067
1068 // UInt_t uFebIdx = (fuGet4Id / fUnpackPar->GetNrOfGet4PerFee());
1069 // UInt_t uFullFebIdx = (fuGdpbNr * fUnpackPar->GetNrOfFeePerGdpb()) + uFebIdx;
1070
1071 UInt_t uChanInGdpb = fuGet4Id * fuNrOfChannelsPerGet4 + uChannel;
1072 UInt_t uChanInSyst = fuGdpbNr * fuNrOfChannelsPerGdpb + uChanInGdpb;
1073 if (fUnpackPar->GetNumberOfChannels() < uChanInSyst) {
1074 LOG(error) << "Invalid mapping index " << uChanInSyst << " VS " << fUnpackPar->GetNumberOfChannels() << ", from "
1075 << fuGdpbNr << ", " << fuGet4Id << ", " << uChannel;
1076 return;
1077 } // if( fUnpackPar->GetNumberOfChannels() < uChanUId )
1078}
1079
1080void CbmStar2019EventBuilderEtofAlgo::ProcessSlCtrl(gdpbv100::Message /*mess*/, uint64_t /*ulCurEpochGdpbGet4*/) {}
1081
1083{
1084 switch (mess.getGdpbSysSubType()) {
1086 uint32_t uData = mess.getGdpbSysErrData();
1091 LOG(debug) << " +++++++ > gDPB: " << std::hex << std::setw(4) << fuGdpbId << std::dec
1092 << ", Chip = " << std::setw(2) << mess.getGdpbGenChipId() << ", Chan = " << std::setw(1)
1093 << mess.getGdpbSysErrChanId() << ", Edge = " << std::setw(1) << mess.getGdpbSysErrEdge()
1094 << ", Empt = " << std::setw(1) << mess.getGdpbSysErrUnused() << ", Data = " << std::hex
1095 << std::setw(2) << uData << std::dec << " -- GET4 V1 Error Event";
1096 else
1097 LOG(debug) << " +++++++ >gDPB: " << std::hex << std::setw(4) << fuGdpbId << std::dec
1098 << ", Chip = " << std::setw(2) << mess.getGdpbGenChipId() << ", Chan = " << std::setw(1)
1099 << mess.getGdpbSysErrChanId() << ", Edge = " << std::setw(1) << mess.getGdpbSysErrEdge()
1100 << ", Empt = " << std::setw(1) << mess.getGdpbSysErrUnused() << ", Data = " << std::hex
1101 << std::setw(2) << uData << std::dec << " -- GET4 V1 Error Event ";
1102 break;
1103 } // case gdpbv100::SYSMSG_GET4_EVENT
1105 LOG(debug) << "Unknown GET4 message, data: " << std::hex << std::setw(8) << mess.getGdpbSysUnkwData() << std::dec
1106 << " Full message: " << std::hex << std::setw(16) << mess.getData() << std::dec;
1107 break;
1108 } // case gdpbv100::SYS_GDPB_UNKWN:
1110 LOG(debug) << "GET4 synchronization pulse missing";
1111 break;
1112 } // case gdpbv100::SYS_GET4_SYNC_MISS:
1113 case gdpbv100::SYS_PATTERN: {
1114 LOG(debug) << "ASIC pattern for missmatch, disable or resync";
1115 break;
1116 } // case gdpbv100::SYS_PATTERN:
1117 default: {
1118 LOG(debug) << "Crazy system message, subtype " << mess.getGdpbSysSubType();
1119 break;
1120 } // default
1121
1122 } // switch( getGdpbSysSubType() )
1123}
1124// -------------------------------------------------------------------------
1126{
1127 uint16_t usType = mess.getGdpbSysPattType();
1128 uint16_t usIndex = mess.getGdpbSysPattIndex();
1129 uint32_t uPattern = mess.getGdpbSysPattPattern();
1130 // std::bitset< 32 > bPattern( mess.getGdpbSysPattPattern() );
1131 UInt_t uNbBits = (kuNbMsgPerPattern - 1 == usIndex ? 16 : 32);
1132
1133 switch (usType) {
1135 // LOG(debug) << Form( "Missmatch pattern message => Type %d, Index %2d, Pattern 0x%08X", usType, usIndex, uPattern );
1136 for (UInt_t uBit = 0; uBit < uNbBits; ++uBit) {
1137 fvbGdpbLastMissmatchPattern[fuGdpbNr][32 * usIndex + uBit] = (uPattern >> uBit) & 0x1;
1138 // fvbGdpbLastMissmatchPattern[ fuGdpbNr ][ 32 * usIndex + uBit ] = bPattern[ uBit ];
1139 } // for( UInt_t uBit = 0; uBit < uNbBits; ++uBit )
1140 break;
1141 } // case gdpbv100::PATT_MISSMATCH:
1142 case gdpbv100::PATT_ENABLE: {
1143 // LOG(debug2) << Form( "Enable pattern message => Type %d, Index %2d, Pattern 0x%08X", usType, usIndex, uPattern );
1144 for (UInt_t uBit = 0; uBit < uNbBits; ++uBit) {
1145 fvbGdpbLastEnablePattern[fuGdpbNr][32 * usIndex + uBit] = (uPattern >> uBit) & 0x1;
1146 // fvbGdpbLastEnablePattern[ fuGdpbNr ][ 32 * usIndex + uBit ] = bPattern[ uBit ];
1147 } // for( UInt_t uBit = 0; uBit < uNbBits; ++uBit )
1148 break;
1149 } // case gdpbv100::PATT_ENABLE:
1150 case gdpbv100::PATT_RESYNC: {
1151 // LOG(debug) << Form( "RESYNC pattern message => Type %d, Index %2d, Pattern 0x%08X", usType, usIndex, uPattern );
1152 for (UInt_t uBit = 0; uBit < uNbBits; ++uBit) {
1153 fvbGdpbLastResyncPattern[fuGdpbNr][32 * usIndex + uBit] = (uPattern >> uBit) & 0x1;
1154 // fvbGdpbLastResyncPattern[ fuGdpbNr ][ 32 * usIndex + uBit ] = bPattern[ uBit ];
1155 } // for( UInt_t uBit = 0; uBit < uNbBits; ++uBit )
1156 break;
1157 } // case gdpbv100::PATT_RESYNC:
1158 default: {
1159 LOG(debug) << "Crazy pattern message, subtype " << usType;
1160 break;
1161 } // default
1162 } // switch( usType )
1163
1164 if (kuNbMsgPerPattern - 1 == usIndex) UpdateStatusPatternCurrGdpb();
1165
1166 return;
1167}
1181 std::pair<uint64_t, std::bitset<kuNbAsicPerGdpb>> statusIn)
1182{
1184 static constexpr std::bitset<kuNbAsicPerGdpb> kbMask32 {0xFFFFFFFF};
1185
1187 gdpbv100::FullMessage outMess(0, statusIn.first);
1188
1190 const std::bitset<kuNbAsicPerGdpb> bPatterBlock {(statusIn.second >> (32 * uIndex)) & kbMask32};
1191
1193 outMess.setGdpbGenGdpbId(uGdpbId);
1197 outMess.setGdpbSysPattIndex(uIndex);
1198 outMess.setGdpbSysPattPattern(bPatterBlock.to_ulong()); // Safe as we masked with 32b
1199
1200 return outMess;
1201}
1202// -------------------------------------------------------------------------
1203
1205{
1210 std::vector<std::vector<CbmTofStarTrigger2019>::iterator> itTrigger;
1211 std::vector<std::vector<gdpbv100::FullMessage>::iterator> itErrorMessStart;
1212 std::vector<std::vector<gdpbv100::FullMessage>::iterator> itMessStart;
1213 std::vector<std::vector<std::pair<uint64_t, std::bitset<kuNbAsicPerGdpb>>>::iterator> itStatUpdtStart;
1214 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
1215 itTrigger.push_back(fvvBufferTriggers[uGdpb].begin());
1216 itErrorMessStart.push_back(fvvBufferMajorAsicErrors[uGdpb].begin());
1217 itMessStart.push_back(fvvBufferMessages[uGdpb].begin());
1218 itStatUpdtStart.push_back(fvSectorStatusPattern[uGdpb].begin());
1219 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1220
1222 Bool_t bAllSectAllTriggDone = kTRUE;
1223 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
1225 if (fvvBufferTriggers[uGdpb].end() != itTrigger[uGdpb]) {
1226 bAllSectAllTriggDone = kFALSE;
1227 break;
1228 } // if( fvvBufferTriggers[ uGdpb ].end() != (*itTrigger[ uGdpb ]) )
1229 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1230
1234 while (kFALSE == bAllSectAllTriggDone) {
1236 Bool_t bAllSectorsMatch = kTRUE;
1237 UInt_t uFirstStarToken = 0;
1238 UShort_t usFirstStarDaqCmd = 0;
1239 UShort_t usFirstStarTrigCmd = 0;
1240 ULong64_t ulFirstFullGdpbTs = 0;
1241
1244 if (0 == fvvBufferTriggers[0].size()) { bAllSectorsMatch = kFALSE; } // if( 0 == fvvBufferTriggers[ 0 ].size() )
1245 else {
1248 uFirstStarToken = (*itTrigger[0]).GetStarToken();
1249 usFirstStarDaqCmd = (*itTrigger[0]).GetStarDaqCmd();
1250 usFirstStarTrigCmd = (*itTrigger[0]).GetStarTrigCmd();
1251 ulFirstFullGdpbTs = (*itTrigger[0]).GetFullGdpbTs();
1252 for (UInt_t uGdpb = 1; uGdpb < fuNrOfGdpbs; ++uGdpb) {
1253 if (itTrigger[uGdpb] == fvvBufferTriggers[uGdpb].end() || (*itTrigger[uGdpb]).GetStarToken() != uFirstStarToken
1254 || (*itTrigger[uGdpb]).GetStarDaqCmd() != usFirstStarDaqCmd
1255 || (*itTrigger[uGdpb]).GetStarTrigCmd() != usFirstStarTrigCmd) {
1256 bAllSectorsMatch = kFALSE;
1257 break;
1258 } // If end of buffer or any field differs for any gDPB/sector current trigger, not all matched!
1259 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1260 } // else of if( 0 == itTrigger[ uGdpb ].size() )
1261
1262 if (kTRUE == bAllSectorsMatch) {
1265 CbmTofStarSubevent2019 starSubEvent;
1266
1268 Double_t dMeanTriggerGdpbTs = 0;
1269 Double_t dMeanTriggerStarTs = 0;
1270 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
1272 dMeanTriggerGdpbTs += (*itTrigger[uGdpb]).GetFullGdpbTs();
1273 dMeanTriggerStarTs += (*itTrigger[uGdpb]).GetFullStarTs();
1274
1276 Long64_t ilTriggerDt =
1277 static_cast<Long64_t>((*itTrigger[uGdpb]).GetFullGdpbTs()) - static_cast<Long64_t>(ulFirstFullGdpbTs);
1278 fvhTriggerDt[uGdpb]->Fill(ilTriggerDt);
1279 } // if( fbMonitorMode && fbDebugMonitorMode )
1280
1282 std::vector<gdpbv100::FullMessage> vTrigMess = (*itTrigger[uGdpb]).GetGdpbMessages();
1283 for (std::vector<gdpbv100::FullMessage>::iterator itMess = vTrigMess.begin(); itMess != vTrigMess.end();
1284 ++itMess) {
1285 starSubEvent.AddMsg((*itMess));
1286 } // for( std::vector< gdpbv100::FullMessage >::iterator itMess = vTrigMess.begin(); itMess != vTrigMess.end(); ++ itMess )
1287
1288 if (fbAddStatusToEvent) {
1289 /*
1290 LOG(info) << Form( "Trying to add status to event, status buffer size is %d, first test status time is %llu, trigger time is %llu Evt size is %u",
1291 fvSectorStatusPattern[ uGdpb ].size(),
1292 (*itStatUpdtStart[ uGdpb ]).first * gdpbv100::kuCoarseCounterSize,
1293 (*itTrigger[ uGdpb ]).GetFullGdpbTs(),
1294 starSubEvent.GetMsgBuffSize() );
1295*/
1300 while (fvSectorStatusPattern[uGdpb].end() != itStatUpdtStart[uGdpb]
1301 && ((*itStatUpdtStart[uGdpb]).first) * gdpbv100::kuCoarseCounterSize
1302 <= (*itTrigger[uGdpb]).GetFullGdpbTs()) {
1303 ++itStatUpdtStart[uGdpb];
1304 } // while( not at the end of updates vector AND Update time under trigger time )
1306 if (fvSectorStatusPattern[uGdpb].begin() != itStatUpdtStart[uGdpb]) {
1308 --itStatUpdtStart[uGdpb];
1310 uint16_t usGdpbId = fUnpackPar->GetGdpbId(uGdpb);
1311 for (uint32_t uIdx = 0; uIdx < kuNbMsgPerPattern; ++uIdx)
1312 starSubEvent.AddMsg(CreateStatusMessage(usGdpbId, uIdx, (*itStatUpdtStart[uGdpb])));
1313 /*
1314 LOG(info) << Form( "found status time is %llu, diff to trigger time is %llu Evt size is %u",
1315 ( (*itStatUpdtStart[ uGdpb ]).first ) * gdpbv100::kuCoarseCounterSize,
1316 (*itTrigger[ uGdpb ]).GetFullGdpbTs() - ( (*itStatUpdtStart[ uGdpb ]).first ) * gdpbv100::kuCoarseCounterSize,
1317 starSubEvent.GetMsgBuffSize() );
1318*/
1319 } // if( fvSectorStatusPattern[ uGdpb ].begin() != itStatUpdtStart[ uGdpb ] )
1320 } // if( fbAddStatusToEvent )
1321
1323 Double_t dWinBeg =
1324 gdpbv100::kdClockCycleSizeNs * (*itTrigger[uGdpb]).GetFullGdpbTs() + fdStarTriggerDelay[uGdpb];
1325 Double_t dWinEnd = gdpbv100::kdClockCycleSizeNs * (*itTrigger[uGdpb]).GetFullGdpbTs()
1326 + fdStarTriggerDelay[uGdpb] + fdStarTriggerWinSize[uGdpb];
1327 Double_t dDeadEnd =
1328 gdpbv100::kdClockCycleSizeNs * (*itTrigger[uGdpb]).GetFullGdpbTs() + fdStarTriggerDeadtime[uGdpb];
1329 Double_t dWinBegErrors = dWinBeg - 2 * gdpbv100::kdEpochInNs;
1330
1332 UInt_t uNbErrorsInEventGdpb = 0;
1333 while (itErrorMessStart[uGdpb] != fvvBufferMajorAsicErrors[uGdpb].end()
1334 && (*itErrorMessStart[uGdpb]).GetFullTimeNs() < dWinBegErrors) {
1335 ++itErrorMessStart[uGdpb];
1336 } // while( not at buffer end && (*itErrorMessStart[ uGdpb ]).GetFullTimeNs() < dWinBeg - 2 * gdpbv100::kdEpochInNs )
1337 while (itErrorMessStart[uGdpb] != fvvBufferMajorAsicErrors[uGdpb].end()
1338 && (*itErrorMessStart[uGdpb]).GetFullTimeNs() < dWinEnd) {
1340 if (uNbErrorsInEventGdpb < kuMaxNbErrorsPerGdpbPerEvent) starSubEvent.AddMsg((*itErrorMessStart[uGdpb]));
1341
1342 ++itErrorMessStart[uGdpb];
1343 ++uNbErrorsInEventGdpb;
1344 } // while( not at buffer end && (*itErrorMessStart[ uGdpb ]).GetFullTimeNs() < dWinEnd )
1345
1348 std::vector<gdpbv100::FullMessage>::iterator itFirstMessOutOfDeadtime = itMessStart[uGdpb];
1349 while (
1350 itMessStart[uGdpb] != fvvBufferMessages[uGdpb].end()
1351 && ((*itMessStart[uGdpb]).GetFullTimeNs() < dDeadEnd || (*itMessStart[uGdpb]).GetFullTimeNs() < dWinEnd)) {
1352 Double_t dMessTime = (*itMessStart[uGdpb]).GetFullTimeNs();
1353 /*
1354 if( 5916 == fulCurrentTsIndex )
1355 LOG(info) << Form( "gDPB %2u Match Mess %12f Start %12f Stop %12f dStart %6f dStop %6f",
1356 uGdpb, dMessTime,
1357 dWinBeg, dWinEnd,
1358 dMessTime - dWinBeg, dWinEnd - dMessTime )
1359 << ( dWinBeg < dMessTime && dMessTime < dWinEnd ? " IN" : "" );
1360*/
1361
1363 if (fbMonitorMode) {
1364 fvhHitsTimeToTriggerRaw[uGdpb]->Fill(dMessTime
1365 - gdpbv100::kdClockCycleSizeNs * (*itTrigger[uGdpb]).GetFullGdpbTs());
1366 if (fbDebugMonitorMode && gdpbv100::MSG_HIT == (*itMessStart[uGdpb]).getMessageType()) {
1367 UInt_t uTot = (*itMessStart[uGdpb]).getGdpbHit32Tot();
1368 Double_t dTimeDiff = dMessTime - gdpbv100::kdClockCycleSizeNs * (*itTrigger[uGdpb]).GetFullGdpbTs();
1369 if (fuMinTotPulser < uTot && uTot < fuMaxTotPulser && TMath::Abs(dTimeDiff) < 20e3)
1370 fvhHitsTimeToTriggerRawPulser[uGdpb]->Fill(dTimeDiff, uTot);
1371 } // if( fbDebugMonitorMode && gdpbv100::MSG_HIT == (*itMessStart[ uGdpb ]).getMessageType() )
1372 } // if( fbMonitorMode )
1373
1374
1376 if (dWinBeg < dMessTime && dMessTime < dWinEnd) {
1377 starSubEvent.AddMsg((*itMessStart[uGdpb]));
1379 if (fbMonitorMode) {
1380 fvhHitsTimeToTriggerSel[uGdpb]->Fill(
1381 dMessTime - gdpbv100::kdClockCycleSizeNs * (*itTrigger[uGdpb]).GetFullGdpbTs());
1382 if (fbDebugMonitorMode) {
1383 fvhHitsTimeToTriggerSelVsDaq[uGdpb]->Fill(
1384 dMessTime - gdpbv100::kdClockCycleSizeNs * (*itTrigger[uGdpb]).GetFullGdpbTs(), usFirstStarDaqCmd);
1385 fvhHitsTimeToTriggerSelVsTrig[uGdpb]->Fill(
1386 dMessTime - gdpbv100::kdClockCycleSizeNs * (*itTrigger[uGdpb]).GetFullGdpbTs(), usFirstStarTrigCmd);
1387 } // if( fbDebugMonitorMode )
1388 } // if( fbMonitorMode )
1389 } // if( dWinBeg < dMessTime && dMessTime < dWinEnd )
1390
1392 if (dMessTime < dDeadEnd) ++itFirstMessOutOfDeadtime;
1393
1394 ++itMessStart[uGdpb];
1395 } // while( not at buffer end && (not out of deadtime || not out of trigg win ) )
1396 itMessStart[uGdpb] = itFirstMessOutOfDeadtime;
1397
1399 ++itTrigger[uGdpb];
1400 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1401
1403 dMeanTriggerGdpbTs /= fuNrOfGdpbs;
1404 dMeanTriggerStarTs /= fuNrOfGdpbs;
1405 CbmTofStarTrigger2019 meanTrigger(static_cast<ULong64_t>(dMeanTriggerGdpbTs),
1406 static_cast<ULong64_t>(dMeanTriggerStarTs), uFirstStarToken, usFirstStarDaqCmd,
1407 usFirstStarTrigCmd);
1408 starSubEvent.SetTrigger(meanTrigger);
1409
1411 fvEventsBuffer.push_back(starSubEvent);
1412 } // if( kTRUE == bAllSectorsMatch )
1413 else {
1416 std::vector<CbmTofStarTrigger2019>::iterator itEarliestTrigger;
1417 ULong64_t ulEarliestGdpbTs = 0xFFFFFFFFFFFFFFFFUL;
1418 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
1419 if (itTrigger[uGdpb] != fvvBufferTriggers[uGdpb].end()
1420 && ((*itTrigger[uGdpb]).GetFullGdpbTs() < ulEarliestGdpbTs
1421 || 0xFFFFFFFFFFFFFFFFUL == (*itTrigger[uGdpb]).GetFullGdpbTs())) {
1422 itEarliestTrigger = itTrigger[uGdpb];
1423 ulEarliestGdpbTs = (*itTrigger[uGdpb]).GetFullGdpbTs();
1424 } // if( not at end of buffer && (*itTrigger[ uGdpb ]).GetFullGdpbTs() < ulEarliestGdpbTs )
1425 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1426 UInt_t uEarliestStarToken = (*itEarliestTrigger).GetStarToken();
1427 UShort_t usEarliestStarDaqCmd = (*itEarliestTrigger).GetStarDaqCmd();
1428 UShort_t usEarliestStarTrigCmd = (*itEarliestTrigger).GetStarTrigCmd();
1429 ULong64_t ulEarliestFullGdpbTs = (*itEarliestTrigger).GetFullGdpbTs();
1430
1432 UInt_t uNrOfMatchedGdpbs = 0;
1433 std::vector<Bool_t> vbMatchingTrigger(fuNrOfGdpbs, kFALSE);
1434 Double_t dMeanTriggerGdpbTs = 0;
1435 Double_t dMeanTriggerStarTs = 0;
1436 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
1437 if (itTrigger[uGdpb] != fvvBufferTriggers[uGdpb].end()
1438 && (*itTrigger[uGdpb]).GetStarToken() == uEarliestStarToken
1439 && (*itTrigger[uGdpb]).GetStarDaqCmd() == usEarliestStarDaqCmd
1440 && (*itTrigger[uGdpb]).GetStarTrigCmd() == usEarliestStarTrigCmd) {
1441 uNrOfMatchedGdpbs++;
1442 vbMatchingTrigger[uGdpb] = kTRUE;
1443
1445 Long64_t ilTriggerDt =
1446 static_cast<Long64_t>((*itTrigger[uGdpb]).GetFullGdpbTs()) - static_cast<Long64_t>(ulEarliestFullGdpbTs);
1447 fvhTriggerDt[uGdpb]->Fill(ilTriggerDt);
1448 } // if( fbMonitorMode && fbDebugMonitorMode )
1449
1450 dMeanTriggerGdpbTs += (*itTrigger[uGdpb]).GetFullGdpbTs();
1451 dMeanTriggerStarTs += (*itTrigger[uGdpb]).GetFullStarTs();
1452 } // if matching trigger
1453 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1454
1456 dMeanTriggerGdpbTs /= uNrOfMatchedGdpbs;
1457 dMeanTriggerStarTs /= uNrOfMatchedGdpbs;
1458 CbmTofStarTrigger2019 meanTrigger(static_cast<ULong64_t>(dMeanTriggerGdpbTs),
1459 static_cast<ULong64_t>(dMeanTriggerStarTs), (*itEarliestTrigger).GetStarToken(),
1460 (*itEarliestTrigger).GetStarDaqCmd(), (*itEarliestTrigger).GetStarTrigCmd());
1461
1463 CbmTofStarSubevent2019 starSubEvent;
1464 starSubEvent.SetTrigger(meanTrigger);
1465 starSubEvent.SetIncompleteEventFlag();
1466
1468 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
1472 ULong64_t ulTriggerTime = static_cast<ULong64_t>(dMeanTriggerGdpbTs);
1473 if (kTRUE == vbMatchingTrigger[uGdpb]) {
1475 ulTriggerTime = (*itTrigger[uGdpb]).GetFullGdpbTs();
1476
1478 std::vector<gdpbv100::FullMessage> vTrigMess = (*itTrigger[uGdpb]).GetGdpbMessages();
1479 for (std::vector<gdpbv100::FullMessage>::iterator itMess = vTrigMess.begin(); itMess != vTrigMess.end();
1480 ++itMess) {
1481 starSubEvent.AddMsg((*itMess));
1482 } // for( std::vector< gdpbv100::FullMessage >::iterator itMess = vTrigMess.begin(); itMess != vTrigMess.end(); ++ itMess )
1483 } // if( kTRUE == vbMatchingTrigger[ uGdpb ] )
1484 else if (fbMonitorMode && fbDebugMonitorMode) {
1486 fhMissingTriggersEvolution->Fill(dMeanTriggerGdpbTs / 1e9 / 60.0,
1488 } // else if( fbMonitorMode && fbDebugMonitorMode )
1489
1490 if (fbAddStatusToEvent) {
1495 while (fvSectorStatusPattern[uGdpb].end() != itStatUpdtStart[uGdpb]
1496 && ((*itStatUpdtStart[uGdpb]).first) * gdpbv100::kuCoarseCounterSize <= ulTriggerTime) {
1497 ++itStatUpdtStart[uGdpb];
1498 } // while( not at the end of updates vector AND Update time under trigger time )
1500 if (fvSectorStatusPattern[uGdpb].begin() != itStatUpdtStart[uGdpb]) {
1502 --itStatUpdtStart[uGdpb];
1504 uint16_t usGdpbId = fUnpackPar->GetGdpbId(uGdpb);
1505 for (uint32_t uIdx = 0; uIdx < kuNbMsgPerPattern; ++uIdx)
1506 starSubEvent.AddMsg(CreateStatusMessage(usGdpbId, uIdx, (*itStatUpdtStart[uGdpb])));
1507 } // if( fvSectorStatusPattern[ uGdpb ].begin() != itStatUpdtStart[ uGdpb ] )
1508 } // if( fbAddStatusToEvent )
1509
1511 Double_t dWinBeg = gdpbv100::kdClockCycleSizeNs * ulTriggerTime + fdStarTriggerDelay[uGdpb];
1512 Double_t dWinEnd =
1514 Double_t dDeadEnd = gdpbv100::kdClockCycleSizeNs * ulTriggerTime + fdStarTriggerDeadtime[uGdpb];
1515 Double_t dWinBegErrors = dWinBeg - 2 * gdpbv100::kdEpochInNs;
1516
1518 UInt_t uNbErrorsInEventGdpb = 0;
1519 while (itErrorMessStart[uGdpb] != fvvBufferMajorAsicErrors[uGdpb].end()
1520 && (*itErrorMessStart[uGdpb]).GetFullTimeNs() < dWinBegErrors) {
1521 ++itErrorMessStart[uGdpb];
1522 } // while( not at buffer end && (*itErrorMessStart[ uGdpb ]).GetFullTimeNs() < dWinBeg - 2 * gdpbv100::kdEpochInNs )
1523 while (itErrorMessStart[uGdpb] != fvvBufferMajorAsicErrors[uGdpb].end()
1524 && (*itErrorMessStart[uGdpb]).GetFullTimeNs() < dWinEnd) {
1526 if (uNbErrorsInEventGdpb < kuMaxNbErrorsPerGdpbPerEvent) starSubEvent.AddMsg((*itErrorMessStart[uGdpb]));
1527
1528 ++itErrorMessStart[uGdpb];
1529 ++uNbErrorsInEventGdpb;
1530 } // while( not at buffer end && (*itErrorMessStart[ uGdpb ]).GetFullTimeNs() < dWinEnd )
1531
1534 std::vector<gdpbv100::FullMessage>::iterator itFirstMessOutOfDeadtime = itMessStart[uGdpb];
1535 while (
1536 itMessStart[uGdpb] != fvvBufferMessages[uGdpb].end()
1537 && ((*itMessStart[uGdpb]).GetFullTimeNs() < dDeadEnd || (*itMessStart[uGdpb]).GetFullTimeNs() < dWinEnd)) {
1538 Double_t dMessTime = (*itMessStart[uGdpb]).GetFullTimeNs();
1539 /*
1540 if( 5916 == fulCurrentTsIndex )
1541 LOG(info) << Form( "gDPB %2u Miss Mess %12f Start %12f Stop %12f",
1542 uGdpb, dMessTime,
1543 dWinBeg,
1544 dWinEnd )
1545 << ( dWinBeg < dMessTime && dMessTime < dWinEnd ? " IN" : "" );
1546*/
1548 if (fbMonitorMode) {
1549 fvhHitsTimeToTriggerRaw[uGdpb]->Fill(dMessTime - gdpbv100::kdClockCycleSizeNs * ulTriggerTime);
1550 if (fbDebugMonitorMode && gdpbv100::MSG_HIT == (*itMessStart[uGdpb]).getMessageType()) {
1551 UInt_t uTot = (*itMessStart[uGdpb]).getGdpbHit32Tot();
1552 Double_t dTimeDiff = dMessTime - gdpbv100::kdClockCycleSizeNs * ulTriggerTime;
1553 if (fuMinTotPulser < uTot && uTot < fuMaxTotPulser && TMath::Abs(dTimeDiff) < 20e3)
1554 fvhHitsTimeToTriggerRawPulser[uGdpb]->Fill(dTimeDiff, uTot);
1555 } // if( fbDebugMonitorMode && gdpbv100::MSG_HIT == (*itMessStart[ uGdpb ]).getMessageType() )
1556 } // if( fbMonitorMode )
1557
1559 if (dWinBeg < dMessTime && dMessTime < dWinEnd) {
1560 starSubEvent.AddMsg((*itMessStart[uGdpb]));
1562 if (fbMonitorMode) {
1563 fvhHitsTimeToTriggerSel[uGdpb]->Fill(dMessTime - gdpbv100::kdClockCycleSizeNs * ulTriggerTime);
1564 if (fbDebugMonitorMode) {
1565 fvhHitsTimeToTriggerSelVsDaq[uGdpb]->Fill(dMessTime - gdpbv100::kdClockCycleSizeNs * ulTriggerTime,
1566 usEarliestStarDaqCmd);
1567 fvhHitsTimeToTriggerSelVsTrig[uGdpb]->Fill(dMessTime - gdpbv100::kdClockCycleSizeNs * ulTriggerTime,
1568 usEarliestStarTrigCmd);
1569 } // if( fbDebugMonitorMode )
1570 } // if( fbMonitorMode )
1571 } // if( dWinBeg < dMessTime && dMessTime < dWinEnd )
1572
1574 if (dMessTime < dDeadEnd) ++itFirstMessOutOfDeadtime;
1575
1576 ++itMessStart[uGdpb];
1577 } // while( not at buffer end && (not out of deadtime || not out of trigg win ) )
1578 itMessStart[uGdpb] = itFirstMessOutOfDeadtime;
1579
1581 if (kTRUE == vbMatchingTrigger[uGdpb]) ++itTrigger[uGdpb];
1582 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1583
1585 fvEventsBuffer.push_back(starSubEvent);
1586 } // else of if( kTRUE == bAllSectorsMatch )
1587
1589 bAllSectAllTriggDone = kTRUE;
1590 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
1592 if (fvvBufferTriggers[uGdpb].end() != itTrigger[uGdpb]) {
1593 bAllSectAllTriggDone = kFALSE;
1594 break;
1595 } // if( fvvBufferTriggers[ uGdpb ].end() != (*itTrigger[ uGdpb ]) )
1596 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1597 } // while( kFALSE == bAllSectAllTriggDone )
1598
1599 return kTRUE;
1600}
1601
1602// -------------------------------------------------------------------------
1603
1605{
1607 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
1608 UInt_t uSector = fUnpackPar->GetGdpbToSectorOffset() + uGdpb;
1609 std::string sFolder = Form("sector%2u", uSector);
1610
1611 LOG(info) << "gDPB " << uGdpb << " is " << sFolder;
1612
1613 fvhHitsTimeToTriggerRaw.push_back(new TH1D(Form("hHitsTimeToTriggerRawSect%2u", uSector),
1614 Form("Time to trigger for all neighboring hits in sector %2u; t "
1615 "- Ttrigg [ns]; Hits []",
1616 uSector),
1617 2000, -5000, 5000));
1618
1619 UInt_t uNbBinsDtSel = fdStarTriggerWinSize[uGdpb];
1620 Double_t dMaxDtSel = fdStarTriggerDelay[uGdpb] + fdStarTriggerWinSize[uGdpb];
1621 fvhHitsTimeToTriggerSel.push_back(new TH1D(Form("hHitsTimeToTriggerSelSect%2u", uSector),
1622 Form("Time to trigger for all selected hits in sector %2u; t - "
1623 "Ttrigg [ns]; Hits []",
1624 uSector),
1625 uNbBinsDtSel, fdStarTriggerDelay[uGdpb], dMaxDtSel));
1626
1630
1631 if (kTRUE == fbDebugMonitorMode) {
1633 new TH2D(Form("hHitsTimeToTriggerRawPulserSect%2u", uSector),
1634 Form("Time to trigger for all neighboring hits within pulser TOT range "
1635 "in sector %2u; t - Ttrigg [ns]; TOT [bins]; Hits []",
1636 uSector),
1637 2000, -5000, 5000, 256, 0, 256));
1638
1639 fvhHitsTimeToTriggerSelVsDaq.push_back(new TH2D(Form("hHitsTimeToTriggerSelVsDaqSect%2u", uSector),
1640 Form("Time to trigger for all selected hits vs DAQ CMD in "
1641 "sector %2u; t - Ttrigg [ns]; DAQ CMD []; Hits []",
1642 uSector),
1643 uNbBinsDtSel, fdStarTriggerDelay[uGdpb], dMaxDtSel, 16, 0., 16.));
1644
1645 fvhHitsTimeToTriggerSelVsTrig.push_back(new TH2D(Form("hHitsTimeToTriggerSelVsTrigSect%2u", uSector),
1646 Form("Time to trigger for all selected hits vs TRIG CMD in "
1647 "sector %2u; t - Ttrigg [ns]; TRIG CMD []; Hits []",
1648 uSector),
1649 uNbBinsDtSel, fdStarTriggerDelay[uGdpb], dMaxDtSel, 16, 0.,
1650 16.));
1651
1652 fvhTriggerDt.push_back(new TH1I(Form("hTriggerDtSect%2u", uSector),
1653 Form("Trigger time difference between sector %2u and the first sector, "
1654 "full events only; Ttrigg%2u - TtriggRef [Clk]; events []",
1655 uSector, uSector),
1656 200, -100, 100));
1657
1660 UInt_t uNbBinsInTs = fdMsSizeInNs * 111 / 1000. / 10.;
1661 UInt_t uNbBinsInMs = fdMsSizeInNs * 20 / 1000. / 10.;
1662
1663 fvhTriggerDistributionInTs.push_back(new TH1I(Form("hTriggerDistInTsSect%2u", uSector),
1664 Form("Trigger distribution inside TS in sector %2u; Time in "
1665 "TS [us]; Trigger [];",
1666 uSector),
1667 uNbBinsInTs, -0.5 - fdMsSizeInNs * 10 / 1000.,
1668 fdMsSizeInNs * 101 / 1000. - 0.5));
1669
1670 fvhTriggerDistributionInMs.push_back(new TH1I(Form("hTriggerDistInMsSect%2u", uSector),
1671 Form("Trigger distribution inside MS in sector %2u; Time in "
1672 "MS [us]; Trigger [];",
1673 uSector),
1674 uNbBinsInMs, -0.5 - fdMsSizeInNs * 10 / 1000.,
1675 fdMsSizeInNs * 10 / 1000. - 0.5));
1676
1677 fvhMessDistributionInMs.push_back(new TH1I(Form("hMessDistInMsSect%2u", uSector),
1678 Form("Messages distribution inside MS in sector %2u; Time in "
1679 "MS [us]; Trigger [];",
1680 uSector),
1681 uNbBinsInMs, -0.5 - fdMsSizeInNs * 10 / 1000.,
1682 fdMsSizeInNs * 10 / 1000. - 0.5));
1683
1688 AddHistoToVector(fvhTriggerDt[uGdpb], sFolder);
1692 } // if( kTRUE == fbDebugMonitorMode )
1693 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1694
1696 fhEventNbPerTs = new TH1I("hEventNbPerTs", "Number of Events per TS; Events []; TS []", 1000, 0, 1000);
1697
1699 new TH1I("hEventSizeDistribution", "Event size distribution; Event size [byte]; Events []",
1701
1702 fhEventSizeEvolution = new TProfile(
1703 "hEventSizeEvolution", "Event size evolution; Time in run [min]; mean Event size [byte];", 14400, 0, 14400);
1704
1706 new TH1I("hEventNbEvolution", "Event number evolution; Time in run [min]; Events [];", 14400, 0, 14400);
1707
1709 AddHistoToVector(fhEventNbPerTs, "eventbuilder");
1711 AddHistoToVector(fhEventSizeEvolution, "eventbuilder");
1712 AddHistoToVector(fhEventNbEvolution, "eventbuilder");
1713
1714 if (kTRUE == fbDebugMonitorMode) {
1717 UInt_t uNbBinsInTs = fdMsSizeInNs * 101 / 1000. / 10.;
1718
1720 new TH1I("hEventNbDistributionInTs", "Event number distribution inside TS; Time in TS [us]; Events [];",
1721 uNbBinsInTs, -0.5, fdMsSizeInNs * 101 / 1000. - 0.5);
1722
1723 fhEventSizeDistributionInTs = new TProfile("hEventSizeDistributionInTs",
1724 "Event size distribution inside TS; Time in TS [us]; mean "
1725 "Event size [Byte];",
1726 uNbBinsInTs, -0.5, fdMsSizeInNs * 101 / 1000. - 0.5);
1727
1728 fhRawTriggersStats = new TH2I("hRawTriggersStats", "Raw triggers statistics per sector; ; Sector []; Messages []",
1729 5, 0, 5, 12, 13, 25);
1730 fhRawTriggersStats->GetXaxis()->SetBinLabel(1, "A");
1731 fhRawTriggersStats->GetXaxis()->SetBinLabel(2, "B");
1732 fhRawTriggersStats->GetXaxis()->SetBinLabel(3, "C");
1733 fhRawTriggersStats->GetXaxis()->SetBinLabel(4, "D");
1734 fhRawTriggersStats->GetXaxis()->SetBinLabel(5, "F");
1735
1737 new TH2I("hRawTriggersStatsCore", "Raw triggers in Core MS statistics per sector; ; Sector []; Messages []", 5, 0,
1738 5, 12, 13, 25);
1739 fhRawTriggersStatsCore->GetXaxis()->SetBinLabel(1, "A");
1740 fhRawTriggersStatsCore->GetXaxis()->SetBinLabel(2, "B");
1741 fhRawTriggersStatsCore->GetXaxis()->SetBinLabel(3, "C");
1742 fhRawTriggersStatsCore->GetXaxis()->SetBinLabel(4, "D");
1743 fhRawTriggersStatsCore->GetXaxis()->SetBinLabel(5, "F");
1744
1745 fhRawTriggersStatsOver = new TH2I("hRawTriggersStatsOver",
1746 "Raw triggers in Overlap MS statistics "
1747 "per sector; ; Sector []; Messages []",
1748 5, 0, 5, 12, 13, 25);
1749 fhRawTriggersStatsOver->GetXaxis()->SetBinLabel(1, "A");
1750 fhRawTriggersStatsOver->GetXaxis()->SetBinLabel(2, "B");
1751 fhRawTriggersStatsOver->GetXaxis()->SetBinLabel(3, "C");
1752 fhRawTriggersStatsOver->GetXaxis()->SetBinLabel(4, "D");
1753 fhRawTriggersStatsOver->GetXaxis()->SetBinLabel(5, "F");
1754
1755 fhRawTriggersStatsSel = new TH2I(
1756 "hRawTriggersStatsSel", "Selected triggers statistics per sector; ; Sector []; Messages []", 3, 0, 3, 12, 13, 25);
1757 fhRawTriggersStatsSel->GetXaxis()->SetBinLabel(1, "All");
1758 fhRawTriggersStatsSel->GetXaxis()->SetBinLabel(2, "Core");
1759 fhRawTriggersStatsSel->GetXaxis()->SetBinLabel(3, "Over");
1760
1761 fhMissingTriggersEvolution = new TH2I("hMissingTriggersEvolution",
1762 "Missing trigger counts per sector vs time in run; Time in run "
1763 "[min]; Sector []; Missing triggers []",
1764 14400, 0, 14400, 12, 13, 25);
1765
1769 AddHistoToVector(fhRawTriggersStats, "eventbuilder");
1772 AddHistoToVector(fhRawTriggersStatsSel, "eventbuilder");
1774 } // if( kTRUE == fbDebugMonitorMode )
1775
1777 Double_t w = 10;
1778 Double_t h = 10;
1779
1781 fcTimeToTrigRaw = new TCanvas("cTimeToTrigRaw", "Raw Time to trig for all sectors", w, h);
1782 fcTimeToTrigRaw->Divide(2, fuNrOfGdpbs / 2);
1783 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
1784 fcTimeToTrigRaw->cd(1 + uGdpb);
1785 gPad->SetGridx();
1786 gPad->SetGridy();
1787 gPad->SetLogy();
1788 fvhHitsTimeToTriggerRaw[uGdpb]->Draw();
1789 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1791
1793 fcTimeToTrigSel = new TCanvas("cTimeToTrigSel", "Selected Time to trig for all sectors", w, h);
1794 fcTimeToTrigSel->Divide(2, fuNrOfGdpbs / 2);
1795 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
1796 fcTimeToTrigSel->cd(1 + uGdpb);
1797 gPad->SetGridx();
1798 gPad->SetGridy();
1799 gPad->SetLogy();
1800 fvhHitsTimeToTriggerSel[uGdpb]->Draw();
1801 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1803
1804 if (kTRUE == fbDebugMonitorMode) {
1806 fcTrigDistMs = new TCanvas("cTrigDistMs", "Trigger time to MS start for all sectors", w, h);
1807 fcTrigDistMs->Divide(2, fuNrOfGdpbs / 2);
1808 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
1809 fcTrigDistMs->cd(1 + uGdpb);
1810 gPad->SetGridx();
1811 gPad->SetGridy();
1812 gPad->SetLogy();
1813 fvhTriggerDistributionInMs[uGdpb]->Draw();
1814 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1815 AddCanvasToVector(fcTrigDistMs, "canvases");
1816
1818 fcMessDistMs = new TCanvas("cMessDistMs", "Message time to MS start for all sectors", w, h);
1819 fcMessDistMs->Divide(2, fuNrOfGdpbs / 2);
1820 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
1821 fcMessDistMs->cd(1 + uGdpb);
1822 gPad->SetGridx();
1823 gPad->SetGridy();
1824 gPad->SetLogy();
1825 fvhMessDistributionInMs[uGdpb]->Draw();
1826 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1827 AddCanvasToVector(fcMessDistMs, "canvases");
1828
1829 fcTriggerStats = new TCanvas("cTriggerStats", "Trigger statistics per sector", w, h);
1830 fcTriggerStats->Divide(2, 2);
1831
1832 fcTriggerStats->cd(1);
1833 gPad->SetGridx();
1834 gPad->SetGridy();
1835 gPad->SetLogz();
1836 fhRawTriggersStats->Draw("colz text");
1837
1838 fcTriggerStats->cd(2);
1839 gPad->SetGridx();
1840 gPad->SetGridy();
1841 gPad->SetLogz();
1842 fhRawTriggersStatsCore->Draw("colz text");
1843
1844 fcTriggerStats->cd(3);
1845 gPad->SetGridx();
1846 gPad->SetGridy();
1847 gPad->SetLogz();
1848 fhRawTriggersStatsOver->Draw("colz text");
1849
1850 fcTriggerStats->cd(4);
1851 gPad->SetGridx();
1852 gPad->SetGridy();
1853 gPad->SetLogz();
1854 fhRawTriggersStatsSel->Draw("colz text");
1855
1856 AddCanvasToVector(fcTriggerStats, "canvases");
1857 } // if( kTRUE == fbDebugMonitorMode )
1858
1860 fcEventBuildStats = new TCanvas("cEvtBuildStats", "Event building statistics", w, h);
1861 if (kTRUE == fbDebugMonitorMode) fcEventBuildStats->Divide(2, 3);
1862 else
1863 fcEventBuildStats->Divide(2, 2);
1864
1865 fcEventBuildStats->cd(1);
1866 gPad->SetGridx();
1867 gPad->SetGridy();
1868 gPad->SetLogy();
1869 fhEventNbPerTs->Draw();
1870
1871 fcEventBuildStats->cd(2);
1872 gPad->SetGridx();
1873 gPad->SetGridy();
1874 gPad->SetLogy();
1876
1877 fcEventBuildStats->cd(3);
1878 gPad->SetGridx();
1879 gPad->SetGridy();
1880 gPad->SetLogy();
1881 fhEventSizeEvolution->Draw();
1882
1883 fcEventBuildStats->cd(4);
1884 gPad->SetGridx();
1885 gPad->SetGridy();
1886 gPad->SetLogy();
1887 fhEventNbEvolution->Draw();
1888
1889 if (kTRUE == fbDebugMonitorMode) {
1890 fcEventBuildStats->cd(5);
1891 gPad->SetGridx();
1892 gPad->SetGridy();
1893 gPad->SetLogy();
1895
1896 fcEventBuildStats->cd(6);
1897 gPad->SetGridx();
1898 gPad->SetGridy();
1899 gPad->SetLogy();
1901 } // if( kTRUE == fbDebugMonitorMode )
1902
1904
1905 return kTRUE;
1906}
1908{
1909 UInt_t uNbEvents = fvEventsBuffer.size();
1910 fhEventNbPerTs->Fill(uNbEvents);
1911
1912 for (UInt_t uEvent = 0; uEvent < uNbEvents; ++uEvent) {
1913 UInt_t uEventSize = fvEventsBuffer[uEvent].GetEventSize();
1914 Double_t dEventTimeSec = fvEventsBuffer[uEvent].GetEventTimeSec();
1915 Double_t dEventTimeMin = dEventTimeSec / 60.0;
1916
1917 fhEventSizeDistribution->Fill(uEventSize);
1918 if (fbDebugMonitorMode) fhEventSizeEvolution->Fill(dEventTimeSec, uEventSize);
1919 else
1920 fhEventSizeEvolution->Fill(dEventTimeMin, uEventSize);
1921 if (fbDebugMonitorMode) fhEventNbEvolution->Fill(dEventTimeSec);
1922 else
1923 fhEventNbEvolution->Fill(dEventTimeMin);
1924
1925 if (kTRUE == fbDebugMonitorMode) {
1926 Double_t dEventTimeInTs =
1927 (fvEventsBuffer[uEvent].GetTrigger().GetFullGdpbTs() * gdpbv100::kdClockCycleSizeNs - fdTsStartTime) / 1000.0;
1928
1929 fhEventNbDistributionInTs->Fill(dEventTimeInTs);
1930 fhEventSizeDistributionInTs->Fill(dEventTimeInTs, uEventSize);
1931 } // if( kTRUE == fbDebugMonitorMode )
1932 } // for( UInt_t uEvent = 0; uEvent < uNbEvents; ++uEvent )
1933
1934 return kTRUE;
1935}
1937{
1938 for (UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb) {
1939 fvhHitsTimeToTriggerRaw[uGdpb]->Reset();
1940 fvhHitsTimeToTriggerSel[uGdpb]->Reset();
1941
1942 if (kTRUE == fbDebugMonitorMode) {
1943 fvhHitsTimeToTriggerRawPulser[uGdpb]->Reset();
1944 fvhHitsTimeToTriggerSelVsDaq[uGdpb]->Reset();
1945 fvhHitsTimeToTriggerSelVsTrig[uGdpb]->Reset();
1946 fvhTriggerDt[uGdpb]->Reset();
1947 fvhTriggerDistributionInTs[uGdpb]->Reset();
1948 fvhTriggerDistributionInMs[uGdpb]->Reset();
1949 fvhMessDistributionInMs[uGdpb]->Reset();
1950 } // if( kTRUE == fbDebugMonitorMode )
1951 } // for( UInt_t uGdpb = 0; uGdpb < fuNrOfGdpbs; ++uGdpb )
1952
1954 fhEventNbPerTs->Reset();
1955 fhEventSizeDistribution->Reset();
1956 fhEventSizeEvolution->Reset();
1957 fhEventNbEvolution->Reset();
1958
1959 if (kTRUE == fbDebugMonitorMode) {
1962 fhRawTriggersStats->Reset();
1963 fhRawTriggersStatsCore->Reset();
1964 fhRawTriggersStatsOver->Reset();
1965 fhRawTriggersStatsSel->Reset();
1967 } // if( kTRUE == fbDebugMonitorMode )
1968
1969 return kTRUE;
1970}
1971// -------------------------------------------------------------------------
Histogram manager.
static Int_t iWarnMess
static constexpr size_t size()
Definition KfSimdPseudo.h:2
bool first
void AddHistoToVector(TNamed *pointer, std::string sFolder="")
std::vector< size_t > fvMsComponentsList
void AddCanvasToVector(TCanvas *pointer, std::string sFolder="")
std::vector< ULong64_t > fvulStarTsMsb
[sector]
std::vector< std::vector< std::pair< uint64_t, std::bitset< kuNbAsicPerGdpb > > > > fvSectorStatusPattern
[sector][asic] Exclude from dictionnary
std::vector< TH2 * > fvhHitsTimeToTriggerSelVsDaq
[sector]
std::vector< std::bitset< kuNbAsicPerGdpb > > fvbGdpbLastResyncPattern
[sector][asic] Exclude from dictionnary
std::vector< UInt_t > fvuStarDaqCmdLast
[sector]
Bool_t ProcessMs(const fles::Timeslice &ts, size_t uMsCompIdx, size_t uMsIdx)
std::vector< UInt_t > fvuStarTokenLastCore
[sector]
std::vector< ULong64_t > fvulGdpbTsFullLastCore
[sector]
std::vector< UInt_t > fvuStarTrigCmdLastCore
[sector]
std::vector< TH1 * > fvhHitsTimeToTriggerSel
[sector], extra monitor for debug
void ProcessSlCtrl(gdpbv100::Message mess, uint64_t ulCurEpochGdpbGet4)
std::vector< TH1 * > fvhTriggerDistributionInMs
[sector], extra monitor for debug
std::vector< std::vector< gdpbv100::FullMessage > > fvvBufferMessages
[sector], buffer to make sure GET4 errors 0-12 are always transmitted
std::vector< TH1 * > fvhTriggerDt
[sector], extra monitor for debug
std::vector< TH1 * > fvhHitsTimeToTriggerRaw
[sector][Update]<time, [asic]> Exclude from dictionnary
std::vector< std::vector< gdpbv100::Message > > fvvmEpSupprBuffer
Epoch + Epoch Cycle.
Bool_t fbAddStatusToEvent
Switch ON the insertion of the LostEvent messages from GET4s with the critical errors (default is fal...
std::vector< ULong64_t > fvulStarTsFullLastCore
[sector]
gdpbv100::FullMessage CreateStatusMessage(uint16_t uGdpbId, uint32_t uIndex, std::pair< uint64_t, std::bitset< kuNbAsicPerGdpb > > statusIn)
std::vector< Double_t > fdStarTriggerDelay
[sector]
Bool_t fbDebugMonitorMode
Switch ON the filling of a minimal set of histograms.
TCanvas * fcMessDistMs
All sectors, extra monitor for debug.
TH2 * fhRawTriggersStatsCore
extra monitor for debug
std::vector< TH1 * > fvhTriggerDistributionInTs
[sector], extra monitor for debug
std::vector< std::vector< gdpbv100::FullMessage > > fvvBufferMajorAsicErrors
[sector]
std::vector< std::bitset< kuNbAsicPerGdpb > > fvbGdpbLastEnablePattern
[sector][asic] Exclude from dictionnary
std::vector< ULong64_t > fvulStarTsFullLast
[sector]
std::vector< TH2 * > fvhHitsTimeToTriggerSelVsTrig
[sector], extra monitor for debug
std::vector< Double_t > fvdTrigCandidateTimeStart
[sector]
TCanvas * fcEventBuildStats
All sectors, extra monitor for debug.
std::vector< Double_t > fvdMessCandidateTimeStart
[sector]
TH2 * fhRawTriggersStatsOver
extra monitor for debug
Bool_t fbStoreLostEventMsg
Switch ON the filling of a additional set of histograms.
void AddMsComponentToList(size_t component, UShort_t usDetectorId)
void StoreMessageInBuffer(gdpbv100::FullMessage fullMess, uint32_t uGdpbNr)
std::vector< std::vector< CbmTofStarTrigger2019 > > fvvBufferTriggers
[sector]
void ProcessEpochCycle(uint64_t ulCycleData)
All sectors, extra monitor for debug.
std::vector< std::bitset< kuNbAsicPerGdpb > > fvbGdpbLastMissmatchPattern
std::vector< ULong64_t > fvulCurrentEpochFull
Epoch cycle from the Ms Start message and Epoch counter flip.
void ProcessSysMess(gdpbv100::Message mess, uint64_t ulCurEpochGdpbGet4)
void ProcessHit(gdpbv100::Message mess, uint64_t ulCurEpochGdpbGet4)
std::vector< ULong64_t > fvulGdpbTsFullLast
[sector]
void ProcessPattern(gdpbv100::Message mess, uint64_t ulCurEpochGdpbGet4)
Bool_t ProcessTs(const fles::Timeslice &ts)
std::vector< TH2 * > fvhHitsTimeToTriggerRawPulser
[sector]
TProfile * fhEventSizeDistributionInTs
extra monitor for debug
std::vector< Double_t > fvdMessCandidateTimeStop
[sector]
std::vector< ULong64_t > fvulGdpbTsLsb
[sector]
TH1 * fhEventNbPerTs
[sector], extra monitor for debug
std::vector< Double_t > fvdTrigCandidateTimeStop
[sector]
std::vector< UInt_t > fvuStarTrigCmdLast
[sector]
std::vector< TH1 * > fvhMessDistributionInMs
[sector], extra monitor for debug
TH2 * fhMissingTriggersEvolution
extra monitor for debug
Double_t fdAllowedTriggersSpread
Event window limits.
std::vector< Double_t > fdStarTriggerWinSize
[sector]
static const UInt_t kuMaxNbErrorsPerGdpbPerEvent
Switch ON the readout and insertion of STATUS pattern message (default is true)
std::vector< CbmTofStarSubevent2019 > fvEventsBuffer
[sector]
std::vector< ULong64_t > fvulCurrentEpoch
Current time references for each GDPB: merged epoch marker, epoch cycle, full epoch [fuNrOfGdpbs].
TH2 * fhRawTriggersStatsSel
extra monitor for debug
TCanvas * fcTimeToTrigRaw
extra monitor for debug
CbmStar2019TofPar * fUnpackPar
Correspond to ~6000 error messages max per event, leaving 2000 for hits and epoch.
static const uint32_t kuNbMsgPerPattern
ASIC status monitoring.
std::vector< ULong64_t > fvulStarTsMid
[sector]
std::vector< UInt_t > fvuStarDaqCmdLastCore
[sector]
Int_t GetRpcType(UInt_t uGbtx)
static constexpr UInt_t GetNrOfGet4PerGdpb()
static constexpr UInt_t GetNrOfFeePerGdpb()
Int_t GetModuleId(UInt_t uGbtx)
Int_t GetNrOfRpc(UInt_t uGbtx)
Double_t GetStarTriggDeadtime(UInt_t uGdpb)
Int_t GetGdpbId(Int_t i)
Double_t GetStarTriggWinSize(UInt_t uGdpb)
static constexpr UInt_t GetNbByteMessage()
Int_t ElinkIdxToGet4Idx(UInt_t uElink)
static constexpr UInt_t GetGdpbToSectorOffset()
Double_t GetStarTriggAllowedSpread()
static constexpr UInt_t GetNrOfChannelsPerGet4()
static constexpr UInt_t GetNrOfGet4PerFee()
Int_t GetRpcSide(UInt_t uGbtx)
Double_t GetStarTriggDelay(UInt_t uGdpb)
Data class with information on a STS local track.
void SetTrigger(CbmTofStarTrigger2019 triggerIn)
void SetIncompleteEventFlag(Bool_t bFlagState=kTRUE)
void AddMsg(const gdpbv100::FullMessage &msgIn)
static uint32_t GetMaxOutputSize()
ULong64_t GetFullGdpbTs() const
double GetFullTimeNs() const
uint16_t getGdpbHitIs24b() const
uint16_t getStarTrigMsgIndex() const
uint32_t getGdpbEpEpochNb() const
uint16_t getGdpbSysPattType() const
uint64_t getStarTsMidStarC() const
uint16_t getGdpbSysSubType() const
uint64_t getGdpbTsLsbStarB() const
void setGdpbSysPattType(uint16_t v)
uint16_t getGdpbSysPattIndex() const
void setGdpbEpEpochNb(uint32_t v)
uint64_t getStarTsMsbStarB() const
uint32_t getStarTrigCmdStarD() const
uint32_t getGdpbSysPattPattern() const
bool isStarTrigger() const
Returns true is message type is MSG_STAR_TRI_A, _B, _C, _D (STAR Trigger message)
uint32_t getStarDaqCmdStarD() const
void setGdpbSysPattPattern(uint32_t v)
void setMessageType(uint8_t v)
Sets the message type field in the current message.
uint64_t getGdpbTsMsbStarA() const
uint16_t getGdpbGenChipId() const
uint16_t getGdpbSysErrUnused() const
void setGdpbSysPattIndex(uint16_t v)
uint32_t getStarTokenStarD() const
void setGdpbGenGdpbId(uint32_t v)
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
bool getGdpbSysErrEdge() const
uint16_t getGdpbSysErrChanId() const
void setGdpbSysSubType(uint16_t v)
uint16_t getGdpbHitChanId() const
const double kdClockCycleSizeNs
@ GET4_V2X_ERR_LOST_EVT
@ GET4_V2X_ERR_TOT_RANGE
@ GET4_V2X_ERR_ADD_RIS_EDG
@ GET4_V2X_ERR_TOT_OVERWRT
@ GET4_V2X_ERR_EVT_DISCARD
@ GET4_V2X_ERR_UNPAIR_FALL
@ GET4_V2X_ERR_SEQUENCE_ER
const uint32_t kuChipIdMergedEpoch
const uint32_t kuEpochCounterSz
const double kdEpochInNs
const uint64_t kulEpochCycleFieldSz
const uint32_t kuCoarseCounterSize
@ SYS_GET4_SYNC_MISS