CbmRoot
Loading...
Searching...
No Matches
CbmAlgoBuildRawEvents.cxx
Go to the documentation of this file.
1/* Copyright (C) 2020-2024 Facility for Antiproton and Ion Research in Europe, Darmstadt
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Pierre-Alain Loizeau [committer], Dominik Smith, Alexandru Bercuci */
4
6
8#include "CbmBmonDigi.h"
9#include "CbmEvent.h"
10#include "CbmFsdDigi.h"
11#include "CbmMuchBeamTimeDigi.h"
12#include "CbmMuchDigi.h"
13#include "CbmPsdDigi.h"
14#include "CbmRichDigi.h"
15#include "CbmStsDigi.h"
16#include "CbmTofDigi.h"
17#include "CbmTrdAddress.h"
18#include "CbmTrdDigi.h"
19#include "TimesliceMetaData.h"
20
22#include <FairRootManager.h>
23#include <FairRunOnline.h>
24#include <Logger.h>
25
27#include <TCanvas.h>
28#include <TClonesArray.h>
29#include <TDirectoryFile.h>
30#include <TH1.h>
31#include <TH2.h>
32#include <THttpServer.h>
33#include <TProfile.h>
34#include <TStopwatch.h>
35
37#include <algorithm>
38
39#define VERBOSE 0
40
41template<>
43
45{
46 LOG(info) << "CbmAlgoBuildRawEvents::InitAlgo => Starting sequence";
47
48 if (fbGetTimings) {
49 fTimer = new TStopwatch;
50 fTimer->Start();
51 }
52
56 if (fSeedTimes == nullptr) {
57 LOG(fatal) << "No reference detector set and no seed times supplied, stopping there!";
58 }
59 }
60 else {
61 if (fSeedTimes != nullptr) {
62 LOG(fatal) << "Cannot have explicit seed times and reference detector, stopping there!";
63 }
64 if (kFALSE == CheckDataAvailable(fRefDet)) {
65 LOG(fatal) << "Reference detector set but no digi input found, stopping there!";
66 }
67 }
68
70 for (std::vector<RawEventBuilderDetector>::iterator det = fvDets.begin(); det != fvDets.end(); ++det) {
71 if (kFALSE == CheckDataAvailable(*det)) {
72 LOG(fatal) << "No digi input for one of selection detector, stopping there!";
73 }
74 }
75
77 if (fbUseTsMetaData) {
79 LOG(fatal) << "No TS metadata input found"
80 << " => Please check in the unpacking macro if the following line was "
81 "present!"
82 << std::endl
83 << "source->SetWriteOutputFlag(kTRUE); // For writing TS metadata";
84 }
85 }
86 if (fbFillHistos) {
88 }
89 if (fTimer != nullptr) {
90 fTimer->Stop();
91 Double_t rtime = fTimer->RealTime();
92 Double_t ctime = fTimer->CpuTime();
93 LOG(info) << "CbmAlgoBuildRawEvents::Init(): Real time " << rtime << " s, CPU time " << ctime << " s";
94 }
95
96 LOG(info) << "CbmAlgoBuildRawEvents::InitAlgo => Done";
97 return kTRUE;
98}
99
101{
102 if (fbGetTimings) {
103 PrintTimings();
104 }
105}
106
108{
109 if (fTimer == nullptr) {
110 LOG(fatal) << "Trying to print timings but timer not set";
111 }
112 else {
113 Double_t rtime = fTimer->RealTime();
114 Double_t ctime = fTimer->CpuTime();
115 LOG(info) << "CbmAlgoBuildRawEvents: Real time " << rtime << " s, CPU time " << ctime << " s";
116 }
117}
118
120{
122 int counter = 0;
123 for (CbmEvent* event : fEventVector) {
124 LOG(debug) << "Event " << counter << " has " << event->GetNofData() << " digis";
125 delete event;
126 counter++;
127 }
128 fEventVector.clear();
129}
130
132{
133 LOG_IF(info, fuNrTs % 1000 == 0) << "Begin of TS " << fuNrTs;
134 TStopwatch timerTs;
135 timerTs.Start();
136
137 if (fTimer != nullptr) {
138 fTimer->Start(kFALSE);
139 }
140 InitTs();
142 BuildEvents();
143
145 if (nullptr != fCurrentEvent) {
147 // fCurrentEvent->SetStartTime( fPrevTime ); // Replace Seed time with time of first digi in event?
149 if (fbBmonInUse) {
152 }
153 fEventVector.push_back(fCurrentEvent);
154 fuCurEv++;
155
157 fCurrentEvent = nullptr;
158 }
159
160 if (fbFillHistos) {
161 timerTs.Stop();
162 fhCpuTimePerTs->Fill(fuNrTs, timerTs.CpuTime() * 1000.);
163 fhRealTimePerTs->Fill(fuNrTs, timerTs.RealTime() * 1000.);
164 timerTs.Start();
165 }
166
167 LOG(debug) << "Found " << fEventVector.size() << " triggered events";
168 if (fbFillHistos) {
169 FillHistos();
170 }
171 if (fTimer != nullptr) {
172 fTimer->Stop();
173 }
174
175 if (fbFillHistos) {
176 timerTs.Stop();
177 fhCpuTimePerTsHist->Fill(fuNrTs, timerTs.CpuTime() * 1000.);
178 fhRealTimePerTsHist->Fill(fuNrTs, timerTs.RealTime() * 1000.);
179 }
180
181 fuNrTs++;
182}
183
185{
187 fuCurEv = 0;
188
193 }
195 for (std::vector<RawEventBuilderDetector>::iterator det = fvDets.begin(); det != fvDets.end(); ++det) {
196 (*det).fuStartIndex = 0;
197 (*det).fuEndIndex = 0;
198 }
199}
200
202{
204 Double_t dTsStartTime = fdTsStartTime;
205 Double_t dOverlapStart = fdTsStartTime + fdTsLength;
206 Double_t dOverlapSize = fdTsOverLength;
207
208 if (fbUseTsMetaData) {
209 const TimesliceMetaData* pTsMetaData = dynamic_cast<TimesliceMetaData*>(fTimeSliceMetaDataArray->At(0));
210 if (nullptr == pTsMetaData)
211 LOG(fatal) << Form("CbmAlgoBuildRawEvents::LoopOnSeeds => "
212 "No TS metadata found for TS %6u.",
213 fuNrTs);
214 dTsStartTime = pTsMetaData->GetStartTime();
215 dOverlapStart = pTsMetaData->GetOverlapStartTime();
216 dOverlapSize = pTsMetaData->GetOverlapDuration();
217 }
218
220 if ((0.0 < fdEarliestTimeWinBeg && dOverlapSize < fdLatestTimeWinEnd) || (dOverlapSize < fdWidestTimeWinRange)) {
221 LOG(warning) << "CbmAlgoBuildRawEvents::LoopOnSeeds => "
222 << Form("Event window not fitting in TS overlap, risk of "
223 "incomplete events: %f %f %f %f",
225 } // if end of event window does not fit in overlap for a seed at edge of TS core
226
229 if (fbIgnoreTsOverlap) {
231 fdSeedWindowEnd = dOverlapStart;
232 }
233 else {
235 fdSeedWindowEnd = dOverlapStart + (0.0 < fdEarliestTimeWinBeg ? 0.0 : -fdEarliestTimeWinBeg);
236 }
237}
238
240{
242 switch (fRefDet.detId) {
243 case ECbmModuleId::kSts: {
245 break;
246 }
247 case ECbmModuleId::kMuch: {
250 }
251 else {
253 }
254 break;
255 }
256 case ECbmModuleId::kTrd2d: // Same data storage as trd 1d
257 case ECbmModuleId::kTrd: {
259 break;
260 }
261 case ECbmModuleId::kTof: {
263 break;
264 }
265 case ECbmModuleId::kRich: {
267 break;
268 }
269 case ECbmModuleId::kPsd: {
271 break;
272 }
273 case ECbmModuleId::kFsd: {
275 break;
276 }
277 case ECbmModuleId::kBmon: {
279 break;
280 }
281 case ECbmModuleId::kNotExist: { //explicit seed times
283 break;
284 }
285 default: {
286 LOG(fatal) << "CbmAlgoBuildRawEvents::BuildEvents => "
287 << "Trying to search event seeds with unsupported det: " << fRefDet.sName;
288 break;
289 }
290 }
291}
292
293template<>
295{
297 const UInt_t uNbSeeds = fSeedTimes->size();
299 for (UInt_t uSeed = 0; uSeed < uNbSeeds; ++uSeed) {
300 LOG(debug) << Form("Checking seed %6u / %6u", uSeed, uNbSeeds);
301 Double_t dTime = fSeedTimes->at(uSeed);
302
304 if (dTime < fdSeedWindowBeg) {
305 LOG(debug4) << Form("CbmAlgoBuildRawEvents :: reject for time %f < %f\n", dTime, fdSeedWindowBeg);
306 continue;
307 }
308 else if (fdSeedWindowEnd < dTime) {
309 break;
310 }
312 CheckSeed(dTime, uSeed);
313 }
314 }
315 else {
316 LOG(fatal) << "Trying to read explicit seeds while reference detector is set.";
317 }
318}
319
320template<class DigiSeed>
322{
323 const UInt_t uNbRefDigis = GetNofDigis(fRefDet.detId);
325 for (UInt_t uDigi = 0; uDigi < uNbRefDigis; ++uDigi) {
326 LOG(debug) << Form("Checking seed %6u / %6u", uDigi, uNbRefDigis);
327 const DigiSeed* pDigi = GetDigi<DigiSeed>(uDigi);
328 // hack for mCBM2024 data and Bmon station selection
329 if (fRefDet.detId == ECbmModuleId::kBmon && !filterBmon(pDigi->GetAddress())) continue;
330
331 const Double_t dTime = pDigi->GetTime();
332 //printf("time = %f %d %d\n", dTime, CbmTofAddress::GetChannelSide(add), CbmTofAddress::GetChannelId(add));
333
335 if (dTime < fdSeedWindowBeg) {
336 LOG(debug4) << Form("CbmAlgoBuildRawEvents :: reject for time %f < %f\n", dTime, fdSeedWindowBeg);
337 continue;
338 }
339 else if (fdSeedWindowEnd < dTime) {
340 break;
341 }
343 CheckSeed(dTime, uDigi);
344 }
345}
346
348{
350 return fRefDet.GetTimeWinRange();
351 }
352 else {
354 }
355}
356
357void CbmAlgoBuildRawEvents::CheckSeed(Double_t dSeedTime, UInt_t uSeedDigiIdx)
358{
361
362 LOG(debug4) << Form("CbmAlgoBuildRawEvents :: CheckSeed(%f, %d)\n", dSeedTime, uSeedDigiIdx);
363 if (nullptr != fCurrentEvent
365 && dSeedTime - fdPrevEvtTime < fdWidestTimeWinRange) {
367 switch (fOverMode) {
370 LOG(debug1) << "Reject seed due to overlap";
371 return;
372 break;
373 }
376 break;
377 }
381 LOG(debug1) << "Reject seed because part of cluster of previous one";
382 return;
383 break;
384 }
385 default: break;
386 }
387 } // if( prev Event exists and mode forbiden overlap present )
388 else {
391 if (nullptr != fCurrentEvent) {
393 // fCurrentEvent->SetStartTime( fPrevTime ); // Replace Seed time with time of first digi in event?
395 if (fbBmonInUse) {
398 }
399 fEventVector.push_back(fCurrentEvent);
400
401 fuCurEv++;
402 }
403 if (fbBmonInUse) {
405 fCurrentEvent = new CbmEvent(fuCurEv, -1, 0.);
406 }
407 else {
408 fCurrentEvent = new CbmEvent(fuCurEv, dSeedTime, 0.);
409 }
410 } // else of if( prev Event exists and mode forbiden overlap present )
411
416 SearchMatches(dSeedTime, fRefDet);
418 if (0 < fRefDet.fdTimeWinBeg) {
419 AddDigiToEvent(fRefDet, uSeedDigiIdx);
420 }
421 }
422 else {
423 AddDigiToEvent(fRefDet, uSeedDigiIdx);
424 }
425 }
426
427
432
434 bool bAllTriggersOk = true;
435 for (std::vector<RawEventBuilderDetector>::iterator det = fvDets.begin(); det != fvDets.end(); ++det) {
436 SearchMatches(dSeedTime, *det);
437
439 if (kFALSE == CheckTriggerConditions(fCurrentEvent, *det)) { //
440 bAllTriggersOk = false;
441 break;
442 }
443 }
444 if (bAllTriggersOk) {
445 fdPrevEvtTime = dSeedTime;
446
453 }
455 for (std::vector<RawEventBuilderDetector>::iterator det = fvDets.begin(); det != fvDets.end(); ++det) {
456 (*det).fuStartIndex = (*det).fuEndIndex;
457 }
458 }
459 // LOG(info) << Form("Accept seed %9.0f due to Selection Trigger requirements", dSeedTime);
460 }
461 else {
462 // LOG(info) << Form("Reject seed %9.0f due to Selection Trigger requirements", dSeedTime);
463 LOG(debug1) << "Reject seed due to Trigger requirements";
464 delete fCurrentEvent;
465 fCurrentEvent = nullptr;
466 }
467 }
468 else {
469 // LOG(info) << Form("Reject seed %9.0f due to Reference Trigger requirements", dSeedTime);
470 LOG(debug1) << "Reject seed due to Trigger requirements";
471 delete fCurrentEvent;
472 fCurrentEvent = nullptr;
473 }
474 /*
476 for (std::vector<RawEventBuilderDetector>::iterator det = fvDets.begin(); det != fvDets.end(); ++det) {
477 SearchMatches(dSeedTime, *det);
478 }
479
480 CheckTriggerCondition(dSeedTime);
481 */
482}
483
484//----------------------------------------------------------------------
486
487template<>
489{
490 return &((*fStsDigis)[uDigi]);
491}
492template<>
494{
495 return &((*fMuchBeamTimeDigis)[uDigi]);
496}
497template<>
499{
500 return &((*fMuchDigis)[uDigi]);
501}
502template<>
504{
505 return &((*fTrdDigis)[uDigi]);
506}
507template<>
509{
510 return &((*fTofDigis)[uDigi]);
511}
512template<>
514{
515 return &((*fRichDigis)[uDigi]);
516}
517template<>
519{
520 return &((*fPsdDigis)[uDigi]);
521}
522template<>
524{
525 return &((*fFsdDigis)[uDigi]);
526}
527template<>
529{
530 return &((*fBmonDigis)[uDigi]);
531}
532
533
534//----------------------------------------------------------------------
535
537{
538 switch (detMatch.detId) {
539 case ECbmModuleId::kSts: {
540 SearchMatches<CbmStsDigi>(dSeedTime, detMatch);
541 break;
542 }
543 case ECbmModuleId::kMuch: {
545 SearchMatches<CbmMuchBeamTimeDigi>(dSeedTime, detMatch);
546 }
547 else {
548 SearchMatches<CbmMuchDigi>(dSeedTime, detMatch);
549 }
550 break;
551 }
552 case ECbmModuleId::kTrd2d: // Same data storage as trd 1d
553 case ECbmModuleId::kTrd: {
554 SearchMatches<CbmTrdDigi>(dSeedTime, detMatch);
555 break;
556 }
557 case ECbmModuleId::kTof: {
558 SearchMatches<CbmTofDigi>(dSeedTime, detMatch);
559 break;
560 }
561 case ECbmModuleId::kRich: {
562 SearchMatches<CbmRichDigi>(dSeedTime, detMatch);
563 break;
564 }
565 case ECbmModuleId::kPsd: {
566 SearchMatches<CbmPsdDigi>(dSeedTime, detMatch);
567 break;
568 }
569 case ECbmModuleId::kFsd: {
570 SearchMatches<CbmFsdDigi>(dSeedTime, detMatch);
571 break;
572 }
573 case ECbmModuleId::kBmon: {
574 SearchMatches<CbmBmonDigi>(dSeedTime, detMatch);
575 break;
576 }
577 default: {
578 LOG(fatal) << "CbmAlgoBuildRawEvents::LoopOnSeeds => "
579 << "Trying to search matches with unsupported det: " << detMatch.sName << std::endl
580 << "You may want to add support for it in the method.";
581 break;
582 }
583 }
584}
585
586template<class DigiCheck>
588{
590 UInt_t uLocalIndexStart = detMatch.fuStartIndex;
591 UInt_t uLocalIndexEnd = detMatch.fuStartIndex;
592 LOG(debug4) << Form("CbmAlgoBuildRawEvents :: SearchMatches(%f, %s)\n", dSeedTime, detMatch.sName.data());
594 const UInt_t uNbSelDigis = GetNofDigis(detMatch.detId);
596 for (UInt_t uDigi = detMatch.fuStartIndex; uDigi < uNbSelDigis; ++uDigi) {
597 const DigiCheck* pDigi = GetDigi<DigiCheck>(uDigi);
598 const Double_t dTime = pDigi->GetTime();
599 const Double_t dTimeDiff = dTime - dSeedTime;
600 LOG(debug4) << detMatch.sName << Form(" => Checking match %6u / %6u, dt %f", uDigi, uNbSelDigis, dTimeDiff);
601 // int32_t add = pDigi->GetAddress(), dId[3] = {-1, -1, -1};
602 // switch (detMatch.detId) {
603 // case ECbmModuleId::kBmon:
604 // dId[0] = CbmTofAddress::GetChannelSide(add);
605 // dId[1] = CbmTofAddress::GetChannelId(add);
606 // break;
607 // case ECbmModuleId::kSts:
608 // dId[0] = CbmStsAddress::GetElementId(add, EStsElementLevel::kStsUnit);
609 // dId[1] = CbmStsAddress::GetElementId(add, EStsElementLevel::kStsLadder);
610 // dId[2] = CbmStsAddress::GetElementId(add, EStsElementLevel::kStsModule);
611 // break;
612 // case ECbmModuleId::kTrd:
613 // case ECbmModuleId::kTrd2d:
614 // dId[0] = CbmTrdAddress::GetLayerId(add);
615 // dId[1] = CbmTrdAddress::GetModuleId(add);
616 // dId[2] = CbmTrdAddress::GetModuleAddress(add);
617 // break;
618 // case ECbmModuleId::kTof:
619 // dId[0] = CbmTofAddress::GetSmId(add);
620 // dId[1] = CbmTofAddress::GetSmType(add);
621 // dId[2] = CbmTofAddress::GetRpcId(add);
622 // break;
623 // default: break;
624 // }
625 // LOG(debug4) << Form("CbmAlgoBuildRawEvents :: Checking match %6u / %6u, dt=%f %s[%d %d %d]\n", uDigi, uNbSelDigis, dTimeDiff, detMatch.sName.data(), dId[0], dId[1], dId[2]);
626
628 if (dTimeDiff < detMatch.fdTimeWinBeg) {
629 ++uLocalIndexStart;
630 continue;
631 }
632 else if (detMatch.fdTimeWinEnd < dTimeDiff) {
635 uLocalIndexEnd = uDigi;
636 break;
637 }
638
639 // Filter TRD2D digis if 1D and reverse
640 if (detMatch.detId == ECbmModuleId::kTrd) {
641 const CbmTrdDigi* pTrdDigi = GetDigi<CbmTrdDigi>(uDigi);
642 if (pTrdDigi->GetType() == CbmTrdDigi::eCbmTrdAsicType::kFASP) { //
643 continue;
644 }
645 }
646 else if (detMatch.detId == ECbmModuleId::kTrd2d) {
647 const CbmTrdDigi* pTrdDigi = GetDigi<CbmTrdDigi>(uDigi);
648 if (pTrdDigi->GetType() == CbmTrdDigi::eCbmTrdAsicType::kSPADIC) { //
649 continue;
650 }
651 }
652
653 AddDigiToEvent(detMatch, uDigi);
654 if (fdPrevEvtEndTime < dTime) fdPrevEvtEndTime = dTime;
655 }
657 if (uLocalIndexEnd < uLocalIndexStart) uLocalIndexEnd = uNbSelDigis;
658
660 detMatch.fuStartIndex = uLocalIndexStart;
661 detMatch.fuEndIndex = uLocalIndexEnd;
662 LOG(debug4) << Form("CbmAlgoBuildRawEvents :: SearchMatches(%d, %d)\n", uLocalIndexStart, uLocalIndexEnd);
663}
664
666{
667 LOG(debug4) << Form("CbmAlgoBuildRawEvents :: :: AddDigiToEvent(%s)\n", det.sName.data());
668 fCurrentEvent->AddData(det.dataType, _entry);
669}
670
672{
673 LOG(debug4) << Form("CbmAlgoBuildRawEvents :: :: CheckTriggerCondition(%f)\n", dSeedTime);
676 fdPrevEvtTime = dSeedTime;
677
684 }
686 for (std::vector<RawEventBuilderDetector>::iterator det = fvDets.begin(); det != fvDets.end(); ++det) {
687 (*det).fuStartIndex = (*det).fuEndIndex;
688 }
689 }
690 }
691 else {
692 LOG(debug1) << "Reject seed due to Trigger requirements";
693 delete fCurrentEvent;
694 fCurrentEvent = nullptr;
695 }
696}
697
699{
702 if (kFALSE == CheckTriggerConditions(event, fRefDet)) {
703 return kFALSE;
704 }
705 }
707 for (std::vector<RawEventBuilderDetector>::iterator det = fvDets.begin(); det != fvDets.end(); ++det) {
708 if (kFALSE == CheckTriggerConditions(event, *det)) return kFALSE;
709 }
711 return kTRUE;
712}
713
715{
716 const int32_t iNbDigis = event->GetNofData(ECbmDataType::kBmonDigi);
717 double eventTime(0.);
718 bool timeSet(false);
719
720 LOG(debug4) << Form("CbmAlgoBuildRawEvents :: SetBmonEventTime(%p, %d)\n", (void*) event, iNbDigis);
721 for (int idigi = 0; idigi < iNbDigis; ++idigi) {
722 uint idx = event->GetIndex(ECbmDataType::kBmonDigi, idigi);
723 const CbmBmonDigi* pDigi = GetDigi<CbmBmonDigi>(idx);
724 if (nullptr == pDigi) continue;
725 if (!filterBmon(pDigi->GetAddress())) continue;
726
727 if (!timeSet) {
728 eventTime = pDigi->GetTime();
729 timeSet = true;
730 }
731 else
732 eventTime = std::min(pDigi->GetTime(), eventTime);
733 }
734
735 if (timeSet) {
736 event->SetStartTime(eventTime);
737 }
738 else
739 LOG(warning) << "CbmAlgoBuildRawEvents::SetBmonEventTime : failed";
740
741 return timeSet;
742}
743
745{
746 LOG(debug4) << Form("CbmAlgoBuildRawEvents :: :: CheckTriggerConditions(%p, %s)\n", (void*) event, det.sName.data());
748 if (0 == det.fuTriggerMinDigis && det.fiTriggerMaxDigis < 0) {
749 return kTRUE;
750 }
751
753 if (!CheckDataAvailable(det.detId)) {
754 LOG(debug2) << "Event does not have digis storage for " << det.sName
755 << " while the following trigger min/max are defined: " << det.fuTriggerMinDigis << " "
756 << det.fiTriggerMaxDigis;
757 return kFALSE;
758 }
759
761 int32_t iNbDigis = event->GetNofData(det.dataType);
762 int32_t iNbFilteredDigis = (det.detId == ECbmModuleId::kBmon ? getNofFilteredBmonDigis(event) : iNbDigis);
763
765 if (0 < det.fuTriggerMinDigis
766 && ((-1 == iNbFilteredDigis) || (static_cast<UInt_t>(iNbFilteredDigis) < det.fuTriggerMinDigis))) {
767 LOG(debug2) << "Event does not have enough digis: " << iNbFilteredDigis << " vs " << det.fuTriggerMinDigis
768 << " for " << det.sName;
769 return kFALSE;
770 }
771
773 if (0 <= det.fiTriggerMaxDigis && det.fiTriggerMaxDigis < iNbFilteredDigis) {
774 LOG(debug2) << "Event Has too many digis: " << iNbFilteredDigis << " vs " << det.fiTriggerMaxDigis << " for "
775 << det.sName;
776 return kFALSE;
777 }
778
780 if (0 < det.fuTriggerMinLayers) {
781 switch (det.detId) {
782 case ECbmModuleId::kSts: {
785 std::set<uint32_t> setStations; // Use set instead of vector as search by value later
786 std::map<uint32_t, int> mModules;
787
788 for (int idigi = 0; idigi < iNbDigis; ++idigi) {
789 uint idx = event->GetIndex(det.dataType, idigi);
790 const CbmStsDigi* pDigi = GetDigi<CbmStsDigi>(idx);
791 if (nullptr == pDigi) continue;
792
793 int iAddr = pDigi->GetAddress();
798 //int iStationAddr = CbmStsAddress::GetElementId(iAddr, EStsElementLevel::kStsUnit) / 2;
799
800 std::map<uint32_t, int>::iterator itModule = mModules.find(iModuleAddr);
801 if (itModule == mModules.end()) {
802 // LOG(info) << Form("Found new module 0x%08x, side %u", iModuleAddr,
803 // static_cast<uint32_t>(pDigi->GetChannel() / 1024));
804 mModules[iModuleAddr] = static_cast<int32_t>(pDigi->GetChannel() / 1024); // extend map
805 }
806 else {
807 // LOG(info) << Form("Check side %u of module 0x%08x: %d ?",
808 // static_cast<int32_t>(pDigi->GetChannel() / 1024),
809 // iModuleAddr, itModule->second);
810 if (static_cast<int32_t>(pDigi->GetChannel() / 1024) == (1 - itModule->second)) {
812 auto itStation = setStations.find(iStationAddr);
813 if (itStation == setStations.end()) {
814 // LOG(info) << Form("Add station 0x%08x ", iStationAddr);
815 setStations.insert(iStationAddr);
816 }
817 }
818 }
819 }
820 // LOG(info) << "Found " << setStations.size() << " Sts stations, " << " in " << iNbDigis << " Sts digis";
821 if (setStations.size() < det.fuTriggerMinLayers) {
822 LOG(debug2) << "Event does not have enough layers fired: " << setStations.size() << " vs "
823 << det.fuTriggerMinLayers << " for " << det.sName;
824 return kFALSE;
825 }
826 break;
827 }
828 case ECbmModuleId::kMuch: {
829 LOG(fatal) << "CbmAlgoBuildRawEvents::CheckTriggerConditions => Fired layers check not implemented yet for "
830 << det.sName;
831 return kFALSE;
832 break;
833 }
834 case ECbmModuleId::kTrd2d: // Same data storage as trd 1d
835 case ECbmModuleId::kTrd: {
836 LOG(fatal) << "CbmAlgoBuildRawEvents::CheckTriggerConditions => Fired layers check not implemented yet for "
837 << det.sName;
838 return kFALSE;
839 break;
840 }
841 case ECbmModuleId::kTof: {
844 std::set<uint32_t> setRpcs; // Use set instead of vector as search by value later
845 std::map<uint32_t, int> mStrips;
846
847 for (int idigi = 0; idigi < iNbDigis; ++idigi) {
848 uint idx = event->GetIndex(det.dataType, idigi);
849 const CbmTofDigi* pDigi = GetDigi<CbmTofDigi>(idx);
850 if (nullptr == pDigi) continue;
851
852 int iAddr = pDigi->GetAddress();
853 int iStripAddr = CbmTofAddress::GetStripFullId(iAddr);
854 int iRpcAddr = CbmTofAddress::GetRpcFullId(iAddr);
855
856 std::map<uint32_t, int>::iterator itStrip = mStrips.find(iStripAddr);
857 if (itStrip == mStrips.end()) {
858 // LOG(info) << Form("Found new strip 0x%08x, side %u", iStripAddr, pDigi->GetSide());
859 mStrips[iStripAddr] = (int) pDigi->GetSide(); // extend map
860 }
861 else {
862 // LOG(info) << Form("Check side %u of strip 0x%08x: %d ?", pDigi->GetSide(), iStripAddr, itStrip->second);
863 if ((int) pDigi->GetSide() == (1 - itStrip->second)) {
865 auto itRpc = setRpcs.find(iRpcAddr);
866 if (itRpc == setRpcs.end()) {
867 // LOG(info) << Form("Add counter 0x%08x ", iRpcAddr);
868 setRpcs.insert(iRpcAddr);
869 }
870 }
871 }
872 }
873 // LOG(info) << "Found " << setRpcs.size() << " Tof RPCs, " << " in " << iNbDigis << " Tof digis";
874 if (setRpcs.size() < det.fuTriggerMinLayers) {
875 LOG(debug2) << "Event does not have enough RPCs fired: " << setRpcs.size() << " vs " << det.fuTriggerMinLayers
876 << " for " << det.sName;
877 return kFALSE;
878 }
879 break;
880 }
881 case ECbmModuleId::kRich: {
882 LOG(fatal) << "CbmAlgoBuildRawEvents::CheckTriggerConditions => Fired layers check not implemented yet for "
883 << det.sName;
884 return kFALSE;
885 break;
886 }
887 case ECbmModuleId::kPsd: {
888 LOG(fatal) << "CbmAlgoBuildRawEvents::CheckTriggerConditions => Fired layers check not implemented yet for "
889 << det.sName;
890 return kFALSE;
891 break;
892 }
893 case ECbmModuleId::kFsd: {
894 LOG(fatal) << "CbmAlgoBuildRawEvents::CheckTriggerConditions => Fired layers check not implemented yet for "
895 << det.sName;
896 return kFALSE;
897 break;
898 }
899 case ECbmModuleId::kBmon: {
900 LOG(fatal) << "CbmAlgoBuildRawEvents::CheckTriggerConditions => Fired layers check not implemented yet for "
901 << det.sName;
902 return kFALSE;
903 break;
904 }
905 default: {
906 LOG(fatal) << "CbmAlgoBuildRawEvents::CheckTriggerConditions => Fired layers check not implemented yet for "
907 << det.sName;
908 return kFALSE;
909 break;
910 }
911 }
912 }
913
915 return kTRUE;
916}
917
918//----------------------------------------------------------------------
919
921{
922 if (!CheckDataAvailable(det.detId)) {
923 LOG(info) << "No " << det.sName << " digi input found.";
924 return kFALSE;
925 }
926 return kTRUE;
927}
928
930{
931 switch (detId) {
932 case ECbmModuleId::kSts: {
933 return fStsDigis != nullptr;
934 }
935 case ECbmModuleId::kMuch: {
937 return fMuchBeamTimeDigis != nullptr;
938 }
939 else {
940 return fMuchDigis != nullptr;
941 }
942 }
943 case ECbmModuleId::kTrd2d: // Same data storage as trd 1d
944 case ECbmModuleId::kTrd: {
945 return fTrdDigis != nullptr;
946 }
947 case ECbmModuleId::kTof: {
948 return fTofDigis != nullptr;
949 }
950 case ECbmModuleId::kRich: {
951 return fRichDigis != nullptr;
952 }
953 case ECbmModuleId::kPsd: {
954 return fPsdDigis != nullptr;
955 }
956 case ECbmModuleId::kFsd: {
957 return fFsdDigis != nullptr;
958 }
959 case ECbmModuleId::kBmon: {
960 return fBmonDigis != nullptr;
961 }
962 default: {
963 LOG(fatal) << "CbmAlgoBuildRawEvents::CheckDataAvailable => "
964 << "Unsupported detector.";
965 return -1;
966 }
967 }
968}
969
971{
972 switch (detId) {
973 case ECbmModuleId::kSts: {
974 return fStsDigis->size();
975 }
976 case ECbmModuleId::kMuch: {
978 return fMuchBeamTimeDigis->size();
979 }
980 else {
981 return fMuchDigis->size();
982 }
983 }
984 case ECbmModuleId::kTrd2d: // Same data storage as trd 1d
985 case ECbmModuleId::kTrd: {
986 return fTrdDigis->size();
987 }
988 case ECbmModuleId::kTof: {
989 return fTofDigis->size();
990 }
991 case ECbmModuleId::kRich: {
992 return fRichDigis->size();
993 }
994 case ECbmModuleId::kPsd: {
995 return fPsdDigis->size();
996 }
997 case ECbmModuleId::kFsd: {
998 return fFsdDigis->size();
999 }
1000 case ECbmModuleId::kBmon: {
1001 return fBmonDigis->size();
1002 }
1003 default: {
1004 LOG(fatal) << "CbmAlgoBuildRawEvents::GetNofDigis => "
1005 << "Trying to get digi number with unsupported detector.";
1006 return -1;
1007 }
1008 }
1009}
1011{
1012 switch (detId) {
1013 case ECbmModuleId::kSts: {
1014 return ulNbDigis * sizeof(CbmStsDigi);
1015 }
1016 case ECbmModuleId::kMuch: {
1018 return ulNbDigis * sizeof(CbmMuchBeamTimeDigi);
1019 }
1020 else {
1021 return ulNbDigis * sizeof(CbmMuchDigi);
1022 }
1023 }
1024 case ECbmModuleId::kTrd2d: // Same data storage as trd 1d
1025 case ECbmModuleId::kTrd: {
1026 return ulNbDigis * sizeof(CbmTrdDigi);
1027 }
1028 case ECbmModuleId::kTof: {
1029 return ulNbDigis * sizeof(CbmTofDigi);
1030 }
1031 case ECbmModuleId::kRich: {
1032 return ulNbDigis * sizeof(CbmRichDigi);
1033 }
1034 case ECbmModuleId::kPsd: {
1035 return ulNbDigis * sizeof(CbmPsdDigi);
1036 }
1037 case ECbmModuleId::kFsd: {
1038 return ulNbDigis * sizeof(CbmFsdDigi);
1039 }
1040 case ECbmModuleId::kBmon: {
1041 return ulNbDigis * sizeof(CbmBmonDigi);
1042 }
1043 default: {
1044 LOG(fatal) << "CbmAlgoBuildRawEvents::GetSizeFromDigisNb => "
1045 << "Trying to get digi number with unsupported detector.";
1046 return -1;
1047 }
1048 }
1049}
1050
1051//----------------------------------------------------------------------
1053{
1054 outFolder = new TDirectoryFile("AlgoBuildRawEventsHist", " AlgoBuildRawEvents Histos");
1055 outFolder->Clear();
1056
1057 fhEventTime = new TH1F("hEventTime", "seed time of the events; Seed time within TS [s]; Events", 100000, 0, 0.2);
1058 // fhEventTime->SetCanExtend(TH1::kAllAxes); // Breaks the MQ histogram server as cannot be merged!
1059
1060 fhEventDt =
1061 new TH1F("fhEventDt", "interval in seed time of consecutive events; Seed time dt [ns]; Events", 10000, 0, 100000);
1062 // fhEventDt->SetCanExtend(TH1::kAllAxes); // Breaks the MQ histogram server as cannot be merged!
1063
1064 double_t dHistMaxTotDigis = fRefDet.fdHistMaxDigiNb;
1065 for (std::vector<RawEventBuilderDetector>::iterator det = fvDets.begin(); det != fvDets.end(); ++det) {
1066 dHistMaxTotDigis += (*det).fdHistMaxDigiNb;
1067 }
1068 fhEventSize = new TH1F("hEventSize", "nb of all digis in the event; Nb Digis []; Events []", dHistMaxTotDigis, 0,
1069 dHistMaxTotDigis);
1070 // fhEventSize->SetCanExtend(TH1::kAllAxes); // Breaks the MQ histogram server as cannot be merged!
1071
1072 fhNbDigiPerEvtTime = new TH2I("hNbDigiPerEvtTime",
1073 "nb of all digis per event vs seed time of the events; Seed time "
1074 "[s]; Nb Digis []; Events []",
1075 1000, 0, 0.2, dHistMaxTotDigis, 0, dHistMaxTotDigis);
1076 // fhNbDigiPerEvtTime->SetCanExtend(TH2::kAllAxes); // Breaks he MQ histogram server as cannot be merged!
1077
1078 fhCpuTimePerTs = new TH1D("hCpuTimePerTs", "CPU Processing time of TS vs TS; Ts; CPU time [ms]", 6000, 0, 6000);
1079 fhRealTimePerTs = new TH1D("hRealTimePerTs", "Real Processing time of TS vs TS; Ts; Real time [ms]", 6000, 0, 6000);
1080
1082 new TH1D("hCpuTimePerTsHist", "CPU Histo filling time of TS vs TS; Ts; CPU time [ms]", 6000, 0, 6000);
1084 new TH1D("hRealTimePerTsHist", "Real Histo filling time of TS vs TS; Ts; Real time [ms]", 6000, 0, 6000);
1085
1086 AddHistoToVector(fhEventTime, "evtbuild");
1087 AddHistoToVector(fhEventDt, "evtbuild");
1088 AddHistoToVector(fhEventSize, "evtbuild");
1090 AddHistoToVector(fhCpuTimePerTs, "evtbuild-eff");
1091 AddHistoToVector(fhRealTimePerTs, "evtbuild-eff");
1092 AddHistoToVector(fhCpuTimePerTsHist, "evtbuild-eff");
1093 AddHistoToVector(fhRealTimePerTsHist, "evtbuild-eff");
1094 outFolder->Add(fhEventTime);
1095 outFolder->Add(fhEventDt);
1096 outFolder->Add(fhEventSize);
1102
1104 fhOverEventShare = new TH1I("fhOverEventShare", "Share of overlap evt; Overlap? []; Events", 2, -0.5, 1.5);
1105 fhOverEventShareTs = new TProfile(
1106 "fhOverEventShareTs", "Share of overlap evt per TS; TS index []; Overlap Events prop. []", 2500, 0, 2500);
1108 new TH2I("fhOverEventSizeTs", "Size of overlap of evt per TS; TS index []; Size of overlap between events [ns]",
1109 2500, 0, 2500, 200, 0, 1000);
1116 }
1117
1119 for (std::vector<RawEventBuilderDetector>::iterator det = fvDets.begin(); det != fvDets.end(); ++det) {
1121 if ("Invalid" == (*det).sName) {
1122 fvhNbDigiPerEvtTimeDet.push_back(nullptr);
1123 fvhNbDigiPerEvtDet.push_back(nullptr);
1124 continue;
1125 }
1126 TH2I* hNbDigiPerEvtTimeDet = new TH2I(Form("hNbDigiPerEvtTime%s", (*det).sName.data()),
1127 Form("nb of %s digis per event vs seed time of the events; Seed time in TS "
1128 "[s]; Nb Digis []; Events []",
1129 (*det).sName.data()),
1130 1000, 0, 0.2, (*det).fdHistMaxDigiNb, 0, (*det).fdHistMaxDigiNb);
1131 // hNbDigiPerEvtTimeDet->SetCanExtend(TH2::kAllAxes); // Breaks he MQ histogram server as cannot be merged!
1132 fvhNbDigiPerEvtTimeDet.push_back(hNbDigiPerEvtTimeDet);
1133
1134 TH1* hNbDigiPerEvtDet = new TH1I(Form("hNbDigiPerEvt%s", (*det).sName.data()),
1135 Form("nb of %s digis per event; Nb Digis []", (*det).sName.data()),
1136 (*det).fdHistMaxDigiNb, 0, (*det).fdHistMaxDigiNb);
1137 fvhNbDigiPerEvtDet.push_back(hNbDigiPerEvtDet);
1138
1139 TH1* hTDiff =
1140 new TH1I(Form("hTDiff%s", (*det).sName.data()),
1141 Form("#DeltaT of %s digis to seed time of event;#DeltaT (ns); Counts []", (*det).sName.data()), 200,
1142 (*det).fdTimeWinBeg, (*det).fdTimeWinEnd);
1143 fvhTDiff.push_back(hTDiff);
1144
1145 // clang-format off
1146 TH1* hSelRatioPerTsNb = new TH1D(Form("hSelRatioPerTsNb%s", (*det).sName.data()),
1147 Form("ratio of sel digis per TS vs TS for %s; TS; Sel Digis Ratio []",
1148 (*det).sName.data()),
1149 6000, 0, 6000);
1150 TH1* hInpRatioPerTsSz = new TH1D(Form("hInpRatioPerTsSz%s", (*det).sName.data()),
1151 Form("ratio of input digi size in total input size vs TS for %s; TS; Size Ratio []",
1152 (*det).sName.data()),
1153 6000, 0, 6000);
1154 TH1* hOutRatioPerTsSz = new TH1D(Form("hOutRatioPerTsSz%s", (*det).sName.data()),
1155 Form("ratio of selected digi size in event size vs TS for %s; TS; Size Ratio []",
1156 (*det).sName.data()),
1157 6000, 0, 6000);
1158 // clang-format on
1159
1160 fvhSelRatioPerTsNb.push_back(hSelRatioPerTsNb);
1161 fvhInpRatioPerTsSz.push_back(hInpRatioPerTsSz);
1162 fvhOutRatioPerTsSz.push_back(hOutRatioPerTsSz);
1163
1164 AddHistoToVector(hSelRatioPerTsNb, "evtbuild-eff");
1165 AddHistoToVector(hInpRatioPerTsSz, "evtbuild-eff");
1166 AddHistoToVector(hOutRatioPerTsSz, "evtbuild-eff");
1167
1168 outFolder->Add(hSelRatioPerTsNb);
1169 outFolder->Add(hInpRatioPerTsSz);
1170 outFolder->Add(hOutRatioPerTsSz);
1171 }
1172
1175 TH2I* hNbDigiPerEvtTimeDet = new TH2I(Form("hNbDigiPerEvtTime%s", fRefDet.sName.data()),
1176 Form("nb of %s digis per event vs seed time of the events; Seed time in TS "
1177 "[s]; Nb Digis []; Events []",
1178 fRefDet.sName.data()),
1180 fvhNbDigiPerEvtTimeDet.push_back(hNbDigiPerEvtTimeDet);
1181
1182 TH1I* hNbDigiPerEvtDet = new TH1I(Form("hNbDigiPerEvt%s", fRefDet.sName.data()),
1183 Form("nb of %s digis per event; Nb Digis []", fRefDet.sName.data()),
1185 fvhNbDigiPerEvtDet.push_back(hNbDigiPerEvtDet);
1186
1187 TH1I* hTDiff =
1188 new TH1I(Form("hTDiff%s", fRefDet.sName.data()),
1189 Form("#DeltaT of %s digis to seed time of event;#DeltaT (ns); Counts []", fRefDet.sName.data()), 200,
1190 fRefDet.fdTimeWinBeg, fRefDet.fdTimeWinEnd); // FIXME, adjust to configured window
1191 fvhTDiff.push_back(hTDiff);
1192
1193 // clang-format off
1194 TH1* hSelRatioPerTsNb = new TH1D(Form("hSelRatioPerTsNb%s", fRefDet.sName.data()),
1195 Form("ratio of sel digis per TS vs TS for %s; TS; Sel Digis Ratio []",
1196 fRefDet.sName.data()),
1197 6000, 0, 6000);
1198 TH1* hInpRatioPerTsSz = new TH1D(Form("hInpRatioPerTsSz%s", fRefDet.sName.data()),
1199 Form("ratio of input digi size in total input size vs TS for %s; TS; Size Ratio []",
1200 fRefDet.sName.data()),
1201 6000, 0, 6000);
1202 TH1* hOutRatioPerTsSz = new TH1D(Form("hOutRatioPerTsSz%s", fRefDet.sName.data()),
1203 Form("ratio of selected digi size in event size vs TS for %s; TS; Size Ratio []",
1204 fRefDet.sName.data()),
1205 6000, 0, 6000);
1206 // clang-format on
1207
1208 fvhSelRatioPerTsNb.push_back(hSelRatioPerTsNb);
1209 fvhInpRatioPerTsSz.push_back(hInpRatioPerTsSz);
1210 fvhOutRatioPerTsSz.push_back(hOutRatioPerTsSz);
1211
1212 AddHistoToVector(hSelRatioPerTsNb, "evtbuild-eff");
1213 AddHistoToVector(hInpRatioPerTsSz, "evtbuild-eff");
1214 AddHistoToVector(hOutRatioPerTsSz, "evtbuild-eff");
1215
1216 outFolder->Add(hSelRatioPerTsNb);
1217 outFolder->Add(hInpRatioPerTsSz);
1218 outFolder->Add(hOutRatioPerTsSz);
1219 }
1220
1222 new TH1D("hSizeReductionPerTs", "ratio of tot. sel. digi size to tot. input digi size vs TS; TS; Size Ratio []",
1223 6000, 0, 6000);
1224 AddHistoToVector(fhSizeReductionPerTs, "evtbuild-eff");
1226
1227 for (std::vector<TH2*>::iterator itHist = fvhNbDigiPerEvtTimeDet.begin(); itHist != fvhNbDigiPerEvtTimeDet.end();
1228 ++itHist) {
1229 if (nullptr != (*itHist)) {
1230 AddHistoToVector((*itHist), "evtbuild");
1231 outFolder->Add((*itHist));
1232 }
1233 }
1234
1235 for (std::vector<TH1*>::iterator itHist = fvhNbDigiPerEvtDet.begin(); itHist != fvhNbDigiPerEvtDet.end(); ++itHist) {
1236 if (nullptr != (*itHist)) {
1237 AddHistoToVector((*itHist), "evtbuild");
1238 outFolder->Add((*itHist));
1239 }
1240 }
1241 for (std::vector<TH1*>::iterator itHist = fvhTDiff.begin(); itHist != fvhTDiff.end(); ++itHist) {
1242 if (nullptr != (*itHist)) {
1243 AddHistoToVector((*itHist), "evtbuild");
1244 outFolder->Add((*itHist));
1245 }
1246 }
1247
1249 // std::vector<std::pair<TCanvas*, std::string>> vCanvases = {};
1250
1251 TCanvas* fcSummary = new TCanvas("cEvBSummary", "EvB monitoring plots");
1252 fcSummary->Divide(2, 2);
1253
1254 fcSummary->cd(1);
1255 gPad->SetGridx();
1256 gPad->SetGridy();
1257 fhEventTime->Draw("hist");
1258
1259 fcSummary->cd(2);
1260 gPad->SetGridx();
1261 gPad->SetGridy();
1262 gPad->SetLogx();
1263 gPad->SetLogy();
1264 fhEventDt->Draw("hist");
1265
1266 fcSummary->cd(3);
1267 gPad->SetGridx();
1268 gPad->SetGridy();
1269 gPad->SetLogy();
1270 fhEventSize->Draw("hist");
1271
1272 fcSummary->cd(4);
1273 gPad->SetGridx();
1274 gPad->SetGridy();
1275 fhNbDigiPerEvtTime->Draw("colz");
1276
1277
1279 AddCanvasToVector(fcSummary, "canvases");
1280
1281 // ------------------------ //
1282 TCanvas* fcNbDigi = new TCanvas("cEvBNbDigi", "EvB NbDigi evolution ");
1283 if (fvhNbDigiPerEvtDet.size() <= 6) { //
1284 fcNbDigi->Divide(2, 3);
1285 }
1286 else { //
1287 fcNbDigi->Divide(3, 3);
1288 }
1289 int iPad = 1;
1290 for (std::vector<TH1*>::iterator itHist = fvhNbDigiPerEvtDet.begin(); itHist != fvhNbDigiPerEvtDet.end(); ++itHist) {
1291 if (nullptr != (*itHist)) {
1292 fcNbDigi->cd(iPad++);
1293 gPad->SetGridx();
1294 gPad->SetGridy();
1295 gPad->SetLogy();
1296 (*itHist)->Draw(); //"colz");
1297 }
1298 }
1299 AddCanvasToVector(fcNbDigi, "canvases");
1300
1301 // ------------------------ //
1302 TCanvas* fcTdif = new TCanvas("cEvBTdif", "EvB Time Difference plots");
1303 if (fvhNbDigiPerEvtDet.size() <= 6) { //
1304 fcTdif->Divide(2, 3);
1305 }
1306 else { //
1307 fcTdif->Divide(3, 3);
1308 }
1309 iPad = 1;
1310 for (std::vector<TH1*>::iterator itHist = fvhTDiff.begin(); itHist != fvhTDiff.end(); ++itHist) {
1311 if (nullptr != (*itHist)) {
1312 fcTdif->cd(iPad++);
1313 gPad->SetGridx();
1314 gPad->SetGridy();
1315 gPad->SetLogy();
1316 (*itHist)->Draw();
1317 }
1318 }
1319 AddCanvasToVector(fcTdif, "canvases");
1320}
1321
1323{
1325 uint32_t uRefDetIdx = fvDets.size();
1326 uint64_t ulTotalInputSize = 0;
1327 uint64_t ulTotalOutputSize = 0;
1328 std::vector<uint64_t> vulTotalInputSizeDet(fvDets.size() + 1, 0);
1329 std::vector<uint64_t> vulTotalOutputSizeDet(fvDets.size() + 1, 0);
1330
1332 Double_t dPreEvtTime = -1.0;
1333 for (CbmEvent* evt : fEventVector) {
1334 fhEventTime->Fill(evt->GetStartTime() * 1e-9);
1335 if (0.0 <= dPreEvtTime) {
1336 fhEventDt->Fill((evt->GetStartTime() - dPreEvtTime) * 1e-9);
1337
1339 if (evt->GetStartTime() - dPreEvtTime < fdWidestTimeWinRange) {
1340 fhOverEventShare->Fill(1);
1341 fhOverEventShareTs->Fill(fuNrTs, 1);
1342 fhOverEventSizeTs->Fill(fuNrTs, fdWidestTimeWinRange - (evt->GetStartTime() - dPreEvtTime));
1343 }
1344 else {
1345 fhOverEventShare->Fill(0);
1346 fhOverEventShareTs->Fill(fuNrTs, 0);
1347 }
1348 }
1349 }
1352 fhOverEventShare->Fill(0);
1353 fhOverEventShareTs->Fill(fuNrTs, 0);
1354 }
1355
1356 fhEventSize->Fill(evt->GetNofData());
1357 fhNbDigiPerEvtTime->Fill(evt->GetStartTime() * 1e-9, evt->GetNofData());
1358
1360 uint32_t uNbDataTrd1d = 0;
1361 uint32_t uNbDataTrd2d = 0;
1362 for (UInt_t uDetIdx = 0; uDetIdx < fvDets.size(); ++uDetIdx) {
1363 if (nullptr == fvhNbDigiPerEvtDet[uDetIdx]) continue;
1364
1365 for (size_t idigi = 0; idigi < evt->GetNofData(fvDets[uDetIdx].dataType); ++idigi) {
1366 double dTimeDiff = 1.E30;
1367 uint idx = evt->GetIndex(fvDets[uDetIdx].dataType, idigi);
1368 switch (fvDets[uDetIdx].dataType) {
1370 auto pDigi = GetDigi<CbmBmonDigi>(idx);
1371 if (nullptr == pDigi) continue;
1372 dTimeDiff = pDigi->GetTime() - evt->GetStartTime();
1373 break;
1374 }
1376 auto pDigi = GetDigi<CbmStsDigi>(idx);
1377 if (nullptr == pDigi) continue;
1378 dTimeDiff = pDigi->GetTime() - evt->GetStartTime();
1379 break;
1380 }
1383 auto pDigi = GetDigi<CbmMuchBeamTimeDigi>(idx);
1384 if (nullptr == pDigi) continue;
1385 dTimeDiff = pDigi->GetTime() - evt->GetStartTime();
1386 }
1387 else {
1388 auto pDigi = GetDigi<CbmMuchDigi>(idx);
1389 if (nullptr == pDigi) continue;
1390 dTimeDiff = pDigi->GetTime() - evt->GetStartTime();
1391 }
1392 break;
1393 }
1395 auto pDigi = GetDigi<CbmTofDigi>(idx);
1396 if (nullptr == pDigi) continue;
1397 dTimeDiff = pDigi->GetTime() - evt->GetStartTime();
1398 break;
1399 }
1401 auto pDigi = GetDigi<CbmTrdDigi>(idx);
1402 if (nullptr == pDigi) continue;
1403 dTimeDiff = pDigi->GetTime() - evt->GetStartTime();
1404 if (pDigi->GetType() == CbmTrdDigi::eCbmTrdAsicType::kSPADIC) {
1405 if (fvDets[uDetIdx].sName == "Trd2D") continue;
1406 ++uNbDataTrd1d;
1407 }
1408 else if (pDigi->GetType() == CbmTrdDigi::eCbmTrdAsicType::kFASP) {
1409 if (fvDets[uDetIdx].sName == "Trd1D") continue;
1410 ++uNbDataTrd2d;
1411 }
1412 break;
1413 }
1415 auto pDigi = GetDigi<CbmRichDigi>(idx); // FIXME, need to find the proper digi template
1416 if (nullptr == pDigi) continue;
1417 dTimeDiff = pDigi->GetTime() - evt->GetStartTime();
1418 break;
1419 }
1421 auto pDigi = GetDigi<CbmPsdDigi>(idx); // FIXME, need to find the proper digi template
1422 if (nullptr == pDigi) continue;
1423 dTimeDiff = pDigi->GetTime() - evt->GetStartTime();
1424 break;
1425 }
1427 auto pDigi = GetDigi<CbmFsdDigi>(idx); // FIXME, need to find the proper digi template
1428 if (nullptr == pDigi) continue;
1429 dTimeDiff = pDigi->GetTime() - evt->GetStartTime();
1430 break;
1431 }
1432 default: LOG(error) << "Unkown dataType " << fvDets[uDetIdx].dataType;
1433 }
1434
1435 if (dTimeDiff < 1.E30) fvhTDiff[uDetIdx]->Fill(dTimeDiff);
1436 }
1437 }
1438
1441 if (nullptr != fvhNbDigiPerEvtDet[uRefDetIdx]) {
1442 for (size_t idigi = 0; idigi < evt->GetNofData(fRefDet.dataType); ++idigi) {
1443 double dTimeDiff = 1.E30;
1444 uint idx = evt->GetIndex(fRefDet.dataType, idigi);
1445 switch (fRefDet.dataType) {
1447 auto pDigi = GetDigi<CbmBmonDigi>(idx);
1448 if (nullptr == pDigi) continue;
1449 dTimeDiff = pDigi->GetTime() - evt->GetStartTime();
1450 break;
1451 }
1453 auto pDigi = GetDigi<CbmStsDigi>(idx);
1454 if (nullptr == pDigi) continue;
1455 dTimeDiff = pDigi->GetTime() - evt->GetStartTime();
1456 break;
1457 }
1460 auto pDigi = GetDigi<CbmMuchBeamTimeDigi>(idx);
1461 if (nullptr == pDigi) continue;
1462 dTimeDiff = pDigi->GetTime() - evt->GetStartTime();
1463 }
1464 else {
1465 auto pDigi = GetDigi<CbmMuchDigi>(idx);
1466 if (nullptr == pDigi) continue;
1467 dTimeDiff = pDigi->GetTime() - evt->GetStartTime();
1468 }
1469 break;
1470 }
1472 auto pDigi = GetDigi<CbmTofDigi>(idx);
1473 if (nullptr == pDigi) continue;
1474 dTimeDiff = pDigi->GetTime() - evt->GetStartTime();
1475 break;
1476 }
1478 auto pDigi = GetDigi<CbmTrdDigi>(idx);
1479 if (nullptr == pDigi) continue;
1480 dTimeDiff = pDigi->GetTime() - evt->GetStartTime();
1481 if (pDigi->GetType() == CbmTrdDigi::eCbmTrdAsicType::kSPADIC) {
1482 if (fRefDet.sName == "Trd2D") continue;
1483 ++uNbDataTrd1d;
1484 }
1485 else if (pDigi->GetType() == CbmTrdDigi::eCbmTrdAsicType::kFASP) {
1486 if (fRefDet.sName == "Trd1D") continue;
1487 ++uNbDataTrd2d;
1488 }
1489 break;
1490 }
1492 auto pDigi = GetDigi<CbmRichDigi>(idx); // FIXME, need to find the proper digi template
1493 if (nullptr == pDigi) continue;
1494 dTimeDiff = pDigi->GetTime() - evt->GetStartTime();
1495 break;
1496 }
1498 auto pDigi = GetDigi<CbmPsdDigi>(idx); // FIXME, need to find the proper digi template
1499 if (nullptr == pDigi) continue;
1500 dTimeDiff = pDigi->GetTime() - evt->GetStartTime();
1501 break;
1502 }
1504 auto pDigi = GetDigi<CbmFsdDigi>(idx); // FIXME, need to find the proper digi template
1505 if (nullptr == pDigi) continue;
1506 dTimeDiff = pDigi->GetTime() - evt->GetStartTime();
1507 break;
1508 }
1509 default: LOG(error) << "Unkown dataType " << fRefDet.dataType;
1510 }
1511
1512 if (dTimeDiff < 1.E30) fvhTDiff[uRefDetIdx]->Fill(dTimeDiff);
1513 }
1514 }
1515 }
1516
1518 for (UInt_t uDetIdx = 0; uDetIdx < fvDets.size(); ++uDetIdx) {
1519 if (nullptr == fvhNbDigiPerEvtTimeDet[uDetIdx]) continue;
1520
1521 if (fvDets[uDetIdx].sName == "Trd1D") {
1522 fvhNbDigiPerEvtDet[uDetIdx]->Fill(uNbDataTrd1d);
1523 fvhNbDigiPerEvtTimeDet[uDetIdx]->Fill(evt->GetStartTime() * 1e-9, uNbDataTrd1d);
1524
1525 if (0 < GetNofDigis(fvDets[uDetIdx].detId)) {
1527 uint64_t ulDigiSizeOut = GetSizeFromDigisNb(fvDets[uDetIdx].detId, uNbDataTrd1d + uNbDataTrd2d);
1528
1529 ulTotalOutputSize += ulDigiSizeOut;
1530 vulTotalOutputSizeDet[uDetIdx] += ulDigiSizeOut;
1531 }
1532 }
1533 else if (fvDets[uDetIdx].sName == "Trd2D") {
1534 fvhNbDigiPerEvtDet[uDetIdx]->Fill(uNbDataTrd2d);
1535 fvhNbDigiPerEvtTimeDet[uDetIdx]->Fill(evt->GetStartTime() * 1e-9, uNbDataTrd2d);
1536 }
1537 else {
1538 fvhNbDigiPerEvtDet[uDetIdx]->Fill(evt->GetNofData(fvDets[uDetIdx].dataType));
1539 fvhNbDigiPerEvtTimeDet[uDetIdx]->Fill(evt->GetStartTime() * 1e-9, evt->GetNofData(fvDets[uDetIdx].dataType));
1540
1541 if (0 < GetNofDigis(fvDets[uDetIdx].detId)) {
1543 uint64_t ulDigiSizeOut = GetSizeFromDigisNb(fvDets[uDetIdx].detId, evt->GetNofData(fvDets[uDetIdx].dataType));
1544
1545 ulTotalOutputSize += ulDigiSizeOut;
1546 vulTotalOutputSizeDet[uDetIdx] += ulDigiSizeOut;
1547 }
1548 }
1549 }
1550
1553 if (fRefDet.sName == "Trd1D") {
1554 fvhNbDigiPerEvtDet[uRefDetIdx]->Fill(uNbDataTrd1d);
1555 fvhNbDigiPerEvtTimeDet[uRefDetIdx]->Fill(evt->GetStartTime() * 1e-9, uNbDataTrd1d);
1556
1557 if (0 < GetNofDigis(fRefDet.detId)) {
1559 uint64_t ulDigiSizeOut = GetSizeFromDigisNb(fRefDet.detId, uNbDataTrd1d + uNbDataTrd2d);
1560
1561 ulTotalOutputSize += ulDigiSizeOut;
1562 vulTotalOutputSizeDet[uRefDetIdx] += ulDigiSizeOut;
1563 }
1564 }
1565 else if (fRefDet.sName == "Trd2D") {
1566 fvhNbDigiPerEvtDet[uRefDetIdx]->Fill(uNbDataTrd2d);
1567 fvhNbDigiPerEvtTimeDet[uRefDetIdx]->Fill(evt->GetStartTime() * 1e-9, uNbDataTrd2d);
1568 }
1569 else {
1570 fvhNbDigiPerEvtDet[uRefDetIdx]->Fill(evt->GetNofData(fRefDet.dataType));
1571 fvhNbDigiPerEvtTimeDet[uRefDetIdx]->Fill(evt->GetStartTime() * 1e-9, evt->GetNofData(fRefDet.dataType));
1572
1573 if (0 < GetNofDigis(fRefDet.detId)) {
1575 uint64_t ulDigiSizeOut = GetSizeFromDigisNb(fRefDet.detId, evt->GetNofData(fRefDet.dataType));
1576
1577 ulTotalOutputSize += ulDigiSizeOut;
1578 vulTotalOutputSizeDet[uRefDetIdx] += ulDigiSizeOut;
1579 }
1580 }
1581 }
1582
1583 dPreEvtTime = evt->GetStartTime();
1584 }
1585
1587 for (UInt_t uDetIdx = 0; uDetIdx < fvDets.size(); ++uDetIdx) {
1588 if (fvDets[uDetIdx].sName == "Trd2D") {
1591 continue;
1592 }
1593 uint64_t ulDigiSizeIn = GetSizeFromDigisNb(fvDets[uDetIdx].detId, GetNofDigis(fvDets[uDetIdx].detId));
1594 ulTotalInputSize += ulDigiSizeIn;
1595 vulTotalInputSizeDet[uDetIdx] += ulDigiSizeIn;
1596 }
1598 uint64_t ulDigiSizeIn = GetSizeFromDigisNb(fRefDet.detId, GetNofDigis(fRefDet.detId));
1599 ulTotalInputSize += ulDigiSizeIn;
1600 vulTotalInputSizeDet[uRefDetIdx] += ulDigiSizeIn;
1601 }
1602
1604 for (UInt_t uDetIdx = 0; uDetIdx < fvDets.size(); ++uDetIdx) {
1605 if (0 != vulTotalInputSizeDet[uDetIdx]) { //
1606 fvhSelRatioPerTsNb[uDetIdx]->Fill(fuNrTs, vulTotalOutputSizeDet[uDetIdx] * 1.0 / vulTotalInputSizeDet[uDetIdx]);
1607 }
1608 if (0 != ulTotalInputSize) { //
1609 fvhInpRatioPerTsSz[uDetIdx]->Fill(fuNrTs, vulTotalInputSizeDet[uDetIdx] * 1.0 / ulTotalInputSize);
1610 }
1611 if (0 != ulTotalOutputSize) { //
1612 fvhOutRatioPerTsSz[uDetIdx]->Fill(fuNrTs, vulTotalOutputSizeDet[uDetIdx] * 1.0 / ulTotalOutputSize);
1613 }
1614 }
1617 if (0 != vulTotalInputSizeDet[uRefDetIdx]) { //
1618 fvhSelRatioPerTsNb[uRefDetIdx]->Fill(fuNrTs,
1619 vulTotalOutputSizeDet[uRefDetIdx] * 1.0 / vulTotalInputSizeDet[uRefDetIdx]);
1620 }
1621 if (0 != ulTotalInputSize) { //
1622 fvhInpRatioPerTsSz[uRefDetIdx]->Fill(fuNrTs, vulTotalInputSizeDet[uRefDetIdx] * 1.0 / ulTotalInputSize);
1623 }
1624 if (0 != ulTotalOutputSize) { //
1625 fvhOutRatioPerTsSz[uRefDetIdx]->Fill(fuNrTs, vulTotalOutputSizeDet[uRefDetIdx] * 1.0 / ulTotalOutputSize);
1626 }
1627 }
1629 if (0 != ulTotalInputSize) { //
1630 fhSizeReductionPerTs->Fill(fuNrTs, ulTotalOutputSize * 1.0 / ulTotalInputSize);
1631 }
1632 LOG(debug) << "I/O Size ratio: " << (ulTotalOutputSize * 1.0 / ulTotalInputSize);
1633}
1634
1636{
1637 fhEventTime->Reset();
1638 fhEventDt->Reset();
1639 fhEventSize->Reset();
1640
1641 fhNbDigiPerEvtTime->Reset();
1643 for (std::vector<TH2*>::iterator itHist = fvhNbDigiPerEvtTimeDet.begin(); itHist != fvhNbDigiPerEvtTimeDet.end();
1644 ++itHist) {
1645 (*itHist)->Reset();
1646 }
1647
1648 for (std::vector<TH1*>::iterator itHist = fvhNbDigiPerEvtDet.begin(); itHist != fvhNbDigiPerEvtDet.end(); ++itHist) {
1649 (*itHist)->Reset();
1650 }
1651
1652 for (std::vector<TH1*>::iterator itHist = fvhTDiff.begin(); itHist != fvhTDiff.end(); ++itHist) {
1653 (*itHist)->Reset();
1654 }
1655
1656 /*
1657 if( kTRUE == bResetTime )
1658 {
1660 fdStartTime = -1.0;
1661 }
1662 */
1663}
1664
1665void CbmAlgoBuildRawEvents::SetReferenceDetector(ECbmModuleId refDet, ECbmDataType dataTypeIn, std::string sNameIn,
1666 UInt_t uTriggerMinDigisIn, Int_t iTriggerMaxDigisIn,
1667 Double_t fdTimeWinBegIn, Double_t fdTimeWinEndIn)
1668{
1670 SetReferenceDetector(RawEventBuilderDetector(refDet, dataTypeIn, sNameIn, uTriggerMinDigisIn, iTriggerMaxDigisIn,
1671 fdTimeWinBegIn, fdTimeWinEndIn));
1672}
1673void CbmAlgoBuildRawEvents::AddDetector(ECbmModuleId selDet, ECbmDataType dataTypeIn, std::string sNameIn,
1674 UInt_t uTriggerMinDigisIn, Int_t iTriggerMaxDigisIn, Double_t fdTimeWinBegIn,
1675 Double_t fdTimeWinEndIn)
1676{
1678 AddDetector(RawEventBuilderDetector(selDet, dataTypeIn, sNameIn, uTriggerMinDigisIn, iTriggerMaxDigisIn,
1679 fdTimeWinBegIn, fdTimeWinEndIn));
1680}
1681
1683{
1685 for (std::vector<RawEventBuilderDetector>::iterator det = fvDets.begin(); det != fvDets.end(); ++det) {
1686 if ((*det) == refDetIn) {
1687 LOG(warning) << "CbmAlgoBuildRawEvents::SetReferenceDetector => "
1688 "Reference detector already in selection detector list! "
1689 << refDetIn.sName;
1690 LOG(warning) << " => "
1691 "It will be automatically removed from selection detector list!";
1692 LOG(warning) << " => "
1693 "Please also remember to update the selection windows to store "
1694 "clusters!";
1695 RemoveDetector(refDetIn);
1696 break;
1697 }
1698 }
1699
1700 if (fRefDet == refDetIn) {
1701 LOG(warning) << "CbmAlgoBuildRawEvents::SetReferenceDetector => "
1702 "Doing nothing, identical reference detector already in use";
1703 }
1704 else {
1705 LOG(info) << "CbmAlgoBuildRawEvents::SetReferenceDetector => "
1706 << "Replacing " << fRefDet.sName << " with " << refDetIn.sName << " as reference detector";
1707 LOG(warning) << " => "
1708 "You may want to use AddDetector after this command to add in "
1709 "selection "
1710 << fRefDet.sName;
1711 LOG(warning) << " => "
1712 "Please also remember to update the selection windows!";
1713 }
1714 fRefDet = refDetIn;
1715
1722
1723 if (fbBmonInUse) {
1724 if (select.size()) {
1725 LOG(info) << "CbmAlgoBuildRawEvents::SetReferenceDetector => Detected Bmon station selection.";
1726 fUseBmonMap.assign(select.size(), false);
1727 int it0(0);
1728 for (auto t0 : select)
1729 SwitchBmonStation(it0++, (t0 > 0));
1730 }
1731 }
1732 else {
1733 if (select.size())
1734 LOG(warning) << "CbmAlgoBuildRawEvents::SetReferenceDetector => Detected use of selector\nfor a reference "
1735 "detector which does not support this option. Skip selection for the moment.";
1736 }
1737}
1738
1740{
1741 if (fRefDet == selDet) {
1742 LOG(fatal) << "CbmAlgoBuildRawEvents::AddDetector => Cannot "
1743 "add the reference detector as selection detector!"
1744 << std::endl
1745 << "=> Maybe first change the reference detector with "
1746 "SetReferenceDetector?";
1747 }
1748
1750 for (std::vector<RawEventBuilderDetector>::iterator det = fvDets.begin(); det != fvDets.end(); ++det) {
1751 if ((*det) == selDet) {
1752 LOG(warning) << "CbmAlgoBuildRawEvents::AddDetector => "
1753 "Doing nothing, selection detector already in list! "
1754 << selDet.sName;
1755 return;
1756 }
1757 }
1758 fvDets.push_back(selDet);
1759
1766}
1767
1769{
1771 for (std::vector<RawEventBuilderDetector>::iterator det = fvDets.begin(); det != fvDets.end(); ++det) {
1772 if ((*det) == selDet) {
1773 fvDets.erase(det);
1776 return;
1777 }
1778 }
1779 LOG(warning) << "CbmAlgoBuildRawEvents::RemoveDetector => Doing "
1780 "nothing, selection detector not in list! "
1781 << selDet.sName;
1782}
1783
1785{
1786 if (ECbmModuleId::kBmon == fRefDet.detId) { //
1787 fbBmonInUse = true;
1788 }
1789 else {
1790 fbBmonInUse = std::any_of(fvDets.cbegin(), fvDets.cend(),
1791 [](RawEventBuilderDetector det) { return ECbmModuleId::kBmon == det.detId; });
1792 }
1793}
1794
1796{
1798 if (fRefDet.detId == selDet) {
1800 LOG(debug) << "Set Trigger min limit for " << fRefDet.sName << " to " << uVal;
1801 return;
1802 }
1803
1805 for (std::vector<RawEventBuilderDetector>::iterator det = fvDets.begin(); det != fvDets.end(); ++det) {
1806 if ((*det).detId == selDet) {
1807 (*det).fuTriggerMinDigis = uVal;
1808 LOG(debug) << "Set Trigger min limit for " << (*det).sName << " to " << uVal;
1809 return;
1810 }
1811 }
1812 LOG(warning) << "CbmAlgoBuildRawEvents::SetTriggerMinNumber => "
1813 "Doing nothing, detector neither reference nor in selection list! "
1814 << selDet;
1815}
1816
1818{
1820 if (fRefDet.detId == selDet) {
1822 LOG(debug) << "Set Trigger min limit for " << fRefDet.sName << " to " << iVal;
1823 return;
1824 }
1825
1827 for (std::vector<RawEventBuilderDetector>::iterator det = fvDets.begin(); det != fvDets.end(); ++det) {
1828 if ((*det).detId == selDet) {
1829 (*det).fiTriggerMaxDigis = iVal;
1830 LOG(debug) << "Set Trigger min limit for " << (*det).sName << " to " << iVal;
1831 return;
1832 }
1833 }
1834 LOG(warning) << "CbmAlgoBuildRawEvents::SetTriggerMaxNumber => "
1835 "Doing nothing, detector neither reference nor in selection list! "
1836 << selDet;
1837}
1838
1840{
1842 if (fRefDet.detId == selDet) {
1844 LOG(debug) << "Set Trigger min fired layers limit for " << fRefDet.sName << " to " << uVal;
1845 return;
1846 }
1847
1849 for (std::vector<RawEventBuilderDetector>::iterator det = fvDets.begin(); det != fvDets.end(); ++det) {
1850 if ((*det).detId == selDet) {
1851 (*det).fuTriggerMinLayers = uVal;
1852 LOG(debug) << "Set Trigger min fired layers limit for " << (*det).sName << " to " << uVal;
1853 return;
1854 }
1855 }
1856 LOG(warning) << "CbmAlgoBuildRawEvents::SetTriggerMinLayersNumber => "
1857 "Doing nothing, detector neither reference nor in selection list! "
1858 << selDet;
1859}
1860
1861void CbmAlgoBuildRawEvents::SetTriggerWindow(ECbmModuleId selDet, Double_t dWinBeg, Double_t dWinEnd)
1862{
1864 if (dWinEnd <= dWinBeg) {
1865 LOG(fatal) << "CbmAlgoBuildRawEvents::SetTriggerWindow => "
1866 "Invalid time window: [ "
1867 << dWinBeg << ", " << dWinEnd << " ]";
1868 }
1869
1870 Bool_t bFound = kFALSE;
1872 if (fRefDet.detId == selDet) {
1873 fRefDet.fdTimeWinBeg = dWinBeg;
1874 fRefDet.fdTimeWinEnd = dWinEnd;
1875 bFound = kTRUE;
1876 }
1877
1879 for (std::vector<RawEventBuilderDetector>::iterator det = fvDets.begin(); det != fvDets.end(); ++det) {
1880 if ((*det).detId == selDet) {
1881 (*det).fdTimeWinBeg = dWinBeg;
1882 (*det).fdTimeWinEnd = dWinEnd;
1883 bFound = kTRUE;
1884 }
1885 }
1886
1887 if (kFALSE == bFound) {
1888 LOG(warning) << "CbmAlgoBuildRawEvents::SetTriggerWindow => "
1889 "Doing nothing, detector neither reference nor in selection list! "
1890 << selDet;
1891 }
1892
1897}
1898
1900{
1902 if (fRefDet.detId == selDet) {
1903 fRefDet.fdHistMaxDigiNb = dVal;
1904 LOG(debug) << "Set histogram max digi nb for " << fRefDet.sName << " to " << dVal;
1905 return;
1906 }
1907
1909 for (std::vector<RawEventBuilderDetector>::iterator det = fvDets.begin(); det != fvDets.end(); ++det) {
1910 if ((*det).detId == selDet) {
1911 (*det).fdHistMaxDigiNb = dVal;
1912 LOG(debug) << "Set histogram max digi nb " << (*det).sName << " to " << dVal;
1913 return;
1914 }
1915 }
1916 LOG(warning) << "CbmAlgoBuildRawEvents::SetHistogramMaxDigiNb => "
1917 "Doing nothing, detector neither reference nor in selection list! "
1918 << selDet;
1919}
1920
1922{
1927 }
1928 else {
1931 }
1932
1934 for (std::vector<RawEventBuilderDetector>::iterator det = fvDets.begin(); det != fvDets.end(); ++det) {
1935 fdEarliestTimeWinBeg = std::min(fdEarliestTimeWinBeg, (*det).fdTimeWinBeg);
1936 fdLatestTimeWinEnd = std::max(fdLatestTimeWinEnd, (*det).fdTimeWinEnd);
1937 }
1938}
1939
1941{
1944
1946 for (std::vector<RawEventBuilderDetector>::iterator det = fvDets.begin(); det != fvDets.end(); ++det) {
1947 fdWidestTimeWinRange = std::max(fdWidestTimeWinRange, (*det).fdTimeWinEnd - (*det).fdTimeWinBeg);
1948 }
1949}
1950
1952{
1953 if (id < 0 || id >= (int) fUseBmonMap.size()) {
1954 LOG(warning) << "CbmAlgoBuildRawEvents::SwitchBmonStation: Bmon station id outside range. Skip.";
1955 return;
1956 }
1957 LOG(info) << "CbmAlgoBuildRawEvents::SwitchBmonStation: Bmon" << id << " station switched " << (on ? "ON" : "OFF")
1958 << " for trigger.";
1959 fUseBmonMap[id] = on;
1960}
1961
1963{
1964 // skip any filtering if not explicitly requested by user.
1965 if (!fUseBmonMap.size()) return true;
1966
1967 // consistency check on the Bmon detector
1968 if (CbmTofAddress::GetSmType(add) != 5) {
1969 LOG(warning) << "CbmAlgoBuildRawEvents::filterBmon: Bmon digi with wrong address [GetSmType() != 5]. Skip.";
1970 return false;
1971 }
1972 size_t mod = (size_t) CbmTofAddress::GetChannelSide(add);
1973 // select Bmon detector
1974 if (!fUseBmonMap[mod]) {
1975 LOG(debug2) << "CbmAlgoBuildRawEvents::filterBmon : reject seed from Bmon" << mod;
1976 return false;
1977 }
1978 return true;
1979}
1980
1982{
1983 // skip any filtering if not explicitly requested by user.
1984 if (!fUseBmonMap.size()) return event->GetNofData(ECbmDataType::kBmonDigi);
1985
1986 int32_t ndigi(0);
1987 for (size_t idigi = 0; idigi < event->GetNofData(ECbmDataType::kBmonDigi); ++idigi) {
1988 uint idx = event->GetIndex(ECbmDataType::kBmonDigi, idigi);
1989 const CbmBmonDigi* pDigi = GetDigi<CbmBmonDigi>(idx);
1990 if (!filterBmon(pDigi->GetAddress())) continue;
1991 ndigi++;
1992 }
1993 return ndigi;
1994}
1995
ClassImp(CbmConverterManager)
ECbmDataType
Definition CbmDefs.h:90
ECbmModuleId
Definition CbmDefs.h:39
@ kTrd
Transition Radiation Detector.
@ kTof
Time-of-flight Detector.
@ kNotExist
If not found.
@ kPsd
Projectile spectator detector.
@ kSts
Silicon Tracking System.
@ kTrd2d
TRD-FASP Detector (FIXME)
@ kMuch
Muon detection system.
@ kFsd
Forward spectator detector.
@ kRich
Ring-Imaging Cherenkov Detector.
@ kStsModule
@ kStsUnit
static double dTsStartTime
Helper class to convert unique channel ID back and forth.
Bool_t CheckTriggerConditions(CbmEvent *event, const RawEventBuilderDetector &det)
UInt_t fuCurEv
histogram with size of overlap between evt vs TS index, AllowOverlap only
const std::vector< CbmTrdDigi > * fTrdDigis
const std::vector< CbmBmonDigi > * fBmonDigis
void AddDigiToEvent(const RawEventBuilderDetector &det, Int_t uIdx)
TH1 * fhCpuTimePerTsHist
Processing time per TS.
std::vector< TH1 * > fvhOutRatioPerTsSz
ratio of input digi size in total input size vs TS in run
void SetHistogramMaxDigiNb(ECbmModuleId selDet, Double_t dDigiNbMax)
int32_t getNofFilteredBmonDigis(CbmEvent *ev)
TH2 * fhOverEventSizeTs
histogram with proportion of overlap evt vs TS index, AllowOverlap only
uint64_t GetSizeFromDigisNb(ECbmModuleId detId, uint64_t ulNbDigis)
const std::vector< Double_t > * fSeedTimes
EOverlapModeRaw fOverMode
bit map for Bmon trigger. Defined by user
TH2 * fhNbDigiPerEvtTime
histogram with the nb of all digis in the event
std::vector< TH2 * > fvhNbDigiPerEvtTimeDet
Plotting time per TS.
TH1 * fhEventSize
histogram with the interval in seed time of consecutive events
void CheckTriggerCondition(Double_t dSeedTime)
Double_t fdPrevEvtEndTime
Save previous time information.
std::vector< TH1 * > fvhSelRatioPerTsNb
void SwitchBmonStation(int id, bool on=true)
Bool_t fbUseMuchBeamtimeDigi
Switch ON/OFF filling of histograms.
const std::vector< CbmTofDigi > * fTofDigis
CbmEvent * fCurrentEvent
Data ouptut.
std::vector< RawEventBuilderDetector > fvDets
TH1 * fhCpuTimePerTs
histogram with the nb of all digis per event vs seed time of the events
UInt_t GetNofDigis(ECbmModuleId detId)
UInt_t fuNrTs
Event Counter.
Double_t fdSeedWindowBeg
Seed window.
void RemoveDetector(RawEventBuilderDetector selDet)
void AddDetector(ECbmModuleId selDet, ECbmDataType dataTypeIn, std::string sNameIn, UInt_t uTriggerMinDigisIn=0, Int_t iTriggerMaxDigisIn=-1, Double_t fdTimeWinBegIn=-100, Double_t fdTimeWinEndIn=100)
std::vector< CbmEvent * > fEventVector
pointer to the event which is currently build
Bool_t fbFillHistos
Ignore data in Overlap part of the TS.
void AddHistoToVector(TNamed *pointer, std::string sFolder="")
For monitor algos.
void SearchMatches(Double_t dSeedTime, RawEventBuilderDetector &detMatch)
TH1 * fhSizeReductionPerTs
ratio of selected digi size in total event size vs TS in run
void SetTriggerMaxNumber(ECbmModuleId selDet, Int_t iVal)
const std::vector< CbmRichDigi > * fRichDigis
TH1 * fhRealTimePerTsHist
Plotting time per TS.
const std::vector< CbmStsDigi > * fStsDigis
TH1 * fhEventTime
Vector of pointers to canvases + optional folder name.
Bool_t fbGetTimings
Switch between MUCH digi classes.
const std::vector< CbmMuchBeamTimeDigi > * fMuchBeamTimeDigis
std::vector< TH1 * > fvhTDiff
histograms with the nb of digis in each detector per event
const std::vector< CbmMuchDigi > * fMuchDigis
void ResetHistograms(Bool_t bResetTime=kTRUE)
TH1 * fhRealTimePerTs
Processing time per TS.
TProfile * fhOverEventShareTs
histogram with proportion of overlap evt, AllowOverlap only
bool filterBmon(int32_t add)
Filter Bmon stations. Hack added for the mCBM2024 data (AB)
void AddCanvasToVector(TCanvas *pointer, std::string sFolder="")
Bool_t fbUseTsMetaData
Measure CPU time using stopwatch.
std::vector< TH1 * > fvhNbDigiPerEvtDet
histograms with the nb of digis in each detector per event vs seed time of the events
void CheckSeed(Double_t dSeedTime, UInt_t uSeedDigiIdx)
void SetTriggerMinLayersNumber(ECbmModuleId selDet, UInt_t uVal)
Double_t fdPrevEvtTime
Timeslice Counter.
void SetTriggerMinNumber(ECbmModuleId selDet, UInt_t uVal)
TH1 * fhOverEventShare
ratio of total selected size to input size selected vs TS in run
TClonesArray * fTimeSliceMetaDataArray
Data input.
Bool_t CheckDataAvailable(const RawEventBuilderDetector &det)
Internal methods.
void SetReferenceDetector(ECbmModuleId refDet, ECbmDataType dataTypeIn, std::string sNameIn, UInt_t uTriggerMinDigisIn=0, Int_t iTriggerMaxDigisIn=-1, Double_t fdTimeWinBegIn=-100, Double_t fdTimeWinEndIn=100)
bool SetBmonEventTime(CbmEvent *event)
const std::vector< CbmPsdDigi > * fPsdDigis
RawEventBuilderDetector fRefDet
is create when fbGetTimings is set before init
void SetTriggerWindow(ECbmModuleId selDet, Double_t dWinBeg, Double_t dWinEnd)
const Digi * GetDigi(UInt_t uDigi)
const std::vector< CbmFsdDigi > * fFsdDigis
TH1 * fhEventDt
histogram with the seed time of the events
std::vector< bool > fUseBmonMap
Read Ts Parameters from input tree.
std::vector< TH1 * > fvhInpRatioPerTsSz
ratio of selected/input digi vs TS in run
Data class for a signal in the t-zero detector.
Definition CbmBmonDigi.h:30
double GetTime() const
Time.
Definition CbmBmonDigi.h:86
int32_t GetAddress() const
Address.
Definition CbmBmonDigi.h:80
Class characterising one event by a collection of links (indices) to data objects,...
Definition CbmEvent.h:34
void AddData(ECbmDataType type, uint32_t index)
Definition CbmEvent.cxx:33
void SetEndTime(double endTime)
Definition CbmEvent.h:157
Data class for FSD digital information.
Definition CbmFsdDigi.h:36
Data class for PSD digital information.
Definition CbmPsdDigi.h:36
Data class for a single-channel message in the STS.
Definition CbmStsDigi.h:40
XPU_D uint16_t GetChannel() const
Channel number in module @value Channel number.
Definition CbmStsDigi.h:93
XPU_D int32_t GetAddress() const
Definition CbmStsDigi.h:74
static int32_t GetStripFullId(uint32_t address)
static int32_t GetSmType(uint32_t address)
static int32_t GetRpcFullId(uint32_t address)
static int32_t GetChannelSide(uint32_t address)
Data class for expanded digital TOF information.
Definition CbmTofDigi.h:47
double GetSide() const
Channel Side.
Definition CbmTofDigi.h:160
int32_t GetAddress() const
Inherited from CbmDigi.
Definition CbmTofDigi.h:112
eCbmTrdAsicType GetType() const
Channel FEE SPADIC/FASP according to CbmTrdAsicType.
Definition CbmTrdDigi.h:173
Int_t fiTriggerMaxDigis
Maximum number of digis per detector needed to generate an event, -1 means no cut,...
UInt_t fuTriggerMinLayers
Minimum number of fired layers needed to generate an event, 0 means do not require for event selectio...
ECbmModuleId detId
Settings.
Double_t fdTimeWinBeg
Selection Window.
UInt_t fuTriggerMinDigis
Minimum number of digis per detector needed to generate an event, 0 means do not use for event select...
Double_t fdHistMaxDigiNb
Histo configuration.
UInt_t fuStartIndex
Book-keeping variables.
uint64_t GetOverlapDuration() const
uint64_t GetOverlapStartTime() const
uint64_t GetStartTime() const
int32_t GetMotherAddress(int32_t address, int32_t level)
Construct the address of an element from the address of a descendant element.
uint32_t GetElementId(int32_t address, int32_t level)
Get the index of an element.