CbmRoot
Loading...
Searching...
No Matches
CbmLitConverter.h
Go to the documentation of this file.
1/* Copyright (C) 2008-2020 GSI/JINR-LIT, Darmstadt/Dubna
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Andrey Lebedev [committer], Timur Ablyazimov */
4
5#ifndef CBMLITCONVERTER_H_
6#define CBMLITCONVERTER_H_
7
8#include "CbmEvent.h"
9#include "CbmGlobalTrack.h"
10#include "CbmHit.h"
12#include "CbmMuchGeoScheme.h"
13#include "CbmMuchTrack.h"
14#include "CbmMvdHit.h"
15#include "CbmPixelHit.h"
16#include "CbmStripHit.h"
17#include "CbmStsTrack.h"
18#include "CbmTofTrack.h"
19#include "CbmTrack.h"
20#include "CbmTrdAddress.h"
21#include "CbmTrdTrack.h"
22#include "CbmVertex.h"
23#include "FairRootManager.h"
24#include "FairTrackParam.h"
25#include "TClonesArray.h"
26#include "base/CbmLitEnums.h"
27#include "base/CbmLitTypes.h"
28#include "data/CbmLitFitNode.h"
29#include "data/CbmLitHit.h"
30#include "data/CbmLitPixelHit.h"
31#include "data/CbmLitStripHit.h"
32#include "data/CbmLitTofTrack.h"
33#include "data/CbmLitTrack.h"
35
36#include <cassert>
37#include <cmath>
38#include <iostream>
39#include <set>
40
41
43 public:
44 static void CbmPixelHitToCbmLitPixelHit(const CbmPixelHit* hit, Int_t index, CbmLitPixelHit* litHit)
45 {
46 assert(hit->GetType() == kTRDHIT || hit->GetType() == kMUCHPIXELHIT || hit->GetType() == kTOFHIT
47 || hit->GetType() == kMVDHIT || hit->GetType() == kSTSHIT || hit->GetType() == kPIXELHIT);
48
49 litHit->SetX(hit->GetX());
50 litHit->SetY(hit->GetY());
51 litHit->SetZ(hit->GetZ());
52
53 litHit->SetT(
54 hit
55 ->GetTime());
56
57 litHit->SetDx(hit->GetDx());
58 litHit->SetDy(hit->GetDy());
59 litHit->SetDz(hit->GetDz());
60 litHit->SetDxy(hit->GetDxy());
61
62 if (hit->GetTimeError() > 0)
63 litHit->SetDt(hit->GetTimeError());
64 else
65 litHit->SetDt(100);
66
67 litHit->SetRefId(index);
68
69 if (hit->GetType() == kTRDHIT) {
70 litHit->SetDetectorId(kLITTRD, hit->GetPlaneId());
71 }
72 else if (hit->GetType() == kMUCHPIXELHIT) {
73 litHit->SetDetectorId(kLITMUCH, (hit->GetPlaneId() - 1) / 2);
74 }
75 else if (hit->GetType() == kTOFHIT) {
76 litHit->SetDetectorId(kLITTOF, 0);
77 }
78 }
79
80 static void CbmMvdHitToCbmLitPixelHit(const CbmMvdHit* hit, Int_t index, CbmLitPixelHit* litHit)
81 {
82 litHit->SetX(hit->GetX());
83 litHit->SetY(hit->GetY());
84 litHit->SetZ(hit->GetZ());
85 litHit->SetT(hit->GetTime());
86 litHit->SetDx(hit->GetDx());
87 litHit->SetDy(hit->GetDy());
88 litHit->SetDz(hit->GetDz());
89 litHit->SetDt(hit->GetTimeError());
90 litHit->SetDxy(0.);
91 litHit->SetRefId(index);
92
93 litHit->SetDetectorId(kLITMVD, hit->GetStationNr());
94 }
95
96 static void CbmStsTrackToCbmLitTrack(const CbmStsTrack* stsTrack, CbmLitTrack* litTrack)
97 {
98 // TODO note that hits are not copied now
99
100 litTrack->SetQuality(kLITGOOD);
101 litTrack->SetChi2(stsTrack->GetChiSq());
102 litTrack->SetNDF(stsTrack->GetNDF());
103 litTrack->SetPreviousTrackId(-1);
104 CbmLitTrackParam paramFirst, paramLast;
105 //TODO remove this const typecasting
106 CbmTrackParam cbmParamFirst;
107 cbmParamFirst.Set(*stsTrack->GetParamFirst(), stsTrack->GetFirstHitTime(), stsTrack->GetFirstHitTimeError());
109 CbmTrackParam cbmParamLast;
110 cbmParamLast.Set(*stsTrack->GetParamLast(), stsTrack->GetLastHitTime(), stsTrack->GetLastHitTimeError());
112 Double_t firstTime;
113 Double_t lastTime;
114 GetStsTrackTimes(stsTrack, firstTime, lastTime);
115 paramFirst.SetTime(firstTime);
116 paramLast.SetTime(lastTime);
117 litTrack->SetParamFirst(&paramFirst);
118 litTrack->SetParamLast(&paramLast);
119 }
120
121 static void CbmTrackToCbmLitTrack(const CbmTrack* track, const HitPtrVector& lhits, CbmLitTrack* ltrack)
122 {
123 for (Int_t iHit = 0; iHit < track->GetNofHits(); iHit++) {
124 // Now we convert only pixel hits
125 if (track->GetHitType(iHit) != kPIXELHIT && track->GetHitType(iHit) != kTRDHIT
126 && track->GetHitType(iHit) != kMUCHPIXELHIT)
127 continue;
128 Int_t hitId = track->GetHitIndex(iHit);
129 ltrack->AddHit(lhits[hitId]);
130 }
131
132 ltrack->SetQuality(kLITGOOD);
133 ltrack->SetChi2(track->GetChiSq());
134 ltrack->SetNDF(track->GetNDF());
135 ltrack->SetPreviousTrackId(track->GetPreviousTrackId());
136 ltrack->SetLastStationId(track->GetFlag());
137 ltrack->SetPDG(track->GetPidHypo());
138 CbmLitTrackParam paramFirst, paramLast;
139 CbmTrackParam cbmParamFirst;
140 cbmParamFirst.Set(*track->GetParamFirst(), track->GetFirstHitTime(), track->GetFirstHitTimeError());
142 CbmTrackParam cbmParamLast;
143 cbmParamLast.Set(*track->GetParamLast(), track->GetLastHitTime(), track->GetLastHitTimeError());
145 ltrack->SetParamFirst(&paramFirst);
146 ltrack->SetParamLast(&paramLast);
147 }
148
149 static void CbmTrackArrayToCbmLitTrackArray(const TClonesArray* tracks, const HitPtrVector& lhits,
150 TrackPtrVector& ltracks)
151 {
152 Int_t nofTracks = tracks->GetEntriesFast();
153 for (Int_t iTrack = 0; iTrack < nofTracks; iTrack++) {
154 const CbmTrack* track = static_cast<const CbmTrack*>(tracks->At(iTrack));
155 CbmLitTrack* ltrack = new CbmLitTrack();
156 CbmTrackToCbmLitTrack(track, lhits, ltrack);
157 ltracks.push_back(ltrack);
158 }
159 }
160
161 static void CbmLitTrackToCbmTrack(const CbmLitTrack* litTrack, CbmTrack* track, LitSystemId systemId)
162 {
163 Double_t chiSq = 0.;
164 Int_t ndf = 0;
165 Int_t firstHit = -1;
166 Int_t lastHit = -1;
167 for (Int_t iHit = 0; iHit < litTrack->GetNofHits(); iHit++) {
168 const CbmLitHit* hit = litTrack->GetHit(iHit);
169 LitHitType type = hit->GetType();
170 LitSystemId det = hit->GetSystem();
171 if (det != systemId) continue;
172 if (firstHit < 0) firstHit = iHit;
173 lastHit = iHit;
174 if (det == kLITMUCH && type == kLITPIXELHIT) {
175 ndf += 2;
176 track->AddHit(hit->GetRefId(), kMUCHPIXELHIT);
177 }
178 else if (det == kLITTRD) {
179 ndf += 2;
180 track->AddHit(hit->GetRefId(), kTRDHIT);
181 }
182 chiSq += litTrack->GetFitNode(iHit)->GetChiSqFiltered();
183 }
184 ndf -= 5;
185 if (ndf <= 0) ndf = 1;
186
187 track->SetChiSq(chiSq);
188 track->SetNDF(ndf);
189 track->SetPreviousTrackId(litTrack->GetPreviousTrackId());
190 track->SetFlag(litTrack->GetQuality());
191 FairTrackParam parLast, parFirst;
193 &parLast);
195 &parFirst);
196 track->SetParamLast(&parLast);
197 track->SetParamFirst(&parFirst);
198 }
199
200 static void LitTrackVectorToGlobalTrackArray(CbmEvent* event, const TrackPtrVector& litTracks,
201 const TofTrackPtrVector& litTofTracks, TClonesArray* globalTracks,
202 TClonesArray* stsTracks, TClonesArray* trdTracks,
203 TClonesArray* muchTracks, TClonesArray* tofTracks)
204 {
205 // Loop over STS tracks and create GlobalTrack for each StsTrack
206 Int_t globalTrackNo = globalTracks->GetEntriesFast();
207 Int_t nofStsTracks = event ? event->GetNofData(ECbmDataType::kStsTrack) : stsTracks->GetEntriesFast();
208 for (Int_t i = 0; i < nofStsTracks; i++) {
209 Int_t iTrack = event ? event->GetIndex(ECbmDataType::kStsTrack, i) : i;
210 CbmGlobalTrack* globalTrack = new ((*globalTracks)[globalTrackNo++]) CbmGlobalTrack();
211 globalTrack->SetStsTrackIndex(iTrack);
212
213 if (event) event->AddData(ECbmDataType::kGlobalTrack, iTrack);
214 }
215
216 Int_t trdTrackNo = (trdTracks != NULL) ? trdTracks->GetEntriesFast() : 0;
217 Int_t muchTrackNo = (muchTracks != NULL) ? muchTracks->GetEntriesFast() : 0;
218 Int_t tofTrackNo = (tofTracks != NULL) ? tofTracks->GetEntriesFast() : 0;
219
220 for (size_t iTrack = 0; iTrack < litTracks.size(); iTrack++) {
221 const CbmLitTrack* litTrack = litTracks[iTrack];
222
223 if (litTrack->GetQuality() == kLITBAD) {
224 continue;
225 }
226 if (litTrack->GetNofHits() < 1) {
227 continue;
228 }
229 if (static_cast<size_t>(litTrack->GetNofHits()) != litTrack->GetFitNodes().size()) {
230 LOG(error) << "CbmLitConverter::LitTrackVectorToGlobalTrackArray: "
231 "unequal number of hits and fit nodes"
232 << std::endl
233 << litTrack->ToString();
234 continue;
235 }
236
237 CbmGlobalTrack* globalTrack = static_cast<CbmGlobalTrack*>(globalTracks->At(litTrack->GetPreviousTrackId()));
238
239 // Set last parameter of the CbmGlobal track to be last parameter of CbmLitTrack
240 FairTrackParam parLast;
242 globalTrack->SetParamLast(&parLast);
243
244 Bool_t isCreateMuchTrack = false, isCreateTrdTrack = false;
245 for (Int_t iHit = 0; iHit < litTrack->GetNofHits(); iHit++) {
246 const CbmLitHit* thisHit = litTrack->GetHit(iHit);
247 LitSystemId thisDetId = thisHit->GetSystem();
248 if (thisDetId == kLITMUCH && muchTracks != NULL) {
249 isCreateMuchTrack = true;
250 }
251 if (thisDetId == kLITTRD && trdTracks != NULL) {
252 isCreateTrdTrack = true;
253 }
254 }
255 if (isCreateTrdTrack) {
256 CbmTrdTrack* track = new ((*trdTracks)[trdTrackNo]) CbmTrdTrack();
257 CbmLitTrackToCbmTrack(litTrack, track, kLITTRD);
258 globalTrack->SetTrdTrackIndex(trdTrackNo);
259
260 if (event) event->AddData(ECbmDataType::kTrdTrack, trdTrackNo);
261
262 ++trdTrackNo;
263 }
264 if (isCreateMuchTrack) {
265 CbmMuchTrack* track = new ((*muchTracks)[muchTrackNo]) CbmMuchTrack();
266 CbmLitTrackToCbmTrack(litTrack, track, kLITMUCH);
267 globalTrack->SetMuchTrackIndex(muchTrackNo);
268
269 if (event) event->AddData(ECbmDataType::kMuchTrack, muchTrackNo);
270
271 ++muchTrackNo;
272 }
273 }
274
275 for (size_t iTrack = 0; iTrack < litTofTracks.size(); iTrack++) {
276 const CbmLitTofTrack* litTofTrack = litTofTracks[iTrack];
277 CbmTofTrack* track = new ((*tofTracks)[tofTrackNo]) CbmTofTrack();
278 if (event) event->AddData(ECbmDataType::kTofTrack, tofTrackNo);
279 tofTrackNo++;
280 Int_t globalTrackId = litTofTrack->GetTrack()->GetPreviousTrackId();
281 Int_t tofHitId = litTofTrack->GetHit()->GetRefId();
282 track->SetTofHitIndex(tofHitId);
283 track->SetTrackIndex(globalTrackId);
284 track->SetDistance(litTofTrack->GetDistance());
285 FairTrackParam par;
287 track->SetTrackParameter(&par);
288
289 CbmGlobalTrack* globalTrack = static_cast<CbmGlobalTrack*>(globalTracks->At(globalTrackId));
290 globalTrack->SetTofHitIndex(tofHitId);
291
292 if (event) event->AddData(ECbmDataType::kTofHit, tofHitId);
293 }
294 }
295
296 static void HitArrayToHitVector(CbmEvent* event, ECbmDataType hitDataType, const TClonesArray* hits,
297 HitPtrVector& litHits)
298 {
299 Int_t nofHits = event ? event->GetNofData(hitDataType) : hits->GetEntriesFast();
300 for (Int_t i = 0; i < nofHits; ++i) {
301 Int_t iHit = event ? event->GetIndex(hitDataType, i) : i;
302 CbmHit* hit = (CbmHit*) hits->At(iHit);
303 if (NULL == hit) {
304 continue;
305 }
306 CbmLitPixelHit* litHit = new CbmLitPixelHit();
307 CbmPixelHit* pixelHit = static_cast<CbmPixelHit*>(hit);
308 CbmPixelHitToCbmLitPixelHit(pixelHit, iHit, litHit);
309 litHits.push_back(litHit);
310 }
311 }
312
313 static void MvdHitArrayToHitVector(const TClonesArray* hits, HitPtrVector& litHits)
314 {
315 Int_t nofHits = hits->GetEntriesFast();
316 for (Int_t iHit = 0; iHit < nofHits; iHit++) {
317 CbmMvdHit* hit = static_cast<CbmMvdHit*>(hits->At(iHit));
318 if (NULL == hit) {
319 continue;
320 }
321 CbmLitPixelHit* litHit = new CbmLitPixelHit();
322 CbmMvdHitToCbmLitPixelHit(hit, iHit, litHit);
323 litHits.push_back(litHit);
324 }
325 }
326
327 static void StsTrackArrayToTrackVector(CbmEvent* event, const TClonesArray* tracks, TrackPtrVector& litTracks)
328 {
329 Int_t nofTracks = event ? event->GetNofData(ECbmDataType::kStsTrack) : tracks->GetEntriesFast();
330 for (Int_t i = 0; i < nofTracks; ++i) {
331 Int_t iTrack = event ? event->GetIndex(ECbmDataType::kStsTrack, i) : i;
332 CbmStsTrack* track = (CbmStsTrack*) tracks->At(iTrack);
333 if (track == NULL) {
334 continue;
335 }
336 if (track->GetParamLast()->GetQp() == 0) {
337 continue;
338 }
339 CbmLitTrack* litTrack = new CbmLitTrack;
340 CbmStsTrackToCbmLitTrack(track, litTrack);
341 litTrack->SetPreviousTrackId(iTrack);
342 litTrack->SetRefId(iTrack);
343 litTracks.push_back(litTrack);
344 }
345 }
346
347 static void GetStsTrackTimes(const CbmStsTrack* track, Double_t& firstTime, Double_t& lastTime)
348 {
349 firstTime = track->GetFirstHitTime();
350 lastTime = track->GetLastHitTime();
351 return;
352
353 static FairRootManager* ioman = 0;
354 static CbmVertex* primVertex = 0;
355 static TClonesArray* stsHits = 0;
356 static TClonesArray* mvdHits = 0;
357 static TrackPropagatorPtr propagator;
358 static TrackUpdatePtr filter;
359 static bool init = false;
360 static Int_t pdg = 211;
361
362 if (!init) {
363 init = true;
364 ioman = FairRootManager::Instance();
365
366 if (0 != ioman) {
367 primVertex = static_cast<CbmVertex*>(ioman->GetObject("PrimaryVertex"));
368 stsHits = static_cast<TClonesArray*>(ioman->GetObject("StsHit"));
369 mvdHits = static_cast<TClonesArray*>(ioman->GetObject("MvdHit"));
370 }
371
373 filter = CbmLitToolFactory::CreateTrackUpdate("kalman");
374 }
375
376 Double_t xVert = primVertex ? primVertex->GetX() : 0;
377 Double_t yVert = primVertex ? primVertex->GetY() : 0;
378 Double_t zVert = primVertex ? primVertex->GetZ() : 0;
379 CbmLitTrackParam paramFirst;
380 CbmLitTrackParam paramLast;
381 CbmTrackParam cbmParamFirst;
382 cbmParamFirst.Set(*track->GetParamFirst(), track->GetFirstHitTime(), track->GetFirstHitTimeError());
384 CbmTrackParam cbmParamLast;
385 cbmParamLast.Set(*track->GetParamLast(), track->GetLastHitTime(), track->GetLastHitTimeError());
387
388 Double_t x = paramFirst.GetX();
389 Double_t y = paramFirst.GetY();
390 Double_t z = paramFirst.GetZ();
391 //Double_t p = paramFirst.GetQp() ? TMath::Abs(1 / paramFirst.GetQp()) : 1; (VF) unused
392 CbmLitTrackParam par = paramFirst;
393 Double_t deltaTFirst = 0;
394
395 if (propagator->Propagate(&par, zVert, pdg) == kLITERROR)
396 deltaTFirst = -TMath::Sqrt(TMath::Power(x - xVert, 2) + TMath::Power(y - yVert, 2) + TMath::Power(z - zVert, 2))
398 else
399 deltaTFirst = par.GetTime() - paramFirst.GetTime();
400
401 if (deltaTFirst > 0) deltaTFirst = -deltaTFirst;
402
403 paramFirst.SetTime(paramFirst.GetTime() - deltaTFirst);
404 firstTime = paramFirst.GetTime();
405 par = paramFirst;
406 //int nofHits = track->GetNofHits();
407 int nofMvdHits = track->GetNofMvdHits();
408 int nofStsHits = track->GetNofStsHits();
409 int nofHits = nofMvdHits + nofStsHits;
410 Double_t deltaTLast = 0;
411
412 for (int i = 1; i < nofHits; ++i) {
413 //HitType hitType = track->GetHitType(i);
414 HitType hitType = (i < nofMvdHits) ? kMVDHIT : kSTSHIT;
415 //Int_t hitInd = track->GetHitIndex(i);
416 Int_t hitInd = (i < nofMvdHits) ? track->GetMvdHitIndex(i) : track->GetStsHitIndex(i - nofMvdHits);
417 CbmPixelHit* hit = static_cast<CbmPixelHit*>(kMVDHIT == hitType ? mvdHits->At(hitInd) : stsHits->At(hitInd));
418
419 if (i == nofHits - 1)
420 z = paramLast.GetZ();
421 else
422 z = hit->GetZ();
423
424 if (propagator->Propagate(&par, z, pdg) == kLITERROR) {
425 deltaTLast = TMath::Sqrt(TMath::Power(paramLast.GetX() - paramFirst.GetX(), 2)
426 + TMath::Power(paramLast.GetY() - paramFirst.GetY(), 2)
427 + TMath::Power(paramLast.GetZ() - paramFirst.GetZ(), 2))
429 break;
430 }
431
432 CbmLitPixelHit litHit;
433 CbmPixelHitToCbmLitPixelHit(hit, hitInd, &litHit);
434 litfloat chi = 0;
435 filter->Update(&par, &litHit, chi);
436 }
437
438 if (0 == deltaTLast)
439 lastTime = par.GetTime();
440 else
441 lastTime = firstTime + deltaTLast;
442 }
443};
444
445#endif /*CBMLITCONVERTER_H_*/
TClonesArray * tracks
ECbmDataType
Definition CbmDefs.h:90
HitType
Definition CbmHit.h:21
@ kSTSHIT
Definition CbmHit.h:25
@ kTOFHIT
Definition CbmHit.h:31
@ kPIXELHIT
Definition CbmHit.h:23
@ kMVDHIT
Definition CbmHit.h:26
@ kTRDHIT
Definition CbmHit.h:30
@ kMUCHPIXELHIT
Definition CbmHit.h:28
Define enumerations used in littrack.
@ kLITBAD
Definition CbmLitEnums.h:40
@ kLITGOOD
Definition CbmLitEnums.h:39
LitHitType
Definition CbmLitEnums.h:19
@ kLITPIXELHIT
Definition CbmLitEnums.h:21
@ kLITERROR
Definition CbmLitEnums.h:31
LitSystemId
Definition CbmLitEnums.h:48
@ kLITTRD
Definition CbmLitEnums.h:50
@ kLITMVD
Definition CbmLitEnums.h:52
@ kLITMUCH
Definition CbmLitEnums.h:49
@ kLITTOF
Definition CbmLitEnums.h:51
Data class for storage of fitted track parameters, transport matrix and chi-square on each detector s...
double litfloat
Definition CbmLitFloat.h:19
Base data class for hits.
Base data class for pixel hits.
Base data class for strip hits.
Data class for track parameters.
Base data class for track.
Typedefs for data structures used in littrack.
vector< CbmLitTofTrack * > TofTrackPtrVector
Definition CbmLitTypes.h:37
Data class for STS tracks.
static vector< vector< QAHit > > hits
boost::shared_ptr< CbmLitTrackPropagator > TrackPropagatorPtr
boost::shared_ptr< CbmLitTrackUpdate > TrackUpdatePtr
std::vector< CbmTofTrack * > TrackPtrVector
Definition CbmTofTypes.h:26
std::vector< CbmTofHit * > HitPtrVector
Definition CbmTofTypes.h:20
Helper class to convert unique channel ID back and forth.
Class characterising one event by a collection of links (indices) to data objects,...
Definition CbmEvent.h:34
void SetParamLast(const FairTrackParam *parLast)
void SetMuchTrackIndex(int32_t iMuch)
void SetTofHitIndex(int32_t iTofHit)
void SetTrdTrackIndex(int32_t iTrd)
void SetStsTrackIndex(int32_t iSts)
double GetTimeError() const
Definition CbmHit.h:77
double GetDz() const
Definition CbmHit.h:72
double GetTime() const
Definition CbmHit.h:76
HitType GetType() const
Definition CbmHit.h:70
virtual int32_t GetPlaneId() const
Definition CbmHit.h:99
double GetZ() const
Definition CbmHit.h:71
static void FairTrackParamToCbmLitTrackParam(const FairTrackParam *par, CbmLitTrackParam *litPar)
static void CbmLitTrackParamToFairTrackParam(const CbmLitTrackParam *litPar, FairTrackParam *par)
static void LitTrackVectorToGlobalTrackArray(CbmEvent *event, const TrackPtrVector &litTracks, const TofTrackPtrVector &litTofTracks, TClonesArray *globalTracks, TClonesArray *stsTracks, TClonesArray *trdTracks, TClonesArray *muchTracks, TClonesArray *tofTracks)
static void MvdHitArrayToHitVector(const TClonesArray *hits, HitPtrVector &litHits)
static void CbmStsTrackToCbmLitTrack(const CbmStsTrack *stsTrack, CbmLitTrack *litTrack)
static void GetStsTrackTimes(const CbmStsTrack *track, Double_t &firstTime, Double_t &lastTime)
static void CbmLitTrackToCbmTrack(const CbmLitTrack *litTrack, CbmTrack *track, LitSystemId systemId)
static void CbmMvdHitToCbmLitPixelHit(const CbmMvdHit *hit, Int_t index, CbmLitPixelHit *litHit)
static void HitArrayToHitVector(CbmEvent *event, ECbmDataType hitDataType, const TClonesArray *hits, HitPtrVector &litHits)
static void CbmTrackArrayToCbmLitTrackArray(const TClonesArray *tracks, const HitPtrVector &lhits, TrackPtrVector &ltracks)
static void StsTrackArrayToTrackVector(CbmEvent *event, const TClonesArray *tracks, TrackPtrVector &litTracks)
static void CbmPixelHitToCbmLitPixelHit(const CbmPixelHit *hit, Int_t index, CbmLitPixelHit *litHit)
static void CbmTrackToCbmLitTrack(const CbmTrack *track, const HitPtrVector &lhits, CbmLitTrack *ltrack)
const CbmLitTrackParam * GetUpdatedParam() const
litfloat GetChiSqFiltered() const
Base data class for hits.
Definition CbmLitHit.h:29
void SetRefId(Int_t refId)
Definition CbmLitHit.h:52
LitSystemId GetSystem() const
Definition CbmLitHit.h:48
LitHitType GetType() const
Definition CbmLitHit.h:43
void SetT(litfloat t)
Definition CbmLitHit.h:56
void SetZ(litfloat z)
Definition CbmLitHit.h:54
void SetDetectorId(LitSystemId sysId, Int_t station)
Definition CbmLitHit.h:58
void SetDt(litfloat dt)
Definition CbmLitHit.h:57
Int_t GetRefId() const
Definition CbmLitHit.h:42
void SetDz(litfloat dz)
Definition CbmLitHit.h:55
Base data class for pixel hits.
void SetY(litfloat y)
void SetDx(litfloat dx)
void SetDxy(litfloat dxy)
void SetDy(litfloat dy)
void SetX(litfloat x)
litfloat GetDistance() const
const CbmLitTrackParam * GetTrackParam() const
const CbmLitHit * GetHit() const
const CbmLitTrack * GetTrack() const
static TrackPropagatorPtr CreateTrackPropagator(const string &name)
Create track propagation tool by name.
static TrackUpdatePtr CreateTrackUpdate(const string &name)
Create track update tool by name.
Data class for track parameters.
litfloat GetZ() const
static litfloat fSpeedOfLight
litfloat GetX() const
litfloat GetY() const
void SetTime(litfloat t)
litfloat GetTime() const
Base data class for track.
Definition CbmLitTrack.h:34
const vector< CbmLitFitNode > & GetFitNodes() const
Definition CbmLitTrack.h:74
void SetRefId(Int_t refId)
Definition CbmLitTrack.h:92
string ToString() const
Return string representation of class.
void SetParamLast(const CbmLitTrackParam *par)
Definition CbmLitTrack.h:86
Int_t GetPreviousTrackId() const
Definition CbmLitTrack.h:66
const CbmLitFitNode * GetFitNode(Int_t index) const
Definition CbmLitTrack.h:73
void SetParamFirst(const CbmLitTrackParam *par)
Definition CbmLitTrack.h:85
void SetNDF(Int_t ndf)
Definition CbmLitTrack.h:82
const CbmLitHit * GetHit(Int_t index) const
Definition CbmLitTrack.h:71
const CbmLitTrackParam * GetParamLast() const
Definition CbmLitTrack.h:69
void SetChi2(litfloat chi2)
Definition CbmLitTrack.h:81
void SetPDG(Int_t pdg)
Definition CbmLitTrack.h:84
void SetPreviousTrackId(Int_t id)
Definition CbmLitTrack.h:83
Int_t GetNofHits() const
Definition CbmLitTrack.h:62
void AddHit(const CbmLitHit *hit)
Add hit to track. No additional memory is allocated for hit.
Definition CbmLitTrack.h:98
LitTrackQa GetQuality() const
Definition CbmLitTrack.h:63
void SetQuality(LitTrackQa quality)
Definition CbmLitTrack.h:80
void SetLastStationId(Int_t lastPlaneId)
Definition CbmLitTrack.h:89
virtual int32_t GetStationNr() const
Definition CbmMvdHit.h:61
double GetDy() const
Definition CbmPixelHit.h:76
double GetDx() const
Definition CbmPixelHit.h:75
double GetY() const
Definition CbmPixelHit.h:74
double GetX() const
Definition CbmPixelHit.h:73
double GetDxy() const
Definition CbmPixelHit.h:77
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
void SetTrackIndex(int32_t trackIndex)
Definition CbmTofTrack.h:82
void SetTofHitIndex(int32_t tofHitIndex)
Definition CbmTofTrack.h:85
void SetTrackParameter(const FairTrackParam *par)
Definition CbmTofTrack.h:93
void SetDistance(double distance)
void Set(const FairTrackParam &ftp, double time=0., double timeError=0.)
double GetFirstHitTime() const
Definition CbmTrack.h:73
double GetLastHitTime() const
Definition CbmTrack.h:75
double GetFirstHitTimeError() const
Definition CbmTrack.h:74
int32_t GetPidHypo() const
Definition CbmTrack.h:61
void AddHit(int32_t index, HitType type)
Definition CbmTrack.cxx:97
int32_t GetFlag() const
Definition CbmTrack.h:62
const FairTrackParam * GetParamLast() const
Definition CbmTrack.h:69
double GetLastHitTimeError() const
Definition CbmTrack.h:76
int32_t GetNDF() const
Definition CbmTrack.h:64
void SetFlag(int32_t flag)
Definition CbmTrack.h:80
virtual int32_t GetNofHits() const
Definition CbmTrack.h:58
void SetParamFirst(const FairTrackParam *par)
Definition CbmTrack.h:86
int32_t GetPreviousTrackId() const
Definition CbmTrack.h:67
const FairTrackParam * GetParamFirst() const
Definition CbmTrack.h:68
void SetChiSq(double chiSq)
Definition CbmTrack.h:81
int32_t GetHitIndex(int32_t iHit) const
Definition CbmTrack.h:59
HitType GetHitType(int32_t iHit) const
Definition CbmTrack.h:60
void SetPreviousTrackId(int32_t previousTrackId)
Definition CbmTrack.h:85
void SetNDF(int32_t ndf)
Definition CbmTrack.h:82
void SetParamLast(const FairTrackParam *par)
Definition CbmTrack.h:87
double GetChiSq() const
Definition CbmTrack.h:63
double GetZ() const
Definition CbmVertex.h:69
double GetY() const
Definition CbmVertex.h:68
double GetX() const
Definition CbmVertex.h:67