CbmRoot
Loading...
Searching...
No Matches
CbmTaskMakeRecoEvents.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
5
7
8#include "CbmDefs.h"
9#include "CbmEvent.h"
10
11#include <FairRootManager.h>
12#include <Logger.h>
13
14#include <TClonesArray.h>
15#include <TStopwatch.h>
16
17#include <algorithm>
18#include <cassert>
19#include <iomanip>
20#include <vector>
21
22
23using namespace std;
24
25
26// ----- Constructor -----------------------------------------------------
27CbmTaskMakeRecoEvents::CbmTaskMakeRecoEvents() : FairTask("MakeRecoEvents") {}
28// ---------------------------------------------------------------------------
29
30
31// ----- Destructor ------------------------------------------------------
33// ---------------------------------------------------------------------------
34
35
36// ----- Execution -------------------------------------------------------
38{
39
40 // --- Timer and counters
41 TStopwatch timer;
42 timer.Start();
43
44 // --- No action if no DigiEvent branch is present
45 if (!fDigiEvents) return;
46
47 // --- Clear output arrays
48 fBmonDigis->clear();
49 fMvdDigis->clear();
50 fStsDigis->clear();
51 fRichDigis->clear();
52 fMuchDigis->clear();
53 fTrdDigis->clear();
54 fTofDigis->clear();
55 fFsdDigis->clear();
56
57 //fRecoEvents->Clear(); //causes memory leak
58 fRecoEvents->Delete();
59
60 // --- Event loop
61 Int_t eventNr = 0;
62 for (auto& digiEvent : *fDigiEvents) {
63
64 // --- Create CbmEvent object
65 CbmEvent* recoEvent = new ((*fRecoEvents)[eventNr]) CbmEvent(eventNr);
66
67 // --- Copy BMON digis
68 FillTree<CbmBmonDigi>(digiEvent.fData.fBmon.fDigis, fBmonDigis, recoEvent, ECbmDataType::kBmonDigi);
69
70 // --- Copy MVD digis
71 FillTreeMvd(digiEvent.fData.fMvd.fDigis, fMvdDigis, recoEvent);
72
73 // --- Copy STS digis
74 FillTree<CbmStsDigi>(digiEvent.fData.fSts.fDigis, fStsDigis, recoEvent, ECbmDataType::kStsDigi);
75
76 // --- Copy RICH digis
77 FillTree<CbmRichDigi>(digiEvent.fData.fRich.fDigis, fRichDigis, recoEvent, ECbmDataType::kRichDigi);
78
79 // --- Copy MUCH digis
80 FillTree<CbmMuchDigi>(digiEvent.fData.fMuch.fDigis, fMuchDigis, recoEvent, ECbmDataType::kMuchDigi);
81
82 // --- Copy TRD digis
83 FillTree<CbmTrdDigi>(digiEvent.fData.fTrd.fDigis, fTrdDigis, recoEvent, ECbmDataType::kTrdDigi);
84
85 // --- Copy TOF digis
86 FillTree<CbmTofDigi>(digiEvent.fData.fTof.fDigis, fTofDigis, recoEvent, ECbmDataType::kTofDigi);
87
88 // --- Copy FSD digis
89 FillTree<CbmFsdDigi>(digiEvent.fData.fFsd.fDigis, fFsdDigis, recoEvent, ECbmDataType::kFsdDigi);
90
91 eventNr++;
92 } //# events
93
94 // --- Timeslice log
95 timer.Stop();
96 stringstream logOut;
97 logOut << setw(20) << left << GetName() << " [";
98 logOut << fixed << setw(8) << setprecision(1) << right << timer.RealTime() * 1000. << " ms] ";
99 logOut << "TS " << fNumTs << ", events " << fDigiEvents->size() << ", Digis: BMON " << fBmonDigis->size() << " MVD "
100 << fMvdDigis->size() << " STS " << fStsDigis->size() << " RICH " << fRichDigis->size() << " MUCH "
101 << fMuchDigis->size() << " TRD " << fTrdDigis->size() << " TOF " << fTofDigis->size() << " FSD "
102 << fFsdDigis->size();
103 LOG(info) << logOut.str();
104
105 // --- Run statistics
106 fNumTs++;
107 fTimeTot += timer.RealTime();
108 assert(fDigiEvents->size() == static_cast<size_t>(fRecoEvents->GetEntriesFast()));
109 fNumEvents += fDigiEvents->size();
110}
111// ----------------------------------------------------------------------------
112
113
114// ----- End-of-timeslice action ------------------------------------------
116{
117 LOG(info) << "=====================================";
118 LOG(info) << GetName() << ": Run summary";
119 LOG(info) << "Timeslices : " << fNumTs;
120 LOG(info) << "Events : " << fNumEvents;
121 LOG(info) << "Time / TS : " << fixed << setprecision(2) << 1000. * fTimeTot / double(fNumTs) << " ms";
122 LOG(info) << "=====================================";
123}
124// ----------------------------------------------------------------------------
125
126
127// ----- Initialisation ---------------------------------------------------
129{
130
131 LOG(info) << "==================================================";
132 LOG(info) << GetName() << ": Initialising ";
133
134 // --- Get FairRootManager instance
135 FairRootManager* frm = FairRootManager::Instance();
136 assert(frm);
137
138 // --- Try to get input vector (CbmDigiEvent)
139 fDigiEvents = frm->InitObjectAs<const std::vector<CbmDigiEvent>*>("DigiEvent");
140
141 // --- If DigiEvents are present, create Digi and CbmEvent branches
142 if (fDigiEvents) {
143 LOG(info) << GetName() << ": Found branch DigiEvent";
144
145 // --- Event
146 if (frm->GetObject("CbmEvent")) {
147 LOG(error) << GetName() << ": Found branch CbmEvent! Aborting...";
148 return kFATAL;
149 }
150 fRecoEvents = new TClonesArray("CbmEvent", 1);
151 frm->Register("CbmEvent", "Reco events", fRecoEvents, IsOutputBranchPersistent("CbmEvent"));
152 if (frm->GetObject("StsDigi")) {
153 LOG(error) << GetName() << ": Found branch StsDigi! Aborting...";
154 return kFATAL;
155 }
156
157 // --- BMON digis
158 fBmonDigis = new std::vector<CbmBmonDigi>;
159 frm->RegisterAny("BmonDigi", fBmonDigis, kFALSE);
160
161 // --- MVD digis
162 fMvdDigis = new std::vector<CbmMvdDigi>;
163 frm->RegisterAny("MvdDigi", fMvdDigis, kFALSE);
164
165 // --- STS digis
166 fStsDigis = new std::vector<CbmStsDigi>;
167 frm->RegisterAny("StsDigi", fStsDigis, kFALSE);
168
169 // --- RICH digis
170 fRichDigis = new std::vector<CbmRichDigi>;
171 frm->RegisterAny("RichDigi", fRichDigis, kFALSE);
172
173 // --- MUCH digis
174 fMuchDigis = new std::vector<CbmMuchDigi>;
175 frm->RegisterAny("MuchDigi", fMuchDigis, kFALSE);
176
177 // --- TRD digis
178 fTrdDigis = new std::vector<CbmTrdDigi>;
179 frm->RegisterAny("TrdDigi", fTrdDigis, kFALSE);
180
181 // --- TOF digis
182 fTofDigis = new std::vector<CbmTofDigi>;
183 frm->RegisterAny("TofDigi", fTofDigis, kFALSE);
184
185 // --- FSD digis
186 fFsdDigis = new std::vector<CbmFsdDigi>;
187 frm->RegisterAny("FsdDigi", fFsdDigis, kFALSE);
188 }
189
190 // --- If no DigiEvent branch is present, there must be a CbmEvent branch
191 else {
192 fRecoEvents = dynamic_cast<TClonesArray*>(frm->GetObject("CbmEvent"));
193 if (fRecoEvents == nullptr) {
194 LOG(error) << GetName() << ": Neither DigiEvent nor CbmEvent branch found! Aborting...";
195 return kFATAL;
196 }
197 }
198
199 LOG(info) << "==================================================";
200
201 return kSUCCESS;
202}
203// ----------------------------------------------------------------------------
204
ClassImp(CbmConverterManager)
int Int_t
Class characterising one event by a collection of links (indices) to data objects,...
Definition CbmEvent.h:34
Task class for interfacing storable raw events in the CbmDigiEvent format to the current offline reco...
std::vector< CbmTrdDigi > * fTrdDigis
std::vector< CbmBmonDigi > * fBmonDigis
size_t fNumTs
Number of processed timeslices.
size_t fNumEvents
Number of events.
std::vector< CbmMvdDigi > * fMvdDigis
virtual void Exec(Option_t *opt)
Task execution.
std::vector< CbmRichDigi > * fRichDigis
const std::vector< CbmDigiEvent > * fDigiEvents
virtual ~CbmTaskMakeRecoEvents()
Destructor.
double fTimeTot
Execution time.
std::vector< CbmFsdDigi > * fFsdDigis
std::vector< CbmTofDigi > * fTofDigis
std::vector< CbmMuchDigi > * fMuchDigis
void FillTree(gsl::span< const Digi > inVec, std::vector< Digi > *outVec, CbmEvent *event, ECbmDataType digiType)
Fill the tree structure with digis from CbmDigiEvent.
virtual void Finish()
Finish timeslice.
virtual InitStatus Init()
Task initialisation.
std::vector< CbmStsDigi > * fStsDigis
void FillTreeMvd(gsl::span< const CbmMvdRawDigi > inVec, std::vector< CbmMvdDigi > *outVec, CbmEvent *event)
Fill the tree with MVD digis.
Hash for CbmL1LinkKey.