CbmRoot
Loading...
Searching...
No Matches
CbmBuildEventsIdealNew.cxx
Go to the documentation of this file.
1/* Copyright (C) 2020 GSI Helmholtzzentrum fuer Schwerionenforschung, Darmstadt
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Volker Friese [committer] */
4
11
12#include "CbmDigiManager.h"
13#include "CbmEvent.h"
14#include "CbmEventStore.h"
15#include "CbmFsdDigi.h"
16#include "CbmLink.h"
17#include "CbmMatch.h"
18#include "CbmModuleList.h"
19#include "CbmMuchDigi.h"
20#include "CbmMvdDigi.h"
21#include "CbmPsdDigi.h"
22#include "CbmRichDigi.h"
23#include "CbmStsDigi.h"
24#include "CbmTofDigi.h"
25#include "CbmTrdDigi.h"
26
27#include <FairRootManager.h>
28#include <Logger.h>
29
30#include <TClonesArray.h>
31#include <TStopwatch.h>
32
33#include <cassert>
34#include <iomanip>
35#include <iostream>
36
37using namespace std;
38
39
40// ===== Constructor =====================================================
41CbmBuildEventsIdealNew::CbmBuildEventsIdealNew() : FairTask("BuildEventsIdealNew") {}
42// ===========================================================================
43
44
45// ===== Destructor ======================================================
47// ===========================================================================
48
49
50// ===== Task execution ==================================================
52{
53
54 TStopwatch timer;
55 timer.Start();
56 fEvents->Delete();
57 std::map<Int_t, CbmEventStore*> eventMap;
58 Int_t nEvents = 0;
59 UInt_t nDigisTot = 0;
60 UInt_t nDigisNoise = 0;
61
62 for (ECbmModuleId& system : fSystems) {
63 if (!fDigiMan->IsPresent(system)) continue;
64 if (!fDigiMan->IsMatchPresent(system)) continue;
65 Int_t nDigis = fDigiMan->GetNofDigis(system);
66 UInt_t nNoise = 0;
67
68 for (Int_t iDigi = 0; iDigi < nDigis; iDigi++) {
69 const CbmMatch* match = fDigiMan->GetMatch(system, iDigi);
70 assert(match);
71
72 // This implementation uses only MC event number from
73 // the matched link, i.e. that with the largest weight.
74 // Can be refined later on.
75 Int_t mcEventNr = match->GetMatchedLink().GetEntry();
76
77 // Ignore digis with missing event number (noise)
78 if (mcEventNr < 0) {
79 nNoise++;
80 continue;
81 }
82
83 // Get event pointer. If event is not yet present, create it.
84 auto mapIt = eventMap.find(mcEventNr);
85 if (mapIt == eventMap.end()) {
86 eventMap[mcEventNr] = new CbmEventStore(nEvents++, kTRUE);
87 mapIt = eventMap.find(mcEventNr);
88 }
89 assert(mapIt != eventMap.end());
90 CbmEventStore* event = mapIt->second;
91
92 // Fill digi into event
93 switch (system) {
94 case ECbmModuleId::kMvd: {
95 const CbmMvdDigi* digi = fDigiMan->Get<CbmMvdDigi>(iDigi);
96 event->AddDigi<CbmMvdDigi>(digi, match);
97 break;
98 }
99 case ECbmModuleId::kSts: {
100 const CbmStsDigi* digi = fDigiMan->Get<CbmStsDigi>(iDigi);
101 event->AddDigi<CbmStsDigi>(digi, match);
102 break;
103 }
104 case ECbmModuleId::kRich: {
105 const CbmRichDigi* digi = fDigiMan->Get<CbmRichDigi>(iDigi);
106 event->AddDigi<CbmRichDigi>(digi, match);
107 break;
108 }
109 case ECbmModuleId::kMuch: {
110 const CbmMuchDigi* digi = fDigiMan->Get<CbmMuchDigi>(iDigi);
111 event->AddDigi<CbmMuchDigi>(digi, match);
112 break;
113 }
114 case ECbmModuleId::kTrd: {
115 const CbmTrdDigi* digi = fDigiMan->Get<CbmTrdDigi>(iDigi);
116 event->AddDigi<CbmTrdDigi>(digi, match);
117 break;
118 }
119 case ECbmModuleId::kTof: {
120 const CbmTofDigi* digi = fDigiMan->Get<CbmTofDigi>(iDigi);
121 event->AddDigi<CbmTofDigi>(digi, match);
122 break;
123 }
124 case ECbmModuleId::kPsd: {
125 const CbmPsdDigi* digi = fDigiMan->Get<CbmPsdDigi>(iDigi);
126 event->AddDigi<CbmPsdDigi>(digi, match);
127 break;
128 }
129 case ECbmModuleId::kFsd: {
130 const CbmFsdDigi* digi = fDigiMan->Get<CbmFsdDigi>(iDigi);
131 event->AddDigi<CbmFsdDigi>(digi, match);
132 break;
133 }
134 default: break;
135 } //? detector
136
137 } //# digis
138
139 nDigisTot += nDigis;
140 nDigisNoise += nNoise;
141 } //# detectors
142
143 // --- Fill output array
144 for (auto& entry : eventMap) {
145 assert(fEvents);
146 UInt_t nEntry = fEvents->GetEntriesFast();
147 new ((*fEvents)[nEntry]) CbmEventStore(*(entry.second));
148 delete entry.second;
149 CbmEventStore* event = (CbmEventStore*) fEvents->At(nEntry);
150 assert(event);
151 LOG(debug) << event->ToString();
152 }
153 LOG(debug) << "Created " << fEvents->GetEntriesFast() << " events";
154
155 timer.Stop();
156
157 // --- Execution log
158 std::cout << std::endl;
159 LOG(info) << "+ " << setw(15) << GetName() << ": Time-slice " << setw(3) << right << fNofEntries
160 << ", events: " << setw(6) << nEvents << ", digis: " << nDigisTot << ", noise: " << nDigisNoise
161 << ". Exec time " << fixed << setprecision(6) << timer.RealTime() << " s.";
162
163 fNofEntries++;
164}
165// ===========================================================================
166
167
168// ===== Task initialisation =============================================
170{
171
172 // --- Get FairRootManager instance
173 FairRootManager* ioman = FairRootManager::Instance();
174 assert(ioman);
175
176 // --- DigiManager instance
178 fDigiMan->Init();
179
180 std::cout << std::endl;
181 LOG(info) << "==================================================";
182 LOG(info) << GetName() << ": Initialising...";
183
184
185 // --- Check input data
186 for (ECbmModuleId system = ECbmModuleId::kMvd; system < ECbmModuleId::kNofSystems; ++system) {
187 if (fDigiMan->IsMatchPresent(system)) {
188 LOG(info) << GetName() << ": Found match branch for " << CbmModuleList::GetModuleNameCaps(system);
189 fSystems.push_back(system);
190 }
191 }
192 if (fSystems.empty()) {
193 LOG(fatal) << GetName() << ": No match branch found!";
194 return kFATAL;
195 }
196
197 // Register output array (CbmEvent)
198 if (ioman->GetObject("CbmEventStore")) {
199 LOG(fatal) << GetName() << ": Branch CbmEventStore already exists!";
200 return kFATAL;
201 }
202 /*
203 fEvents = new std::vector<CbmEventStore>();
204 ioman->RegisterAny("CbmEventStore", fEvents, kTRUE);
205 if ( ! fEvents ) {
206 LOG(fatal) << GetName() << ": Output branch could not be created!";
207 return kFATAL;
208 }
209 */
210 fEvents = new TClonesArray("CbmEventStore", 100);
211 ioman->Register("CbmEventStore", "Events", fEvents, kTRUE);
212
213
214 LOG(info) << "==================================================";
215 std::cout << std::endl;
216
217 return kSUCCESS;
218}
219// ===========================================================================
220
221
ClassImp(CbmConverterManager)
ECbmModuleId
Definition CbmDefs.h:39
@ kMvd
Micro-Vertex Detector.
@ kTrd
Transition Radiation Detector.
@ kTof
Time-of-flight Detector.
@ kPsd
Projectile spectator detector.
@ kSts
Silicon Tracking System.
@ kMuch
Muon detection system.
@ kFsd
Forward spectator detector.
@ kNofSystems
For loops over active systems.
@ kRich
Ring-Imaging Cherenkov Detector.
virtual void Exec(Option_t *opt)
std::vector< ECbmModuleId > fSystems
static Int_t GetNofDigis(ECbmModuleId systemId)
static Bool_t IsPresent(ECbmModuleId systemId)
Presence of a digi branch.
static Bool_t IsMatchPresent(ECbmModuleId systemId)
Presence of a digi match branch.
InitStatus Init()
Initialisation.
const Digi * Get(Int_t index) const
Get a digi object.
static CbmDigiManager * Instance()
Static instance.
const CbmMatch * GetMatch(ECbmModuleId systemId, UInt_t index) const
Get a match object.
Storable event class for CBM.
Data class for FSD digital information.
Definition CbmFsdDigi.h:36
const CbmLink & GetMatchedLink() const
Definition CbmMatch.h:41
static TString GetModuleNameCaps(ECbmModuleId moduleId)
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
Hash for CbmL1LinkKey.