CbmRoot
Loading...
Searching...
No Matches
CbmSourceDigiTimeslice.cxx
Go to the documentation of this file.
1/* Copyright (C) 2024 GSI Helmholtzzentrum fuer Schwerionenforschung, Darmstadt
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Volker Friese[committer] */
4
6
7#include "CbmTimeSlice.h"
8#include "CbmTsEventHeader.h"
9
10#include <FairRootManager.h>
11#include <FairRun.h>
12#include <Logger.h>
13
14#include <utility>
15
16
17// ----- Constructor ------------------------------------------------------
19// ----------------------------------------------------------------------------
20
21
22// ----- Close ------------------------------------------------------------
24{
25 LOG(info) << "Source: Closing after " << fNumTs << " timeslices";
26 // Clear output vectors
27 fBmonDigis->clear();
28 fMvdDigis->clear();
29 fStsDigis->clear();
30 fMuchDigis->clear();
31 fTrdDigis->clear();
32 fTofDigis->clear();
33 fRichDigis->clear();
34 fFsdDigis->clear();
35}
36// ----------------------------------------------------------------------------
37
38
39// ----- Initialisation ---------------------------------------------------
40template<typename TVecobj>
41Bool_t CbmSourceDigiTimeslice::RegisterVector(FairRootManager* ioman, std::vector<TVecobj>*& vec)
42{
43 if (ioman->GetObject(TVecobj::GetBranchName())) {
44 LOG(fatal) << "Source: Branch " << TVecobj::GetBranchName() << " already exists!";
45 return kFALSE;
46 }
47
48 ioman->RegisterAny(TVecobj::GetBranchName(), vec, kTRUE);
49 LOG(info) << "Source: Registered branch " << TVecobj::GetBranchName() << " at " << vec;
50
51 return kTRUE;
52}
53
55{
56 // Create input archive
57 fArchive = std::make_unique<cbm::algo::RecoResultsInputArchive>(fInputFileName);
58 LOG(info) << "Source: Reading from input archive " << fInputFileName;
59 auto desc = fArchive->descriptor();
60 LOG(info) << " - Time created: " << desc.time_created();
61 LOG(info) << " - Host name : " << desc.hostname();
62 LOG(info) << " - User name : " << desc.username();
63
64 // Create and register the Digi tree branches
65 FairRootManager* ioman = FairRootManager::Instance();
66 assert(ioman);
67
68 if (!(fTsEventHeader = dynamic_cast<CbmTsEventHeader*>(FairRun::Instance()->GetEventHeader()))) {
69 LOG(fatal) << "CbmSourceDigiTimeslice::Init() no CbmTsEventHeader was added to the run. Without it, we can not "
70 "store the UTC of the "
71 "Timeslices correctly. Hence, this causes a fatal. Please add it in the steering macro to the Run.";
72 return kFALSE;
73 }
74
75 // TimeSlice. branch initialization
76 if (ioman->GetObject("TimeSlice.")) {
77 LOG(fatal) << "Source: Branch TimeSlice. already exists!";
78 return kFALSE;
79 }
80 else {
81 // NOTE: the max time of timeslice is 1.28e8, taken from CbmRecoUnpack.cxx
82 fTimeslice = new CbmTimeSlice(0., 1.28e8 + 1.28e6);
83 ioman->Register("TimeSlice.", "DAQ", fTimeslice, kTRUE);
84 }
85
86 fBmonDigis = new std::vector<CbmBmonDigi>();
87 if (kFALSE == RegisterVector<CbmBmonDigi>(ioman, fBmonDigis)) {
88 return kFALSE;
89 }
90
91
92 fMvdDigis = new std::vector<CbmMvdDigi>();
93 if (kFALSE == RegisterVector<CbmMvdDigi>(ioman, fMvdDigis)) {
94 return kFALSE;
95 }
96
97
98 fStsDigis = new std::vector<CbmStsDigi>();
99 if (kFALSE == RegisterVector<CbmStsDigi>(ioman, fStsDigis)) {
100 return kFALSE;
101 }
102
103
104 fMuchDigis = new std::vector<CbmMuchDigi>();
105 if (kFALSE == RegisterVector<CbmMuchDigi>(ioman, fMuchDigis)) {
106 return kFALSE;
107 }
108
109
110 fTrdDigis = new std::vector<CbmTrdDigi>();
111 if (kFALSE == RegisterVector<CbmTrdDigi>(ioman, fTrdDigis)) {
112 return kFALSE;
113 }
114
115
116 fTofDigis = new std::vector<CbmTofDigi>();
117 if (kFALSE == RegisterVector<CbmTofDigi>(ioman, fTofDigis)) {
118 return kFALSE;
119 }
120
121
122 fRichDigis = new std::vector<CbmRichDigi>();
123 if (kFALSE == RegisterVector<CbmRichDigi>(ioman, fRichDigis)) {
124 return kFALSE;
125 }
126
127
128 fFsdDigis = new std::vector<CbmFsdDigi>();
129 if (kFALSE == RegisterVector<CbmFsdDigi>(ioman, fFsdDigis)) {
130 return kFALSE;
131 }
132
133
134 return kTRUE;
135}
136// ----------------------------------------------------------------------------
137
138
139// ----- Read one time slice from archive ---------------------------------
141{
142
143 LOG(debug) << "CbmSourceTimeslice::ReadEvent";
144 // Clear output vectors
145 fBmonDigis->clear();
146 fMvdDigis->clear();
147 fStsDigis->clear();
148 fMuchDigis->clear();
149 fTrdDigis->clear();
150 fTofDigis->clear();
151 fRichDigis->clear();
152 fFsdDigis->clear();
153 LOG(debug) << "CbmSourceTimeslice: Cleared digi vectors";
154
155
156 // Get next timeslice data from archive. Stop run if end of archive is reached.
157 auto results = fArchive->get();
158 if (fArchive->eos()) {
159 LOG(info) << "Source: End of input archive; terminating run";
160 return 1;
161 }
162 LOG(debug) << "CbmSourceTimeslice: read fles timeslice";
163
164 // Move event data from input archive to ROOT tree
165 if (results == nullptr) {
166 LOG(error) << "Source: Failed to read RecoResults from archive";
167 return 1;
168 }
169 LOG(info) << "Source: Reading TS " << fNumTs << ", index " << results->TsIndex() << ", start "
170 << results->TsStartTime() << ", contains Digis: \n"
171 << "MVD = " << results->MvdDigis().size() //
172 << " STS = " << results->StsDigis().size() << " MUCH = " << results->MuchDigis().size() //
173 << " TOF = " << results->TofDigis().size() << " BMON = " << results->BmonDigis().size() //
174 << " TRD = " << results->TrdDigis().size() << " TRD2D=" << results->Trd2dDigis().size() //
175 << " RICH=" << results->RichDigis().size() << " FSD=" << results->FsdDigis().size();
176
177 fTsEventHeader->SetTsIndex(results->TsIndex());
178 fTsEventHeader->SetTsStartTime(results->TsStartTime());
179
180 fTimeslice->SetStartTime(results->TsStartTime());
181
182 std::move(results->BmonDigis().begin(), results->BmonDigis().end(), std::back_inserter(*fBmonDigis));
183 std::move(results->StsDigis().begin(), results->StsDigis().end(), std::back_inserter(*fStsDigis));
184 std::move(results->MuchDigis().begin(), results->MuchDigis().end(), std::back_inserter(*fMuchDigis));
185 std::move(results->Trd2dDigis().begin(), results->Trd2dDigis().end(), std::back_inserter(*fTrdDigis));
186 std::move(results->TrdDigis().begin(), results->TrdDigis().end(), std::back_inserter(*fTrdDigis));
187 std::move(results->TofDigis().begin(), results->TofDigis().end(), std::back_inserter(*fTofDigis));
188 std::move(results->RichDigis().begin(), results->RichDigis().end(), std::back_inserter(*fRichDigis));
189 std::move(results->FsdDigis().begin(), results->FsdDigis().end(), std::back_inserter(*fFsdDigis));
190
191 // CbmMvdRawDigi from the source ts needs to be converted into CbmMvdDigi (temporary fix)
192 for (const auto& mvdRawDigi : results->MvdDigis()) {
193 fMvdDigis->emplace_back(mvdRawDigi);
194 }
195
196
197 // Time-sort the TRD vector as we merged TRD1D and TRD2D (needed for compatibility with legacy unpackers)
199
200 // Update counters
201 fNumTs++;
202
203 return 0;
204}
205// ----------------------------------------------------------------------------
206
207
ClassImp(CbmConverterManager)
int Int_t
bool Bool_t
Source class for reading from files resulting from online processing (containing raw Digis)
std::vector< CbmStsDigi > * fStsDigis
std::vector< CbmBmonDigi > * fBmonDigis
virtual void Close()
Close source after end of run.
CbmSourceDigiTimeslice(const char *fileName="")
Constructor.
Bool_t RegisterVector(FairRootManager *ioman, std::vector< TVecobj > *&vec)
std::vector< CbmFsdDigi > * fFsdDigis
std::vector< CbmTofDigi > * fTofDigis
Int_t ReadEvent(UInt_t=0)
Read one time slice from file.
std::vector< CbmMvdDigi > * fMvdDigis
std::enable_if< std::is_member_function_pointer< decltype(&TVecobj::GetTime)>::value, void >::type Timesort(std::vector< TVecobj > *vec=nullptr)
std::vector< CbmMuchDigi > * fMuchDigis
std::vector< CbmRichDigi > * fRichDigis
std::vector< CbmTrdDigi > * fTrdDigis
CbmTsEventHeader * fTsEventHeader
virtual Bool_t Init()
Initialisation.
std::unique_ptr< cbm::algo::RecoResultsInputArchive > fArchive
Bookkeeping of time-slice content.