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 }
61}
62
64{
65 SetSlidingWindowSeedFinder(1, 0.0, 0.0, 0.0);
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 }
79
80 fSeedFinderSlidingWindow = new CbmSeedFinderSlidingWindow(fSeedTimes, minDigis, dWindDur, dDeadT);
82}
83
85{
86 if (bFlagIn == kTRUE) {
87 if (fSeedFinderSlidingWindow == nullptr) {
88 std::cout << "SetSeedFinderQa(): Cannot activate Qa when seed finder not active. Exiting." << std::endl;
89 exit(1);
90 }
91 if (fvDigiMatchQa == nullptr) {
92 fvDigiMatchQa = new std::vector<CbmMatch>;
93 }
94 }
95 else //kFALSE
96 {
97 if (fvDigiMatchQa != nullptr) {
98 delete fvDigiMatchQa;
99 }
100 fvDigiMatchQa = nullptr;
101 }
103}
104
105template<class TDigi>
106void CbmTaskBuildRawEvents::InitDigis(ECbmModuleId detId, std::vector<TDigi>** vDigi)
107{
108 TString detName = CbmModuleList::GetModuleNameCaps(detId);
109 if (!fDigiMan->IsPresent(detId)) {
110 LOG(info) << "No " << detName << " digi input.";
111 }
112 else {
113 LOG(info) << detName << " digi input.";
114 *vDigi = new std::vector<TDigi>;
115 fpAlgo->SetDigis(*vDigi);
116 }
117}
118
120{
121 if (fbGetTimings) {
122 fTimer = new TStopwatch;
123 fTimer->Start();
124 fCopyTimer = new TStopwatch;
125 fCopyTimer->Reset();
126 }
127
129 FairRootManager* ioman = FairRootManager::Instance();
130
131 // Get a pointer to the previous already existing data level
135 }
136 fDigiMan->Init();
137
138 //Init digis
141 }
142 else {
144 }
152
154 if (fbDigiEvtOut) {
155 if (fbUseMuchBeamtimeDigi) LOG(fatal) << "DigiEvent output branch not compatible with MuchBeamtimeDigi";
156
157 fDigiEvents = new std::vector<CbmDigiEvent>();
158 ioman->RegisterAny("DigiEvent", fDigiEvents, kTRUE);
159 if (!fDigiEvents) LOG(fatal) << "Output branch was not created";
160 LOG(info) << "DigiEvent out instead of CbmEvent, you will need an instance of CbmTaskMakeRecoEvents in reco macro";
162 LOG(info) << "Exclusive TRD extraction, DigiEvents will be comparable to CbmEvents but slower";
163 }
164 else {
165 LOG(info) << "Inclusive TRD extraction, faster but DigiEvents will noy be comparable to CbmEvents (extra digis)";
166 }
167 }
168 else {
169 fEvents = new TClonesArray("CbmEvent", 100);
170 ioman->Register("CbmEvent", "Cbm_Event", fEvents, IsOutputBranchPersistent("CbmEvent"));
171 if (!fEvents) LOG(fatal) << "Output branch was not created";
172 LOG(info) << "CbmEvent oupput, you will need an instance of CbmTaskEventsCloneInToOut in reco macro to update them";
173 }
174
175 // Set timeslice meta data
176 fpAlgo->SetTimeSliceMetaDataArray(dynamic_cast<TClonesArray*>(ioman->GetObject("TimesliceMetaData")));
177
178 if (fTimer != nullptr) {
179 fTimer->Stop();
180 Double_t rtime = fTimer->RealTime();
181 Double_t ctime = fTimer->CpuTime();
182 LOG(info) << "CbmTaskBuildRawEvents::Init(): Real time " << rtime << " s, CPU time " << ctime << " s";
183 }
184
185 // Init seed finder
188 }
189
191 if (kTRUE == fpAlgo->InitAlgo())
192 return kSUCCESS;
193 else
194 return kFATAL;
195}
196
197InitStatus CbmTaskBuildRawEvents::ReInit() { return kSUCCESS; }
198
199template<class TDigi>
200void CbmTaskBuildRawEvents::ReadDigis(ECbmModuleId detId, std::vector<TDigi>* vDigis)
201{
202 //Warning: Int_t must be used for the loop counters instead of UInt_t,
203 //as the digi manager can return -1, which would be casted to +1
204 //during comparison, leading to an error.
205 if (fCopyTimer != nullptr) {
206 fCopyTimer->Start(kFALSE);
207 }
208
209 const TString detName = CbmModuleList::GetModuleNameCaps(detId);
210
211 if (fDigiMan->IsPresent(detId)) {
212 vDigis->clear();
213 const Int_t nDigis = fDigiMan->GetNofDigis(detId);
214
215 for (Int_t i = 0; i < nDigis; i++) {
216 const TDigi* Digi = fDigiMan->Get<TDigi>(i);
217 vDigis->push_back(*Digi);
218 }
219 LOG(debug) << "Read: " << vDigis->size() << " " << detName << " digis.";
220 LOG(debug) << "In DigiManager: " << nDigis << " " << detName << " digis.";
221 }
222
223 if (fCopyTimer != nullptr) {
224 fCopyTimer->Stop();
225 }
226}
227
228void CbmTaskBuildRawEvents::Exec(Option_t* /*option*/)
229{
230 if (fTimer != nullptr) {
231 fTimer->Start(kFALSE);
232 }
233 TStopwatch timer;
234 timer.Start();
235 LOG(debug2) << "CbmTaskBuildRawEvents::Exec => Starting sequence";
236
237 // Process Timeslice
238 BuildEvents();
239 if (fTimer != nullptr) {
240 fTimer->Stop();
241 }
242
243 // --- Timeslice log and statistics
244 timer.Stop();
245 std::stringstream logOut;
246 logOut << std::setw(20) << std::left << GetName() << " [";
247 logOut << std::fixed << std::setw(8) << std::setprecision(1) << std::right << timer.RealTime() * 1000. << " ms] ";
248 logOut << "TS " << fNofTs;
249 if (fbDigiEvtOut) {
250 logOut << ", events " << fDigiEvents->size();
251 fNofEvents += fDigiEvents->size();
252 }
253 else {
254 logOut << ", events " << fEvents->GetEntriesFast();
255 fNofEvents += fEvents->GetEntriesFast();
256 }
257 LOG(info) << logOut.str();
259 const size_t seedCount = fSeedFinderSlidingWindow->GetNofSeeds();
260 LOG(info) << seedCount << " trigger(s) for this TS.";
261 fTotalSeedCount += seedCount;
262 }
263 fNofTs++;
264 fTime += timer.RealTime();
265
266 LOG(debug2) << "CbmTaskBuildRawEvents::Exec => Done";
267}
268
270{
271 //Reset explicit seed times if set
272 if (fSeedTimes != nullptr) {
273 fSeedTimes->clear();
274 }
275
276 //Read digis
279 }
280 else {
282 }
290
291 //Fill seeds
292 if (fSeedFinderSlidingWindow != nullptr) {
294 }
295 else if (fSeedTimeDetList.size() > 0) {
297 }
298 //DumpSeedTimesFromDetList();
299
301 fpAlgo->ProcessTs();
302
304 FillOutput();
305}
306
307void CbmTaskBuildRawEvents::FillSeedTimesFromDetList(std::vector<Double_t>* vdSeedTimes,
308 std::vector<CbmMatch>* vDigiMatch)
309{
310 std::map<ECbmModuleId, UInt_t> DigiNumbers;
311 std::map<ECbmModuleId, UInt_t> DigiCounters;
312 vdSeedTimes->clear();
313
314 if (vDigiMatch != nullptr) vDigiMatch->clear();
315
317 DigiNumbers[system.detId] = GetNofDigis(system.detId);
318 DigiCounters[system.detId] = 0;
319 }
320
321 do {
322 ECbmModuleId nextAddedSystem;
323 Double_t earliestTime = -1;
324
326 if (DigiCounters[system.detId] < DigiNumbers[system.detId]) {
327
328 Double_t thisTime = GetDigiTime(system.detId, DigiCounters[system.detId]);
329 if (thisTime < earliestTime || earliestTime == -1) {
330 nextAddedSystem = system.detId;
331 earliestTime = thisTime;
332 }
333 }
334 }
335 if (earliestTime != -1) {
336
337 if (vDigiMatch != nullptr) {
338 const CbmMatch* digiMatch = fDigiMan->GetMatch(nextAddedSystem, DigiCounters[nextAddedSystem]);
339 vDigiMatch->push_back(*digiMatch);
340 }
341 vdSeedTimes->push_back(earliestTime);
342 DigiCounters[nextAddedSystem]++;
343 }
344 else {
345 break;
346 }
347 } while (true);
348}
349
351{
352 if (fSeedTimeDetList.size() == 0) {
355 return;
356 }
357 else {
358 std::cout << "FillSeedTimesFromSlidingWindow(): Error, seed detector list empty." << std::endl;
359 exit(1);
360 }
361 }
362 if (fSeedTimeDetList.size() == 1) {
363 const RawEventBuilderDetector seedDet = fSeedTimeDetList.at(0);
365 }
366 else { // more than one seed detector
367 if (!fTempDigiTimes) {
368 fTempDigiTimes = new std::vector<Double_t>;
369 }
372 }
373}
374
376{
377 if (fvDigiMatchQa != nullptr) {
378 if (!fDigiMan->IsMatchPresent(seedDet->detId)) {
379 std::cout << "FillSeedTimesFromSlidingWindow(): Error, QA set but no CbmMatch found for seed detector."
380 << std::endl;
381 exit(1);
382 }
383 fvDigiMatchQa->clear();
384 for (Int_t i = 0; i < fDigiMan->GetNofDigis(seedDet->detId); i++) {
385 const CbmMatch* digiMatch = fDigiMan->GetMatch(seedDet->detId, i);
386 fvDigiMatchQa->push_back(*digiMatch);
387 }
388 }
389
390 switch (seedDet->detId) {
394 break;
395 }
396 else {
398 break;
399 }
407 default: break;
408 }
409}
410
411Double_t CbmTaskBuildRawEvents::GetDigiTime(ECbmModuleId _system, UInt_t _entry)
412{
413 switch (_system) {
416 return (fMuchBeamTimeDigis->at(_entry)).GetTime();
417 }
418 else {
419 return (fMuchDigis->at(_entry)).GetTime();
420 }
421 case ECbmModuleId::kSts: return (fStsDigis->at(_entry)).GetTime();
422 case ECbmModuleId::kTrd: return (fTrdDigis->at(_entry)).GetTime();
423 case ECbmModuleId::kTof: return (fTofDigis->at(_entry)).GetTime();
424 case ECbmModuleId::kRich: return (fRichDigis->at(_entry)).GetTime();
425 case ECbmModuleId::kPsd: return (fPsdDigis->at(_entry)).GetTime();
426 case ECbmModuleId::kFsd: return (fFsdDigis->at(_entry)).GetTime();
427 case ECbmModuleId::kBmon: return (fBmonDigis->at(_entry)).GetTime();
428 default: break;
429 }
430 return -1;
431}
432
434{
435 switch (_system) {
438 return fMuchBeamTimeDigis->size();
439 }
440 else {
441 return fMuchDigis->size();
442 }
443 case ECbmModuleId::kSts: return fStsDigis->size();
444 case ECbmModuleId::kTrd: return fTrdDigis->size();
445 case ECbmModuleId::kTof: return fTofDigis->size();
446 case ECbmModuleId::kRich: return fRichDigis->size();
447 case ECbmModuleId::kPsd: return fPsdDigis->size();
448 case ECbmModuleId::kFsd: return fFsdDigis->size();
449 case ECbmModuleId::kBmon: return fBmonDigis->size();
450 default: break;
451 }
452 return 0;
453}
454
456{
457 if (fTimer == nullptr) {
458 LOG(fatal) << "Trying to print timings but timer not set";
459 }
460 else {
461 Double_t rtime = fTimer->RealTime();
462 Double_t ctime = fTimer->CpuTime();
463 LOG(info) << "CbmTaskBuildRawEvents: Real time " << rtime << " s, CPU time " << ctime << " s";
464 }
465 if (fCopyTimer == nullptr) {
466 LOG(fatal) << "Trying to print timings but timer not set";
467 }
468 else {
469 Double_t rtime = fCopyTimer->RealTime();
470 Double_t ctime = fCopyTimer->CpuTime();
471 LOG(info) << "CbmTaskBuildRawEvents (digi copy only): Real time " << rtime << " s, CPU time " << ctime << " s";
472 }
473}
474
476{
477 if ((fvDigiMatchQa != nullptr) && (fSeedFinderSlidingWindow != nullptr)) {
479 }
480
482 fpAlgo->Finish();
483 if (fbFillHistos) {
484 SaveHistos();
485 }
486 if (fbGetTimings) {
487 PrintTimings();
488 }
489
490 std::cout << std::endl;
491 LOG(info) << "=====================================";
492 LOG(info) << GetName() << ": Run summary";
493 LOG(info) << "Time slices : " << fNofTs;
494 LOG(info) << "Events : " << fNofEvents;
496 LOG(info) << "Triggers : " << fTotalSeedCount;
497 }
498 LOG(info) << "Time / TS : " << std::fixed << std::setprecision(2) << 1000. * fTime / Double_t(fNofTs)
499 << " ms";
500 LOG(info) << "=====================================";
501}
502
504{
506 std::vector<CbmEvent*> vEvents = fpAlgo->GetEventVector();
507
508 if (fbDigiEvtOut) {
510 fDigiEvents->clear();
511
513 ExtractSelectedData(vEvents);
514 }
515 else {
517 fEvents->Delete();
518
520 for (CbmEvent* event : vEvents) {
521 LOG(debug) << "Vector: " << event->ToString();
522 new ((*fEvents)[fEvents->GetEntriesFast()]) CbmEvent(std::move(*event));
523 LOG(debug) << "TClonesArray: " << static_cast<CbmEvent*>(fEvents->At(fEvents->GetEntriesFast() - 1))->ToString();
524 }
525 }
528}
529
531{
533 if (!FairRootManager::Instance() || !FairRootManager::Instance()->GetSink()) {
534 LOG(error) << "No sink found";
535 return;
536 }
537 FairSink* sink = FairRootManager::Instance()->GetSink();
538 sink->WriteObject(dynamic_cast<TObject*>(fpAlgo->GetOutFolder()), nullptr);
539 }
540 else {
541
543 std::vector<std::pair<TNamed*, std::string>> vHistos = fpAlgo->GetHistoVector();
544
546 TDirectory* oldDir = NULL;
547 TFile* histoFile = NULL;
549 oldDir = gDirectory;
551 histoFile = new TFile(fsOutFileName, "RECREATE");
552 histoFile->cd();
553
555 for (UInt_t uHisto = 0; uHisto < vHistos.size(); ++uHisto) {
557 const TString sFolder = vHistos[uHisto].second.data();
558 if (nullptr == gDirectory->Get(sFolder)) gDirectory->mkdir(sFolder);
559 gDirectory->cd(sFolder);
560
562 vHistos[uHisto].first->Write();
563 histoFile->cd();
564 }
565
567 std::vector<std::pair<TCanvas*, std::string>> vCanvases = fpAlgo->GetCanvasVector();
568
570 for (UInt_t uCanvas = 0; uCanvas < vCanvases.size(); ++uCanvas) {
572 TString sFolder = vCanvases[uCanvas].second.data();
573 if (nullptr == gDirectory->Get(sFolder)) gDirectory->mkdir(sFolder);
574 gDirectory->cd(sFolder);
575
577 vCanvases[uCanvas].first->Write();
578
579 histoFile->cd();
580 }
581
583 oldDir->cd();
584 histoFile->Close();
585 }
586}
587
589{
590 std::ofstream timesUnsorted("digiTimesUnsorted.dat", std::ofstream::out);
591 timesUnsorted << std::setprecision(16);
592
594 for (UInt_t i = 0; i < GetNofDigis(system.detId); i++) {
595 timesUnsorted << GetDigiTime(system.detId, i) << std::endl;
596 }
597 }
598 timesUnsorted.close();
599 LOG(info) << "Completed write of unsorted digi list.";
600
601 std::ofstream timesSorted("digiTimesSorted.dat", std::ofstream::out);
602 timesSorted << std::setprecision(16);
603
604 for (UInt_t i = 0; i < fSeedTimes->size(); i++) {
605 timesSorted << fSeedTimes->at(i) << std::endl;
606 }
607 timesSorted.close();
608 LOG(info) << "Completed DumpSeedTimesFromDetList(). Closing.";
609 exit(0); //terminate as this method should only be used for diagnostics
610}
611
612void CbmTaskBuildRawEvents::ExtractSelectedData(std::vector<CbmEvent*> vEvents)
613{
615 LOG(debug) << "In Vector size: " << vEvents.size();
616
617 fDigiEvents->reserve(vEvents.size());
618 for (CbmEvent* event : vEvents) {
619 CbmDigiEvent selEvent;
620 selEvent.fTime = event->GetStartTime();
621 selEvent.fNumber = event->GetNumber();
622
625 event->SortIndices();
626
630 uint32_t uNbDigis =
631 (0 < event->GetNofData(ECbmDataType::kBmonDigi) ? event->GetNofData(ECbmDataType::kBmonDigi) : 0);
632 if (0 < uNbDigis) {
633 auto startIt = fBmonDigis->begin() + event->GetIndex(ECbmDataType::kBmonDigi, 0);
634 auto stopIt = fBmonDigis->begin() + event->GetIndex(ECbmDataType::kBmonDigi, uNbDigis - 1);
635 ++stopIt;
636 selEvent.fData.fBmon.fDigis.assign(startIt, stopIt);
637 }
638
640 uNbDigis = (0 < event->GetNofData(ECbmDataType::kStsDigi) ? event->GetNofData(ECbmDataType::kStsDigi) : 0);
641 if (0 < uNbDigis) {
642 auto startIt = fStsDigis->begin() + event->GetIndex(ECbmDataType::kStsDigi, 0);
643 auto stopIt = fStsDigis->begin() + event->GetIndex(ECbmDataType::kStsDigi, uNbDigis - 1);
644 ++stopIt;
645 selEvent.fData.fSts.fDigis.assign(startIt, stopIt);
646 }
647
649 uNbDigis = (0 < event->GetNofData(ECbmDataType::kMuchDigi) ? event->GetNofData(ECbmDataType::kMuchDigi) : 0);
650 if (0 < uNbDigis) {
651 auto startIt = fMuchDigis->begin() + event->GetIndex(ECbmDataType::kMuchDigi, 0);
652 auto stopIt = fMuchDigis->begin() + event->GetIndex(ECbmDataType::kMuchDigi, uNbDigis - 1);
653 ++stopIt;
654 selEvent.fData.fMuch.fDigis.assign(startIt, stopIt);
655 }
656
658 uNbDigis = (0 < event->GetNofData(ECbmDataType::kTrdDigi) ? event->GetNofData(ECbmDataType::kTrdDigi) : 0);
659 if (0 < uNbDigis) {
661 for (uint32_t uDigiInEvt = 0; uDigiInEvt < uNbDigis; ++uDigiInEvt) {
666 selEvent.fData.fTrd.fDigis.push_back(fTrdDigis->at(event->GetIndex(ECbmDataType::kTrdDigi, uDigiInEvt)));
667 }
668 }
669 else {
674 auto startIt = fTrdDigis->begin() + event->GetIndex(ECbmDataType::kTrdDigi, 0);
675 auto stopIt = fTrdDigis->begin() + event->GetIndex(ECbmDataType::kTrdDigi, uNbDigis - 1);
676 ++stopIt;
677 selEvent.fData.fTrd.fDigis.assign(startIt, stopIt);
678 }
679 }
680
682 uNbDigis = (0 < event->GetNofData(ECbmDataType::kTofDigi) ? event->GetNofData(ECbmDataType::kTofDigi) : 0);
683 if (0 < uNbDigis) {
684 auto startIt = fTofDigis->begin() + event->GetIndex(ECbmDataType::kTofDigi, 0);
685 auto stopIt = fTofDigis->begin() + event->GetIndex(ECbmDataType::kTofDigi, uNbDigis - 1);
686 ++stopIt;
687 selEvent.fData.fTof.fDigis.assign(startIt, stopIt);
688 }
689
691 uNbDigis = (0 < event->GetNofData(ECbmDataType::kRichDigi) ? event->GetNofData(ECbmDataType::kRichDigi) : 0);
692 if (0 < uNbDigis) {
693 auto startIt = fRichDigis->begin() + event->GetIndex(ECbmDataType::kRichDigi, 0);
694 auto stopIt = fRichDigis->begin() + event->GetIndex(ECbmDataType::kRichDigi, uNbDigis - 1);
695 ++stopIt;
696 selEvent.fData.fRich.fDigis.assign(startIt, stopIt);
697 }
698
700 uNbDigis = (0 < event->GetNofData(ECbmDataType::kPsdDigi) ? event->GetNofData(ECbmDataType::kPsdDigi) : 0);
701 if (0 < uNbDigis) {
702 auto startIt = fPsdDigis->begin() + event->GetIndex(ECbmDataType::kPsdDigi, 0);
703 auto stopIt = fPsdDigis->begin() + event->GetIndex(ECbmDataType::kPsdDigi, uNbDigis - 1);
704 ++stopIt;
705 selEvent.fData.fPsd.fDigis.assign(startIt, stopIt);
706 }
707
708 fDigiEvents->push_back(std::move(selEvent));
709 }
710}
711
ClassImp(CbmConverterManager)
ECbmModuleId
Definition CbmDefs.h:39
@ kTrd
Transition Radiation Detector.
@ kTof
Time-of-flight Detector.
@ kPsd
Projectile spectator detector.
@ kSts
Silicon Tracking System.
@ kMuch
Muon detection system.
@ kFsd
Forward spectator detector.
@ kRich
Ring-Imaging Cherenkov Detector.
Class for sliding window seed finder.
void SetTimeSliceMetaDataArray(TClonesArray *TimeSliceMetaDataArray)
std::vector< std::pair< TNamed *, std::string > > GetHistoVector()
std::vector< std::pair< TCanvas *, std::string > > GetCanvasVector()
void SetSeedTimes(std::vector< Double_t > *SeedTimes)
TDirectoryFile * GetOutFolder()
std::vector< CbmEvent * > & GetEventVector()
Data output access.
void SetDigis(std::vector< CbmBmonDigi > *BmonDigis)
Set digi containers.
std::vector< CbmBmonDigi > fDigis
Data vector.
CbmPsdDigiData fPsd
PSD data.
Definition CbmDigiData.h:42
CbmTrdDigiData fTrd
TRD data.
Definition CbmDigiData.h:39
CbmTofDigiData fTof
TOF data.
Definition CbmDigiData.h:41
CbmStsDigiData fSts
STS data.
Definition CbmDigiData.h:36
CbmRichDigiData fRich
RICH data.
Definition CbmDigiData.h:38
CbmMuchDigiData fMuch
MUCH data.
Definition CbmDigiData.h:37
CbmBmonDigiData fBmon
Beam monitor data.
Definition CbmDigiData.h:35
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 Int_t GetNofDigis(ECbmModuleId systemId)
static Bool_t IsPresent(ECbmModuleId systemId)
Presence of a digi branch.
static Bool_t IsMatchPresent(ECbmModuleId systemId)
Presence of a digi match branch.
void UseMuchBeamTimeDigi(Bool_t)
Use CbmMuchBeamTimeDigi instead of CbmMuchDigi for MUCH.
InitStatus Init()
Initialisation.
const Digi * Get(Int_t index) const
Get a digi object.
static CbmDigiManager * Instance()
Static instance.
const CbmMatch * GetMatch(ECbmModuleId systemId, UInt_t index) const
Get a match object.
Class characterising one event by a collection of links (indices) to data objects,...
Definition CbmEvent.h:34
std::string ToString() const
Definition CbmEvent.cxx:96
static TString GetModuleNameCaps(ECbmModuleId moduleId)
std::vector< CbmMuchDigi > fDigis
Data vector.
std::vector< CbmPsdDigi > fDigis
Data vector.
std::vector< CbmRichDigi > fDigis
Data vector.
void SetIdealMode(const int32_t fileId=-1)
Switches to `‘ideal mode’' in which event times from MC data are used as triggers (no algorithm is ru...
bool IsIdealMode()
Is `‘ideal mode’' switched on?
void OutputQa()
Output QA Information.
void SetOffset(double offset)
Sets a global constant offset which is applied to each trigger time @params Value of offset.
size_t GetNofSeeds()
Returns number of seed times currently stored in buffer.
void SetQa(bool doQA=true)
Enable or disable the generation of QA information.
void FillSeedTimes(const std::vector< inType > *vIn, const std::vector< CbmMatch > *vDigiMatch=nullptr)
Function which builds event seeds @params Vector of input data (either digis or digi times)....
void Init()
Initializes QA object if set.
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.