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"
27 class AliKFParticle :public AliKFParticleBase
36 //* Set magnetic field for all particles
38 static void SetField( Double_t Bz );
40 //* Constructor (empty)
42 AliKFParticle():AliKFParticleBase(){ ; }
44 //* Destructor (empty)
48 //* Construction of mother particle by its 2-3-4 daughters
50 AliKFParticle( const AliKFParticle &d1, const AliKFParticle &d2 );
52 AliKFParticle( const AliKFParticle &d1, const AliKFParticle &d2,
53 const AliKFParticle &d3 );
55 AliKFParticle( const AliKFParticle &d1, const AliKFParticle &d2,
56 const AliKFParticle &d3, const AliKFParticle &d4 );
58 //* Initialisation from "cartesian" coordinates ( X Y Z Px Py Pz )
59 //* Parameters, covariance matrix, charge and PID hypothesis should be provided
61 void Create( const Double_t Param[], const Double_t Cov[], Int_t Charge, Int_t PID );
63 //* Initialisation from ALICE track, PID hypothesis shoould be provided
65 AliKFParticle( const AliVTrack &track, Int_t PID );
67 //* Initialisation from VVertex
69 AliKFParticle( const AliVVertex &vertex );
71 //* Initialise covariance matrix and set current parameters to 0.0
75 //* Set decay vertex parameters for linearisation
77 void SetVtxGuess( Double_t x, Double_t y, Double_t z );
85 Double_t GetX () const ; //* x of current position
86 Double_t GetY () const ; //* y of current position
87 Double_t GetZ () const ; //* z of current position
88 Double_t GetPx () const ; //* x-compoment of 3-momentum
89 Double_t GetPy () const ; //* y-compoment of 3-momentum
90 Double_t GetPz () const ; //* z-compoment of 3-momentum
91 Double_t GetE () const ; //* energy
92 Double_t GetS () const ; //* decay length / momentum
93 Int_t GetQ () const ; //* charge
94 Double_t GetChi2 () const ; //* chi^2
95 Int_t GetNDF () const ; //* Number of Degrees of Freedom
97 Double_t GetParameter ( int i ) const ;
98 Double_t GetCovariance( int i ) const ;
99 Double_t GetCovariance( int i, int j ) const ;
101 //* Accessors with calculations, value returned w/o error flag
103 Double_t GetP () const; //* momentum
104 Double_t GetPt () const; //* transverse momentum
105 Double_t GetEta () const; //* pseudorapidity
106 Double_t GetPhi () const; //* phi
107 Double_t GetMomentum () const; //* momentum (same as GetP() )
108 Double_t GetMass () const; //* mass
109 Double_t GetDecayLength () const; //* decay length
110 Double_t GetDecayLengthXY () const; //* decay length in XY
111 Double_t GetLifeTime () const; //* life time
112 Double_t GetR () const; //* distance to the origin
114 //* Accessors to estimated errors
116 Double_t GetErrX () const ; //* x of current position
117 Double_t GetErrY () const ; //* y of current position
118 Double_t GetErrZ () const ; //* z of current position
119 Double_t GetErrPx () const ; //* x-compoment of 3-momentum
120 Double_t GetErrPy () const ; //* y-compoment of 3-momentum
121 Double_t GetErrPz () const ; //* z-compoment of 3-momentum
122 Double_t GetErrE () const ; //* energy
123 Double_t GetErrS () const ; //* decay length / momentum
124 Double_t GetErrP () const ; //* momentum
125 Double_t GetErrPt () const ; //* transverse momentum
126 Double_t GetErrEta () const ; //* pseudorapidity
127 Double_t GetErrPhi () const ; //* phi
128 Double_t GetErrMomentum () const ; //* momentum
129 Double_t GetErrMass () const ; //* mass
130 Double_t GetErrDecayLength () const ; //* decay length
131 Double_t GetErrDecayLengthXY () const ; //* decay length in XY
132 Double_t GetErrLifeTime () const ; //* life time
133 Double_t GetErrR () const ; //* distance to the origin
135 //* Accessors with calculations( &value, &estimated sigma )
136 //* error flag returned (0 means no error during calculations)
138 int GetP ( Double_t &P, Double_t &SigmaP ) const ; //* momentum
139 int GetPt ( Double_t &Pt, Double_t &SigmaPt ) const ; //* transverse momentum
140 int GetEta ( Double_t &Eta, Double_t &SigmaEta ) const ; //* pseudorapidity
141 int GetPhi ( Double_t &Phi, Double_t &SigmaPhi ) const ; //* phi
142 int GetMomentum ( Double_t &P, Double_t &SigmaP ) const ; //* momentum
143 int GetMass ( Double_t &M, Double_t &SigmaM ) const ; //* mass
144 int GetDecayLength ( Double_t &L, Double_t &SigmaL ) const ; //* decay length
145 int GetDecayLengthXY ( Double_t &L, Double_t &SigmaL ) const ; //* decay length in XY
146 int GetLifeTime ( Double_t &T, Double_t &SigmaT ) const ; //* life time
147 int GetR ( Double_t &R, Double_t &SigmaR ) const ; //* R
166 Double_t & Parameter ( int i ) ;
167 Double_t & Covariance( int i ) ;
168 Double_t & Covariance( int i, int j ) ;
169 Double_t * Parameters () ;
170 Double_t * CovarianceMatrix() ;
173 //* CONSTRUCTION OF THE PARTICLE BY ITS DAUGHTERS AND MOTHER
174 //* USING THE KALMAN FILTER METHOD
178 //* Add daughter to the particle
180 void AddDaughter( const AliKFParticle &Daughter );
182 //* Add daughter via += operator: ex.{ D0; D0+=Pion; D0+= Kaon; }
184 void operator +=( const AliKFParticle &Daughter );
186 //* Set production vertex
188 void SetProductionVertex( const AliKFParticle &Vtx );
190 //* Set mass constraint
192 void SetMassConstraint( Double_t Mass, Double_t SigmaMass = 0 );
194 //* Set no decay length for resonances
196 void SetNoDecayLength();
198 //* Everything in one go
200 void Construct( const AliKFParticle *vDaughters[], int NDaughters,
201 const AliKFParticle *ProdVtx=0, Double_t Mass=-1, Bool_t IsConstrained=0 );
206 //* ( main transportation parameter is S = SignedPath/Momentum )
207 //* ( parameters of decay & production vertices are stored locally )
210 //* Transport the particle to its decay vertex
212 void TransportToDecayVertex();
214 //* Transport the particle to its production vertex
216 void TransportToProductionVertex();
218 //* Transport the particle close to xyz[] point
220 void TransportToPoint( const Double_t xyz[] );
222 //* Transport the particle close to VVertex
224 void TransportToVertex( const AliVVertex &v );
226 //* Transport the particle close to another particle p
228 void TransportToParticle( const AliKFParticle &p );
230 //* Transport the particle on dS parameter (SignedPath/Momentum)
232 void TransportToDS( Double_t dS );
234 //* Get dS to a certain space point
236 Double_t GetDStoPoint( const Double_t xyz[] ) const ;
238 //* Get dS to other particle p (dSp for particle p also returned)
240 void GetDStoParticle( const AliKFParticle &p,
241 Double_t &DS, Double_t &DSp ) const ;
243 //* Get dS to other particle p in XY-plane
245 void GetDStoParticleXY( const AliKFParticleBase &p,
246 Double_t &DS, Double_t &DSp ) const ;
253 //* Calculate distance from another object [cm]
255 Double_t GetDistanceFromVertex( const Double_t vtx[] ) const ;
256 Double_t GetDistanceFromVertex( const AliKFParticle &Vtx ) const ;
257 Double_t GetDistanceFromVertex( const AliVVertex &Vtx ) const ;
258 Double_t GetDistanceFromParticle( const AliKFParticle &p ) const ;
260 //* Calculate sqrt(Chi2/ndf) deviation from another object
261 //* ( v = [xyz]-vertex, Cv=[Cxx,Cxy,Cyy,Cxz,Cyz,Czz]-covariance matrix )
263 Double_t GetDeviationFromVertex( const Double_t v[], const Double_t Cv[]=0 ) const ;
264 Double_t GetDeviationFromVertex( const AliKFParticle &Vtx ) const ;
265 Double_t GetDeviationFromVertex( const AliVVertex &Vtx ) const ;
266 Double_t GetDeviationFromParticle( const AliKFParticle &p ) const ;
268 //* Calculate distance from another object [cm] in XY-plane
270 Bool_t GetDistanceFromVertexXY( const Double_t vtx[], Double_t &val, Double_t &err ) const ;
271 Bool_t GetDistanceFromVertexXY( const Double_t vtx[], const Double_t Cv[], Double_t &val, Double_t &err ) const ;
272 Bool_t GetDistanceFromVertexXY( const AliKFParticle &Vtx, Double_t &val, Double_t &err ) const ;
273 Bool_t GetDistanceFromVertexXY( const AliVVertex &Vtx, Double_t &val, Double_t &err ) const ;
275 Double_t GetDistanceFromVertexXY( const Double_t vtx[] ) const ;
276 Double_t GetDistanceFromVertexXY( const AliKFParticle &Vtx ) const ;
277 Double_t GetDistanceFromVertexXY( const AliVVertex &Vtx ) const ;
278 Double_t GetDistanceFromParticleXY( const AliKFParticle &p ) const ;
280 //* Calculate sqrt(Chi2/ndf) deviation from another object in XY plane
281 //* ( v = [xyz]-vertex, Cv=[Cxx,Cxy,Cyy,Cxz,Cyz,Czz]-covariance matrix )
283 Double_t GetDeviationFromVertexXY( const Double_t v[], const Double_t Cv[]=0 ) const ;
284 Double_t GetDeviationFromVertexXY( const AliKFParticle &Vtx ) const ;
285 Double_t GetDeviationFromVertexXY( const AliVVertex &Vtx ) const ;
286 Double_t GetDeviationFromParticleXY( const AliKFParticle &p ) const ;
288 //* Calculate opennig angle between two particles
290 Double_t GetAngle ( const AliKFParticle &p ) const ;
291 Double_t GetAngleXY( const AliKFParticle &p ) const ;
292 Double_t GetAngleRZ( const AliKFParticle &p ) const ;
294 //* Subtract the particle from the vertex
296 void SubtractFromVertex( AliKFParticle &v ) const ;
298 //* Special method for creating gammas
300 void ConstructGamma( const AliKFParticle &daughter1,
301 const AliKFParticle &daughter2 );
309 //* Method to access ALICE field
311 static Double_t GetFieldAlice();
313 //* Other methods required by the abstract AliKFParticleBase class
315 void GetFieldValue( const Double_t xyz[], Double_t B[] ) const ;
316 void GetDStoParticle( const AliKFParticleBase &p, Double_t &DS, Double_t &DSp )const ;
317 void Transport( Double_t dS, Double_t P[], Double_t C[] ) const ;
318 static void GetExternalTrackParam( const AliKFParticleBase &p, Double_t &X, Double_t &Alpha, Double_t P[5] ) ;
320 //void GetDStoParticleALICE( const AliKFParticleBase &p, Double_t &DS, Double_t &DS1 ) const;
325 static Double_t fgBz; //* Bz compoment of the magnetic field
327 ClassDef( AliKFParticle, 1 );
333 //---------------------------------------------------------------------
335 // Inline implementation of the AliKFParticle methods
337 //---------------------------------------------------------------------
340 inline void AliKFParticle::SetField( Double_t Bz )
346 inline AliKFParticle::AliKFParticle( const AliKFParticle &d1,
347 const AliKFParticle &d2 )
349 AliKFParticle mother;
355 inline AliKFParticle::AliKFParticle( const AliKFParticle &d1,
356 const AliKFParticle &d2,
357 const AliKFParticle &d3 )
359 AliKFParticle mother;
366 inline AliKFParticle::AliKFParticle( const AliKFParticle &d1,
367 const AliKFParticle &d2,
368 const AliKFParticle &d3,
369 const AliKFParticle &d4 )
371 AliKFParticle mother;
380 inline void AliKFParticle::Initialize()
382 AliKFParticleBase::Initialize();
385 inline void AliKFParticle::SetVtxGuess( Double_t x, Double_t y, Double_t z )
387 AliKFParticleBase::SetVtxGuess(x,y,z);
390 inline Double_t AliKFParticle::GetX () const
392 return AliKFParticleBase::GetX();
395 inline Double_t AliKFParticle::GetY () const
397 return AliKFParticleBase::GetY();
400 inline Double_t AliKFParticle::GetZ () const
402 return AliKFParticleBase::GetZ();
405 inline Double_t AliKFParticle::GetPx () const
407 return AliKFParticleBase::GetPx();
410 inline Double_t AliKFParticle::GetPy () const
412 return AliKFParticleBase::GetPy();
415 inline Double_t AliKFParticle::GetPz () const
417 return AliKFParticleBase::GetPz();
420 inline Double_t AliKFParticle::GetE () const
422 return AliKFParticleBase::GetE();
425 inline Double_t AliKFParticle::GetS () const
427 return AliKFParticleBase::GetS();
430 inline Int_t AliKFParticle::GetQ () const
432 return AliKFParticleBase::GetQ();
435 inline Double_t AliKFParticle::GetChi2 () const
437 return AliKFParticleBase::GetChi2();
440 inline Int_t AliKFParticle::GetNDF () const
442 return AliKFParticleBase::GetNDF();
445 inline Double_t AliKFParticle::GetParameter ( int i ) const
447 return AliKFParticleBase::GetParameter(i);
450 inline Double_t AliKFParticle::GetCovariance( int i ) const
452 return AliKFParticleBase::GetCovariance(i);
455 inline Double_t AliKFParticle::GetCovariance( int i, int j ) const
457 return AliKFParticleBase::GetCovariance(i,j);
461 inline Double_t AliKFParticle::GetP () const
464 if( AliKFParticleBase::GetMomentum( par, err ) ) return 0;
468 inline Double_t AliKFParticle::GetPt () const
471 if( AliKFParticleBase::GetPt( par, err ) ) return 0;
475 inline Double_t AliKFParticle::GetEta () const
478 if( AliKFParticleBase::GetEta( par, err ) ) return 0;
482 inline Double_t AliKFParticle::GetPhi () const
485 if( AliKFParticleBase::GetPhi( par, err ) ) return 0;
489 inline Double_t AliKFParticle::GetMomentum () const
492 if( AliKFParticleBase::GetMomentum( par, err ) ) return 0;
496 inline Double_t AliKFParticle::GetMass () const
499 if( AliKFParticleBase::GetMass( par, err ) ) return 0;
503 inline Double_t AliKFParticle::GetDecayLength () const
506 if( AliKFParticleBase::GetDecayLength( par, err ) ) return 0;
510 inline Double_t AliKFParticle::GetDecayLengthXY () const
513 if( AliKFParticleBase::GetDecayLengthXY( par, err ) ) return 0;
517 inline Double_t AliKFParticle::GetLifeTime () const
520 if( AliKFParticleBase::GetLifeTime( par, err ) ) return 0;
524 inline Double_t AliKFParticle::GetR () const
527 if( AliKFParticleBase::GetR( par, err ) ) return 0;
531 inline Double_t AliKFParticle::GetErrX () const
533 return TMath::Sqrt(TMath::Abs( GetCovariance(0,0) ));
536 inline Double_t AliKFParticle::GetErrY () const
538 return TMath::Sqrt(TMath::Abs( GetCovariance(1,1) ));
541 inline Double_t AliKFParticle::GetErrZ () const
543 return TMath::Sqrt(TMath::Abs( GetCovariance(2,2) ));
546 inline Double_t AliKFParticle::GetErrPx () const
548 return TMath::Sqrt(TMath::Abs( GetCovariance(3,3) ));
551 inline Double_t AliKFParticle::GetErrPy () const
553 return TMath::Sqrt(TMath::Abs( GetCovariance(4,4) ));
556 inline Double_t AliKFParticle::GetErrPz () const
558 return TMath::Sqrt(TMath::Abs( GetCovariance(5,5) ));
561 inline Double_t AliKFParticle::GetErrE () const
563 return TMath::Sqrt(TMath::Abs( GetCovariance(6,6) ));
566 inline Double_t AliKFParticle::GetErrS () const
568 return TMath::Sqrt(TMath::Abs( GetCovariance(7,7) ));
571 inline Double_t AliKFParticle::GetErrP () const
574 if( AliKFParticleBase::GetMomentum( par, err ) ) return 1.e10;
578 inline Double_t AliKFParticle::GetErrPt () const
581 if( AliKFParticleBase::GetPt( par, err ) ) return 1.e10;
585 inline Double_t AliKFParticle::GetErrEta () const
588 if( AliKFParticleBase::GetEta( par, err ) ) return 1.e10;
592 inline Double_t AliKFParticle::GetErrPhi () const
595 if( AliKFParticleBase::GetPhi( par, err ) ) return 1.e10;
599 inline Double_t AliKFParticle::GetErrMomentum () const
602 if( AliKFParticleBase::GetMomentum( par, err ) ) return 1.e10;
606 inline Double_t AliKFParticle::GetErrMass () const
609 if( AliKFParticleBase::GetMass( par, err ) ) return 1.e10;
613 inline Double_t AliKFParticle::GetErrDecayLength () const
616 if( AliKFParticleBase::GetDecayLength( par, err ) ) return 1.e10;
620 inline Double_t AliKFParticle::GetErrDecayLengthXY () const
623 if( AliKFParticleBase::GetDecayLengthXY( par, err ) ) return 1.e10;
627 inline Double_t AliKFParticle::GetErrLifeTime () const
630 if( AliKFParticleBase::GetLifeTime( par, err ) ) return 1.e10;
634 inline Double_t AliKFParticle::GetErrR () const
637 if( AliKFParticleBase::GetR( par, err ) ) return 1.e10;
642 inline int AliKFParticle::GetP( Double_t &P, Double_t &SigmaP ) const
644 return AliKFParticleBase::GetMomentum( P, SigmaP );
647 inline int AliKFParticle::GetPt( Double_t &Pt, Double_t &SigmaPt ) const
649 return AliKFParticleBase::GetPt( Pt, SigmaPt );
652 inline int AliKFParticle::GetEta( Double_t &Eta, Double_t &SigmaEta ) const
654 return AliKFParticleBase::GetEta( Eta, SigmaEta );
657 inline int AliKFParticle::GetPhi( Double_t &Phi, Double_t &SigmaPhi ) const
659 return AliKFParticleBase::GetPhi( Phi, SigmaPhi );
662 inline int AliKFParticle::GetMomentum( Double_t &P, Double_t &SigmaP ) const
664 return AliKFParticleBase::GetMomentum( P, SigmaP );
667 inline int AliKFParticle::GetMass( Double_t &M, Double_t &SigmaM ) const
669 return AliKFParticleBase::GetMass( M, SigmaM );
672 inline int AliKFParticle::GetDecayLength( Double_t &L, Double_t &SigmaL ) const
674 return AliKFParticleBase::GetDecayLength( L, SigmaL );
677 inline int AliKFParticle::GetDecayLengthXY( Double_t &L, Double_t &SigmaL ) const
679 return AliKFParticleBase::GetDecayLengthXY( L, SigmaL );
682 inline int AliKFParticle::GetLifeTime( Double_t &T, Double_t &SigmaT ) const
684 return AliKFParticleBase::GetLifeTime( T, SigmaT );
687 inline int AliKFParticle::GetR( Double_t &R, Double_t &SigmaR ) const
689 return AliKFParticleBase::GetR( R, SigmaR );
692 inline Double_t & AliKFParticle::X()
694 return AliKFParticleBase::X();
697 inline Double_t & AliKFParticle::Y()
699 return AliKFParticleBase::Y();
702 inline Double_t & AliKFParticle::Z()
704 return AliKFParticleBase::Z();
707 inline Double_t & AliKFParticle::Px()
709 return AliKFParticleBase::Px();
712 inline Double_t & AliKFParticle::Py()
714 return AliKFParticleBase::Py();
717 inline Double_t & AliKFParticle::Pz()
719 return AliKFParticleBase::Pz();
722 inline Double_t & AliKFParticle::E()
724 return AliKFParticleBase::E();
727 inline Double_t & AliKFParticle::S()
729 return AliKFParticleBase::S();
732 inline Int_t & AliKFParticle::Q()
734 return AliKFParticleBase::Q();
737 inline Double_t & AliKFParticle::Chi2()
739 return AliKFParticleBase::Chi2();
742 inline Int_t & AliKFParticle::NDF()
744 return AliKFParticleBase::NDF();
747 inline Double_t & AliKFParticle::Parameter ( int i )
749 return AliKFParticleBase::Parameter(i);
752 inline Double_t & AliKFParticle::Covariance( int i )
754 return AliKFParticleBase::Covariance(i);
757 inline Double_t & AliKFParticle::Covariance( int i, int j )
759 return AliKFParticleBase::Covariance(i,j);
762 inline Double_t * AliKFParticle::Parameters ()
767 inline Double_t * AliKFParticle::CovarianceMatrix()
773 inline void AliKFParticle::operator +=( const AliKFParticle &Daughter )
775 AliKFParticleBase::operator +=( Daughter );
779 inline void AliKFParticle::AddDaughter( const AliKFParticle &Daughter )
781 AliKFParticleBase::AddDaughter( Daughter );
784 inline void AliKFParticle::SetProductionVertex( const AliKFParticle &Vtx )
786 AliKFParticleBase::SetProductionVertex( Vtx );
789 inline void AliKFParticle::SetMassConstraint( Double_t Mass, Double_t SigmaMass )
791 AliKFParticleBase::SetMassConstraint( Mass, SigmaMass );
794 inline void AliKFParticle::SetNoDecayLength()
796 AliKFParticleBase::SetNoDecayLength();
799 inline void AliKFParticle::Construct( const AliKFParticle *vDaughters[], int NDaughters,
800 const AliKFParticle *ProdVtx, Double_t Mass, Bool_t IsConstrained )
802 AliKFParticleBase::Construct( ( const AliKFParticleBase**)vDaughters, NDaughters,
803 ( const AliKFParticleBase*)ProdVtx, Mass, IsConstrained );
806 inline void AliKFParticle::TransportToDecayVertex()
808 AliKFParticleBase::TransportToDecayVertex();
811 inline void AliKFParticle::TransportToProductionVertex()
813 AliKFParticleBase::TransportToProductionVertex();
816 inline void AliKFParticle::TransportToPoint( const Double_t xyz[] )
818 TransportToDS( GetDStoPoint(xyz) );
821 inline void AliKFParticle::TransportToVertex( const AliVVertex &v )
823 TransportToPoint( AliKFParticle(v).fP );
826 inline void AliKFParticle::TransportToParticle( const AliKFParticle &p )
829 GetDStoParticle( p, dS, dSp );
833 inline void AliKFParticle::TransportToDS( Double_t dS )
835 AliKFParticleBase::TransportToDS( dS );
838 inline Double_t AliKFParticle::GetDStoPoint( const Double_t xyz[] ) const
840 return AliKFParticleBase::GetDStoPointBz( GetFieldAlice(), xyz );
844 inline void AliKFParticle::GetDStoParticle( const AliKFParticle &p,
845 Double_t &DS, Double_t &DSp ) const
847 GetDStoParticleXY( p, DS, DSp );
851 inline Double_t AliKFParticle::GetDistanceFromVertex( const Double_t vtx[] ) const
853 return AliKFParticleBase::GetDistanceFromVertex( vtx );
856 inline Double_t AliKFParticle::GetDeviationFromVertex( const Double_t v[],
857 const Double_t Cv[] ) const
859 return AliKFParticleBase::GetDeviationFromVertex( v, Cv);
862 inline Double_t AliKFParticle::GetDistanceFromVertex( const AliKFParticle &Vtx ) const
864 return AliKFParticleBase::GetDistanceFromVertex( Vtx );
867 inline Double_t AliKFParticle::GetDeviationFromVertex( const AliKFParticle &Vtx ) const
869 return AliKFParticleBase::GetDeviationFromVertex( Vtx );
872 inline Double_t AliKFParticle::GetDistanceFromVertex( const AliVVertex &Vtx ) const
874 return GetDistanceFromVertex( AliKFParticle(Vtx) );
877 inline Double_t AliKFParticle::GetDeviationFromVertex( const AliVVertex &Vtx ) const
879 return GetDeviationFromVertex( AliKFParticle(Vtx) );
882 inline Double_t AliKFParticle::GetDistanceFromParticle( const AliKFParticle &p ) const
884 return AliKFParticleBase::GetDistanceFromParticle( p );
887 inline Double_t AliKFParticle::GetDeviationFromParticle( const AliKFParticle &p ) const
889 return AliKFParticleBase::GetDeviationFromParticle( p );
892 inline void AliKFParticle::SubtractFromVertex( AliKFParticle &v ) const
894 AliKFParticleBase::SubtractFromVertex( v );
897 inline Double_t AliKFParticle::GetFieldAlice()
902 inline void AliKFParticle::GetFieldValue( const Double_t * /*xyz*/, Double_t B[] ) const
905 B[2] = GetFieldAlice();
908 inline void AliKFParticle::GetDStoParticle( const AliKFParticleBase &p,
909 Double_t &DS, Double_t &DSp )const
911 GetDStoParticleXY( p, DS, DSp );
914 inline void AliKFParticle::GetDStoParticleXY( const AliKFParticleBase &p,
915 Double_t &DS, Double_t &DSp ) const
917 AliKFParticleBase::GetDStoParticleBz( GetFieldAlice(), p, DS, DSp ) ;
918 //GetDStoParticleALICE( p, DS, DSp ) ;
921 inline void AliKFParticle::Transport( Double_t dS, Double_t P[], Double_t C[] ) const
923 AliKFParticleBase::TransportBz( GetFieldAlice(), dS, P, C );
926 inline void AliKFParticle::ConstructGamma( const AliKFParticle &daughter1,
927 const AliKFParticle &daughter2 )
929 AliKFParticleBase::ConstructGammaBz( daughter1, daughter2, GetFieldAlice() );