1 //---------------------------------------------------------------------------------
2 // The AliKFParticle class
4 // @author S.Gorbunov, I.Kisel
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
13 // This class is ALICE interface to general mathematics in AliKFParticleBase
15 // -= Copyright © ALICE HLT Group =-
16 //_________________________________________________________________________________
18 #ifndef ALIKFPARTICLE_H
19 #define ALIKFPARTICLE_H
21 #include "AliKFParticleBase.h"
26 class AliExternalTrackParam;
28 class AliKFParticle :public AliKFParticleBase
37 //* Set magnetic field for all particles
39 static void SetField( Double_t Bz );
41 //* Constructor (empty)
43 AliKFParticle():AliKFParticleBase(){ ; }
45 //* Destructor (empty)
49 //* Construction of mother particle by its 2-3-4 daughters
51 AliKFParticle( const AliKFParticle &d1, const AliKFParticle &d2, Bool_t gamma = kFALSE );
53 AliKFParticle( const AliKFParticle &d1, const AliKFParticle &d2,
54 const AliKFParticle &d3 );
56 AliKFParticle( const AliKFParticle &d1, const AliKFParticle &d2,
57 const AliKFParticle &d3, const AliKFParticle &d4 );
59 //* Initialisation from "cartesian" coordinates ( X Y Z Px Py Pz )
60 //* Parameters, covariance matrix, charge and PID hypothesis should be provided
62 void Create( const Double_t Param[], const Double_t Cov[], Int_t Charge, Int_t PID );
63 void Create( const Double_t Param[], const Double_t Cov[], Int_t Charge, Double_t Mass );
65 //* Initialisation from ALICE track, PID hypothesis shoould be provided
67 AliKFParticle( const AliVTrack &track, Int_t PID );
68 AliKFParticle( const AliExternalTrackParam &track, Double_t Mass, Int_t Charge );
70 //* Initialisation from VVertex
72 AliKFParticle( const AliVVertex &vertex );
74 //* Initialise covariance matrix and set current parameters to 0.0
78 //* Set decay vertex parameters for linearisation
80 void SetVtxGuess( Double_t x, Double_t y, Double_t z );
88 Double_t GetX () const ; //* x of current position
89 Double_t GetY () const ; //* y of current position
90 Double_t GetZ () const ; //* z of current position
91 Double_t GetPx () const ; //* x-compoment of 3-momentum
92 Double_t GetPy () const ; //* y-compoment of 3-momentum
93 Double_t GetPz () const ; //* z-compoment of 3-momentum
94 Double_t GetE () const ; //* energy
95 Double_t GetS () const ; //* decay length / momentum
96 Int_t GetQ () const ; //* charge
97 Double_t GetChi2 () const ; //* chi^2
98 Int_t GetNDF () const ; //* Number of Degrees of Freedom
100 const Double_t& X () const { return fP[0]; }
101 const Double_t& Y () const { return fP[1]; }
102 const Double_t& Z () const { return fP[2]; }
103 const Double_t& Px () const { return fP[3]; }
104 const Double_t& Py () const { return fP[4]; }
105 const Double_t& Pz () const { return fP[5]; }
106 const Double_t& E () const { return fP[6]; }
107 const Double_t& S () const { return fP[7]; }
108 const Int_t & Q () const { return fQ; }
109 const Double_t& Chi2 () const { return fChi2; }
110 const Int_t & NDF () const { return fNDF; }
112 Double_t GetParameter ( int i ) const ;
113 Double_t GetCovariance( int i ) const ;
114 Double_t GetCovariance( int i, int j ) const ;
116 //* Accessors with calculations, value returned w/o error flag
118 Double_t GetP () const; //* momentum
119 Double_t GetPt () const; //* transverse momentum
120 Double_t GetEta () const; //* pseudorapidity
121 Double_t GetPhi () const; //* phi
122 Double_t GetMomentum () const; //* momentum (same as GetP() )
123 Double_t GetMass () const; //* mass
124 Double_t GetDecayLength () const; //* decay length
125 Double_t GetDecayLengthXY () const; //* decay length in XY
126 Double_t GetLifeTime () const; //* life time
127 Double_t GetR () const; //* distance to the origin
129 //* Accessors to estimated errors
131 Double_t GetErrX () const ; //* x of current position
132 Double_t GetErrY () const ; //* y of current position
133 Double_t GetErrZ () const ; //* z of current position
134 Double_t GetErrPx () const ; //* x-compoment of 3-momentum
135 Double_t GetErrPy () const ; //* y-compoment of 3-momentum
136 Double_t GetErrPz () const ; //* z-compoment of 3-momentum
137 Double_t GetErrE () const ; //* energy
138 Double_t GetErrS () const ; //* decay length / momentum
139 Double_t GetErrP () const ; //* momentum
140 Double_t GetErrPt () const ; //* transverse momentum
141 Double_t GetErrEta () const ; //* pseudorapidity
142 Double_t GetErrPhi () const ; //* phi
143 Double_t GetErrMomentum () const ; //* momentum
144 Double_t GetErrMass () const ; //* mass
145 Double_t GetErrDecayLength () const ; //* decay length
146 Double_t GetErrDecayLengthXY () const ; //* decay length in XY
147 Double_t GetErrLifeTime () const ; //* life time
148 Double_t GetErrR () const ; //* distance to the origin
150 //* Accessors with calculations( &value, &estimated sigma )
151 //* error flag returned (0 means no error during calculations)
153 int GetP ( Double_t &P, Double_t &SigmaP ) const ; //* momentum
154 int GetPt ( Double_t &Pt, Double_t &SigmaPt ) const ; //* transverse momentum
155 int GetEta ( Double_t &Eta, Double_t &SigmaEta ) const ; //* pseudorapidity
156 int GetPhi ( Double_t &Phi, Double_t &SigmaPhi ) const ; //* phi
157 int GetMomentum ( Double_t &P, Double_t &SigmaP ) const ; //* momentum
158 int GetMass ( Double_t &M, Double_t &SigmaM ) const ; //* mass
159 int GetDecayLength ( Double_t &L, Double_t &SigmaL ) const ; //* decay length
160 int GetDecayLengthXY ( Double_t &L, Double_t &SigmaL ) const ; //* decay length in XY
161 int GetLifeTime ( Double_t &T, Double_t &SigmaT ) const ; //* life time
162 int GetR ( Double_t &R, Double_t &SigmaR ) const ; //* R
181 Double_t & Parameter ( int i ) ;
182 Double_t & Covariance( int i ) ;
183 Double_t & Covariance( int i, int j ) ;
184 Double_t * Parameters () ;
185 Double_t * CovarianceMatrix() ;
188 //* CONSTRUCTION OF THE PARTICLE BY ITS DAUGHTERS AND MOTHER
189 //* USING THE KALMAN FILTER METHOD
193 //* Add daughter to the particle
195 void AddDaughter( const AliKFParticle &Daughter );
197 //* Add daughter via += operator: ex.{ D0; D0+=Pion; D0+= Kaon; }
199 void operator +=( const AliKFParticle &Daughter );
201 //* Set production vertex
203 void SetProductionVertex( const AliKFParticle &Vtx );
205 //* Set mass constraint
207 void SetMassConstraint( Double_t Mass, Double_t SigmaMass = 0 );
209 //* Set no decay length for resonances
211 void SetNoDecayLength();
213 //* Everything in one go
215 void Construct( const AliKFParticle *vDaughters[], int NDaughters,
216 const AliKFParticle *ProdVtx=0, Double_t Mass=-1, Bool_t IsConstrained=0 );
221 //* ( main transportation parameter is S = SignedPath/Momentum )
222 //* ( parameters of decay & production vertices are stored locally )
225 //* Transport the particle to its decay vertex
227 void TransportToDecayVertex();
229 //* Transport the particle to its production vertex
231 void TransportToProductionVertex();
233 //* Transport the particle close to xyz[] point
235 void TransportToPoint( const Double_t xyz[] );
237 //* Transport the particle close to VVertex
239 void TransportToVertex( const AliVVertex &v );
241 //* Transport the particle close to another particle p
243 void TransportToParticle( const AliKFParticle &p );
245 //* Transport the particle on dS parameter (SignedPath/Momentum)
247 void TransportToDS( Double_t dS );
249 //* Get dS to a certain space point
251 Double_t GetDStoPoint( const Double_t xyz[] ) const ;
253 //* Get dS to other particle p (dSp for particle p also returned)
255 void GetDStoParticle( const AliKFParticle &p,
256 Double_t &DS, Double_t &DSp ) const ;
258 //* Get dS to other particle p in XY-plane
260 void GetDStoParticleXY( const AliKFParticleBase &p,
261 Double_t &DS, Double_t &DSp ) const ;
268 //* Calculate distance from another object [cm]
270 Double_t GetDistanceFromVertex( const Double_t vtx[] ) const ;
271 Double_t GetDistanceFromVertex( const AliKFParticle &Vtx ) const ;
272 Double_t GetDistanceFromVertex( const AliVVertex &Vtx ) const ;
273 Double_t GetDistanceFromParticle( const AliKFParticle &p ) const ;
275 //* Calculate sqrt(Chi2/ndf) deviation from another object
276 //* ( v = [xyz]-vertex, Cv=[Cxx,Cxy,Cyy,Cxz,Cyz,Czz]-covariance matrix )
278 Double_t GetDeviationFromVertex( const Double_t v[], const Double_t Cv[]=0 ) const ;
279 Double_t GetDeviationFromVertex( const AliKFParticle &Vtx ) const ;
280 Double_t GetDeviationFromVertex( const AliVVertex &Vtx ) const ;
281 Double_t GetDeviationFromParticle( const AliKFParticle &p ) const ;
283 //* Calculate distance from another object [cm] in XY-plane
285 Bool_t GetDistanceFromVertexXY( const Double_t vtx[], Double_t &val, Double_t &err ) const ;
286 Bool_t GetDistanceFromVertexXY( const Double_t vtx[], const Double_t Cv[], Double_t &val, Double_t &err ) const ;
287 Bool_t GetDistanceFromVertexXY( const AliKFParticle &Vtx, Double_t &val, Double_t &err ) const ;
288 Bool_t GetDistanceFromVertexXY( const AliVVertex &Vtx, Double_t &val, Double_t &err ) const ;
290 Double_t GetDistanceFromVertexXY( const Double_t vtx[] ) const ;
291 Double_t GetDistanceFromVertexXY( const AliKFParticle &Vtx ) const ;
292 Double_t GetDistanceFromVertexXY( const AliVVertex &Vtx ) const ;
293 Double_t GetDistanceFromParticleXY( const AliKFParticle &p ) const ;
295 //* Calculate sqrt(Chi2/ndf) deviation from another object in XY plane
296 //* ( v = [xyz]-vertex, Cv=[Cxx,Cxy,Cyy,Cxz,Cyz,Czz]-covariance matrix )
298 Double_t GetDeviationFromVertexXY( const Double_t v[], const Double_t Cv[]=0 ) const ;
299 Double_t GetDeviationFromVertexXY( const AliKFParticle &Vtx ) const ;
300 Double_t GetDeviationFromVertexXY( const AliVVertex &Vtx ) const ;
301 Double_t GetDeviationFromParticleXY( const AliKFParticle &p ) const ;
303 //* Calculate opennig angle between two particles
305 Double_t GetAngle ( const AliKFParticle &p ) const ;
306 Double_t GetAngleXY( const AliKFParticle &p ) const ;
307 Double_t GetAngleRZ( const AliKFParticle &p ) const ;
309 //* Subtract the particle from the vertex
311 void SubtractFromVertex( AliKFParticle &v ) const ;
313 //* Special method for creating gammas
315 void ConstructGamma( const AliKFParticle &daughter1,
316 const AliKFParticle &daughter2 );
318 // * Pseudo Proper Time of decay = (r*pt) / |pt| * M/|pt|
319 // @primVertex - primary vertex
320 // @mass - mass of the mother particle (in the case of "Hb -> JPsi" it would be JPsi mass)
321 // @*timeErr2 - squared error of the decay time. If timeErr2 = 0 it isn't calculated
322 Double_t GetPseudoProperDecayTime( const AliKFParticle &primVertex, const Double_t& mass, Double_t* timeErr2 = 0 ) const;
330 //* Method to access ALICE field
332 static Double_t GetFieldAlice();
334 //* Other methods required by the abstract AliKFParticleBase class
336 void GetFieldValue( const Double_t xyz[], Double_t B[] ) const ;
337 void GetDStoParticle( const AliKFParticleBase &p, Double_t &DS, Double_t &DSp )const ;
338 void Transport( Double_t dS, Double_t P[], Double_t C[] ) const ;
339 static void GetExternalTrackParam( const AliKFParticleBase &p, Double_t &X, Double_t &Alpha, Double_t P[5] ) ;
341 //void GetDStoParticleALICE( const AliKFParticleBase &p, Double_t &DS, Double_t &DS1 ) const;
346 static Double_t fgBz; //* Bz compoment of the magnetic field
348 ClassDef( AliKFParticle, 1 );
354 //---------------------------------------------------------------------
356 // Inline implementation of the AliKFParticle methods
358 //---------------------------------------------------------------------
361 inline void AliKFParticle::SetField( Double_t Bz )
366 inline AliKFParticle::AliKFParticle( const AliKFParticle &d1,
367 const AliKFParticle &d2,
368 const AliKFParticle &d3 )
370 AliKFParticle mother;
377 inline AliKFParticle::AliKFParticle( const AliKFParticle &d1,
378 const AliKFParticle &d2,
379 const AliKFParticle &d3,
380 const AliKFParticle &d4 )
382 AliKFParticle mother;
391 inline void AliKFParticle::Initialize()
393 AliKFParticleBase::Initialize();
396 inline void AliKFParticle::SetVtxGuess( Double_t x, Double_t y, Double_t z )
398 AliKFParticleBase::SetVtxGuess(x,y,z);
401 inline Double_t AliKFParticle::GetX () const
403 return AliKFParticleBase::GetX();
406 inline Double_t AliKFParticle::GetY () const
408 return AliKFParticleBase::GetY();
411 inline Double_t AliKFParticle::GetZ () const
413 return AliKFParticleBase::GetZ();
416 inline Double_t AliKFParticle::GetPx () const
418 return AliKFParticleBase::GetPx();
421 inline Double_t AliKFParticle::GetPy () const
423 return AliKFParticleBase::GetPy();
426 inline Double_t AliKFParticle::GetPz () const
428 return AliKFParticleBase::GetPz();
431 inline Double_t AliKFParticle::GetE () const
433 return AliKFParticleBase::GetE();
436 inline Double_t AliKFParticle::GetS () const
438 return AliKFParticleBase::GetS();
441 inline Int_t AliKFParticle::GetQ () const
443 return AliKFParticleBase::GetQ();
446 inline Double_t AliKFParticle::GetChi2 () const
448 return AliKFParticleBase::GetChi2();
451 inline Int_t AliKFParticle::GetNDF () const
453 return AliKFParticleBase::GetNDF();
456 inline Double_t AliKFParticle::GetParameter ( int i ) const
458 return AliKFParticleBase::GetParameter(i);
461 inline Double_t AliKFParticle::GetCovariance( int i ) const
463 return AliKFParticleBase::GetCovariance(i);
466 inline Double_t AliKFParticle::GetCovariance( int i, int j ) const
468 return AliKFParticleBase::GetCovariance(i,j);
472 inline Double_t AliKFParticle::GetP () const
475 if( AliKFParticleBase::GetMomentum( par, err ) ) return 0;
479 inline Double_t AliKFParticle::GetPt () const
482 if( AliKFParticleBase::GetPt( par, err ) ) return 0;
486 inline Double_t AliKFParticle::GetEta () const
489 if( AliKFParticleBase::GetEta( par, err ) ) return 0;
493 inline Double_t AliKFParticle::GetPhi () const
496 if( AliKFParticleBase::GetPhi( par, err ) ) return 0;
500 inline Double_t AliKFParticle::GetMomentum () const
503 if( AliKFParticleBase::GetMomentum( par, err ) ) return 0;
507 inline Double_t AliKFParticle::GetMass () const
510 if( AliKFParticleBase::GetMass( par, err ) ) return 0;
514 inline Double_t AliKFParticle::GetDecayLength () const
517 if( AliKFParticleBase::GetDecayLength( par, err ) ) return 0;
521 inline Double_t AliKFParticle::GetDecayLengthXY () const
524 if( AliKFParticleBase::GetDecayLengthXY( par, err ) ) return 0;
528 inline Double_t AliKFParticle::GetLifeTime () const
531 if( AliKFParticleBase::GetLifeTime( par, err ) ) return 0;
535 inline Double_t AliKFParticle::GetR () const
538 if( AliKFParticleBase::GetR( par, err ) ) return 0;
542 inline Double_t AliKFParticle::GetErrX () const
544 return TMath::Sqrt(TMath::Abs( GetCovariance(0,0) ));
547 inline Double_t AliKFParticle::GetErrY () const
549 return TMath::Sqrt(TMath::Abs( GetCovariance(1,1) ));
552 inline Double_t AliKFParticle::GetErrZ () const
554 return TMath::Sqrt(TMath::Abs( GetCovariance(2,2) ));
557 inline Double_t AliKFParticle::GetErrPx () const
559 return TMath::Sqrt(TMath::Abs( GetCovariance(3,3) ));
562 inline Double_t AliKFParticle::GetErrPy () const
564 return TMath::Sqrt(TMath::Abs( GetCovariance(4,4) ));
567 inline Double_t AliKFParticle::GetErrPz () const
569 return TMath::Sqrt(TMath::Abs( GetCovariance(5,5) ));
572 inline Double_t AliKFParticle::GetErrE () const
574 return TMath::Sqrt(TMath::Abs( GetCovariance(6,6) ));
577 inline Double_t AliKFParticle::GetErrS () const
579 return TMath::Sqrt(TMath::Abs( GetCovariance(7,7) ));
582 inline Double_t AliKFParticle::GetErrP () const
585 if( AliKFParticleBase::GetMomentum( par, err ) ) return 1.e10;
589 inline Double_t AliKFParticle::GetErrPt () const
592 if( AliKFParticleBase::GetPt( par, err ) ) return 1.e10;
596 inline Double_t AliKFParticle::GetErrEta () const
599 if( AliKFParticleBase::GetEta( par, err ) ) return 1.e10;
603 inline Double_t AliKFParticle::GetErrPhi () const
606 if( AliKFParticleBase::GetPhi( par, err ) ) return 1.e10;
610 inline Double_t AliKFParticle::GetErrMomentum () const
613 if( AliKFParticleBase::GetMomentum( par, err ) ) return 1.e10;
617 inline Double_t AliKFParticle::GetErrMass () const
620 if( AliKFParticleBase::GetMass( par, err ) ) return 1.e10;
624 inline Double_t AliKFParticle::GetErrDecayLength () const
627 if( AliKFParticleBase::GetDecayLength( par, err ) ) return 1.e10;
631 inline Double_t AliKFParticle::GetErrDecayLengthXY () const
634 if( AliKFParticleBase::GetDecayLengthXY( par, err ) ) return 1.e10;
638 inline Double_t AliKFParticle::GetErrLifeTime () const
641 if( AliKFParticleBase::GetLifeTime( par, err ) ) return 1.e10;
645 inline Double_t AliKFParticle::GetErrR () const
648 if( AliKFParticleBase::GetR( par, err ) ) return 1.e10;
653 inline int AliKFParticle::GetP( Double_t &P, Double_t &SigmaP ) const
655 return AliKFParticleBase::GetMomentum( P, SigmaP );
658 inline int AliKFParticle::GetPt( Double_t &Pt, Double_t &SigmaPt ) const
660 return AliKFParticleBase::GetPt( Pt, SigmaPt );
663 inline int AliKFParticle::GetEta( Double_t &Eta, Double_t &SigmaEta ) const
665 return AliKFParticleBase::GetEta( Eta, SigmaEta );
668 inline int AliKFParticle::GetPhi( Double_t &Phi, Double_t &SigmaPhi ) const
670 return AliKFParticleBase::GetPhi( Phi, SigmaPhi );
673 inline int AliKFParticle::GetMomentum( Double_t &P, Double_t &SigmaP ) const
675 return AliKFParticleBase::GetMomentum( P, SigmaP );
678 inline int AliKFParticle::GetMass( Double_t &M, Double_t &SigmaM ) const
680 return AliKFParticleBase::GetMass( M, SigmaM );
683 inline int AliKFParticle::GetDecayLength( Double_t &L, Double_t &SigmaL ) const
685 return AliKFParticleBase::GetDecayLength( L, SigmaL );
688 inline int AliKFParticle::GetDecayLengthXY( Double_t &L, Double_t &SigmaL ) const
690 return AliKFParticleBase::GetDecayLengthXY( L, SigmaL );
693 inline int AliKFParticle::GetLifeTime( Double_t &T, Double_t &SigmaT ) const
695 return AliKFParticleBase::GetLifeTime( T, SigmaT );
698 inline int AliKFParticle::GetR( Double_t &R, Double_t &SigmaR ) const
700 return AliKFParticleBase::GetR( R, SigmaR );
703 inline Double_t & AliKFParticle::X()
705 return AliKFParticleBase::X();
708 inline Double_t & AliKFParticle::Y()
710 return AliKFParticleBase::Y();
713 inline Double_t & AliKFParticle::Z()
715 return AliKFParticleBase::Z();
718 inline Double_t & AliKFParticle::Px()
720 return AliKFParticleBase::Px();
723 inline Double_t & AliKFParticle::Py()
725 return AliKFParticleBase::Py();
728 inline Double_t & AliKFParticle::Pz()
730 return AliKFParticleBase::Pz();
733 inline Double_t & AliKFParticle::E()
735 return AliKFParticleBase::E();
738 inline Double_t & AliKFParticle::S()
740 return AliKFParticleBase::S();
743 inline Int_t & AliKFParticle::Q()
745 return AliKFParticleBase::Q();
748 inline Double_t & AliKFParticle::Chi2()
750 return AliKFParticleBase::Chi2();
753 inline Int_t & AliKFParticle::NDF()
755 return AliKFParticleBase::NDF();
758 inline Double_t & AliKFParticle::Parameter ( int i )
760 return AliKFParticleBase::Parameter(i);
763 inline Double_t & AliKFParticle::Covariance( int i )
765 return AliKFParticleBase::Covariance(i);
768 inline Double_t & AliKFParticle::Covariance( int i, int j )
770 return AliKFParticleBase::Covariance(i,j);
773 inline Double_t * AliKFParticle::Parameters ()
778 inline Double_t * AliKFParticle::CovarianceMatrix()
784 inline void AliKFParticle::operator +=( const AliKFParticle &Daughter )
786 AliKFParticleBase::operator +=( Daughter );
790 inline void AliKFParticle::AddDaughter( const AliKFParticle &Daughter )
792 AliKFParticleBase::AddDaughter( Daughter );
795 inline void AliKFParticle::SetProductionVertex( const AliKFParticle &Vtx )
797 AliKFParticleBase::SetProductionVertex( Vtx );
800 inline void AliKFParticle::SetMassConstraint( Double_t Mass, Double_t SigmaMass )
802 AliKFParticleBase::SetMassConstraint( Mass, SigmaMass );
805 inline void AliKFParticle::SetNoDecayLength()
807 AliKFParticleBase::SetNoDecayLength();
810 inline void AliKFParticle::Construct( const AliKFParticle *vDaughters[], int NDaughters,
811 const AliKFParticle *ProdVtx, Double_t Mass, Bool_t IsConstrained )
813 AliKFParticleBase::Construct( ( const AliKFParticleBase**)vDaughters, NDaughters,
814 ( const AliKFParticleBase*)ProdVtx, Mass, IsConstrained );
817 inline void AliKFParticle::TransportToDecayVertex()
819 AliKFParticleBase::TransportToDecayVertex();
822 inline void AliKFParticle::TransportToProductionVertex()
824 AliKFParticleBase::TransportToProductionVertex();
827 inline void AliKFParticle::TransportToPoint( const Double_t xyz[] )
829 TransportToDS( GetDStoPoint(xyz) );
832 inline void AliKFParticle::TransportToVertex( const AliVVertex &v )
834 TransportToPoint( AliKFParticle(v).fP );
837 inline void AliKFParticle::TransportToParticle( const AliKFParticle &p )
840 GetDStoParticle( p, dS, dSp );
844 inline void AliKFParticle::TransportToDS( Double_t dS )
846 AliKFParticleBase::TransportToDS( dS );
849 inline Double_t AliKFParticle::GetDStoPoint( const Double_t xyz[] ) const
851 return AliKFParticleBase::GetDStoPointBz( GetFieldAlice(), xyz );
855 inline void AliKFParticle::GetDStoParticle( const AliKFParticle &p,
856 Double_t &DS, Double_t &DSp ) const
858 GetDStoParticleXY( p, DS, DSp );
862 inline Double_t AliKFParticle::GetDistanceFromVertex( const Double_t vtx[] ) const
864 return AliKFParticleBase::GetDistanceFromVertex( vtx );
867 inline Double_t AliKFParticle::GetDeviationFromVertex( const Double_t v[],
868 const Double_t Cv[] ) const
870 return AliKFParticleBase::GetDeviationFromVertex( v, Cv);
873 inline Double_t AliKFParticle::GetDistanceFromVertex( const AliKFParticle &Vtx ) const
875 return AliKFParticleBase::GetDistanceFromVertex( Vtx );
878 inline Double_t AliKFParticle::GetDeviationFromVertex( const AliKFParticle &Vtx ) const
880 return AliKFParticleBase::GetDeviationFromVertex( Vtx );
883 inline Double_t AliKFParticle::GetDistanceFromVertex( const AliVVertex &Vtx ) const
885 return GetDistanceFromVertex( AliKFParticle(Vtx) );
888 inline Double_t AliKFParticle::GetDeviationFromVertex( const AliVVertex &Vtx ) const
890 return GetDeviationFromVertex( AliKFParticle(Vtx) );
893 inline Double_t AliKFParticle::GetDistanceFromParticle( const AliKFParticle &p ) const
895 return AliKFParticleBase::GetDistanceFromParticle( p );
898 inline Double_t AliKFParticle::GetDeviationFromParticle( const AliKFParticle &p ) const
900 return AliKFParticleBase::GetDeviationFromParticle( p );
903 inline void AliKFParticle::SubtractFromVertex( AliKFParticle &v ) const
905 AliKFParticleBase::SubtractFromVertex( v );
908 inline Double_t AliKFParticle::GetFieldAlice()
913 inline void AliKFParticle::GetFieldValue( const Double_t * /*xyz*/, Double_t B[] ) const
916 B[2] = GetFieldAlice();
919 inline void AliKFParticle::GetDStoParticle( const AliKFParticleBase &p,
920 Double_t &DS, Double_t &DSp )const
922 GetDStoParticleXY( p, DS, DSp );
925 inline void AliKFParticle::GetDStoParticleXY( const AliKFParticleBase &p,
926 Double_t &DS, Double_t &DSp ) const
928 AliKFParticleBase::GetDStoParticleBz( GetFieldAlice(), p, DS, DSp ) ;
929 //GetDStoParticleALICE( p, DS, DSp ) ;
932 inline void AliKFParticle::Transport( Double_t dS, Double_t P[], Double_t C[] ) const
934 AliKFParticleBase::TransportBz( GetFieldAlice(), dS, P, C );
937 inline void AliKFParticle::ConstructGamma( const AliKFParticle &daughter1,
938 const AliKFParticle &daughter2 )
940 AliKFParticleBase::ConstructGammaBz( daughter1, daughter2, GetFieldAlice() );