CbmRoot
Loading...
Searching...
No Matches
CbmAlgoBuildRawEvents.h
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
5#ifndef CBMALGOBUILDRAWEVENTS_H
6#define CBMALGOBUILDRAWEVENTS_H
7
9#include "CbmDefs.h"
10
12
14#include "TFolder.h"
15
17#include "CbmMvdDigi.h"
18
19#include <boost/any.hpp>
20
21#include <array>
22#include <map>
23#include <set>
24#include <tuple>
25#include <vector>
26
28class CbmEvent;
29class CbmMuchDigi;
31class CbmPsdDigi;
32class CbmFsdDigi;
33class CbmRichDigi;
34class CbmStsDigi;
35class CbmTofDigi;
36class CbmTrdDigi;
37class CbmBmonDigi;
38class TClonesArray;
39class TH1;
40class TH2;
41class TProfile;
42class TNamed;
43class TStopwatch;
44class TCanvas;
45class TDirectoryFile;
46
54
56 public:
58
59 RawEventBuilderDetector(ECbmModuleId detIdIn, ECbmDataType dataTypeIn, std::string sNameIn)
60 : detId{detIdIn}
61 , dataType{dataTypeIn}
62 , sName{sNameIn}
63 {
64 }
65
66 RawEventBuilderDetector(ECbmModuleId detIdIn, ECbmDataType dataTypeIn, std::string sNameIn, UInt_t uTriggerMinDigisIn,
67 Int_t iTriggerMaxDigisIn, Double_t fdTimeWinBegIn, Double_t fdTimeWinEndIn,
68 UInt_t uTriggerMinLayersIn = 0, Double_t fdHistMaxDigiNbIn = 1000)
69 : RawEventBuilderDetector(detIdIn, dataTypeIn, sNameIn)
70 {
71 fuTriggerMinDigis = uTriggerMinDigisIn;
72 fiTriggerMaxDigis = iTriggerMaxDigisIn;
73 fuTriggerMinLayers = uTriggerMinLayersIn;
74 fdTimeWinBeg = fdTimeWinBegIn;
75 fdTimeWinEnd = fdTimeWinEndIn;
76 fdHistMaxDigiNb = fdHistMaxDigiNbIn;
77 }
78
79 bool operator==(const RawEventBuilderDetector& other) const { return (other.detId == this->detId); }
80 bool operator!=(const RawEventBuilderDetector& other) const { return (other.detId != this->detId); }
81
82 Double_t GetTimeWinRange() { return fdTimeWinEnd - fdTimeWinBeg; }
83
87 std::string sName = "Invalid";
95 Double_t fdTimeWinBeg = -100;
96 Double_t fdTimeWinEnd = 100;
98 Double_t fdHistMaxDigiNb = 1000;
100 UInt_t fuStartIndex = 0;
101 UInt_t fuEndIndex = 0;
102};
103
126
128 public:
131
134
137
140
142 void ProcessTs();
143
145 void Finish();
146
147 void SetFillHistos(Bool_t var) { fbFillHistos = var; }
148 void ResetHistograms(Bool_t bResetTime = kTRUE);
149
151 void SetTimings(Bool_t var) { fbGetTimings = var; }
152 void PrintTimings();
153
154 void SetReferenceDetector(ECbmModuleId refDet, ECbmDataType dataTypeIn, std::string sNameIn,
155 UInt_t uTriggerMinDigisIn = 0, Int_t iTriggerMaxDigisIn = -1,
156 Double_t fdTimeWinBegIn = -100, Double_t fdTimeWinEndIn = 100);
157 void AddDetector(ECbmModuleId selDet, ECbmDataType dataTypeIn, std::string sNameIn, UInt_t uTriggerMinDigisIn = 0,
158 Int_t iTriggerMaxDigisIn = -1, Double_t fdTimeWinBegIn = -100, Double_t fdTimeWinEndIn = 100);
159
160 void SetReferenceDetector(RawEventBuilderDetector refDetIn, std::vector<bool> select = {});
161 void AddDetector(RawEventBuilderDetector selDet);
162 void RemoveDetector(RawEventBuilderDetector selDet);
163
164 void SetTriggerMinNumber(ECbmModuleId selDet, UInt_t uVal);
165 void SetTriggerMaxNumber(ECbmModuleId selDet, Int_t iVal);
166 void SetTriggerMinLayersNumber(ECbmModuleId selDet, UInt_t uVal);
167 void SetTriggerWindow(ECbmModuleId selDet, Double_t dWinBeg, Double_t dWinEnd);
168 void SetHistogramMaxDigiNb(ECbmModuleId selDet, Double_t dDigiNbMax);
169 void SetTsParameters(Double_t dTsStartTime, Double_t dTsLength, Double_t dTsOverLength)
170 {
172 fdTsLength = dTsLength;
173 fdTsOverLength = dTsOverLength;
174 fbUseTsMetaData = kFALSE;
175 }
176
177 void SetSeedTimeWindow(Double_t timeWinBeg, Double_t timeWinEnd)
178 {
179 fdSeedTimeWinBeg = timeWinBeg;
180 fdSeedTimeWinEnd = timeWinEnd;
183 }
184
187 void SetIgnoreTsOverlap(Bool_t bFlagIn = kTRUE) { fbIgnoreTsOverlap = bFlagIn; }
188
189 void ChangeMuchBeamtimeDigiFlag(Bool_t bFlagIn = kFALSE) { fbUseMuchBeamtimeDigi = bFlagIn; }
190
192 void AddHistoToVector(TNamed* pointer, std::string sFolder = "")
193 {
194 fvpAllHistoPointers.push_back(std::pair<TNamed*, std::string>(pointer, sFolder));
195 }
196 std::vector<std::pair<TNamed*, std::string>> GetHistoVector() { return fvpAllHistoPointers; }
197 void AddCanvasToVector(TCanvas* pointer, std::string sFolder = "")
198 {
199 fvpAllCanvasPointers.push_back(std::pair<TCanvas*, std::string>(pointer, sFolder));
200 }
201 std::vector<std::pair<TCanvas*, std::string>> GetCanvasVector() { return fvpAllCanvasPointers; }
202
204 void SetDigis(std::vector<CbmBmonDigi>* BmonDigis) { fBmonDigis = BmonDigis; }
205 void SetDigis(std::vector<CbmMvdDigi>* MvdDigis) { fMvdDigis = MvdDigis; }
206 void SetDigis(std::vector<CbmStsDigi>* StsDigis) { fStsDigis = StsDigis; }
207 void SetDigis(std::vector<CbmMuchDigi>* MuchDigis)
208 {
209 fMuchDigis = MuchDigis;
210 fbUseMuchBeamtimeDigi = kFALSE;
211 }
212 void SetDigis(std::vector<CbmTrdDigi>* TrdDigis) { fTrdDigis = TrdDigis; }
213 void SetDigis(std::vector<CbmTofDigi>* TofDigis) { fTofDigis = TofDigis; }
214 void SetDigis(std::vector<CbmRichDigi>* RichDigis) { fRichDigis = RichDigis; }
215 void SetDigis(std::vector<CbmPsdDigi>* PsdDigis) { fPsdDigis = PsdDigis; }
216 void SetDigis(std::vector<CbmFsdDigi>* FsdDigis) { fFsdDigis = FsdDigis; }
217 void SetDigis(std::vector<CbmMuchBeamTimeDigi>* MuchBeamTimeDigis)
218 {
219 fMuchBeamTimeDigis = MuchBeamTimeDigis;
220 fbUseMuchBeamtimeDigi = kTRUE;
221 }
222
223 void SetSeedTimes(std::vector<Double_t>* SeedTimes) { fSeedTimes = SeedTimes; }
224
225 // TS metadata
226 void SetTimeSliceMetaDataArray(TClonesArray* TimeSliceMetaDataArray)
227 {
228 fTimeSliceMetaDataArray = TimeSliceMetaDataArray;
229 }
230
231 // Output folder for histograms
232 TDirectoryFile* GetOutFolder() { return outFolder; }
233
235 std::vector<CbmEvent*>& GetEventVector() { return fEventVector; }
236 void ClearEventVector();
237
238 private:
241 void InitTs();
242 void InitSeedWindow();
243 void BuildEvents();
244
245 void CreateHistograms();
246 void FillHistos();
247
248 template<class DigiSeed>
249 void LoopOnSeeds();
250
251 void CheckSeed(Double_t dSeedTime, UInt_t uSeedDigiIdx);
252 void CheckTriggerCondition(Double_t dSeedTime);
253
254 template<class DigiCheck>
255 void SearchMatches(Double_t dSeedTime, RawEventBuilderDetector& detMatch);
256 void SearchMatches(Double_t dSeedTime, RawEventBuilderDetector& detMatch);
257 void AddDigiToEvent(const RawEventBuilderDetector& det, Int_t uIdx);
260
263
264 void SwitchBmonStation(int id, bool on = true);
265 bool SetBmonEventTime(CbmEvent* event);
266
267 void CheckBmonInUse();
268
273 bool filterBmon(int32_t add);
275
276 TDirectoryFile* outFolder; // oputput folder to store histograms
277
279 static constexpr Double_t kdDefaultTimeWinBeg = -100.0;
280 static constexpr Double_t kdDefaultTimeWinEnd = 100.0;
281
289 std::vector<bool> fUseBmonMap = {};
290
293
294 TStopwatch* fTimer = nullptr;
295
300 bool fbBmonInUse = false;
301
306 Double_t fdSeedWindowBeg = 0;
307 Double_t fdSeedWindowEnd = 0;
308
309 Double_t fdTsStartTime = -1;
310 Double_t fdTsLength = -1;
311 Double_t fdTsOverLength = -1;
312
314 TClonesArray* fTimeSliceMetaDataArray = nullptr;
315
316 const std::vector<CbmBmonDigi>* fBmonDigis = nullptr;
317 const std::vector<CbmMuchDigi>* fMuchDigis = nullptr;
318 const std::vector<CbmMuchBeamTimeDigi>* fMuchBeamTimeDigis = nullptr;
319 const std::vector<CbmMvdDigi>* fMvdDigis = nullptr;
320 const std::vector<CbmStsDigi>* fStsDigis = nullptr;
321 const std::vector<CbmTrdDigi>* fTrdDigis = nullptr;
322 const std::vector<CbmTofDigi>* fTofDigis = nullptr;
323 const std::vector<CbmRichDigi>* fRichDigis = nullptr;
324 const std::vector<CbmPsdDigi>* fPsdDigis = nullptr;
325 const std::vector<CbmFsdDigi>* fFsdDigis = nullptr;
326
327 // If explicit seed times are supplied.
328 const std::vector<Double_t>* fSeedTimes = nullptr;
329 Double_t fdSeedTimeWinBeg = -100.0;
330 Double_t fdSeedTimeWinEnd = 100.0;
331
333 UInt_t GetNofDigis(ECbmModuleId detId);
334 template<class Digi>
335 const Digi* GetDigi(UInt_t uDigi);
336 uint64_t GetSizeFromDigisNb(ECbmModuleId detId, uint64_t ulNbDigis);
337
338 Double_t GetSeedTimeWinRange();
339
342 std::vector<CbmEvent*> fEventVector = {};
343
355 std::vector<std::pair<TNamed*, std::string>>
357 std::vector<std::pair<TCanvas*, std::string>>
359
360 TH1* fhEventTime = nullptr;
361 TH1* fhEventDt = nullptr;
362 TH1* fhEventSize = nullptr;
363 TH2* fhNbDigiPerEvtTime = nullptr;
364 TH1* fhCpuTimePerTs = nullptr;
365 TH1* fhRealTimePerTs = nullptr;
366
367 TH1* fhCpuTimePerTsHist = nullptr;
368 TH1* fhRealTimePerTsHist = nullptr;
369
370 std::vector<TH2*> fvhNbDigiPerEvtTimeDet =
371 {};
372 std::vector<TH1*> fvhNbDigiPerEvtDet = {};
373 std::vector<TH1*> fvhTDiff = {}; // digi time difference to seed
374
375 std::vector<TH1*> fvhSelRatioPerTsNb = {};
376 std::vector<TH1*> fvhInpRatioPerTsSz = {};
377 std::vector<TH1*> fvhOutRatioPerTsSz = {};
378 TH1* fhSizeReductionPerTs = nullptr;
379
380 TH1* fhOverEventShare = nullptr;
381 TProfile* fhOverEventShareTs = nullptr;
382 TH2* fhOverEventSizeTs = nullptr;
383
385 UInt_t fuCurEv = 0;
386 UInt_t fuNrTs = 0;
387 Double_t fdPrevEvtTime = 0.;
388 Double_t fdPrevEvtEndTime = 0.;
389
391};
392
393#endif // CBMALGOBUILDRAWEVENTS_H
static const RawEventBuilderDetector kRawEventBuilderDetSts
static const RawEventBuilderDetector kRawEventBuilderDetTrd2D
static const RawEventBuilderDetector kRawEventBuilderDetMvd
Pre-defined detector types.
static const RawEventBuilderDetector kRawEventBuilderDetBmon
static const RawEventBuilderDetector kRawEventBuilderDetFsd
static const RawEventBuilderDetector kRawEventBuilderDetTrd
static const RawEventBuilderDetector kRawEventBuilderDetPsd
static const RawEventBuilderDetector kRawEventBuilderDetRich
static const RawEventBuilderDetector kRawEventBuilderDetTof
static const RawEventBuilderDetector kRawEventBuilderDetMuch
static const RawEventBuilderDetector kRawEventBuilderDetUndef
ECbmDataType
Enumerator for CBM data types.
Definition CbmDefs.h:122
ECbmModuleId
Enumerator for module Identifiers.
Definition CbmDefs.h:45
@ kMvd
Micro-Vertex Detector.
Definition CbmDefs.h:47
@ kTrd
Transition Radiation Detector.
Definition CbmDefs.h:51
@ kTof
Time-of-flight Detector.
Definition CbmDefs.h:52
@ kNotExist
If not found.
Definition CbmDefs.h:68
@ kPsd
Projectile spectator detector.
Definition CbmDefs.h:54
@ kSts
Silicon Tracking System.
Definition CbmDefs.h:48
@ kTrd2d
TRD-FASP Detector (FIXME)
Definition CbmDefs.h:58
@ 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
static double dTsStartTime
@ Undefined
Definition LmvmDef.h:29
int Int_t
bool Bool_t
CbmAlgoBuildRawEvents()=default
Bool_t CheckTriggerConditions(CbmEvent *event, const RawEventBuilderDetector &det)
void SetEventOverlapMode(EOverlapModeRaw mode)
Control flags.
static constexpr Double_t kdDefaultTimeWinBeg
Constants.
void SetTimeSliceMetaDataArray(TClonesArray *TimeSliceMetaDataArray)
void SetDigis(std::vector< CbmStsDigi > *StsDigis)
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
std::vector< std::pair< TNamed *, std::string > > GetHistoVector()
void AddDigiToEvent(const RawEventBuilderDetector &det, Int_t uIdx)
TH1 * fhCpuTimePerTsHist
Processing time per TS.
void SetDigis(std::vector< CbmRichDigi > *RichDigis)
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
void SetIgnoreTsOverlap(Bool_t bFlagIn=kTRUE)
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< std::pair< TNamed *, std::string > > fvpAllHistoPointers
vector with all created events
void SetDigis(std::vector< CbmTofDigi > *TofDigis)
std::vector< TH2 * > fvhNbDigiPerEvtTimeDet
Plotting time per TS.
static constexpr Double_t kdDefaultTimeWinEnd
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 SetDigis(std::vector< CbmTrdDigi > *TrdDigis)
std::vector< std::pair< TCanvas *, std::string > > GetCanvasVector()
void SwitchBmonStation(int id, bool on=true)
void ChangeMuchBeamtimeDigiFlag(Bool_t bFlagIn=kFALSE)
Bool_t fbUseMuchBeamtimeDigi
Switch ON/OFF filling of histograms.
const std::vector< CbmTofDigi > * fTofDigis
CbmEvent * fCurrentEvent
Data ouptut.
void SetTsParameters(Double_t dTsStartTime, Double_t dTsLength, Double_t dTsOverLength)
void SetDigis(std::vector< CbmMuchBeamTimeDigi > *MuchBeamTimeDigis)
ClassDefNV(CbmAlgoBuildRawEvents, 2)
Save previous event last digi time information.
void SetSeedTimeWindow(Double_t timeWinBeg, Double_t timeWinEnd)
CbmAlgoBuildRawEvents operator=(const CbmAlgoBuildRawEvents &)=delete
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)
void SetDigis(std::vector< CbmFsdDigi > *FsdDigis)
void SetSeedTimes(std::vector< Double_t > *SeedTimes)
UInt_t fuNrTs
Event Counter.
std::vector< std::pair< TCanvas *, std::string > > fvpAllCanvasPointers
Vector of pointers to histograms + optional folder name.
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
TDirectoryFile * GetOutFolder()
Bool_t fbFillHistos
Ignore data in Overlap part of the TS.
std::vector< CbmEvent * > & GetEventVector()
Data output access.
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
const std::vector< CbmMvdDigi > * fMvdDigis
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 SetDigis(std::vector< CbmPsdDigi > *PsdDigis)
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="")
void SetDigis(std::vector< CbmMvdDigi > *MvdDigis)
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 SetDigis(std::vector< CbmBmonDigi > *BmonDigis)
Set digi containers.
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)
CbmAlgoBuildRawEvents(const CbmAlgoBuildRawEvents &)=delete
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)
void SetDigis(std::vector< CbmMuchDigi > *MuchDigis)
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:31
Class characterising one event by a collection of links (indices) to data objects,...
Definition CbmEvent.h:34
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
Data class for expanded digital TOF information.
Definition CbmTofDigi.h:47
RawEventBuilderDetector(ECbmModuleId detIdIn, ECbmDataType dataTypeIn, std::string sNameIn, UInt_t uTriggerMinDigisIn, Int_t iTriggerMaxDigisIn, Double_t fdTimeWinBegIn, Double_t fdTimeWinEndIn, UInt_t uTriggerMinLayersIn=0, Double_t fdHistMaxDigiNbIn=1000)
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...
RawEventBuilderDetector()=default
RawEventBuilderDetector(ECbmModuleId detIdIn, ECbmDataType dataTypeIn, std::string sNameIn)
ECbmModuleId detId
Settings.
Double_t fdTimeWinBeg
Selection Window.
bool operator!=(const RawEventBuilderDetector &other) const
UInt_t fuTriggerMinDigis
Minimum number of digis per detector needed to generate an event, 0 means do not use for event select...
bool operator==(const RawEventBuilderDetector &other) const
Double_t fdHistMaxDigiNb
Histo configuration.
UInt_t fuStartIndex
Book-keeping variables.