BmnRoot
Loading...
Searching...
No Matches
BmnParticle.cxx
Go to the documentation of this file.
1
6
7#include "BmnParticle.h"
8// #include "MpdHelix.h"
9// #include "MpdKalmanFilter.h"
10// #include "MpdKalmanHit.h"
11// #include "MpdKalmanTrack.h"
12#include "BmnMotherFitterPart.h"
13#include "CbmKF.h"
14#include "CbmKFTrack.h"
15#include "CbmStsTrack.h"
16
17#include <TDatabasePDG.h>
18#include <TMatrixDSym.h>
19
20//__________________________________________________________________________
22 : TObject()
23 , fIndx(-1)
24 , fCharge(0)
25 , fMass(TDatabasePDG::Instance()->GetParticle("pi+")->Mass())
26 , fChi2(0.0)
27 , fChi2ver(-1.0)
28 , fFlag(0)
29 , fPoint00(kTRUE)
30 , fZ(0.0)
31 , fNDF(2)
32{
34
35 // fKF = MpdKalmanFilter::Instance();
36 fKF = CbmKF::Instance();
37 fMeas.ResizeTo(5, 1);
38 fq.ResizeTo(3, 1);
39 fx.ResizeTo(3, 1);
40 fJ.ResizeTo(3, 3);
41 fJinv.ResizeTo(3, 3);
42 fD.ResizeTo(3, 3);
43 fE.ResizeTo(3, 3);
44 // fQ.ResizeTo(3,3);
45 fA.ResizeTo(5, 3);
46 fB.ResizeTo(5, 3);
47 fC.ResizeTo(3, 3);
48 fG.ResizeTo(5, 5);
49 fW.ResizeTo(3, 3);
50
51 fDaughtersInds.clear();
52 /* AZ-020725
53 Double_t pos[3] = {0.0}, magf[3] = {0.0};
54 fKF->GetField(pos, magf);
55 fieldConst = 0.3 * 0.01 * magf[2] / 10;
56 */
57 fXY0[0] = fXY0[1] = 0.0;
58 fG(0, 0) = fG(1, 1) = fG(2, 2) = fG(3, 3) = fG(4, 4) = 0.001; // AZ-139725
59}
60
61//__________________________________________________________________________
62
63BmnParticle::BmnParticle(CbmKFTrack& track, Int_t indx, Double_t mass, Double_t* orig)
64 // : TObject(track), fIndx(indx), fCharge(0),
65 : TObject()
66 , fIndx(indx)
67 , fCharge(0)
68 , fMass(mass)
69 , fChi2(0.0)
70 , fChi2ver(-1.0)
71 , fFlag(0)
72 , fPoint00(kTRUE)
73 , fZ(0.0)
74 , fNDF(2)
75{
77
78 fDaughtersInds.clear();
79 AddDaughter(indx);
80 Track2Part(track, kTRUE, orig);
81}
82
83//__________________________________________________________________________
84
85BmnParticle::BmnParticle(CbmStsTrack& track, Int_t indx, Double_t mass, Double_t* orig)
86 : TObject(track)
87 , fIndx(indx)
88 , fCharge(0)
89 , fMass(mass)
90 , fChi2(0.0)
91 , fChi2ver(-1.0)
92 , fFlag(0)
93 , fPoint00(kTRUE)
94 , fNDF(2)
95{
97
98 fDaughtersInds.clear();
99 AddDaughter(indx);
100 CbmKFTrack kftr(track);
101 Track2Part(kftr, kTRUE, orig);
102}
103
104//__________________________________________________________________________
106 : TObject(part)
107 , fIndx(part.fIndx)
108 , fPdg(part.fPdg)
109 , fCharge(part.fCharge)
110 , fMass(part.fMass)
111 , fieldConst(part.fieldConst)
112 , fMeas(part.fMeas)
113 , fq(part.fq)
114 , fx(part.fx)
115 , fJ(part.fJ)
116 , fJinv(part.fJinv)
117 , fD(part.fD)
118 , fE(part.fE)
119 , fA(part.fA)
120 , fB(part.fB)
121 , fC(part.fC)
122 , fG(part.fG)
123 , fW(part.fW)
124 , fChi2(part.fChi2)
125 , fChi2ver(part.fChi2ver)
126 , fFlag(part.fFlag)
127 , fPoint00(part.fPoint00)
128 , fZ(part.fZ)
129 , fNDF(part.fNDF)
130{
132
133 fXY0[0] = part.fXY0[0];
134 fXY0[1] = part.fXY0[1];
135 fDaughtersInds = part.fDaughtersInds;
136
137 fKF = part.fKF;
138}
139
140//__________________________________________________________________________
142{
144
145 // check assignement to self
146 if (this == &part)
147 return *this;
148
149 // base class assignement
150 TObject::operator=(part);
151
152 fIndx = part.fIndx;
153 fPdg = part.fPdg;
154 fCharge = part.fCharge;
155 fMass = part.fMass;
156 fieldConst = part.fieldConst;
157 fMeas = part.fMeas;
158 fq = part.fq;
159 fx = part.fx;
160 fJ = part.fJ;
161 fJinv = part.fJinv;
162 fD = part.fD;
163 fE = part.fE;
164 fA = part.fA;
165 fB = part.fB;
166 fC = part.fC;
167 fG = part.fG;
168 fW = part.fW;
169 fChi2 = part.fChi2;
170 fChi2ver = part.fChi2ver;
171 fFlag = part.fFlag;
172 fKF = part.fKF;
173 fNDF = part.fNDF;
174 fZ = part.fZ;
175
176 fXY0[0] = part.fXY0[0];
177 fXY0[1] = part.fXY0[1];
178 fDaughtersInds = part.fDaughtersInds;
179
180 return *this;
181}
182
183//__________________________________________________________________________
188
189/*
190//__________________________________________________________________________
191Int_t MpdParticle::Compare(const TObject* part) const
192{
194
195 MpdKalmanTrack *trackKF = (MpdKalmanTrack*) track;
196 Double_t pt = 1. / TMath::Abs(trackKF->GetParam(4));
197 Double_t ptthis = 1. / TMath::Abs((*fParam)(4,0));
198 if (ptthis < pt) return 1;
199 else if (ptthis > pt) return -1;
200 return 0;
201}
202*/
203//__________________________________________________________________________
204void BmnParticle::SetMass(Double_t mass)
205{
207
208 fMass = (mass > -1.0) ? mass : TDatabasePDG::Instance()->GetParticle(fPdg)->Mass();
209}
210
211//__________________________________________________________________________
212void BmnParticle::Track2Part(CbmKFTrack& track, Bool_t setWeight, Double_t* orig)
213{
215
216 // fKF = MpdKalmanFilter::Instance();
217 fKF = CbmKF::Instance();
218 fCharge = TMath::Sign(1.5, track.GetTrack()[4]); // Qp
219 fMeas.ResizeTo(5, 1);
220 fq.ResizeTo(3, 1);
221 fx.ResizeTo(3, 1);
222 fJ.ResizeTo(3, 3);
223 fJinv.ResizeTo(3, 3);
224 fD.ResizeTo(3, 3);
225 fE.ResizeTo(3, 3);
226 // fQ.ResizeTo(3,3);
227 fA.ResizeTo(5, 3);
228 fB.ResizeTo(5, 3);
229 fC.ResizeTo(3, 3);
230 fG.ResizeTo(5, 5);
231 fW.ResizeTo(3, 3);
232
233 /* AZ-020725
234 Double_t pos[3] = {0.0}, magf[3] = {0.0};
235 fKF->GetField(pos, magf);
236 fieldConst = 0.3 * 0.01 * magf[2] / 10;
237 */
238 fXY0[0] = fXY0[1] = 0.0;
239 // AZ if (orig) fPoint00 = kFALSE;
240 fXY0[0] = track.GetTrack()[5]; // AZ-031225 - Z-ccordinate of the first hit
241
242 Double_t vert[3] = {0};
243 if (orig == NULL)
244 orig = vert;
245
246 CbmKFTrack kftr(track);
247 // FindPca (kftr, orig);
248 kftr.Extrapolate(0.0); // extrapolate to Z=0 - AZ-070725
249 // kftr.Propagate(0.0); // extrapolate to Z=0
250 FairTrackParam param;
251 kftr.GetTrackParam(param);
252 fMeas(0, 0) = param.GetX(); // X
253 fMeas(1, 0) = param.GetY(); // Y
254 fMeas(2, 0) = param.GetTx(); // Tx
255 fMeas(3, 0) = param.GetTy(); // Ty
256 fMeas(4, 0) = param.GetQp(); // 1/p
257 fZ = param.GetZ(); // Z
258
259 if (!setWeight)
260 return; // skip all the rest
261
262 fq(0, 0) = fMeas(2, 0);
263 fq(1, 0) = fMeas(3, 0);
264 fq(2, 0) = fMeas(4, 0);
265
266 TMatrixDSym covar(5);
267
268 for (int i = 0; i < 5; ++i) {
269 covar(i, i) = param.GetCovariance(i, i);
270
271 for (int j = i + 1; j < 5; ++j)
272 covar(i, j) = covar(j, i) = param.GetCovariance(i, j);
273 }
274 covar.Invert();
275 TMatrixD weight(5, 5);
276 weight.SetMatrixArray(covar.GetMatrixArray());
277
278 SetG(weight);
279}
280
281//__________________________________________________________________________
282
283// void BmnParticle::FindPca(CbmKFTrackInterface& tr, Double_t* vert)
284void BmnParticle::FindPca(CbmKFTrack& tr, Double_t* vert)
285{
286 // Find track 2-D PCA w.r.t. vert using parabolic approximation
287
288 // Make 2 iterations
289 Double_t zpca = vert[2], dz = 5.0;
290
291 for (int iter = 0; iter < 2; ++iter) {
292 TVector3 pos[3];
293 if (iter > 0)
294 dz = 2.0;
295
296 for (int i = -1; i < 2; ++i) {
297 Double_t z = zpca + dz * i;
298 tr.Extrapolate(z);
299 Double_t* pars = tr.GetTrack();
300 Double_t dx = pars[0] - vert[0], dy = pars[1] - vert[1]; //, dz = pars[5] - vert[2];
301 // Double_t dist = TMath::Sqrt(dx * dx + dy * dy + dz * dz);
302 Double_t dist = TMath::Sqrt(dx * dx + dy * dy);
303 pos[i + 1].SetX(dist);
304 pos[i + 1].SetZ(z);
305 }
306 TVector3 cba = Parabola(pos[0], pos[1], pos[2]);
307 zpca = -cba[1] / cba[2] / 2;
308 }
309 tr.Extrapolate(zpca);
310}
311
312//__________________________________________________________________________
313
314TVector3 BmnParticle::Parabola(TVector3& pos0, TVector3& pos1, TVector3& pos2)
315{
316 // Get parabolic track approximation from 3 points (X vs Z)
317 // y = a*x^2 + bx + c
318
319 Double_t x[3] = {pos0[2], pos1[2], pos2[2]};
320 Double_t y[3] = {pos0[0], pos1[0], pos2[0]};
321
322 Double_t denom = (x[0] - x[1]) * (x[0] - x[2]) * (x[1] - x[2]);
323 Double_t dy10 = y[1] - y[0];
324 Double_t dy02 = y[0] - y[2];
325 Double_t dy21 = y[2] - y[1];
326 Double_t a = x[2] * dy10 + x[1] * dy02 + x[0] * dy21;
327 a /= denom;
328 Double_t b = -x[0] * x[0] * dy21 - x[2] * x[2] * dy10 - x[1] * x[1] * dy02;
329 b /= denom;
330 Double_t c = x[1] * x[1] * (x[2] * y[0] - x[0] * y[2]) + x[1] * (x[0] * x[0] * y[2] - x[2] * x[2] * y[0])
331 + x[0] * x[2] * (x[2] - x[0]) * y[1];
332 c /= denom;
333
334 return TVector3(c, b, a);
335 // 𝐴=𝑥3(𝑦2−𝑦1)+𝑥2(𝑦1−𝑦3)+𝑥1(𝑦3−𝑦2) / (𝑥1−𝑥2)(𝑥1−𝑥3)(𝑥2−𝑥3)
336 // 𝐵=𝑥21(𝑦2−𝑦3)+𝑥23(𝑦1−𝑦2)+𝑥22(𝑦3−𝑦1) / (𝑥1−𝑥2)(𝑥1−𝑥3)(𝑥2−𝑥3)
337 // 𝐶=𝑥22(𝑥3𝑦1−𝑥1𝑦3)+𝑥2(𝑥21𝑦3−𝑥23𝑦1)+𝑥1𝑥3(𝑥3−𝑥1)𝑦2 / (𝑥1−𝑥2)(𝑥1−𝑥3)(𝑥2−𝑥3)
338}
339
340//__________________________________________________________________________
341
342// void BmnParticle::WeightAtDca(CbmStsTrack &tr, Double_t *vert)
343void BmnParticle::WeightAtDca(CbmKFTrack& tr, Double_t* vert)
344{
345 // Obtain BmnParticle weight at DCA from CbmStsTrack weight
346
347 /*
348 FairTrackParam par0;
349 tr.GetTrackParam (par0);
350 Double_t pars[5] = { par0.GetX(), par0.GetY(), par0.GetTx(), par0.GetTy(), par0.GetQp() };
351
352 TMatrixDSym covar;
353 //covar.SetMatrixArray (par0.GetCovMatrix());
354
355 //TMatrixDSym *covar = tr.Weight2Cov(); // get covariance matrix
356 TMatrixD g = *tr.GetWeight(); // track weight matrix
357 TMatrixD meas0 = GetMeas();
358 //TMatrixD param0 = *tr.GetParamNew();
359 TMatrixD param0 (5, 1, pars);
360 TMatrixD jacob (5, 5);
361 // Double_t vert[3] = {0.0};
362
363 Double_t dPar;
364 // Loop over parameters to find change of the propagated vs initial ones
365 Bool_t ok = kTRUE;
366
367 for (Int_t i = 0; i < 5; ++i) {
368 dPar = TMath::Sqrt (covar(i, i));
369 if (i < 4)
370 dPar = TMath::Min(dPar, 0.1);
371 else
372 dPar = TMath::Min(dPar, 0.1 * TMath::Abs(param0(4, 0)));
373
374 tr.SetParam(param0);
375 if (i == 4) dPar *= TMath::Sign(1., -param0(4, 0)); // 1/p
376 // else if (i == 2) dPar *= sign;
377 else if (i == 3)
378 dPar *= TMath::Sign(1., -param0(3, 0)); // dip-angle
379 tr.SetParam(i, param0(i, 0) + dPar);
380 ok = fKF->FindPca(&tr, vert);
381 tr.SetParam(*tr.GetParamNew());
382 BmnParticle tmpPart;
383 tmpPart.Track2Part(tr, kFALSE, vert);
384
385 Double_t meas = 0.0;
386 for (Int_t j = 0; j < 5; ++j) {
387 if (j == 2)
388 meas = fKF->Proxim(meas0(j, 0), tmpPart.GetMeas(j));
389 else
390 meas = tmpPart.GetMeas(j);
391 jacob(j, i) = (meas - meas0(j, 0)) / dPar;
392 // cout << i << " " << j << " " << dPar << " " << track->GetParamNew(j) << " " << paramNew0(j,0)
393 // << " " << track->GetPos() << " " << track->GetPosNew() << " " << jacob(j,i) << endl;
394 }
395 }
396
397 // jacob.Print();
398 jacob.Invert();
399
400 TMatrixD tmp (g, TMatrixD::kMult, jacob); // WD
401 TMatrixD weight1 (jacob, TMatrixD::kTransposeMult, tmp); // DtWD
402
403 SetG (weight1);
404 */
405}
406
407//__________________________________________________________________________
408Double_t BmnParticle::BuildMother(vector<BmnParticle*>& vDaught)
409{
412
413 return BmnMotherFitterPart::Instance()->BuildMother(this, vDaught);
414}
415
416//__________________________________________________________________________
417
418Double_t BmnParticle::BuildMother(vector<CbmStsTrack*>& vTracks, vector<BmnParticle*>& vDaught)
419{
424
425 // AZ-130725 return BmnMotherFitterPart::Instance()->BuildMother(this, vTracks, vDaught);
426 return kTRUE;
427}
428
429//__________________________________________________________________________
430
432{
434
435 Double_t p0 = Momentum();
436 Double_t tx = fq(0, 0);
437 Double_t ty = fq(1, 0);
438 Double_t denom = TMath::Sqrt(tx * tx + ty * ty + 1);
439 Double_t denom2 = denom * denom;
440 Double_t pz = p0 / denom;
441 fJ = 0.0;
442
443 if (fCharge) {
444 // Charged particle
445 // d/d(tx)
446 fJ(0, 0) = pz - pz * tx * tx / denom2;
447 fJ(1, 0) = -pz * tx * ty / denom2;
448 fJ(2, 0) = -pz * tx / denom2;
449 // d/d(ty)
450 fJ(0, 1) = -pz * tx * ty / denom2;
451 fJ(1, 1) = pz - pz * ty * ty / denom2;
452 fJ(2, 1) = -pz * ty / denom2;
453 // d/d(q/p) - CHECK ME !!!
454 // fJ(0, 2) = -tx * pz * p0;
455 /* AZ-050825 fJ(0, 2) = -tx * pz * fq(2, 0);
456 fJ(1, 2) = -ty * pz * fq(2, 0);
457 fJ(2, 2) = -pz * fq(2, 0); */
458 fJ(0, 2) = -tx * pz / fq(2, 0);
459 fJ(1, 2) = -ty * pz / fq(2, 0);
460 fJ(2, 2) = -pz / fq(2, 0);
461 } else {
462 // neutral
463 // d/d(tx)
464 fJ(0, 0) = pz - pz * tx * tx / denom2;
465 fJ(1, 0) = -pz * tx * ty / denom2;
466 fJ(2, 0) = -pz * tx / denom2;
467 // d/d(ty)
468 fJ(0, 1) = -pz * tx * ty / denom2;
469 fJ(1, 1) = pz - pz * ty * ty / denom2;
470 fJ(2, 1) = -pz * ty / denom2;
471 // d/dp
472 fJ(0, 2) = tx / denom;
473 fJ(1, 2) = ty / denom;
474 fJ(2, 2) = 1 / denom;
475 }
476
477 // AZ-120725 TMatrixD tmp(TMatrixD::kTransposed, fJ);
478 // fJ = tmp;
479}
480
481//__________________________________________________________________________
482void BmnParticle::FillJinv(TVector3& mom3)
483{
485
486 Double_t px0 = mom3.X();
487 Double_t py0 = mom3.Y();
488 Double_t pz0 = mom3.Z();
489 Double_t p0 = mom3.Mag();
490 Double_t p03 = p0 * p0 * p0;
491 Double_t pz02 = pz0 * pz0;
492 fJinv = 0.0;
493
494 if (fCharge) {
495 // Charged mother - CHECK ME !!!
496 // d/d(px)
497 fJinv(0, 0) = 1 / pz0;
498 fJinv(1, 0) = 0;
499 fJinv(2, 0) = -fCharge * px0 / p03;
500 // d/d(py)
501 fJinv(0, 1) = 0;
502 fJinv(1, 1) = 1 / pz0;
503 fJinv(2, 1) = -fCharge * py0 / p03;
504 // d/d(pz)
505 fJinv(0, 2) = -px0 / pz02;
506 fJinv(1, 2) = -py0 / pz02;
507 fJinv(2, 2) = -fCharge * pz0 / p03;
508 } else {
509 // neutral mother
510 // d/d(px)
511 fJinv(0, 0) = 1 / pz0;
512 fJinv(1, 0) = 0;
513 fJinv(2, 0) = px0 / p0;
514 // d/d(py)
515 fJinv(0, 1) = 0;
516 fJinv(1, 1) = 1 / pz0;
517 fJinv(2, 1) = py0 / p0;
518 // d/d(pz)
519 fJinv(0, 2) = -px0 / pz02;
520 fJinv(1, 2) = -py0 / pz02;
521 fJinv(2, 2) = pz0 / p0;
522 }
523
524 // AZ-120725 TMatrixD tmp(TMatrixD::kTransposed, fJinv);
525 // fJinv = tmp;
526}
527
528//__________________________________________________________________________
530{
532
533 fChi2ver = BmnMotherFitterPart::Instance()->Chi2Vertex(this, vtx);
534 return fChi2ver;
535}
536
537//__________________________________________________________________________
538Double_t BmnParticle::Energy() const
539{
541
542 Double_t mom = Momentum();
543 return TMath::Sqrt(mom * mom + fMass * fMass);
544}
545
546//__________________________________________________________________________
547Double_t BmnParticle::Rapidity() const
548{
550
551 TVector3 mom3 = Momentum3();
552 Double_t mom = Momentum();
553 Double_t e = TMath::Sqrt(mom * mom + fMass * fMass);
554 Double_t pz = mom3.Z();
555 Double_t y = 0.5 * TMath::Log((e + pz) / (e - pz));
556 return y;
557}
558
559//__________________________________________________________________________
561{
563}
564
565//__________________________________________________________________________
566
568{
570
571 TMatrixD meas(5, 1);
572
573 meas(0, 0) = Getx()(0, 0);
574 meas(1, 0) = Getx()(1, 0);
575 meas(2, 0) = Getq()(0, 0);
576 meas(3, 0) = Getq()(1, 0);
577 meas(4, 0) = Getq()(2, 0);
578 SetMeas(meas);
579 fZ = Getx()(2, 0);
580 //* AZ-130725
581 CbmKFTrack tr = GetKFTrack();
582 tr.Propagate(0.0);
583
584 FairTrackParam par;
585 tr.GetTrackParam(par);
586
587 meas(0, 0) = par.GetX();
588 meas(1, 0) = par.GetY();
589 meas(2, 0) = par.GetTx();
590 meas(3, 0) = par.GetTy();
591 Double_t qp = (GetCharge()) ? par.GetQp() : Getq()(2, 0);
592 meas(4, 0) = qp;
593 SetMeas(meas);
594 fZ = par.GetZ();
595 //*/
596
597 // AZ-090825 - set cov. matrix
598 TMatrixDSym covar(5);
599
600 for (int i = 0; i < 5; ++i) {
601 covar(i, i) = par.GetCovariance(i, i);
602
603 for (int j = i + 1; j < 5; ++j)
604 covar(i, j) = covar(j, i) = par.GetCovariance(i, j);
605 }
606 covar.Invert();
607 TMatrixD weight(5, 5);
608 weight.SetMatrixArray(covar.GetMatrixArray());
609
610 SetG(weight);
611}
612
613//__________________________________________________________________________
614
616{
618
619 TMatrixFSym covMatr(5);
620
621 for (int i = 0; i < 5; ++i) {
622 covMatr(i, i) = fG(i, i);
623
624 for (int j = i + 1; j < 5; ++j) {
625 covMatr(i, j) = covMatr(j, i) = fG(i, j);
626 }
627 }
628 covMatr.Invert();
629 Double_t qp = GetMeas(4);
630 if (GetCharge() == 0)
631 qp = 0.0; // !!! AZ-090725 - this is for CbmKFTrack !!!
632 FairTrackParam par(GetMeas(0), GetMeas(1), fZ, GetMeas(2), GetMeas(3), qp, covMatr);
633 CbmKFTrack tr(par);
634 return tr;
635}
vector< Double_t > dist(vector< Double_t > qp, Double_t mu)
Definition BmnMath.cxx:869
ClassImp(BmnParticle)
int i
Definition P4_F32vec4.h:22
Double_t Chi2Vertex(BmnParticle *part, const CbmVertex *vtx)
compute Chi2 w.r.t. vertex
static BmnMotherFitterPart * Instance()
get singleton instance
Double_t BuildMother(BmnParticle *mother, vector< BmnParticle * > &vDaught)
TMatrixD & Getx()
Double_t Energy() const
void SetG(TMatrixD &matr)
Double_t BuildMother(vector< BmnParticle * > &vDaught)
CbmKFTrack GetKFTrack() const
Double_t Rapidity() const
BmnParticle()
Default ctor.
const Double_t Chi2Vertex()
return Chi2 w.r.t. vertex
Definition BmnParticle.h:95
Int_t GetCharge() const
Definition BmnParticle.h:64
void Track2Part(CbmKFTrack &track, Bool_t setWeight, Double_t *orig)
TMatrixD & GetMeas()
TVector3 Momentum3() const
Definition BmnParticle.h:85
virtual ~BmnParticle()
void SetMass(Double_t mass=-2.0)
TMatrixD & Getq()
void SetMeas(TMatrixD &matr)
void FillJinv(TVector3 &mom3)
void AddDaughter(Int_t indx)
void ParamsAtDca()
Double_t Momentum() const
Definition BmnParticle.h:83
BmnParticle & operator=(const BmnParticle &part)
assignment operator
Int_t Propagate(Double_t z_out, Double_t QP0, Bool_t line=false)
Int_t Extrapolate(Double_t z, Double_t *QP0=0, Bool_t line=false)
Access to i-th hit.
void GetTrackParam(FairTrackParam &track)
Double_t * GetTrack()
Is it electron.
Definition CbmKFTrack.h:51
static CbmKF * Instance()
Definition CbmKF.h:35