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 fStsDigis->clear();
50 fRichDigis->clear();
51 fMuchDigis->clear();
52 fTrdDigis->clear();
53 fTofDigis->clear();
54 fPsdDigis->clear();
55
56 //fRecoEvents->Clear(); //causes memory leak
57 fRecoEvents->Delete();
58
59 // --- Event loop
60 Int_t eventNr = 0;
61 for (auto& digiEvent : *fDigiEvents) {
62
63 // --- Create CbmEvent object
64 CbmEvent* recoEvent = new ((*fRecoEvents)[eventNr]) CbmEvent(eventNr);
65
66 // --- Copy Bmon digis
67 FillTree<CbmBmonDigi>(digiEvent.fData.fBmon.fDigis, fBmonDigis, recoEvent, ECbmDataType::kBmonDigi);
68
69 // --- Copy STS digis
70 FillTree<CbmStsDigi>(digiEvent.fData.fSts.fDigis, fStsDigis, recoEvent, ECbmDataType::kStsDigi);
71
72 // --- Copy RICH digis
73 FillTree<CbmRichDigi>(digiEvent.fData.fRich.fDigis, fRichDigis, recoEvent, ECbmDataType::kRichDigi);
74
75 // --- Copy MUCH digis
76 FillTree<CbmMuchDigi>(digiEvent.fData.fMuch.fDigis, fMuchDigis, recoEvent, ECbmDataType::kMuchDigi);
77
78 // --- Copy TRD digis
79 FillTree<CbmTrdDigi>(digiEvent.fData.fTrd.fDigis, fTrdDigis, recoEvent, ECbmDataType::kTrdDigi);
80
81 // --- Copy TOF digis
82 FillTree<CbmTofDigi>(digiEvent.fData.fTof.fDigis, fTofDigis, recoEvent, ECbmDataType::kTofDigi);
83
84 // --- Copy PSD digis
85 FillTree<CbmPsdDigi>(digiEvent.fData.fPsd.fDigis, fPsdDigis, recoEvent, ECbmDataType::kPsdDigi);
86
87 eventNr++;
88 } //# events
89
90 // --- Timeslice log
91 timer.Stop();
92 stringstream logOut;
93 logOut << setw(20) << left << GetName() << " [";
94 logOut << fixed << setw(8) << setprecision(1) << right << timer.RealTime() * 1000. << " ms] ";
95 logOut << "TS " << fNumTs << ", events " << fDigiEvents->size() << ", Digis: Bmon " << fBmonDigis->size() << " STS "
96 << fStsDigis->size() << " RICH " << fRichDigis->size() << " MUCH " << fMuchDigis->size() << " TRD "
97 << fTrdDigis->size() << " TOF " << fTofDigis->size() << " PSD " << fPsdDigis->size();
98 LOG(info) << logOut.str();
99
100 // --- Run statistics
101 fNumTs++;
102 fTimeTot += timer.RealTime();
103 assert(fDigiEvents->size() == static_cast<size_t>(fRecoEvents->GetEntriesFast()));
104 fNumEvents += fDigiEvents->size();
105}
106// ----------------------------------------------------------------------------
107
108
109// ----- End-of-timeslice action ------------------------------------------
111{
112 LOG(info) << "=====================================";
113 LOG(info) << GetName() << ": Run summary";
114 LOG(info) << "Timeslices : " << fNumTs;
115 LOG(info) << "Events : " << fNumEvents;
116 LOG(info) << "Time / TS : " << fixed << setprecision(2) << 1000. * fTimeTot / double(fNumTs) << " ms";
117 LOG(info) << "=====================================";
118}
119// ----------------------------------------------------------------------------
120
121
122// ----- Initialisation ---------------------------------------------------
124{
125
126 LOG(info) << "==================================================";
127 LOG(info) << GetName() << ": Initialising ";
128
129 // --- Get FairRootManager instance
130 FairRootManager* frm = FairRootManager::Instance();
131 assert(frm);
132
133 // --- Try to get input vector (CbmDigiEvent)
134 fDigiEvents = frm->InitObjectAs<const std::vector<CbmDigiEvent>*>("DigiEvent");
135
136 // --- If DigiEvents are present, create Digi and CbmEvent branches
137 if (fDigiEvents) {
138 LOG(info) << GetName() << ": Found branch DigiEvent";
139
140 // --- Event
141 if (frm->GetObject("CbmEvent")) {
142 LOG(error) << GetName() << ": Found branch CbmEvent! Aborting...";
143 return kFATAL;
144 }
145 fRecoEvents = new TClonesArray("CbmEvent", 1);
146 frm->Register("CbmEvent", "Reco events", fRecoEvents, IsOutputBranchPersistent("CbmEvent"));
147 if (frm->GetObject("StsDigi")) {
148 LOG(error) << GetName() << ": Found branch StsDigi! Aborting...";
149 return kFATAL;
150 }
151
152 // --- Bmon digis
153 fBmonDigis = new std::vector<CbmBmonDigi>;
154 frm->RegisterAny("BmonDigi", fBmonDigis, kFALSE);
155
156 // --- STS digis
157 fStsDigis = new std::vector<CbmStsDigi>;
158 frm->RegisterAny("StsDigi", fStsDigis, kFALSE);
159
160 // --- RICH digis
161 fRichDigis = new std::vector<CbmRichDigi>;
162 frm->RegisterAny("RichDigi", fRichDigis, kFALSE);
163
164 // --- MUCH digis
165 fMuchDigis = new std::vector<CbmMuchDigi>;
166 frm->RegisterAny("MuchDigi", fMuchDigis, kFALSE);
167
168 // --- TRD digis
169 fTrdDigis = new std::vector<CbmTrdDigi>;
170 frm->RegisterAny("TrdDigi", fTrdDigis, kFALSE);
171
172 // --- TOF digis
173 fTofDigis = new std::vector<CbmTofDigi>;
174 frm->RegisterAny("TofDigi", fTofDigis, kFALSE);
175
176 // --- PSD digis
177 fPsdDigis = new std::vector<CbmPsdDigi>;
178 frm->RegisterAny("PsdDigi", fPsdDigis, kFALSE);
179 }
180
181 // --- If no DigiEvent branch is present, there must be a CbmEvent branch
182 else {
183 fRecoEvents = dynamic_cast<TClonesArray*>(frm->GetObject("CbmEvent"));
184 if (fRecoEvents == nullptr) {
185 LOG(error) << GetName() << ": Neither DigiEvent nor CbmEvent branch found! Aborting...";
186 return kFATAL;
187 }
188 }
189
190 LOG(info) << "==================================================";
191
192 return kSUCCESS;
193}
194// ----------------------------------------------------------------------------
195
ClassImp(CbmConverterManager)
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.
virtual void Exec(Option_t *opt)
Task execution.
std::vector< CbmRichDigi > * fRichDigis
std::vector< CbmPsdDigi > * fPsdDigis
const std::vector< CbmDigiEvent > * fDigiEvents
virtual ~CbmTaskMakeRecoEvents()
Destructor.
double fTimeTot
Execution time.
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
Hash for CbmL1LinkKey.