CbmRoot
Loading...
Searching...
No Matches
CbmTaskBuildEvents.cxx
Go to the documentation of this file.
1/* Copyright (C) 2021 GSI Helmholtzzentrum fuer Schwerionenforschung, Darmstadt
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Volker Friese [committer] */
4
6
7#include "CbmDigiBranchBase.h"
8#include "CbmDigiManager.h"
9#include "CbmDigiTimeslice.h"
10
11#include <FairRootManager.h>
12#include <Logger.h>
13
14#include <TStopwatch.h>
15
16#include <cassert>
17#include <iomanip>
18
19using namespace std;
20
21// ----- Constructor -----------------------------------------------------
22CbmTaskBuildEvents::CbmTaskBuildEvents() : FairTask("BuildEvents") {}
23// ---------------------------------------------------------------------------
24
25
26// ----- Destructor ------------------------------------------------------
31// ---------------------------------------------------------------------------
32
33
34// ------ Construct a DigiTimeslice from the data in CbmDigiManager ----------
36{
38 for (const auto& entry : fConfig->fWindows) {
39 auto system = entry.first;
40 CbmDigiBranchBase* digiBranch = fDigiMan->GetBranch(system);
41 if (digiBranch == nullptr) continue;
42 switch (system) {
43 case ECbmModuleId::kSts: {
44 const vector<CbmStsDigi>* digiVec =
45 boost::any_cast<const vector<CbmStsDigi>*>(digiBranch->GetBranchContainer());
46 assert(digiVec);
47 std::copy(digiVec->begin(), digiVec->end(), std::back_inserter(ts.fData.fSts.fDigis));
48 break;
49 }
51 const vector<CbmRichDigi>* digiVec =
52 boost::any_cast<const vector<CbmRichDigi>*>(digiBranch->GetBranchContainer());
53 assert(digiVec);
54 std::copy(digiVec->begin(), digiVec->end(), std::back_inserter(ts.fData.fRich.fDigis));
55 break;
56 }
58 const vector<CbmMuchDigi>* digiVec =
59 boost::any_cast<const vector<CbmMuchDigi>*>(digiBranch->GetBranchContainer());
60 assert(digiVec);
61 std::copy(digiVec->begin(), digiVec->end(), std::back_inserter(ts.fData.fMuch.fDigis));
62 break;
63 }
64 case ECbmModuleId::kTrd: {
65 const vector<CbmTrdDigi>* digiVec =
66 boost::any_cast<const vector<CbmTrdDigi>*>(digiBranch->GetBranchContainer());
67 assert(digiVec);
68 std::copy(digiVec->begin(), digiVec->end(), std::back_inserter(ts.fData.fTrd.fDigis));
69 break;
70 }
72 const vector<CbmTrdDigi>* digiVec =
73 boost::any_cast<const vector<CbmTrdDigi>*>(digiBranch->GetBranchContainer());
74 assert(digiVec);
75 std::copy(digiVec->begin(), digiVec->end(), std::back_inserter(ts.fData.fTrd2d.fDigis));
76 break;
77 }
78 case ECbmModuleId::kTof: {
79 const vector<CbmTofDigi>* digiVec =
80 boost::any_cast<const vector<CbmTofDigi>*>(digiBranch->GetBranchContainer());
81 assert(digiVec);
82 std::copy(digiVec->begin(), digiVec->end(), std::back_inserter(ts.fData.fTof.fDigis));
83 break;
84 }
85 case ECbmModuleId::kFsd: {
86 const vector<CbmFsdDigi>* digiVec =
87 boost::any_cast<const vector<CbmFsdDigi>*>(digiBranch->GetBranchContainer());
88 assert(digiVec);
89 std::copy(digiVec->begin(), digiVec->end(), std::back_inserter(ts.fData.fFsd.fDigis));
90 break;
91 }
92 case ECbmModuleId::kPsd: {
93 const vector<CbmPsdDigi>* digiVec =
94 boost::any_cast<const vector<CbmPsdDigi>*>(digiBranch->GetBranchContainer());
95 assert(digiVec);
96 ts.fData.fPsd.fDigis = *digiVec;
97 break;
98 }
99 case ECbmModuleId::kBmon: { //Bmon has Tof digis
100 const vector<CbmBmonDigi>* digiVec =
101 boost::any_cast<const vector<CbmBmonDigi>*>(digiBranch->GetBranchContainer());
102 assert(digiVec);
103 std::copy(digiVec->begin(), digiVec->end(), std::back_inserter(ts.fData.fTof.fDigis));
104 break;
105 }
106 default: LOG(fatal) << GetName() << ": Unknown detector type!";
107 }
108 }
109 return ts;
110}
111// ---------------------------------------------------------------------------
112
113
114// ----- Execution -------------------------------------------------------
116{
117
118 // --- Timer and counters
119 TStopwatch timerStep;
120 TStopwatch timerTot;
121 timerTot.Start();
122 std::map<ECbmModuleId, size_t> numDigisTs;
123 std::map<ECbmModuleId, size_t> numDigisEv;
124
125 // --- Clear output vector
126 fEvents->clear();
127
128 // --- If the input is already CbmDigiTimeslice (from unpacking), use that directly
129 if (fTimeslice) {
130 timerStep.Start();
132 (*fAlgo)(cbm::algo::DigiData(fTimeslice->fData), *fTriggers, std::nullopt).first);
133 timerStep.Stop();
134 fTimeBuildEvt += timerStep.RealTime();
135 for (const auto& entry : fConfig->fWindows)
136 numDigisTs[entry.first] = GetNumDigis(fTimeslice->fData, entry.first);
137 }
138
139 // --- If input is not DigiTimeslice: construct a transientDigiTimeslice from the data in CbmDigiManager
140 else {
141 timerStep.Start();
143 for (const auto& entry : fConfig->fWindows)
144 numDigisTs[entry.first] = ts.fData.Size(entry.first);
145 timerStep.Stop();
146 fTimeFillTs += timerStep.RealTime();
147 timerStep.Start();
148 *fEvents =
150 timerStep.Stop();
151 fTimeBuildEvt += timerStep.RealTime();
152 }
153
154 // Apply event selector if desired
155 if (fSelector) {
156 timerStep.Start();
157 auto noTrigger = [&](CbmDigiEvent& ev) { return !(*fSelector)(cbm::algo::DigiEvent(ev)); };
158 auto removeIt = std::remove_if(fEvents->begin(), fEvents->end(), noTrigger);
159 fEvents->erase(removeIt, fEvents->end());
160 timerStep.Stop();
161 fTimeSelectorEvt += timerStep.RealTime();
162 }
163
164 // --- Timeslice statistics
165 size_t numTriggers = fTriggers->size();
166 size_t numEvents = fEvents->size();
167 for (const auto& entry : fConfig->fWindows) {
168 for (auto& event : (*fEvents)) {
169 numDigisEv[entry.first] += GetNumDigis(event.fData, entry.first);
170 }
171 }
172
173 // --- Timeslice log
174 timerTot.Stop();
175 fTimeTot += timerTot.RealTime();
176 stringstream logOut;
177 logOut << setw(15) << left << GetName() << " [";
178 logOut << fixed << setw(8) << setprecision(1) << right << timerTot.RealTime() * 1000. << " ms] ";
179 logOut << "TS " << fNumTs << ", triggers " << numTriggers << ", events " << numEvents;
180
181 for (const auto& entry : fConfig->fWindows) {
182 auto system = entry.first;
183 logOut << ", frac " << ToString(system) << " digis "
184 << 100. * double(numDigisEv[system]) / double(numDigisTs[system]) << " %";
185 }
186
187 LOG(info) << logOut.str();
188
189 // --- Run statistics
190 fNumTs++;
191 fNumTriggers += numTriggers;
192 fNumEvents += numEvents;
193
194 for (const auto& entry : fConfig->fWindows) {
195 auto system = entry.first;
196 fNumDigisTs[system] += numDigisTs[system];
197 fNumDigisEv[system] += numDigisEv[system];
198 }
199}
200// ----------------------------------------------------------------------------
201
202
203// ----- End-of-timeslice action ------------------------------------------
205{
206 LOG(info) << "=====================================";
207 LOG(info) << GetName() << ": Run summary";
208 LOG(info) << "Timeslices : " << fNumTs;
209 LOG(info) << "Triggers : " << fNumTriggers;
210 LOG(info) << "Events : " << fNumEvents;
211 for (const auto& entry : fConfig->fWindows) {
212 auto system = entry.first;
213 LOG(info) << setw(4) << left << ToString(system) << " digis in TS : " << fNumDigisTs[system];
214 LOG(info) << setw(4) << left << ToString(system) << " digis in events : " << fNumDigisEv[system] << " = " << fixed
215 << setprecision(2) << 100. * double(fNumDigisEv[system]) / double(fNumDigisTs[system]) << " %";
216 }
217 LOG(info) << "Time / TS : " << fixed << setprecision(2) << 1000. * fTimeTot / double(fNumTs) << " ms";
218 LOG(info) << "Time fill TS : " << fixed << setprecision(2) << 1000. * fTimeFillTs / double(fNumTs)
219 << " ms = " << 100. * fTimeFillTs / fTimeTot << " %";
220 LOG(info) << "Time build events : " << fixed << setprecision(2) << 1000. * fTimeBuildEvt / double(fNumTs)
221 << " ms = " << 100. * fTimeBuildEvt / fTimeTot << " %";
222 LOG(info) << "Time selector events : " << fixed << setprecision(2) << 1000. * fTimeSelectorEvt / double(fNumTs)
223 << " ms = " << 100. * fTimeSelectorEvt / fTimeTot << " %";
224 LOG(info) << "=====================================";
225}
226// ----------------------------------------------------------------------------
227
228
229// ----- Number of digis in the timeslice ----------------------------------
231{
232 size_t result = 0;
233 switch (system) {
234 case ECbmModuleId::kSts: result = data.fSts.fDigis.size(); break;
235 case ECbmModuleId::kRich: result = data.fRich.fDigis.size(); break;
236 case ECbmModuleId::kMuch: result = data.fMuch.fDigis.size(); break;
237 case ECbmModuleId::kTrd: result = data.fTrd.fDigis.size(); break;
238 case ECbmModuleId::kTrd2d: result = data.fTrd2d.fDigis.size(); break;
239 case ECbmModuleId::kTof: result = data.fTof.fDigis.size(); break;
240 case ECbmModuleId::kPsd: result = data.fPsd.fDigis.size(); break;
241 case ECbmModuleId::kBmon: result = data.fBmon.fDigis.size(); break;
242 default: result = 0; break;
243 }
244 return result;
245}
246// ----------------------------------------------------------------------------
247
248
249// ----- Initialisation ---------------------------------------------------
251{
252 // --- Get FairRootManager instance
253 FairRootManager* ioman = FairRootManager::Instance();
254 assert(ioman);
255
256 LOG(info) << "==================================================";
257 LOG(info) << GetName() << ": Initialising...";
258
259 // --- Check input data
260 // --- DigiTimeslice: Unpacked data from FLES
261 fTimeslice = ioman->InitObjectAs<const CbmDigiTimeslice*>("DigiTimeslice.");
262 if (fTimeslice) {
263 LOG(info) << "--- Found branch DigiTimeslice.";
264 }
265 // --- DigiManager: Simulated digi data
266 else {
268 fDigiMan->Init();
269 for (const auto& entry : fConfig->fWindows) {
270 auto system = entry.first;
271 if (!fDigiMan->IsPresent(system)) {
272 LOG(warn) << GetName() << ": No digi branch for " << ToString(system);
273 }
274 LOG(info) << "--- Found digi branch for " << ToString(system);
275 }
276 }
277
278 // --- Initialize diagnostics
279 for (const auto& entry : fConfig->fWindows) {
280 auto system = entry.first;
281 fNumDigisTs[system] = 0;
282 fNumDigisEv[system] = 0;
283 }
284
285 // --- Get input data (triggers)
286 fTriggers = ioman->InitObjectAs<std::vector<double> const*>("Trigger");
287 if (!fTriggers) {
288 LOG(fatal) << GetName() << ": No Trigger branch!" << endl;
289 return kFATAL;
290 }
291 LOG(info) << "--- Found branch Trigger";
292
293 // --- Register output array (CbmDigiEvent)
294 if (ioman->GetObject("DigiEvent")) {
295 LOG(fatal) << GetName() << ": Branch DigiEvent already exists!";
296 return kFATAL;
297 }
298 fEvents = new vector<CbmDigiEvent>;
299 ioman->RegisterAny("DigiEvent", fEvents, IsOutputBranchPersistent("DigiEvent"));
300 if (!fEvents) {
301 LOG(fatal) << GetName() << ": Output branch could not be created!";
302 return kFATAL;
303 }
304 LOG(info) << "--- Registered branch DigiEvent";
305
306 // --- Configure algorithm
308 LOG(info) << fAlgo->ToString();
309
310 // --- Log selector
311 if (fSelector)
312 LOG(info) << fSelector->ToString();
313 else
314 LOG(info) << "--- No event selection";
315
316 LOG(info) << "==================================================";
317
318 return kSUCCESS;
319}
320// ----------------------------------------------------------------------------
321
ClassImp(CbmConverterManager)
std::string ToString(ECbmModuleId modId)
Definition CbmDefs.cxx:70
ECbmModuleId
Definition CbmDefs.h:39
@ kTrd
Transition Radiation Detector.
@ kTof
Time-of-flight Detector.
@ 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.
std::vector< CbmBmonDigi > fDigis
Data vector.
Abstract base class for CBM digi branches.
virtual boost::any GetBranchContainer() const
Get branch pointer.
Collection of digis from all detector systems.
Definition CbmDigiData.h:32
CbmPsdDigiData fPsd
PSD data.
Definition CbmDigiData.h:42
CbmTrdDigiData fTrd2d
TRD2D data.
Definition CbmDigiData.h:40
CbmTrdDigiData fTrd
TRD data.
Definition CbmDigiData.h:39
CbmTofDigiData fTof
TOF data.
Definition CbmDigiData.h:41
CbmStsDigiData fSts
STS data.
Definition CbmDigiData.h:36
CbmFsdDigiData fFsd
FSD data.
Definition CbmDigiData.h:43
size_t Size(ECbmModuleId system) const
Size of detector data.
Definition CbmDigiData.h:82
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.
static Bool_t IsPresent(ECbmModuleId systemId)
Presence of a digi branch.
CbmDigiBranchBase * GetBranch(ECbmModuleId system)
Access to a digi branch.
InitStatus Init()
Initialisation.
static CbmDigiManager * Instance()
Static instance.
Collection of digis from all detector systems within one timeslice.
CbmDigiData fData
Timeslice data.
std::vector< CbmFsdDigi > fDigis
Data vector.
std::vector< CbmMuchDigi > fDigis
Data vector.
std::vector< CbmPsdDigi > fDigis
Data vector.
std::vector< CbmRichDigi > fDigis
Data vector.
std::vector< CbmStsDigi > fDigis
Data vector.
Task class for associating digis to events.
std::map< ECbmModuleId, size_t > fNumDigisTs
Event builder configuration.
virtual InitStatus Init()
Task initialisation.
const CbmDigiTimeslice * fTimeslice
CbmDigiTimeslice FillTimeSlice()
Construct a DigiTimeslice from the data in CbmDigiManager.
CbmTaskBuildEvents()
Constructor.
std::unique_ptr< cbm::algo::evbuild::EventBuilder > fAlgo
Event selector.
size_t GetNumDigis(const CbmDigiData &data, ECbmModuleId system)
Number of digis for a given system.
virtual ~CbmTaskBuildEvents()
Destructor.
const std::vector< double > * fTriggers
Input data (from simulation)
std::map< ECbmModuleId, size_t > fNumDigisEv
std::vector< CbmDigiEvent > * fEvents
Input data (triggers)
virtual void Finish()
Finish timeslice.
std::unique_ptr< cbm::algo::evbuild::DigiEventSelector > fSelector
Output data (events)
virtual void Exec(Option_t *opt)
Task execution.
CbmDigiManager * fDigiMan
Input data (from unpacking)
std::unique_ptr< cbm::algo::evbuild::EventBuilderConfig > fConfig
Algorithm.
std::vector< CbmTofDigi > fDigis
Data vector.
std::vector< CbmTrdDigi > fDigis
Data vector.
Constructs CbmDigiEvents out of CbmDigiTimeslices.
Hash for CbmL1LinkKey.
Collection of digis from all detector systems.
Definition DigiData.h:31
static std::vector< CbmDigiEvent > ToCbmDigiEvents(const std::vector< DigiEvent > &events)
Definition DigiData.cxx:108