PandaRoot
KFParticle.h
Go to the documentation of this file.
1 //****************************************************************************
2 //* This file is part of PandaRoot. *
3 //* *
4 //* PandaRoot is distributed under the terms of the *
5 //* GNU General Public License (GPL) version 3, *
6 //* copied verbatim in the file "LICENSE". *
7 //* *
8 //* Copyright (C) 2006 - 2024 FAIR GmbH and copyright holders of PandaRoot *
9 //* The copyright holders are listed in the file "COPYRIGHTHOLDERS". *
10 //* The authors are listed in the file "AUTHORS". *
11 //****************************************************************************
12 
13 //---------------------------------------------------------------------------------
14 // The KFParticle class
15 // .
16 // @author S.Gorbunov, I.Kisel
17 // @version 1.0
18 // @since 13.05.07
19 //
20 // Class to reconstruct and store the decayed particle parameters.
21 // The method is described in CBM-SOFT note 2007-003,
22 // ``Reconstruction of decayed particles based on the Kalman filter'',
23 // http://www.gsi.de/documents/DOC-2007-May-14-1.pdf
24 //
25 // This class is ALICE interface to general mathematics in KFParticleBase
26 //
27 // -= Copyright &copy ALICE HLT and CBM L1 Groups =-
28 //_________________________________________________________________________________
29 
30 //#define NonhomogeneousField
31 #define HomogeneousField
32 
33 #ifndef ALIKFPARTICLE_H
34 #define ALIKFPARTICLE_H
35 
36 #include "KFParticleBase.h"
37 // #include "TMath.h"
38 #include <cmath>
39 
40 #ifdef HomogeneousField
41 class KFPTrack;
42 class KFPVertex;
43 #endif
44 
45 #ifdef NonhomogeneousField
46 class CbmKFTrackInterface;
47 class CbmKFVertexInterface;
48 #endif
49 
50 class KFParticle : public KFParticleBase {
51 
52  public:
53  //*
54  //* INITIALIZATION
55  //*
56 
57  //* Set magnetic field for all particles
58 #ifdef HomogeneousField
59  static void SetField(Double_t Bz);
60 #endif
61  //* Constructor (empty)
62 
64 
65  //* Destructor (empty)
66 
67  ~KFParticle() { ; }
68 
69  //* Construction of mother particle by its 2-3-4 daughters
70 
71  KFParticle(const KFParticle &d1, const KFParticle &d2, Bool_t gamma = false);
72 
73  KFParticle(const KFParticle &d1, const KFParticle &d2, const KFParticle &d3);
74 
75  KFParticle(const KFParticle &d1, const KFParticle &d2, const KFParticle &d3, const KFParticle &d4);
76 
77  //* Initialisation from "cartesian" coordinates ( X Y Z Px Py Pz )
78  //* Parameters, covariance matrix, charge and PID hypothesis should be provided
79 
80  void Create(const Double_t Param[], const Double_t Cov[], Int_t Charge, Double_t mass /*Int_t PID*/);
81 
82 #ifdef HomogeneousField
83  //* Initialisation from ALICE track, PID hypothesis shoould be provided
84 
85  KFParticle(const KFPTrack &track, Int_t PID);
86 
87  //* Initialisation from VVertex
88 
89  KFParticle(const KFPVertex &vertex);
90 #endif
91 
92 #ifdef NonhomogeneousField
93  KFParticle(CbmKFTrackInterface *Track, Double_t *z0 = 0, Int_t *qHypo = 0, Int_t *PID = 0);
94  KFParticle(CbmKFVertexInterface &vertex);
95 #endif
96 
97  void CleanDaughtersId() { fDaughtersIds.clear(); }
98  void SetNDaughters(int n) { fDaughtersIds.reserve(n); }
99  void AddDaughter(int id) { fDaughtersIds.push_back(id); }
100 
101  //* Initialise covariance matrix and set current parameters to 0.0
102 
103  void Initialize();
104 
105  //* Set decay vertex parameters for linearisation
106 
107  void SetVtxGuess(Double_t x, Double_t y, Double_t z);
108 
109  //*
110  //* ACCESSORS
111  //*
112 
113  //* Simple accessors
114 
115  Double_t GetX() const; //* x of current position
116  Double_t GetY() const; //* y of current position
117  Double_t GetZ() const; //* z of current position
118  Double_t GetPx() const; //* x-compoment of 3-momentum
119  Double_t GetPy() const; //* y-compoment of 3-momentum
120  Double_t GetPz() const; //* z-compoment of 3-momentum
121  Double_t GetE() const; //* energy
122  Double_t GetS() const; //* decay length / momentum
123  Int_t GetQ() const; //* charge
124  Double_t GetChi2() const; //* chi^2
125  Int_t GetNDF() const; //* Number of Degrees of Freedom
126 
127  Bool_t GetAtProductionVertex() const { return fAtProductionVertex; }
129 
130 #ifdef NonhomogeneousField
131  float *GetFieldCoeff() { return fieldRegion; }
132  void SetFieldCoeff(float c, int i) { fieldRegion[i] = c; }
133 #endif
134 
135  const Double_t &X() const { return fP[0]; }
136  const Double_t &Y() const { return fP[1]; }
137  const Double_t &Z() const { return fP[2]; }
138  const Double_t &Px() const { return fP[3]; }
139  const Double_t &Py() const { return fP[4]; }
140  const Double_t &Pz() const { return fP[5]; }
141  const Double_t &E() const { return fP[6]; }
142  const Double_t &S() const { return fP[7]; }
143  const Int_t &Q() const { return fQ; }
144  const Double_t &Chi2() const { return fChi2; }
145  const Int_t &NDF() const { return fNDF; }
146 
147  Double_t GetParameter(int i) const;
148  Double_t GetCovariance(int i) const;
149  Double_t GetCovariance(int i, int j) const;
150 
151  //* Accessors with calculations, value returned w/o error flag
152 
153  Double_t GetP() const; //* momentum
154  Double_t GetPt() const; //* transverse momentum
155  Double_t GetEta() const; //* pseudorapidity
156  Double_t GetPhi() const; //* phi
157  Double_t GetMomentum() const; //* momentum (same as GetP() )
158  Double_t GetMass() const; //* mass
159  Double_t GetDecayLength() const; //* decay length
160  Double_t GetDecayLengthXY() const; //* decay length in XY
161  Double_t GetLifeTime() const; //* life time
162  Double_t GetR() const; //* distance to the origin
163 
164  //* Accessors to estimated errors
165 
166  Double_t GetErrX() const; //* x of current position
167  Double_t GetErrY() const; //* y of current position
168  Double_t GetErrZ() const; //* z of current position
169  Double_t GetErrPx() const; //* x-compoment of 3-momentum
170  Double_t GetErrPy() const; //* y-compoment of 3-momentum
171  Double_t GetErrPz() const; //* z-compoment of 3-momentum
172  Double_t GetErrE() const; //* energy
173  Double_t GetErrS() const; //* decay length / momentum
174  Double_t GetErrP() const; //* momentum
175  Double_t GetErrPt() const; //* transverse momentum
176  Double_t GetErrEta() const; //* pseudorapidity
177  Double_t GetErrPhi() const; //* phi
178  Double_t GetErrMomentum() const; //* momentum
179  Double_t GetErrMass() const; //* mass
180  Double_t GetErrDecayLength() const; //* decay length
181  Double_t GetErrDecayLengthXY() const; //* decay length in XY
182  Double_t GetErrLifeTime() const; //* life time
183  Double_t GetErrR() const; //* distance to the origin
184 
185  //* Accessors with calculations( &value, &estimated sigma )
186  //* error flag returned (0 means no error during calculations)
187 
188  int GetP(Double_t &P, Double_t &SigmaP) const; //* momentum
189  int GetPt(Double_t &Pt, Double_t &SigmaPt) const; //* transverse momentum
190  int GetEta(Double_t &Eta, Double_t &SigmaEta) const; //* pseudorapidity
191  int GetPhi(Double_t &Phi, Double_t &SigmaPhi) const; //* phi
192  int GetMomentum(Double_t &P, Double_t &SigmaP) const; //* momentum
193  int GetMass(Double_t &M, Double_t &SigmaM) const; //* mass
194  int GetDecayLength(Double_t &L, Double_t &SigmaL) const; //* decay length
195  int GetDecayLengthXY(Double_t &L, Double_t &SigmaL) const; //* decay length in XY
196  int GetLifeTime(Double_t &T, Double_t &SigmaT) const; //* life time
197  int GetR(Double_t &R, Double_t &SigmaR) const; //* R
198  Double_t GetRapidity() const { return 0.5 * log((fP[6] + fP[5]) / (fP[6] - fP[5])); }
199  Double_t GetTheta() const { return atan2(GetPt(), fP[5]); }
200 
201  //*
202  //* MODIFIERS
203  //*
204 
205  Double_t &X();
206  Double_t &Y();
207  Double_t &Z();
208  Double_t &Px();
209  Double_t &Py();
210  Double_t &Pz();
211  Double_t &E();
212  Double_t &S();
213  Int_t &Q();
214  Double_t &Chi2();
215  Int_t &NDF();
216 
217  Double_t &Parameter(int i);
218  Double_t &Covariance(int i);
219  Double_t &Covariance(int i, int j);
220  Double_t *Parameters();
221  Double_t *CovarianceMatrix();
222 
223  //*
224  //* CONSTRUCTION OF THE PARTICLE BY ITS DAUGHTERS AND MOTHER
225  //* USING THE KALMAN FILTER METHOD
226  //*
227 
228  //* Add daughter to the particle
229 
230  void AddDaughter(const KFParticle &Daughter);
231 
232  //* Add daughter via += operator: ex.{ D0; D0+=Pion; D0+= Kaon; }
233 
234  void operator+=(const KFParticle &Daughter);
235 
236  //* Set production vertex
237 
238  void SetProductionVertex(const KFParticle &Vtx);
239 
240  //* Set mass constraint
241 
242  void SetMassConstraint(Double_t Mass, Double_t SigmaMass = 0);
243 
244  //* Set no decay length for resonances
245 
246  void SetNoDecayLength();
247 
248  //* Everything in one go
249 
250  void Construct(const KFParticle *vDaughters[], int nDaughters, const KFParticle *ProdVtx = nullptr, Double_t Mass = -1, Bool_t IsConstrained = 0);
251 
252  //*
253  //* TRANSPORT
254  //*
255  //* ( main transportation parameter is S = SignedPath/Momentum )
256  //* ( parameters of decay & production vertices are stored locally )
257  //*
258 
259  //* Transport the particle to its decay vertex
260 
261  void TransportToDecayVertex();
262 
263  //* Transport the particle to its production vertex
264 
266 
267  //* Transport the particle close to xyz[] point
268 
269  void TransportToPoint(const Double_t xyz[]);
270 
271  //* Transport the particle close to VVertex
272 #ifdef HomogeneousField
273  void TransportToVertex(const KFPVertex &v);
274 #endif
275  //* Transport the particle close to another particle p
276 
277  void TransportToParticle(const KFParticle &p);
278 
279  //* Transport the particle on dS parameter (SignedPath/Momentum)
280 
281  void TransportToDS(Double_t dS);
282 
283  //* Get dS to a certain space point
284 
285  Double_t GetDStoPoint(const Double_t xyz[]) const;
286 
287  //* Get dS to other particle p (dSp for particle p also returned)
288 
289  void GetDStoParticle(const KFParticle &p, Double_t &DS, Double_t &DSp) const;
290 
291  //* Get dS to other particle p in XY-plane
292 
293  void GetDStoParticleXY(const KFParticleBase &p, Double_t &DS, Double_t &DSp) const;
294 
295  //*
296  //* OTHER UTILITIES
297  //*
298 
299  //* Calculate distance from another object [cm]
300 
301  Double_t GetDistanceFromVertex(const Double_t vtx[]) const;
302  Double_t GetDistanceFromVertex(const KFParticle &Vtx) const;
303 #ifdef HomogeneousField
304  Double_t GetDistanceFromVertex(const KFPVertex &Vtx) const;
305 #endif
306  Double_t GetDistanceFromParticle(const KFParticle &p) const;
307 
308  //* Calculate sqrt(Chi2/ndf) deviation from another object
309  //* ( v = [xyz]-vertex, Cv=[Cxx,Cxy,Cyy,Cxz,Cyz,Czz]-covariance matrix )
310 
311  Double_t GetDeviationFromVertex(const Double_t v[], const Double_t Cv[] = 0) const;
312  Double_t GetDeviationFromVertex(const KFParticle &Vtx) const;
313 #ifdef HomogeneousField
314  Double_t GetDeviationFromVertex(const KFPVertex &Vtx) const;
315 #endif
316  Double_t GetDeviationFromParticle(const KFParticle &p) const;
317 
318  //* Calculate distance from another object [cm] in XY-plane
319 
320  Bool_t GetDistanceFromVertexXY(const Double_t vtx[], Double_t &val, Double_t &err) const;
321  Bool_t GetDistanceFromVertexXY(const Double_t vtx[], const Double_t Cv[], Double_t &val, Double_t &err) const;
322  Bool_t GetDistanceFromVertexXY(const KFParticle &Vtx, Double_t &val, Double_t &err) const;
323 #ifdef HomogeneousField
324  Bool_t GetDistanceFromVertexXY(const KFPVertex &Vtx, Double_t &val, Double_t &err) const;
325 #endif
326 
327  Double_t GetDistanceFromVertexXY(const Double_t vtx[]) const;
328  Double_t GetDistanceFromVertexXY(const KFParticle &Vtx) const;
329 #ifdef HomogeneousField
330  Double_t GetDistanceFromVertexXY(const KFPVertex &Vtx) const;
331 #endif
332  Double_t GetDistanceFromParticleXY(const KFParticle &p) const;
333 
334  //* Calculate sqrt(Chi2/ndf) deviation from another object in XY plane
335  //* ( v = [xyz]-vertex, Cv=[Cxx,Cxy,Cyy,Cxz,Cyz,Czz]-covariance matrix )
336 
337  Double_t GetDeviationFromVertexXY(const Double_t v[], const Double_t Cv[] = 0) const;
338  Double_t GetDeviationFromVertexXY(const KFParticle &Vtx) const;
339 #ifdef HomogeneousField
340  Double_t GetDeviationFromVertexXY(const KFPVertex &Vtx) const;
341 #endif
342  Double_t GetDeviationFromParticleXY(const KFParticle &p) const;
343 
344  //* Calculate opennig angle between two particles
345 
346  Double_t GetAngle(const KFParticle &p) const;
347  Double_t GetAngleXY(const KFParticle &p) const;
348  Double_t GetAngleRZ(const KFParticle &p) const;
349 
350  //* Subtract the particle from the vertex
351 
352  void SubtractFromVertex(KFParticle &v) const;
353 
354  //* Special method for creating gammas
355 
356  void ConstructGamma(const KFParticle &daughter1, const KFParticle &daughter2);
357 
358  // * Pseudo Proper Time of decay = (r*pt) / |pt| * M/|pt|
359  // @primVertex - primary vertex
360  // @mass - mass of the mother particle (in the case of "Hb -> JPsi" it would be JPsi mass)
361  // @*timeErr2 - squared error of the decay time. If timeErr2 = 0 it isn't calculated
362  Double_t GetPseudoProperDecayTime(const KFParticle &primVertex, const Double_t &mass, Double_t *timeErr2 = 0) const;
363 
364  void GetFieldValue(const Double_t xyz[], Double_t B[]) const;
365 
366  protected:
367  //*
368  //* INTERNAL STUFF
369  //*
370 
371  //* Method to access ALICE field
372 #ifdef HomogeneousField
373  static Double_t GetFieldAlice();
374 #endif
375  //* Other methods required by the abstract KFParticleBase class
376 
377  void GetDStoParticle(const KFParticleBase &p, Double_t &DS, Double_t &DSp) const;
378  void Transport(Double_t dS, Double_t P[], Double_t C[]) const;
379  static void GetExternalTrackParam(const KFParticleBase &p, Double_t &X, Double_t &Alpha, Double_t P[5]);
380 
381  // void GetDStoParticleALICE( const KFParticleBase &p, Double_t &DS, Double_t &DS1 ) const;
382 
383  private:
384 #ifdef HomogeneousField
385  static Double_t fgBz; //* Bz compoment of the magnetic field
386 #endif
387 #ifdef NonhomogeneousField
388  float fieldRegion[10];
389 #endif
390 
391  // ClassDef( KFParticle, 1 );
392 };
393 
394 //---------------------------------------------------------------------
395 //
396 // Inline implementation of the KFParticle methods
397 //
398 //---------------------------------------------------------------------
399 
400 #ifdef HomogeneousField
401 inline void KFParticle::SetField(Double_t Bz)
402 {
403  fgBz = Bz;
404 }
405 #endif
406 
407 inline KFParticle::KFParticle(const KFParticle &d1, const KFParticle &d2, const KFParticle &d3)
408 {
409  KFParticle mother;
410  mother += d1;
411  mother += d2;
412  mother += d3;
413  *this = mother;
414 }
415 
416 inline KFParticle::KFParticle(const KFParticle &d1, const KFParticle &d2, const KFParticle &d3, const KFParticle &d4)
417 {
418  KFParticle mother;
419  mother += d1;
420  mother += d2;
421  mother += d3;
422  mother += d4;
423  *this = mother;
424 }
425 
427 {
429 }
430 
431 inline void KFParticle::SetVtxGuess(Double_t x, Double_t y, Double_t z)
432 {
434 }
435 
436 inline Double_t KFParticle::GetX() const
437 {
438  return KFParticleBase::GetX();
439 }
440 
441 inline Double_t KFParticle::GetY() const
442 {
443  return KFParticleBase::GetY();
444 }
445 
446 inline Double_t KFParticle::GetZ() const
447 {
448  return KFParticleBase::GetZ();
449 }
450 
451 inline Double_t KFParticle::GetPx() const
452 {
453  return KFParticleBase::GetPx();
454 }
455 
456 inline Double_t KFParticle::GetPy() const
457 {
458  return KFParticleBase::GetPy();
459 }
460 
461 inline Double_t KFParticle::GetPz() const
462 {
463  return KFParticleBase::GetPz();
464 }
465 
466 inline Double_t KFParticle::GetE() const
467 {
468  return KFParticleBase::GetE();
469 }
470 
471 inline Double_t KFParticle::GetS() const
472 {
473  return KFParticleBase::GetS();
474 }
475 
476 inline Int_t KFParticle::GetQ() const
477 {
478  return KFParticleBase::GetQ();
479 }
480 
481 inline Double_t KFParticle::GetChi2() const
482 {
483  return KFParticleBase::GetChi2();
484 }
485 
486 inline Int_t KFParticle::GetNDF() const
487 {
488  return KFParticleBase::GetNDF();
489 }
490 
491 inline Double_t KFParticle::GetParameter(int i) const
492 {
494 }
495 
496 inline Double_t KFParticle::GetCovariance(int i) const
497 {
499 }
500 
501 inline Double_t KFParticle::GetCovariance(int i, int j) const
502 {
503  return KFParticleBase::GetCovariance(i, j);
504 }
505 
506 inline Double_t KFParticle::GetP() const
507 {
508  Double_t par, err;
509  if (KFParticleBase::GetMomentum(par, err))
510  return 0;
511  else
512  return par;
513 }
514 
515 inline Double_t KFParticle::GetPt() const
516 {
517  Double_t par, err;
518  if (KFParticleBase::GetPt(par, err))
519  return 0;
520  else
521  return par;
522 }
523 
524 inline Double_t KFParticle::GetEta() const
525 {
526  Double_t par, err;
527  if (KFParticleBase::GetEta(par, err))
528  return 0;
529  else
530  return par;
531 }
532 
533 inline Double_t KFParticle::GetPhi() const
534 {
535  Double_t par, err;
536  if (KFParticleBase::GetPhi(par, err))
537  return 0;
538  else
539  return par;
540 }
541 
542 inline Double_t KFParticle::GetMomentum() const
543 {
544  Double_t par, err;
545  if (KFParticleBase::GetMomentum(par, err))
546  return 0;
547  else
548  return par;
549 }
550 
551 inline Double_t KFParticle::GetMass() const
552 {
553  Double_t par, err;
554  if (KFParticleBase::GetMass(par, err))
555  return 0;
556  else
557  return par;
558 }
559 
560 inline Double_t KFParticle::GetDecayLength() const
561 {
562  Double_t par, err;
563  if (KFParticleBase::GetDecayLength(par, err))
564  return 0;
565  else
566  return par;
567 }
568 
569 inline Double_t KFParticle::GetDecayLengthXY() const
570 {
571  Double_t par, err;
572  if (KFParticleBase::GetDecayLengthXY(par, err))
573  return 0;
574  else
575  return par;
576 }
577 
578 inline Double_t KFParticle::GetLifeTime() const
579 {
580  Double_t par, err;
581  if (KFParticleBase::GetLifeTime(par, err))
582  return 0;
583  else
584  return par;
585 }
586 
587 inline Double_t KFParticle::GetR() const
588 {
589  Double_t par, err;
590  if (KFParticleBase::GetR(par, err))
591  return 0;
592  else
593  return par;
594 }
595 
596 inline Double_t KFParticle::GetErrX() const
597 {
598  return sqrt(fabs(GetCovariance(0, 0)));
599 }
600 
601 inline Double_t KFParticle::GetErrY() const
602 {
603  return sqrt(fabs(GetCovariance(1, 1)));
604 }
605 
606 inline Double_t KFParticle::GetErrZ() const
607 {
608  return sqrt(fabs(GetCovariance(2, 2)));
609 }
610 
611 inline Double_t KFParticle::GetErrPx() const
612 {
613  return sqrt(fabs(GetCovariance(3, 3)));
614 }
615 
616 inline Double_t KFParticle::GetErrPy() const
617 {
618  return sqrt(fabs(GetCovariance(4, 4)));
619 }
620 
621 inline Double_t KFParticle::GetErrPz() const
622 {
623  return sqrt(fabs(GetCovariance(5, 5)));
624 }
625 
626 inline Double_t KFParticle::GetErrE() const
627 {
628  return sqrt(fabs(GetCovariance(6, 6)));
629 }
630 
631 inline Double_t KFParticle::GetErrS() const
632 {
633  return sqrt(fabs(GetCovariance(7, 7)));
634 }
635 
636 inline Double_t KFParticle::GetErrP() const
637 {
638  Double_t par, err;
639  if (KFParticleBase::GetMomentum(par, err))
640  return 1.e10;
641  else
642  return err;
643 }
644 
645 inline Double_t KFParticle::GetErrPt() const
646 {
647  Double_t par, err;
648  if (KFParticleBase::GetPt(par, err))
649  return 1.e10;
650  else
651  return err;
652 }
653 
654 inline Double_t KFParticle::GetErrEta() const
655 {
656  Double_t par, err;
657  if (KFParticleBase::GetEta(par, err))
658  return 1.e10;
659  else
660  return err;
661 }
662 
663 inline Double_t KFParticle::GetErrPhi() const
664 {
665  Double_t par, err;
666  if (KFParticleBase::GetPhi(par, err))
667  return 1.e10;
668  else
669  return err;
670 }
671 
672 inline Double_t KFParticle::GetErrMomentum() const
673 {
674  Double_t par, err;
675  if (KFParticleBase::GetMomentum(par, err))
676  return 1.e10;
677  else
678  return err;
679 }
680 
681 inline Double_t KFParticle::GetErrMass() const
682 {
683  Double_t par, err;
684  if (KFParticleBase::GetMass(par, err))
685  return 1.e10;
686  else
687  return err;
688 }
689 
690 inline Double_t KFParticle::GetErrDecayLength() const
691 {
692  Double_t par, err;
693  if (KFParticleBase::GetDecayLength(par, err))
694  return 1.e10;
695  else
696  return err;
697 }
698 
699 inline Double_t KFParticle::GetErrDecayLengthXY() const
700 {
701  Double_t par, err;
702  if (KFParticleBase::GetDecayLengthXY(par, err))
703  return 1.e10;
704  else
705  return err;
706 }
707 
708 inline Double_t KFParticle::GetErrLifeTime() const
709 {
710  Double_t par, err;
711  if (KFParticleBase::GetLifeTime(par, err))
712  return 1.e10;
713  else
714  return err;
715 }
716 
717 inline Double_t KFParticle::GetErrR() const
718 {
719  Double_t par, err;
720  if (KFParticleBase::GetR(par, err))
721  return 1.e10;
722  else
723  return err;
724 }
725 
726 inline int KFParticle::GetP(Double_t &P, Double_t &SigmaP) const
727 {
728  return KFParticleBase::GetMomentum(P, SigmaP);
729 }
730 
731 inline int KFParticle::GetPt(Double_t &Pt, Double_t &SigmaPt) const
732 {
733  return KFParticleBase::GetPt(Pt, SigmaPt);
734 }
735 
736 inline int KFParticle::GetEta(Double_t &Eta, Double_t &SigmaEta) const
737 {
738  return KFParticleBase::GetEta(Eta, SigmaEta);
739 }
740 
741 inline int KFParticle::GetPhi(Double_t &Phi, Double_t &SigmaPhi) const
742 {
743  return KFParticleBase::GetPhi(Phi, SigmaPhi);
744 }
745 
746 inline int KFParticle::GetMomentum(Double_t &P, Double_t &SigmaP) const
747 {
748  return KFParticleBase::GetMomentum(P, SigmaP);
749 }
750 
751 inline int KFParticle::GetMass(Double_t &M, Double_t &SigmaM) const
752 {
753  return KFParticleBase::GetMass(M, SigmaM);
754 }
755 
756 inline int KFParticle::GetDecayLength(Double_t &L, Double_t &SigmaL) const
757 {
758  return KFParticleBase::GetDecayLength(L, SigmaL);
759 }
760 
761 inline int KFParticle::GetDecayLengthXY(Double_t &L, Double_t &SigmaL) const
762 {
763  return KFParticleBase::GetDecayLengthXY(L, SigmaL);
764 }
765 
766 inline int KFParticle::GetLifeTime(Double_t &T, Double_t &SigmaT) const
767 {
768  return KFParticleBase::GetLifeTime(T, SigmaT);
769 }
770 
771 inline int KFParticle::GetR(Double_t &R, Double_t &SigmaR) const
772 {
773  return KFParticleBase::GetR(R, SigmaR);
774 }
775 
776 inline Double_t &KFParticle::X()
777 {
778  return KFParticleBase::X();
779 }
780 
781 inline Double_t &KFParticle::Y()
782 {
783  return KFParticleBase::Y();
784 }
785 
786 inline Double_t &KFParticle::Z()
787 {
788  return KFParticleBase::Z();
789 }
790 
791 inline Double_t &KFParticle::Px()
792 {
793  return KFParticleBase::Px();
794 }
795 
796 inline Double_t &KFParticle::Py()
797 {
798  return KFParticleBase::Py();
799 }
800 
801 inline Double_t &KFParticle::Pz()
802 {
803  return KFParticleBase::Pz();
804 }
805 
806 inline Double_t &KFParticle::E()
807 {
808  return KFParticleBase::E();
809 }
810 
811 inline Double_t &KFParticle::S()
812 {
813  return KFParticleBase::S();
814 }
815 
816 inline Int_t &KFParticle::Q()
817 {
818  return KFParticleBase::Q();
819 }
820 
821 inline Double_t &KFParticle::Chi2()
822 {
823  return KFParticleBase::Chi2();
824 }
825 
826 inline Int_t &KFParticle::NDF()
827 {
828  return KFParticleBase::NDF();
829 }
830 
831 inline Double_t &KFParticle::Parameter(int i)
832 {
833  return KFParticleBase::Parameter(i);
834 }
835 
836 inline Double_t &KFParticle::Covariance(int i)
837 {
838  return KFParticleBase::Covariance(i);
839 }
840 
841 inline Double_t &KFParticle::Covariance(int i, int j)
842 {
843  return KFParticleBase::Covariance(i, j);
844 }
845 
846 inline Double_t *KFParticle::Parameters()
847 {
848  return fP;
849 }
850 
852 {
853  return fC;
854 }
855 
856 inline void KFParticle::operator+=(const KFParticle &Daughter)
857 {
858  KFParticleBase::operator+=(Daughter);
859 }
860 
861 inline void KFParticle::AddDaughter(const KFParticle &Daughter)
862 {
863  KFParticleBase::AddDaughter(Daughter);
864 }
865 
867 {
869 }
870 
871 inline void KFParticle::SetMassConstraint(Double_t Mass, Double_t SigmaMass)
872 {
873  KFParticleBase::SetMassConstraint(Mass, SigmaMass);
874 }
875 
877 {
879 }
880 
881 inline void KFParticle::Construct(const KFParticle *vDaughters[], int nDaughters, const KFParticle *ProdVtx, Double_t Mass, Bool_t IsConstrained)
882 {
883  KFParticleBase::Construct((const KFParticleBase **)vDaughters, nDaughters, (const KFParticleBase *)ProdVtx, Mass, IsConstrained);
884 }
885 
887 {
889 }
890 
892 {
894 }
895 
896 inline void KFParticle::TransportToPoint(const Double_t xyz[])
897 {
899 }
900 #ifdef HomogeneousField
902 {
904 }
905 #endif
907 {
908  Double_t dS, dSp;
909  GetDStoParticle(p, dS, dSp);
910  TransportToDS(dS);
911 }
912 
913 inline void KFParticle::TransportToDS(Double_t dS)
914 {
916 }
917 
918 inline Double_t KFParticle::GetDStoPoint(const Double_t xyz[]) const
919 {
920 #ifdef HomogeneousField
922 #endif
923 #ifdef NonhomogeneousField
925 #endif
926 }
927 
928 inline void KFParticle::GetDStoParticle(const KFParticle &p, Double_t &DS, Double_t &DSp) const
929 {
930  GetDStoParticleXY(p, DS, DSp);
931 }
932 
933 inline Double_t KFParticle::GetDistanceFromVertex(const Double_t vtx[]) const
934 {
936 }
937 
938 inline Double_t KFParticle::GetDeviationFromVertex(const Double_t v[], const Double_t Cv[]) const
939 {
941 }
942 
943 inline Double_t KFParticle::GetDistanceFromVertex(const KFParticle &Vtx) const
944 {
946 }
947 
948 inline Double_t KFParticle::GetDeviationFromVertex(const KFParticle &Vtx) const
949 {
951 }
952 #ifdef HomogeneousField
953 inline Double_t KFParticle::GetDistanceFromVertex(const KFPVertex &Vtx) const
954 {
955  return GetDistanceFromVertex(KFParticle(Vtx));
956 }
957 
958 inline Double_t KFParticle::GetDeviationFromVertex(const KFPVertex &Vtx) const
959 {
960  return GetDeviationFromVertex(KFParticle(Vtx));
961 }
962 #endif
963 inline Double_t KFParticle::GetDistanceFromParticle(const KFParticle &p) const
964 {
966 }
967 
968 inline Double_t KFParticle::GetDeviationFromParticle(const KFParticle &p) const
969 {
971 }
972 
974 {
976 }
977 
978 #ifdef HomogeneousField
979 inline Double_t KFParticle::GetFieldAlice()
980 {
981  return fgBz;
982 }
983 #endif
984 
985 #ifdef HomogeneousField
986 inline void KFParticle::GetFieldValue(const Double_t * /*xyz*/, Double_t B[]) const
987 {
988  B[0] = B[1] = 0;
989  B[2] = GetFieldAlice();
990 }
991 #endif
992 
993 #ifdef NonhomogeneousField
994 inline void KFParticle::GetFieldValue(const Double_t xyz[], Double_t B[]) const
995 {
996  // FairField *MF = CbmKF::Instance()->GetMagneticField();
997  // MF->GetFieldValue( xyz, B );
998 }
999 #endif
1000 
1001 inline void KFParticle::GetDStoParticle(const KFParticleBase &p, Double_t &DS, Double_t &DSp) const
1002 {
1003  GetDStoParticleXY(p, DS, DSp);
1004 }
1005 
1006 inline void KFParticle::GetDStoParticleXY(const KFParticleBase &p, Double_t &DS, Double_t &DSp) const
1007 {
1008 #ifdef HomogeneousField
1010 #endif
1011 #ifdef NonhomogeneousField
1013 #endif
1014  // GetDStoParticleALICE( p, DS, DSp ) ;
1015 }
1016 
1017 inline void KFParticle::Transport(Double_t dS, Double_t P[], Double_t C[]) const
1018 {
1019 #ifdef HomogeneousField
1021 #endif
1022 #ifdef NonhomogeneousField
1023  KFParticleBase::TransportCBM(dS, P, C);
1024 #endif
1025 }
1026 
1027 inline void KFParticle::ConstructGamma(const KFParticle &daughter1, const KFParticle &daughter2)
1028 {
1029 #ifdef HomogeneousField
1030  KFParticleBase::ConstructGammaBz(daughter1, daughter2, GetFieldAlice());
1031 #endif
1032 }
1033 
1034 #endif
Double_t GetPy() const
void GetDStoParticleXY(const KFParticleBase &p, Double_t &DS, Double_t &DSp) const
Definition: KFParticle.h:1006
Double_t GetDistanceFromVertex(const Double_t vtx[]) const
static Double_t GetFieldAlice()
Definition: KFParticle.h:979
void SetNoDecayLength()
Definition: KFParticle.h:876
Double_t GetZ() const
Double_t GetDistanceFromParticle(const KFParticle &p) const
Definition: KFParticle.h:963
Double_t GetDecayLength() const
Definition: KFParticle.h:560
Double_t GetMomentum() const
Definition: KFParticle.h:542
void TransportToProductionVertex()
Double_t GetPseudoProperDecayTime(const KFParticle &primVertex, const Double_t &mass, Double_t *timeErr2=0) const
void TransportToDecayVertex()
Double_t GetErrPx() const
Definition: KFParticle.h:611
void SetMassConstraint(Double_t Mass, Double_t SigmaMass=0)
Definition: KFParticle.h:871
void SetProductionVertex(const KFParticle &Vtx)
Definition: KFParticle.h:866
const Int_t & NDF() const
Definition: KFParticle.h:145
Double_t GetChi2() const
Definition: KFParticle.h:481
Double_t GetDecayLengthXY() const
Definition: KFParticle.h:569
Int_t GetQ() const
Double_t GetDeviationFromVertex(const Double_t v[], const Double_t Cv[]=0) const
Definition: KFParticle.h:938
Double_t GetPhi() const
Definition: KFParticle.h:533
Double_t GetPx() const
Double_t GetErrPhi() const
Definition: KFParticle.h:663
Int_t GetMomentum(Double_t &P, Double_t &SigmaP) const
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:40
Double_t GetPz() const
Double_t fP[8]
void TransportCBM(Double_t dS, Double_t P[], Double_t C[]) const
Double_t GetE() const
Definition: KFParticle.h:466
Double_t GetErrPz() const
Definition: KFParticle.h:621
Int_t GetPhi(Double_t &Phi, Double_t &SigmaPhi) const
const Int_t & Q() const
void ConstructGamma(const KFParticle &daughter1, const KFParticle &daughter2)
Definition: KFParticle.h:1027
Double_t GetDStoPointBz(Double_t Bz, const Double_t xyz[]) const
const Double_t & Pz() const
void SetVtxGuess(Double_t x, Double_t y, Double_t z)
Definition: KFParticle.h:431
void GetDStoParticleCBM(const KFParticleBase &p, Double_t &dS, Double_t &dS1) const
Double_t * CovarianceMatrix()
Definition: KFParticle.h:851
Int_t GetMass(Double_t &M, Double_t &SigmaM) const
friend F32vec4 log(const F32vec4 &a)
Definition: P4_F32vec4.h:129
Double_t GetY() const
Double_t GetPx() const
Definition: KFParticle.h:451
const Double_t & Chi2() const
Definition: KFParticle.h:144
const Double_t & Z() const
Definition: KFParticle.h:137
const Double_t & Py() const
void Transport(Double_t dS, Double_t P[], Double_t C[]) const
Definition: KFParticle.h:1017
void Initialize()
Definition: KFParticle.h:426
Double_t GetMass() const
Definition: KFParticle.h:551
Double_t GetAngle(const KFParticle &p) const
Double_t GetZ() const
Definition: KFParticle.h:446
void SubtractFromVertex(KFParticleBase &Vtx) const
Double_t & Covariance(int i)
Definition: KFParticle.h:836
Double_t GetY() const
Definition: KFParticle.h:441
Double_t GetErrLifeTime() const
Definition: KFParticle.h:708
void SetVtxGuess(Double_t x, Double_t y, Double_t z)
Double_t GetParameter(Int_t i) const
const Double_t & S() const
Definition: KFParticle.h:142
void operator+=(const KFParticleBase &Daughter)
const Double_t & Chi2() const
Double_t * Parameters()
Definition: KFParticle.h:846
const Double_t & X() const
Definition: KFParticle.h:135
Double_t GetDStoPoint(const Double_t xyz[]) const
Definition: KFParticle.h:918
Double_t GetCovariance(int i) const
Definition: KFParticle.h:496
const Double_t & Y() const
__m128 v
Definition: P4_F32vec4.h:15
unsigned int i
Definition: P4_F32vec4.h:33
Double_t GetErrEta() const
Definition: KFParticle.h:654
Int_t GetNDF() const
void TransportToPoint(const Double_t xyz[])
Definition: KFParticle.h:896
Double_t GetP() const
Definition: KFParticle.h:506
void SetNoDecayLength()
const Double_t & Y() const
Definition: KFParticle.h:136
Bool_t GetAtProductionVertex() const
Definition: KFParticle.h:127
Double_t GetX() const
Definition: KFParticle.h:436
Double_t GetErrMomentum() const
Definition: KFParticle.h:672
Double_t GetEta() const
Definition: KFParticle.h:524
void TransportToProductionVertex()
Definition: KFParticle.h:891
Int_t GetEta(Double_t &Eta, Double_t &SigmaEta) const
Int_t GetPt(Double_t &Pt, Double_t &SigmaPt) const
void SetProductionVertex(const KFParticleBase &Vtx)
Double_t GetTheta() const
Definition: KFParticle.h:199
Int_t GetQ() const
Definition: KFParticle.h:476
Double_t GetErrS() const
Definition: KFParticle.h:631
Double_t GetErrY() const
Definition: KFParticle.h:601
void GetDStoParticleBz(Double_t Bz, const KFParticleBase &p, Double_t &dS, Double_t &dS1) const
Double_t GetErrE() const
Definition: KFParticle.h:626
void SetMassConstraint(Double_t Mass, Double_t SigmaMass=0)
Int_t GetDecayLengthXY(Double_t &L, Double_t &SigmaL) const
void TransportToDS(Double_t dS)
Definition: KFParticle.h:913
void SubtractFromVertex(KFParticle &v) const
Definition: KFParticle.h:973
Double_t GetErrMass() const
Definition: KFParticle.h:681
Double_t GetErrR() const
Definition: KFParticle.h:717
Double_t GetAngleXY(const KFParticle &p) const
const Int_t & NDF() const
const Double_t & Pz() const
Definition: KFParticle.h:140
const Int_t & Q() const
Definition: KFParticle.h:143
void SetAtProductionVertex(Bool_t b)
Definition: KFParticle.h:128
void SetNDaughters(int n)
Definition: KFParticle.h:98
Double_t GetPy() const
Definition: KFParticle.h:456
Double_t GetDeviationFromVertex(const Double_t v[], const Double_t Cv[]=0) const
Double_t GetDeviationFromParticleXY(const KFParticle &p) const
Double_t fC[36]
const Double_t & Py() const
Definition: KFParticle.h:139
Double_t GetErrZ() const
Definition: KFParticle.h:606
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:58
Int_t GetNDF() const
Definition: KFParticle.h:486
Double_t GetCovariance(Int_t i) const
Int_t GetR(Double_t &R, Double_t &SigmaR) const
Double_t GetDeviationFromParticle(const KFParticle &p) const
Definition: KFParticle.h:968
Double_t GetChi2() const
void TransportToDecayVertex()
Definition: KFParticle.h:886
const Double_t & Px() const
Double_t GetDStoPointCBM(const Double_t xyz[]) const
Double_t GetParameter(int i) const
Definition: KFParticle.h:491
friend F32vec4 atan2(const F32vec4 &y, const F32vec4 &x)
Definition: P4_F32vec4.h:136
const Double_t & E() const
Definition: KFParticle.h:141
Double_t GetPz() const
Definition: KFParticle.h:461
Double_t GetErrDecayLengthXY() const
Definition: KFParticle.h:699
Double_t GetX() const
void TransportToVertex(const KFPVertex &v)
Definition: KFParticle.h:901
Double_t GetErrPt() const
Definition: KFParticle.h:645
Double_t GetRapidity() const
Definition: KFParticle.h:198
Double_t GetDistanceFromVertex(const Double_t vtx[]) const
Definition: KFParticle.h:933
Double_t GetDistanceFromParticleXY(const KFParticle &p) const
Bool_t GetDistanceFromVertexXY(const Double_t vtx[], Double_t &val, Double_t &err) const
Double_t GetE() const
void Construct(const KFParticleBase *vDaughters[], Int_t nDaughters, const KFParticleBase *ProdVtx=nullptr, Double_t Mass=-1, Bool_t IsConstrained=0)
void Construct(const KFParticle *vDaughters[], int nDaughters, const KFParticle *ProdVtx=nullptr, Double_t Mass=-1, Bool_t IsConstrained=0)
Definition: KFParticle.h:881
void GetDStoParticle(const KFParticle &p, Double_t &DS, Double_t &DSp) const
Definition: KFParticle.h:928
const Double_t & Px() const
Definition: KFParticle.h:138
void TransportToDS(Double_t dS)
Double_t GetS() const
const Double_t & X() const
Double_t GetR() const
Definition: KFParticle.h:587
Double_t GetAngleRZ(const KFParticle &p) const
Double_t GetDeviationFromParticle(const KFParticleBase &p) const
void operator+=(const KFParticle &Daughter)
Definition: KFParticle.h:856
Double_t GetErrX() const
Definition: KFParticle.h:596
Bool_t fAtProductionVertex
Double_t GetLifeTime() const
Definition: KFParticle.h:578
Double_t GetErrDecayLength() const
Definition: KFParticle.h:690
Double_t GetErrP() const
Definition: KFParticle.h:636
void AddDaughter(const KFParticleBase &Daughter)
static void GetExternalTrackParam(const KFParticleBase &p, Double_t &X, Double_t &Alpha, Double_t P[5])
Double_t & Parameter(Int_t i)
Double_t GetS() const
Definition: KFParticle.h:471
void Create(const Double_t Param[], const Double_t Cov[], Int_t Charge, Double_t mass)
std::vector< int > fDaughtersIds
Double_t GetDistanceFromParticle(const KFParticleBase &p) const
void AddDaughter(int id)
Definition: KFParticle.h:99
static void SetField(Double_t Bz)
Definition: KFParticle.h:401
Double_t & Parameter(int i)
Definition: KFParticle.h:831
void GetFieldValue(const Double_t xyz[], Double_t B[]) const
Definition: KFParticle.h:986
const Double_t & S() const
Double_t GetDeviationFromVertexXY(const Double_t v[], const Double_t Cv[]=0) const
Int_t GetDecayLength(Double_t &L, Double_t &SigmaL) const
const Double_t & E() const
const Double_t & Z() const
void TransportToParticle(const KFParticle &p)
Definition: KFParticle.h:906
void ConstructGammaBz(const KFParticleBase &daughter1, const KFParticleBase &daughter2, double Bz)
void CleanDaughtersId()
Definition: KFParticle.h:97
Double_t GetErrPy() const
Definition: KFParticle.h:616
Double_t & Covariance(Int_t i)
void TransportBz(Double_t Bz, Double_t dS, Double_t P[], Double_t C[]) const
Int_t GetLifeTime(Double_t &T, Double_t &SigmaT) const
Double_t GetPt() const
Definition: KFParticle.h:515