CbmRoot
Loading...
Searching...
No Matches
CbmLitFindGlobalTracksIdeal.cxx
Go to the documentation of this file.
1/* Copyright (C) 2009-2020 GSI/JINR-LIT, Darmstadt/Dubna
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Andrey Lebedev [committer] */
4
11
12#include "CbmGlobalTrack.h"
13#include "CbmMCTrack.h"
14#include "CbmMatch.h"
15#include "CbmTofHit.h"
16#include "CbmTofPoint.h"
17#include "CbmTrackMatchNew.h"
18#include "FairRootManager.h"
19#include "TClonesArray.h"
20
21#include <iostream>
22
24 : fDet()
25 ,
26
27 fMCTracks(NULL)
28 , fStsMatches(NULL)
29 , fMuchMatches(NULL)
30 , fTrdMatches(NULL)
31 , fTofMCPoints(NULL)
32 , fTofHits(NULL)
33 , fTofHitsMatches(NULL)
34 , fGlobalTracks(NULL)
35 ,
36
37 fMcStsMap()
38 , fMcTrdMap()
39 , fMcMuchMap()
40 , fMcTofMap()
41 ,
42
43 fEventNo(0)
44{
45}
46
48
50{
52 std::cout << fDet.ToString();
53
55
56 return kSUCCESS;
57}
58
60{
63 }
66 }
69 }
71 FillMapTof();
72 }
73
75
76 std::cout << "Event: " << fEventNo++ << std::endl;
77}
78
80
82{
83 FairRootManager* ioman = FairRootManager::Instance();
84 if (NULL == ioman) {
85 LOG(fatal) << GetName() << "::Init: CbmRootManager is not instantiated";
86 }
87
88 fMCTracks = (TClonesArray*) ioman->GetObject("MCTrack");
89 if (NULL == fMCTracks) {
90 LOG(fatal) << GetName() << "::Init: No MCTrack array!";
91 }
92
93 //STS data
95 fStsMatches = (TClonesArray*) ioman->GetObject("StsTrackMatch");
96 if (NULL == fStsMatches) {
97 LOG(fatal) << GetName() << "::Init: No StsTrackMatch array!";
98 }
99 }
100
101 //MUCH data
103 fMuchMatches = (TClonesArray*) ioman->GetObject("MuchTrackMatch");
104 if (NULL == fMuchMatches) {
105 LOG(fatal) << GetName() << "::Init: No MuchTrackMatch array!";
106 }
107 }
108
109 //TRD data
111 fTrdMatches = (TClonesArray*) ioman->GetObject("TrdTrackMatch");
112 if (NULL == fTrdMatches) {
113 LOG(fatal) << GetName() << "::Init: No TrdTrackMatch array!";
114 }
115 }
116
117 //TOF data
119 fTofMCPoints = (TClonesArray*) ioman->GetObject("TofPoint");
120 if (NULL == fTofMCPoints) {
121 LOG(fatal) << GetName() << "::Init: No TofPoint array!";
122 }
123 fTofHits = (TClonesArray*) ioman->GetObject("TofHit");
124 if (NULL == fTofHits) {
125 LOG(fatal) << GetName() << "::Init: No TofHit array!";
126 }
127 fTofHitsMatches = (TClonesArray*) ioman->GetObject("TofHitMatch");
128 if (NULL == fTofHitsMatches) {
129 LOG(fatal) << GetName() << "::Init: No TofHitMatch array!";
130 }
131 }
132
133 // Create and register CbmGlobalTrack array
134 fGlobalTracks = new TClonesArray("CbmGlobalTrack", 100);
135 ioman->Register("GlobalTrack", "Global", fGlobalTracks, IsOutputBranchPersistent("GlobalTrack"));
136}
137
138void CbmLitFindGlobalTracksIdeal::FillTrackMap(std::map<Int_t, Int_t>& mcMap, const TClonesArray* matches)
139{
140 mcMap.clear();
141 Int_t nofTracks = matches->GetEntriesFast();
142 for (Int_t iTrack = 0; iTrack < nofTracks; iTrack++) {
143 CbmTrackMatchNew* trackMatch = (CbmTrackMatchNew*) matches->At(iTrack);
144 if (trackMatch == NULL) {
145 continue;
146 }
147 Int_t mcId = trackMatch->GetMatchedLink().GetIndex();
148 if (mcId == -1) {
149 continue;
150 }
151 mcMap.insert(std::pair<Int_t, Int_t>(mcId, iTrack));
152 }
153}
154
156{
157 fMcTofMap.clear();
158 Int_t nofTofHits = fTofHits->GetEntriesFast();
159 for (Int_t iTofHit = 0; iTofHit < nofTofHits; iTofHit++) {
160 CbmTofHit* tofHit = (CbmTofHit*) fTofHits->At(iTofHit);
161 if (tofHit == NULL) {
162 continue;
163 }
164 CbmMatch* tofHitMatch = (CbmMatch*) fTofHitsMatches->At(iTofHit);
165 if (tofHitMatch == NULL) {
166 continue;
167 }
168 Int_t tofPointIndex = tofHitMatch->GetMatchedLink().GetIndex();
169 CbmTofPoint* tofPoint = (CbmTofPoint*) fTofMCPoints->At(tofPointIndex);
170 if (tofPoint == NULL) {
171 continue;
172 }
173 Int_t mcId = tofPoint->GetTrackID();
174 if (mcId == -1) {
175 continue;
176 }
177 fMcTofMap.insert(std::pair<Int_t, Int_t>(mcId, iTofHit));
178 }
179}
180
182{
183 fGlobalTracks->Clear();
184 Int_t nGlobalTracks = 0;
185 Int_t nofMCTracks = fMCTracks->GetEntriesFast();
186 for (Int_t iMCTrack = 0; iMCTrack < nofMCTracks; iMCTrack++) {
187 CbmMCTrack* mcTrack = (CbmMCTrack*) fMCTracks->At(iMCTrack);
188 if (mcTrack == NULL) {
189 continue;
190 }
191 Int_t stsId = -1, trdId = -1, muchId = -1, tofId = -1;
192 if (fDet.GetDet(ECbmModuleId::kSts) && (fMcStsMap.find(iMCTrack) != fMcStsMap.end())) {
193 stsId = fMcStsMap[iMCTrack];
194 }
195 if (fDet.GetDet(ECbmModuleId::kTrd) && (fMcTrdMap.find(iMCTrack) != fMcTrdMap.end())) {
196 trdId = fMcTrdMap[iMCTrack];
197 }
198 if (fDet.GetDet(ECbmModuleId::kMuch) && (fMcMuchMap.find(iMCTrack) != fMcMuchMap.end())) {
199 muchId = fMcMuchMap[iMCTrack];
200 }
201 if (fDet.GetDet(ECbmModuleId::kTof) && (fMcTofMap.find(iMCTrack) != fMcTofMap.end())) {
202 tofId = fMcTofMap[iMCTrack];
203 }
204
205 if (stsId == -1 && trdId == -1 && muchId == -1 && tofId == -1) {
206 continue;
207 }
208
209 CbmGlobalTrack* globalTrack = new ((*fGlobalTracks)[nGlobalTracks++]) CbmGlobalTrack();
210 globalTrack->SetStsTrackIndex(stsId);
211 globalTrack->SetTrdTrackIndex(trdId);
212 globalTrack->SetTofHitIndex(tofId);
213 globalTrack->SetMuchTrackIndex(muchId);
214 }
215}
216
@ kTrd
Transition Radiation Detector.
@ kTof
Time-of-flight Detector.
@ kSts
Silicon Tracking System.
@ kMuch
Muon detection system.
ClassImp(CbmLitFindGlobalTracksIdeal)
FairTask for ideal global track reconstruction.
void SetMuchTrackIndex(int32_t iMuch)
void SetTofHitIndex(int32_t iTofHit)
void SetTrdTrackIndex(int32_t iTrd)
void SetStsTrackIndex(int32_t iSts)
void DetermineSetup()
Determines detector presence using TGeoManager.
string ToString() const
Return string representation of class.
bool GetDet(ECbmModuleId detId) const
Return detector presence in setup.
FairTask for ideal global track reconstruction.
void ReadDataBranches()
Read necessary data branches from the input data files.
virtual void Exec(Option_t *opt)
Derived from FairTask.
virtual InitStatus Init()
Derived from FairTask.
virtual void Finish()
Derived from FairTask.
void CreateGlobalTracks()
Create output CbmGlobalTracks and write them to output array.
void FillMapTof()
Fill map from <MC track index> to <TOF hit index>.
void FillTrackMap(std::map< Int_t, Int_t > &mcMap, const TClonesArray *matches)
Fill map from <MC track index> to <reconstructed track index>.
const CbmLink & GetMatchedLink() const
Definition CbmMatch.h:41
Geometric intersection of a MC track with a TOFb detector.
Definition CbmTofPoint.h:44