CbmRoot
Loading...
Searching...
No Matches
CbmRecoTracks.cxx
Go to the documentation of this file.
1/* Copyright (C) 2016-2021 Laboratory of Information Technologies, Joint Institute for Nuclear Research, Dubna
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Timur Ablyazimov [committer] */
4
5/********************************************************************************
6 * Copyright (C) 2016 GSI Helmholtzzentrum fuer Schwerionenforschung GmbH *
7 * *
8 * This software is distributed under the terms of the *
9 * GNU Lesser General Public Licence version 3 (LGPL) version 3, *
10 * copied verbatim in the file "LICENSE" *
11 ********************************************************************************/
12// -------------------------------------------------------------------------
13// ----- CbmRecoTracks source file -----
14// ----- Created 12/02/16 by T. Ablyazimov -----
15// -------------------------------------------------------------------------
16#include "CbmRecoTracks.h"
17
18#include "CbmGlobalTrack.h" // for CbmGlobalTrack
19#include "CbmHit.h" // for HitType, kMUCHPIXELHIT, kRICHHIT
20#include "CbmPixelHit.h" // for CbmPixelHit
21#include "CbmRichRing.h" // for CbmRichRing
22#include "CbmStsTrack.h" // for CbmStsTrack
23#include "CbmTrack.h" // for CbmTrack
24
25#include <FairEventManager.h> // for FairEventManager
26#include <FairRootManager.h> // for FairRootManager
27#include <FairTask.h> // for FairTask, InitStatus, kERROR, kSUCCESS
28#include <FairTrackParam.h> // for FairTrackParam
29#include <Logger.h> // for LOG, Logger
30
31#include <Rtypes.h> // for ClassImp
32#include <TClonesArray.h> // for TClonesArray
33#include <TEveManager.h> // for TEveManager, gEve
34#include <TEvePathMark.h> // for TEvePathMark
35#include <TEveTrack.h> // for TEveTrackList, TEveTrack
36#include <TEveTrackPropagator.h> // for TEveTrackPropagator
37#include <TEveVector.h> // for TEveVector, TEveVectorT
38#include <TGenericClassInfo.h> // for TGenericClassInfo
39#include <TParticle.h> // for TParticle
40#include <TString.h> // for TString
41#include <TVector3.h> // for TVector3
42
43#include <stdio.h> // for sprintf
44#include <string.h> // for strcmp
45
46// -------------------------------------------------------------------------
48{
49 LOG(debug) << "FairMCTracks::Init()";
50 FairRootManager* fManager = FairRootManager::Instance();
51 fGlobalTracks = (TClonesArray*) fManager->GetObject("GlobalTrack");
52 fStsTracks = (TClonesArray*) fManager->GetObject("StsTrack");
53 fMvdHits = (TClonesArray*) fManager->GetObject("MvdHit");
54 fStsHits = (TClonesArray*) fManager->GetObject("StsHit");
55 fRichRings = (TClonesArray*) fManager->GetObject("RichRing");
56 fRichHits = (TClonesArray*) fManager->GetObject("RichHit");
57 fMuchPixelHits = (TClonesArray*) fManager->GetObject("MuchPixelHit");
58 fMuchTracks = (TClonesArray*) fManager->GetObject("MuchTrack");
59 fTrdHits = (TClonesArray*) fManager->GetObject("TrdHit");
60 fTrdTracks = (TClonesArray*) fManager->GetObject("TrdTrack");
61 fTofHits = (TClonesArray*) fManager->GetObject("TofHit");
62 fTofTracks = (TClonesArray*) fManager->GetObject("TofTrack");
63
64 if (fGlobalTracks == 0) {
65 LOG(error) << "FairMCTracks::Init() branch " << GetName()
66 << "GlobalTrack branch not found! Task will be deactivated";
67 SetActive(kFALSE);
68 }
69
70 LOG(debug1) << "FairMCTracks::Init() get track list" << fStsTracks;
71 LOG(debug1) << "FairMCTracks::Init() create propagator";
72 fEventManager = FairEventManager::Instance();
73 LOG(debug1) << "FairMCTracks::Init() get instance of FairEventManager ";
74 fEvent = "Current Event";
75 MinEnergyLimit = fEventManager->GetEvtMinEnergy();
76 MaxEnergyLimit = fEventManager->GetEvtMaxEnergy();
77 PEnergy = 0;
78 if (IsActive()) { return kSUCCESS; }
79 else {
80 return kERROR;
81 }
82}
83
84void CbmRecoTracks::HandlePixelHit(TEveTrack* eveTrack, Int_t& n, const CbmPixelHit* hit, TEveVector* pMom = 0)
85{
86 eveTrack->SetPoint(n, hit->GetX(), hit->GetY(), hit->GetZ());
87 TEveVector pos = TEveVector(hit->GetX(), hit->GetY(), hit->GetZ());
88 TEvePathMark path;
89 path.fV = pos;
90 path.fTime = 0;
91
92 if (pMom) path.fP = *pMom;
93
94 eveTrack->AddPathMark(path);
95 ++n;
96}
97
98void CbmRecoTracks::HandleTrack(TEveTrack* eveTrack, Int_t& n, const CbmTrack* recoTrack)
99{
100 Int_t nofHits = recoTrack->GetNofHits();
101
102 for (Int_t i = 0; i < nofHits; ++i) {
103 HitType hitType = recoTrack->GetHitType(i);
104 Int_t hitIndex = recoTrack->GetHitIndex(i);
105 const CbmPixelHit* pixelHit = 0;
106
107 switch (hitType) {
108 case kRICHHIT: pixelHit = static_cast<const CbmPixelHit*>(fRichHits->At(hitIndex)); break;
109
110 case kMUCHPIXELHIT: pixelHit = static_cast<const CbmPixelHit*>(fMuchPixelHits->At(hitIndex)); break;
111
112 case kTRDHIT: pixelHit = static_cast<const CbmPixelHit*>(fTrdHits->At(hitIndex)); break;
113
114 case kTOFHIT: pixelHit = static_cast<const CbmPixelHit*>(fTofHits->At(hitIndex)); break;
115 default: LOG(warn) << "Pixel type " << hitType << " not supported.";
116 }
117
118 if (0 != pixelHit) HandlePixelHit(eveTrack, n, pixelHit);
119 }
120}
121
122void CbmRecoTracks::HandleStsTrack(TEveTrack* eveTrack, Int_t& n, const CbmStsTrack* stsTrack)
123{
124 for (Int_t i = 0; i < stsTrack->GetNofMvdHits(); ++i) {
125 const CbmPixelHit* hit = static_cast<const CbmPixelHit*>(fMvdHits->At(stsTrack->GetMvdHitIndex(i)));
126
127 if (0 == n) {
128 TVector3 mom3;
129 stsTrack->GetParamFirst()->Momentum(mom3);
130 TEveVector mom = TEveVector(mom3.X(), mom3.Y(), mom3.Z());
131 HandlePixelHit(eveTrack, n, hit, &mom);
132 }
133 else
134 HandlePixelHit(eveTrack, n, hit);
135 }
136
137 for (Int_t i = 0; i < stsTrack->GetNofStsHits(); ++i) {
138 const CbmPixelHit* hit = static_cast<const CbmPixelHit*>(fStsHits->At(stsTrack->GetStsHitIndex(i)));
139
140 if (0 == n) {
141 TVector3 mom3;
142 stsTrack->GetParamFirst()->Momentum(mom3);
143 TEveVector mom = TEveVector(mom3.X(), mom3.Y(), mom3.Z());
144 HandlePixelHit(eveTrack, n, hit, &mom);
145 }
146 else
147 HandlePixelHit(eveTrack, n, hit);
148 }
149}
150
151void CbmRecoTracks::Exec(Option_t* /*option*/)
152{
153
154 if (IsActive()) {
155
156 LOG(debug1) << " CbmRecoTracks::Exec ";
157
158 Reset();
159
160 Int_t nofGlobalTracks = fGlobalTracks->GetEntriesFast();
161
162 for (Int_t i = 0; i < nofGlobalTracks; ++i) {
163 LOG(debug3) << "CbmRecoTracks::Exec " << i;
164 const CbmGlobalTrack* globalTrack = static_cast<const CbmGlobalTrack*>(fGlobalTracks->At(i));
165 Int_t stsId = globalTrack->GetStsTrackIndex();
166 Int_t richId = globalTrack->GetRichRingIndex();
167 Int_t muchId = globalTrack->GetMuchTrackIndex();
168 Int_t trdId = globalTrack->GetTrdTrackIndex();
169 Int_t tofId = globalTrack->GetTofTrackIndex();
170 Int_t tofHitId = globalTrack->GetTofHitIndex();
171
172 if (0 > stsId) continue;
173
174 const CbmStsTrack* stsTrack = static_cast<const CbmStsTrack*>(fStsTracks->At(stsId));
175
176 if (0 == stsTrack) continue;
177
178 Int_t pdg = stsTrack->GetPidHypo();
179 TParticle P;
180 P.SetPdgCode(pdg);
181 fTrList = GetTrGroup(&P);
182 TEveTrack* eveTrack = new TEveTrack(&P, pdg, fTrPr);
183 eveTrack->SetLineColor(fEventManager->Color(pdg));
184 Int_t n = 0;
185 HandleStsTrack(eveTrack, n, stsTrack);
186 //LOG(info) << "GetPidHypo: " << stsTrack->GetPidHypo();
187
188 if (-1 < richId) {
189 const CbmRichRing* r = static_cast<const CbmRichRing*>(fRichRings->At(richId));
190 const CbmPixelHit* rh = static_cast<const CbmPixelHit*>(fRichHits->At(r->GetHit(0)));
191 CbmPixelHit h(*rh);
192 h.SetX(r->GetCenterX());
193 h.SetY(r->GetCenterY());
194 HandlePixelHit(eveTrack, n, &h);
195 }
196 else if (-1 < muchId)
197 HandleTrack(eveTrack, n, static_cast<const CbmTrack*>(fMuchTracks->At(muchId)));
198
199 if (-1 < trdId) HandleTrack(eveTrack, n, static_cast<const CbmTrack*>(fTrdTracks->At(trdId)));
200
201 if (-1 < tofId) HandleTrack(eveTrack, n, static_cast<const CbmTrack*>(fTofTracks->At(tofId)));
202 else if (-1 < tofHitId)
203 HandlePixelHit(eveTrack, n, static_cast<const CbmPixelHit*>(fTofHits->At(tofHitId)));
204
205 fTrList->AddElement(eveTrack);
206 LOG(debug3) << "track added " << eveTrack->GetName();
207 }
208
209 //fEventManager->SetEvtMaxEnergy(MaxEnergyLimit);
210 //fEventManager->SetEvtMinEnergy(MinEnergyLimit);
211 gEve->Redraw3D(kFALSE);
212 }
213}
214// ----- Destructor ----------------------------------------------------
216// -------------------------------------------------------------------------
218
219// -------------------------------------------------------------------------
221// -------------------------------------------------------------------------
223{
224 for (Int_t i = 0; i < fEveTrList->GetEntriesFast(); i++) {
225 TEveTrackList* ele = (TEveTrackList*) fEveTrList->At(i);
226 gEve->RemoveElement(ele, fEventManager);
227 }
228 fEveTrList->Clear();
229}
230
231TEveTrackList* CbmRecoTracks::GetTrGroup(TParticle* P)
232{
233 size_t buf_size = 128;
234 char name_buf[buf_size];
235 snprintf(name_buf, buf_size - 1, "reco_%s", P->GetName());
236 fTrList = 0;
237 for (Int_t i = 0; i < fEveTrList->GetEntriesFast(); i++) {
238 TEveTrackList* TrListIn = (TEveTrackList*) fEveTrList->At(i);
239 if (strcmp(TrListIn->GetName(), name_buf) == 0) {
240 fTrList = TrListIn;
241 break;
242 }
243 }
244 if (fTrList == 0) {
245 fTrPr = new TEveTrackPropagator();
246 fTrList = new TEveTrackList(name_buf, fTrPr);
247 fTrList->SetMainColor(fEventManager->Color(P->GetPdgCode()));
248 fEveTrList->Add(fTrList);
249 gEve->AddElement(fTrList, fEventManager);
250 fTrList->SetRnrLine(kTRUE);
251 }
252 return fTrList;
253}
254
ClassImp(CbmConverterManager)
HitType
Definition CbmHit.h:21
@ kTOFHIT
Definition CbmHit.h:31
@ kTRDHIT
Definition CbmHit.h:30
@ kRICHHIT
Definition CbmHit.h:27
@ kMUCHPIXELHIT
Definition CbmHit.h:28
Data class for STS tracks.
int32_t GetStsTrackIndex() const
int32_t GetRichRingIndex() const
int32_t GetTofTrackIndex() const
int32_t GetTofHitIndex() const
int32_t GetMuchTrackIndex() const
int32_t GetTrdTrackIndex() const
double GetZ() const
Definition CbmHit.h:71
double GetY() const
Definition CbmPixelHit.h:74
double GetX() const
Definition CbmPixelHit.h:73
TClonesArray * fMuchPixelHits
TClonesArray * fMuchTracks
virtual ~CbmRecoTracks()
TClonesArray * fTofTracks
TEveTrackPropagator * fTrPr
virtual void Exec(Option_t *option)
virtual InitStatus Init()
virtual void Finish()
TClonesArray * fTrdHits
TEveTrackList * fTrList
TClonesArray * fMvdHits
TClonesArray * fStsTracks
FairEventManager * fEventManager
TEveTrackList * GetTrGroup(TParticle *P)
Double_t MinEnergyLimit
TClonesArray * fGlobalTracks
TClonesArray * fRichHits
TClonesArray * fRichRings
virtual void SetParContainers()
TClonesArray * fStsHits
TObjArray * fEveTrList
void HandleTrack(TEveTrack *eveTrack, Int_t &n, const CbmTrack *recoTrack)
void HandlePixelHit(TEveTrack *eveTrack, Int_t &n, const CbmPixelHit *hit, TEveVector *pMom)
Double_t MaxEnergyLimit
TClonesArray * fTrdTracks
void HandleStsTrack(TEveTrack *eveTrack, Int_t &n, const CbmStsTrack *stsTrack)
TClonesArray * fTofHits
uint32_t GetHit(int32_t i) const
Definition CbmRichRing.h:39
float GetCenterX() const
Definition CbmRichRing.h:77
float GetCenterY() const
Definition CbmRichRing.h:78
Data class with information on a STS local track.
int32_t GetNofMvdHits() const
Definition CbmStsTrack.h:87
int32_t GetMvdHitIndex(int32_t iHit) const
Definition CbmStsTrack.h:75
int32_t GetStsHitIndex(int32_t iHit) const
int32_t GetNofStsHits() const
Definition CbmStsTrack.h:93
int32_t GetPidHypo() const
Definition CbmTrack.h:61
virtual int32_t GetNofHits() const
Definition CbmTrack.h:58
const FairTrackParam * GetParamFirst() const
Definition CbmTrack.h:68
int32_t GetHitIndex(int32_t iHit) const
Definition CbmTrack.h:59
HitType GetHitType(int32_t iHit) const
Definition CbmTrack.h:60