CbmRoot
Loading...
Searching...
No Matches
CbmTaskBuildRawEvents.cxx
Go to the documentation of this file.
1/* Copyright (C) 2007-2021 Facility for Antiproton and Ion Research in Europe, Darmstadt
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Dominik Smith [committer] */
4
6
7#include "CbmDigiManager.h"
8#include "CbmEvent.h"
9#include "CbmMatch.h"
10#include "CbmModuleList.h"
12
13#include <FairRootManager.h>
14#include <FairRunOnline.h>
15#include <Logger.h>
16
17#include <TCanvas.h>
18#include <TClonesArray.h>
19#include <TFile.h>
20#include <TFolder.h>
21#include <TH1.h>
22#include <TH2.h>
23#include <THttpServer.h>
24#include <TStopwatch.h>
25
26#include <iomanip>
27#include <iostream>
28
30{
31 if (fpAlgo) delete fpAlgo;
33 if (fSeedTimes) delete fSeedTimes;
35 if (fTimer) delete fTimer;
36 if (fCopyTimer) delete fCopyTimer;
37 if (fDigiEvents) {
38 fDigiEvents->clear();
39 delete fDigiEvents;
40 }
41 else if (fEvents) {
42 fEvents->Delete();
43 delete fEvents;
44 }
45}
46
48{
51}
52
54{
55 fSeedTimeDetList.push_back(seedDet);
56
57 if (fSeedTimes == nullptr) {
58 fSeedTimes = new std::vector<Double_t>;
59 }
60 fpAlgo->SetSeedTimes(fSeedTimes);
61}
62
64{
65 SetSlidingWindowSeedFinder(1, 0.0, 0.0, 0.0);
66 fSeedFinderSlidingWindow->SetIdealMode(fileId);
67}
68
69void CbmTaskBuildRawEvents::SetSlidingWindowSeedFinder(int32_t minDigis, double dWindDur, double dDeadT, double dOffset)
70{
74 }
75 if (fSeedTimes == nullptr) {
76 fSeedTimes = new std::vector<Double_t>;
77 }
78 fpAlgo->SetSeedTimes(fSeedTimes);
79
80 fSeedFinderSlidingWindow = new CbmSeedFinderSlidingWindow(fSeedTimes, minDigis, dWindDur, dDeadT);
81 fSeedFinderSlidingWindow->SetOffset(dOffset);
82 // Set seed time window/deadtime for event building checks
83 fpAlgo->SetSeedTimeWindow(0.0, dDeadT < dWindDur ? dDeadT : dWindDur);
84}
85
87{
88 if (bFlagIn == kTRUE) {
89 if (fSeedFinderSlidingWindow == nullptr) {
90 std::cout << "SetSeedFinderQa(): Cannot activate Qa when seed finder not active. Exiting." << std::endl;
91 exit(1);
92 }
93 if (fvDigiMatchQa == nullptr) {
94 fvDigiMatchQa = new std::vector<CbmMatch>;
95 }
96 }
97 else //kFALSE
98 {
99 if (fvDigiMatchQa != nullptr) {
100 delete fvDigiMatchQa;
101 }
102 fvDigiMatchQa = nullptr;
103 }
104 if (fSeedFinderSlidingWindow != nullptr) fSeedFinderSlidingWindow->SetQa(bFlagIn);
105}
106
107template<class TDigi>
108void CbmTaskBuildRawEvents::InitDigis(ECbmModuleId detId, std::vector<TDigi>** vDigi)
109{
110 TString detName = CbmModuleList::GetModuleNameCaps(detId);
111 if (!fDigiMan->IsPresent(detId)) {
112 LOG(info) << "No " << detName << " digi input.";
113 }
114 else {
115 LOG(info) << detName << " digi input.";
116 *vDigi = new std::vector<TDigi>;
117 fpAlgo->SetDigis(*vDigi);
118 }
119}
120
122{
123 if (fbGetTimings) {
124 fTimer = new TStopwatch;
125 fTimer->Start();
126 fCopyTimer = new TStopwatch;
127 fCopyTimer->Reset();
128 }
129
131 FairRootManager* ioman = FairRootManager::Instance();
132
133 // Get a pointer to the previous already existing data level
136 fDigiMan->UseMuchBeamTimeDigi();
137 }
138 fDigiMan->Init();
139
140 //Init digis
143 }
144 else {
146 }
154
156 if (fbDigiEvtOut) {
157 if (fbUseMuchBeamtimeDigi) LOG(fatal) << "DigiEvent output branch not compatible with MuchBeamtimeDigi";
158
159 fDigiEvents = new std::vector<CbmDigiEvent>();
160 ioman->RegisterAny("DigiEvent", fDigiEvents, kTRUE);
161 if (!fDigiEvents) LOG(fatal) << "Output branch was not created";
162 LOG(info) << "DigiEvent out instead of CbmEvent, you will need an instance of CbmTaskMakeRecoEvents in reco macro";
164 LOG(info) << "Exclusive TRD extraction, DigiEvents will be comparable to CbmEvents but slower";
165 }
166 else {
167 LOG(info) << "Inclusive TRD extraction, faster but DigiEvents will noy be comparable to CbmEvents (extra digis)";
168 }
169 }
170 else {
171 fEvents = new TClonesArray("CbmEvent", 100);
172 ioman->Register("CbmEvent", "Cbm_Event", fEvents, IsOutputBranchPersistent("CbmEvent"));
173 if (!fEvents) LOG(fatal) << "Output branch was not created";
174 LOG(info) << "CbmEvent oupput, you will need an instance of CbmTaskEventsCloneInToOut in reco macro to update them";
175 }
176
177 // Set timeslice meta data
178 fpAlgo->SetTimeSliceMetaDataArray(dynamic_cast<TClonesArray*>(ioman->GetObject("TimesliceMetaData")));
179
180 if (fTimer != nullptr) {
181 fTimer->Stop();
182 Double_t rtime = fTimer->RealTime();
183 Double_t ctime = fTimer->CpuTime();
184 LOG(info) << "CbmTaskBuildRawEvents::Init(): Real time " << rtime << " s, CPU time " << ctime << " s";
185 }
186
187 // Init seed finder
190 }
191
193 if (kTRUE == fpAlgo->InitAlgo())
194 return kSUCCESS;
195 else
196 return kFATAL;
197}
198
199InitStatus CbmTaskBuildRawEvents::ReInit() { return kSUCCESS; }
200
201template<class TDigi>
202void CbmTaskBuildRawEvents::ReadDigis(ECbmModuleId detId, std::vector<TDigi>* vDigis)
203{
204 //Warning: Int_t must be used for the loop counters instead of UInt_t,
205 //as the digi manager can return -1, which would be casted to +1
206 //during comparison, leading to an error.
207 if (fCopyTimer != nullptr) {
208 fCopyTimer->Start(kFALSE);
209 }
210
211 const TString detName = CbmModuleList::GetModuleNameCaps(detId);
212
213 if (fDigiMan->IsPresent(detId)) {
214 vDigis->clear();
215 const Int_t nDigis = fDigiMan->GetNofDigis(detId);
216
217 for (Int_t i = 0; i < nDigis; i++) {
218 const TDigi* Digi = fDigiMan->Get<TDigi>(i);
219 vDigis->push_back(*Digi);
220 }
221 LOG(debug) << "Read: " << vDigis->size() << " " << detName << " digis.";
222 LOG(debug) << "In DigiManager: " << nDigis << " " << detName << " digis.";
223 }
224
225 if (fCopyTimer != nullptr) {
226 fCopyTimer->Stop();
227 }
228}
229
230void CbmTaskBuildRawEvents::Exec(Option_t* /*option*/)
231{
232 if (fTimer != nullptr) {
233 fTimer->Start(kFALSE);
234 }
235 TStopwatch timer;
236 timer.Start();
237 LOG(debug2) << "CbmTaskBuildRawEvents::Exec => Starting sequence";
238
239 // Process Timeslice
240 BuildEvents();
241 if (fTimer != nullptr) {
242 fTimer->Stop();
243 }
244
245 // --- Timeslice log and statistics
246 timer.Stop();
247 std::stringstream logOut;
248 logOut << std::setw(20) << std::left << GetName() << " [";
249 logOut << std::fixed << std::setw(8) << std::setprecision(1) << std::right << timer.RealTime() * 1000. << " ms] ";
250 logOut << "TS " << fNofTs;
251 if (fbDigiEvtOut) {
252 logOut << ", events " << fDigiEvents->size();
253 fNofEvents += fDigiEvents->size();
254 }
255 else {
256 logOut << ", events " << fEvents->GetEntriesFast();
257 fNofEvents += fEvents->GetEntriesFast();
258 }
259 LOG(info) << logOut.str();
261 const size_t seedCount = fSeedFinderSlidingWindow->GetNofSeeds();
262 LOG(info) << seedCount << " trigger(s) for this TS.";
263 fTotalSeedCount += seedCount;
264 }
265 fNofTs++;
266 fTime += timer.RealTime();
267
268 LOG(debug2) << "CbmTaskBuildRawEvents::Exec => Done";
269}
270
272{
273 //Reset explicit seed times if set
274 if (fSeedTimes != nullptr) {
275 fSeedTimes->clear();
276 }
277
278 //Read digis
281 }
282 else {
284 }
292
293 //Fill seeds
294 if (fSeedFinderSlidingWindow != nullptr) {
296 }
297 else if (fSeedTimeDetList.size() > 0) {
299 }
300 //DumpSeedTimesFromDetList();
301
303 fpAlgo->ProcessTs();
304
306 FillOutput();
307}
308
309void CbmTaskBuildRawEvents::FillSeedTimesFromDetList(std::vector<Double_t>* vdSeedTimes,
310 std::vector<CbmMatch>* vDigiMatch)
311{
312 std::map<ECbmModuleId, UInt_t> DigiNumbers;
313 std::map<ECbmModuleId, UInt_t> DigiCounters;
314 vdSeedTimes->clear();
315
316 if (vDigiMatch != nullptr) vDigiMatch->clear();
317
319 DigiNumbers[system.detId] = GetNofDigis(system.detId);
320 DigiCounters[system.detId] = 0;
321 }
322
323 do {
324 ECbmModuleId nextAddedSystem;
325 Double_t earliestTime = -1;
326
328 if (DigiCounters[system.detId] < DigiNumbers[system.detId]) {
329
330 Double_t thisTime = GetDigiTime(system.detId, DigiCounters[system.detId]);
331 if (thisTime < earliestTime || earliestTime == -1) {
332 nextAddedSystem = system.detId;
333 earliestTime = thisTime;
334 }
335 }
336 }
337 if (earliestTime != -1) {
338
339 if (vDigiMatch != nullptr) {
340 const CbmMatch* digiMatch = fDigiMan->GetMatch(nextAddedSystem, DigiCounters[nextAddedSystem]);
341 vDigiMatch->push_back(*digiMatch);
342 }
343 vdSeedTimes->push_back(earliestTime);
344 DigiCounters[nextAddedSystem]++;
345 }
346 else {
347 break;
348 }
349 } while (true);
350}
351
353{
354 if (fSeedTimeDetList.size() == 0) {
355 if (fSeedFinderSlidingWindow->IsIdealMode()) {
356 fSeedFinderSlidingWindow->FillSeedTimes();
357 return;
358 }
359 else {
360 std::cout << "FillSeedTimesFromSlidingWindow(): Error, seed detector list empty." << std::endl;
361 exit(1);
362 }
363 }
364 if (fSeedTimeDetList.size() == 1) {
365 const RawEventBuilderDetector seedDet = fSeedTimeDetList.at(0);
367 }
368 else { // more than one seed detector
369 if (!fTempDigiTimes) {
370 fTempDigiTimes = new std::vector<Double_t>;
371 }
374 }
375}
376
378{
379 if (fvDigiMatchQa != nullptr) {
380 if (!fDigiMan->IsMatchPresent(seedDet->detId)) {
381 std::cout << "FillSeedTimesFromSlidingWindow(): Error, QA set but no CbmMatch found for seed detector."
382 << std::endl;
383 exit(1);
384 }
385 fvDigiMatchQa->clear();
386 for (Int_t i = 0; i < fDigiMan->GetNofDigis(seedDet->detId); i++) {
387 const CbmMatch* digiMatch = fDigiMan->GetMatch(seedDet->detId, i);
388 fvDigiMatchQa->push_back(*digiMatch);
389 }
390 }
391
392 switch (seedDet->detId) {
396 break;
397 }
398 else {
400 break;
401 }
409 default: break;
410 }
411}
412
413Double_t CbmTaskBuildRawEvents::GetDigiTime(ECbmModuleId _system, UInt_t _entry)
414{
415 switch (_system) {
418 return (fMuchBeamTimeDigis->at(_entry)).GetTime();
419 }
420 else {
421 return (fMuchDigis->at(_entry)).GetTime();
422 }
423 case ECbmModuleId::kSts: return (fStsDigis->at(_entry)).GetTime();
424 case ECbmModuleId::kTrd: return (fTrdDigis->at(_entry)).GetTime();
425 case ECbmModuleId::kTof: return (fTofDigis->at(_entry)).GetTime();
426 case ECbmModuleId::kRich: return (fRichDigis->at(_entry)).GetTime();
427 case ECbmModuleId::kPsd: return (fPsdDigis->at(_entry)).GetTime();
428 case ECbmModuleId::kFsd: return (fFsdDigis->at(_entry)).GetTime();
429 case ECbmModuleId::kBmon: return (fBmonDigis->at(_entry)).GetTime();
430 default: break;
431 }
432 return -1;
433}
434
436{
437 switch (_system) {
440 return fMuchBeamTimeDigis->size();
441 }
442 else {
443 return fMuchDigis->size();
444 }
445 case ECbmModuleId::kSts: return fStsDigis->size();
446 case ECbmModuleId::kTrd: return fTrdDigis->size();
447 case ECbmModuleId::kTof: return fTofDigis->size();
448 case ECbmModuleId::kRich: return fRichDigis->size();
449 case ECbmModuleId::kPsd: return fPsdDigis->size();
450 case ECbmModuleId::kFsd: return fFsdDigis->size();
451 case ECbmModuleId::kBmon: return fBmonDigis->size();
452 default: break;
453 }
454 return 0;
455}
456
458{
459 if (fTimer == nullptr) {
460 LOG(fatal) << "Trying to print timings but timer not set";
461 }
462 else {
463 Double_t rtime = fTimer->RealTime();
464 Double_t ctime = fTimer->CpuTime();
465 LOG(info) << "CbmTaskBuildRawEvents: Real time " << rtime << " s, CPU time " << ctime << " s";
466 }
467 if (fCopyTimer == nullptr) {
468 LOG(fatal) << "Trying to print timings but timer not set";
469 }
470 else {
471 Double_t rtime = fCopyTimer->RealTime();
472 Double_t ctime = fCopyTimer->CpuTime();
473 LOG(info) << "CbmTaskBuildRawEvents (digi copy only): Real time " << rtime << " s, CPU time " << ctime << " s";
474 }
475}
476
478{
479 if ((fvDigiMatchQa != nullptr) && (fSeedFinderSlidingWindow != nullptr)) {
480 fSeedFinderSlidingWindow->OutputQa();
481 }
482
484 fpAlgo->Finish();
485 if (fbFillHistos) {
486 SaveHistos();
487 }
488 if (fbGetTimings) {
489 PrintTimings();
490 }
491
492 std::cout << std::endl;
493 LOG(info) << "=====================================";
494 LOG(info) << GetName() << ": Run summary";
495 LOG(info) << "Time slices : " << fNofTs;
496 LOG(info) << "Events : " << fNofEvents;
498 LOG(info) << "Triggers : " << fTotalSeedCount;
499 }
500 LOG(info) << "Time / TS : " << std::fixed << std::setprecision(2) << 1000. * fTime / Double_t(fNofTs)
501 << " ms";
502 LOG(info) << "=====================================";
503}
504
506{
508 std::vector<CbmEvent*> vEvents = fpAlgo->GetEventVector();
509
510 if (fbDigiEvtOut) {
512 fDigiEvents->clear();
513
515 ExtractSelectedData(vEvents);
516 }
517 else {
519 fEvents->Delete();
520
522 for (CbmEvent* event : vEvents) {
523 LOG(debug) << "Vector: " << event->ToString();
524 new ((*fEvents)[fEvents->GetEntriesFast()]) CbmEvent(std::move(*event));
525 LOG(debug) << "TClonesArray: " << static_cast<CbmEvent*>(fEvents->At(fEvents->GetEntriesFast() - 1))->ToString();
526 }
527 }
529 fpAlgo->ClearEventVector();
530}
531
533{
535 if (!FairRootManager::Instance() || !FairRootManager::Instance()->GetSink()) {
536 LOG(error) << "No sink found";
537 return;
538 }
539 FairSink* sink = FairRootManager::Instance()->GetSink();
540 sink->WriteObject(dynamic_cast<TObject*>(fpAlgo->GetOutFolder()), nullptr);
541 }
542 else {
543
545 std::vector<std::pair<TNamed*, std::string>> vHistos = fpAlgo->GetHistoVector();
546
548 TDirectory* oldDir = NULL;
549 TFile* histoFile = NULL;
551 oldDir = gDirectory;
553 histoFile = new TFile(fsOutFileName, "RECREATE");
554 histoFile->cd();
555
557 for (UInt_t uHisto = 0; uHisto < vHistos.size(); ++uHisto) {
559 const TString sFolder = vHistos[uHisto].second.data();
560 if (nullptr == gDirectory->Get(sFolder)) gDirectory->mkdir(sFolder);
561 gDirectory->cd(sFolder);
562
564 vHistos[uHisto].first->Write();
565 histoFile->cd();
566 }
567
569 std::vector<std::pair<TCanvas*, std::string>> vCanvases = fpAlgo->GetCanvasVector();
570
572 for (UInt_t uCanvas = 0; uCanvas < vCanvases.size(); ++uCanvas) {
574 TString sFolder = vCanvases[uCanvas].second.data();
575 if (nullptr == gDirectory->Get(sFolder)) gDirectory->mkdir(sFolder);
576 gDirectory->cd(sFolder);
577
579 vCanvases[uCanvas].first->Write();
580
581 histoFile->cd();
582 }
583
585 oldDir->cd();
586 histoFile->Close();
587 }
588}
589
591{
592 std::ofstream timesUnsorted("digiTimesUnsorted.dat", std::ofstream::out);
593 timesUnsorted << std::setprecision(16);
594
596 for (UInt_t i = 0; i < GetNofDigis(system.detId); i++) {
597 timesUnsorted << GetDigiTime(system.detId, i) << std::endl;
598 }
599 }
600 timesUnsorted.close();
601 LOG(info) << "Completed write of unsorted digi list.";
602
603 std::ofstream timesSorted("digiTimesSorted.dat", std::ofstream::out);
604 timesSorted << std::setprecision(16);
605
606 for (UInt_t i = 0; i < fSeedTimes->size(); i++) {
607 timesSorted << fSeedTimes->at(i) << std::endl;
608 }
609 timesSorted.close();
610 LOG(info) << "Completed DumpSeedTimesFromDetList(). Closing.";
611 exit(0); //terminate as this method should only be used for diagnostics
612}
613
614void CbmTaskBuildRawEvents::ExtractSelectedData(std::vector<CbmEvent*> vEvents)
615{
617 LOG(debug) << "In Vector size: " << vEvents.size();
618
619 fDigiEvents->reserve(vEvents.size());
620 for (CbmEvent* event : vEvents) {
621 CbmDigiEvent selEvent;
622 selEvent.fTime = event->GetStartTime();
623 selEvent.fNumber = event->GetNumber();
624
627 event->SortIndices();
628
632 uint32_t uNbDigis =
633 (0 < event->GetNofData(ECbmDataType::kBmonDigi) ? event->GetNofData(ECbmDataType::kBmonDigi) : 0);
634 if (0 < uNbDigis) {
635 auto startIt = fBmonDigis->begin() + event->GetIndex(ECbmDataType::kBmonDigi, 0);
636 auto stopIt = fBmonDigis->begin() + event->GetIndex(ECbmDataType::kBmonDigi, uNbDigis - 1);
637 ++stopIt;
638 selEvent.fData.fBmon.fDigis.assign(startIt, stopIt);
639 }
640
642 uNbDigis = (0 < event->GetNofData(ECbmDataType::kStsDigi) ? event->GetNofData(ECbmDataType::kStsDigi) : 0);
643 if (0 < uNbDigis) {
644 auto startIt = fStsDigis->begin() + event->GetIndex(ECbmDataType::kStsDigi, 0);
645 auto stopIt = fStsDigis->begin() + event->GetIndex(ECbmDataType::kStsDigi, uNbDigis - 1);
646 ++stopIt;
647 selEvent.fData.fSts.fDigis.assign(startIt, stopIt);
648 }
649
651 uNbDigis = (0 < event->GetNofData(ECbmDataType::kMuchDigi) ? event->GetNofData(ECbmDataType::kMuchDigi) : 0);
652 if (0 < uNbDigis) {
653 auto startIt = fMuchDigis->begin() + event->GetIndex(ECbmDataType::kMuchDigi, 0);
654 auto stopIt = fMuchDigis->begin() + event->GetIndex(ECbmDataType::kMuchDigi, uNbDigis - 1);
655 ++stopIt;
656 selEvent.fData.fMuch.fDigis.assign(startIt, stopIt);
657 }
658
660 uNbDigis = (0 < event->GetNofData(ECbmDataType::kTrdDigi) ? event->GetNofData(ECbmDataType::kTrdDigi) : 0);
661 if (0 < uNbDigis) {
663 for (uint32_t uDigiInEvt = 0; uDigiInEvt < uNbDigis; ++uDigiInEvt) {
668 selEvent.fData.fTrd.fDigis.push_back(fTrdDigis->at(event->GetIndex(ECbmDataType::kTrdDigi, uDigiInEvt)));
669 }
670 }
671 else {
676 auto startIt = fTrdDigis->begin() + event->GetIndex(ECbmDataType::kTrdDigi, 0);
677 auto stopIt = fTrdDigis->begin() + event->GetIndex(ECbmDataType::kTrdDigi, uNbDigis - 1);
678 ++stopIt;
679 selEvent.fData.fTrd.fDigis.assign(startIt, stopIt);
680 }
681 }
682
684 uNbDigis = (0 < event->GetNofData(ECbmDataType::kTofDigi) ? event->GetNofData(ECbmDataType::kTofDigi) : 0);
685 if (0 < uNbDigis) {
686 auto startIt = fTofDigis->begin() + event->GetIndex(ECbmDataType::kTofDigi, 0);
687 auto stopIt = fTofDigis->begin() + event->GetIndex(ECbmDataType::kTofDigi, uNbDigis - 1);
688 ++stopIt;
689 selEvent.fData.fTof.fDigis.assign(startIt, stopIt);
690 }
691
693 uNbDigis = (0 < event->GetNofData(ECbmDataType::kRichDigi) ? event->GetNofData(ECbmDataType::kRichDigi) : 0);
694 if (0 < uNbDigis) {
695 auto startIt = fRichDigis->begin() + event->GetIndex(ECbmDataType::kRichDigi, 0);
696 auto stopIt = fRichDigis->begin() + event->GetIndex(ECbmDataType::kRichDigi, uNbDigis - 1);
697 ++stopIt;
698 selEvent.fData.fRich.fDigis.assign(startIt, stopIt);
699 }
700
702 uNbDigis = (0 < event->GetNofData(ECbmDataType::kPsdDigi) ? event->GetNofData(ECbmDataType::kPsdDigi) : 0);
703 if (0 < uNbDigis) {
704 auto startIt = fPsdDigis->begin() + event->GetIndex(ECbmDataType::kPsdDigi, 0);
705 auto stopIt = fPsdDigis->begin() + event->GetIndex(ECbmDataType::kPsdDigi, uNbDigis - 1);
706 ++stopIt;
707 selEvent.fData.fPsd.fDigis.assign(startIt, stopIt);
708 }
709
710 fDigiEvents->push_back(std::move(selEvent));
711 }
712}
713
ClassImp(CbmConverterManager)
ECbmModuleId
Enumerator for module Identifiers.
Definition CbmDefs.h:45
@ kTrd
Transition Radiation Detector.
Definition CbmDefs.h:51
@ kTof
Time-of-flight Detector.
Definition CbmDefs.h:52
@ kPsd
Projectile spectator detector.
Definition CbmDefs.h:54
@ kSts
Silicon Tracking System.
Definition CbmDefs.h:48
@ kMuch
Muon detection system.
Definition CbmDefs.h:50
@ kFsd
Forward spectator detector.
Definition CbmDefs.h:59
@ kBmon
Bmon Counter.
Definition CbmDefs.h:57
@ kRich
Ring-Imaging Cherenkov Detector.
Definition CbmDefs.h:49
Class for sliding window seed finder.
int Int_t
bool Bool_t
std::vector< CbmBmonDigi > fDigis
Data vector.
CbmPsdDigiData fPsd
PSD data.
Definition CbmDigiData.h:44
CbmTrdDigiData fTrd
TRD data.
Definition CbmDigiData.h:41
CbmTofDigiData fTof
TOF data.
Definition CbmDigiData.h:43
CbmStsDigiData fSts
STS data.
Definition CbmDigiData.h:37
CbmRichDigiData fRich
RICH data.
Definition CbmDigiData.h:40
CbmMuchDigiData fMuch
MUCH data.
Definition CbmDigiData.h:38
CbmBmonDigiData fBmon
Beam monitor data.
Definition CbmDigiData.h:36
Collection of digis from all detector systems within one event.
double fTime
Event trigger time [ns].
CbmDigiData fData
Event data.
uint64_t fNumber
Event identifier.
static CbmDigiManager * Instance()
Static instance.
Class characterising one event by a collection of links (indices) to data objects,...
Definition CbmEvent.h:34
std::string ToString() const
Definition CbmEvent.cxx:101
static TString GetModuleNameCaps(ECbmModuleId moduleId)
std::vector< CbmMuchDigi > fDigis
Data vector.
std::vector< CbmPsdDigi > fDigis
Data vector.
std::vector< CbmRichDigi > fDigis
Data vector.
std::vector< CbmStsDigi > fDigis
Data vector.
void SetSeedFinderQa(Bool_t bFlagIn=kTRUE)
CbmAlgoBuildRawEvents * fpAlgo
timing only for filling of std::vector<Digi> fields
void SetSlidingWindowSeedFinder(int32_t minDigis, double dWindDur, double dDeadT, double dOffset=0.0)
UInt_t GetNofDigis(ECbmModuleId _system)
void ReadDigis(ECbmModuleId detId, std::vector< TDigi > *vDigis)
std::vector< CbmMuchBeamTimeDigi > * fMuchBeamTimeDigis
std::vector< CbmStsDigi > * fStsDigis
std::vector< CbmDigiEvent > * fDigiEvents
output container of CbmEvents
std::vector< CbmPsdDigi > * fPsdDigis
TStopwatch * fCopyTimer
is created when fbGetTimings is set before init
std::vector< Double_t > * fTempDigiTimes
Bool_t fbWriteHistosToFairSink
Switch ON/OFF filling of histograms.
TString fsOutFileName
Measure CPU time using stopwatch.
std::vector< CbmTrdDigi > * fTrdDigis
TClonesArray * fEvents
Enable/disabled loop based extraction of TRD digis due to 1D/2D.
std::vector< CbmTofDigi > * fTofDigis
CbmSeedFinderSlidingWindow * fSeedFinderSlidingWindow
Switch between MUCH digi classes.
std::vector< CbmMatch > * fvDigiMatchQa
std::vector< Double_t > * fSeedTimes
void ExtractSelectedData(std::vector< CbmEvent * > vEvents)
output container of CbmEvents
Bool_t fbGetTimings
Write histos to FairRootManager instead of separate file.
void InitDigis(ECbmModuleId detId, std::vector< TDigi > **vDigi)
std::vector< CbmMuchDigi > * fMuchDigis
virtual void Exec(Option_t *)
std::vector< RawEventBuilderDetector > fSeedTimeDetList
void AddSeedTimeFillerToList(RawEventBuilderDetector seedDet)
void SetIdealSeedFinder(const int32_t fileId=-1)
std::vector< CbmRichDigi > * fRichDigis
Double_t GetDigiTime(ECbmModuleId _system, UInt_t _entry)
void FillSeedTimesFromDetList(std::vector< Double_t > *vdSeedTimes, std::vector< CbmMatch > *vDigiMatch=nullptr)
std::vector< CbmBmonDigi > * fBmonDigis
std::vector< CbmFsdDigi > * fFsdDigis
std::vector< CbmTofDigi > fDigis
Data vector.
std::vector< CbmTrdDigi > fDigis
Data vector.
ECbmModuleId detId
Settings.