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