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"
22 #include "AliESDVertex.h"
26 class AliKFParticle :public AliKFParticleBase
35 //* Set magnetic field for all particles
37 static void SetField( Double_t Bz );
39 //* Constructor (empty)
41 AliKFParticle():AliKFParticleBase(){ ; }
43 //* Destructor (empty)
47 //* Construction of mother particle by its 2-3-4 daughters
49 AliKFParticle( const AliKFParticle &d1, const AliKFParticle &d2 );
51 AliKFParticle( const AliKFParticle &d1, const AliKFParticle &d2,
52 const AliKFParticle &d3 );
54 AliKFParticle( const AliKFParticle &d1, const AliKFParticle &d2,
55 const AliKFParticle &d3, const AliKFParticle &d4 );
57 //* Initialisation from "cartesian" coordinates ( X Y Z Px Py Pz )
58 //* Parameters, covariance matrix, charge and PID hypothesis should be provided
60 void Create( const Double_t Param[], const Double_t Cov[], Int_t Charge, Int_t PID );
62 //* Initialisation from ALICE track, PID hypothesis shoould be provided
64 AliKFParticle( const AliVTrack &track, Int_t PID );
66 //* Initialisation from VVertex
68 AliKFParticle( const AliVVertex &vertex );
70 //* Copy position part to ESD vertex
72 void CopyToESDVertex( AliESDVertex &Vtx ) const ;
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 Double_t GetParameter ( int i ) const ;
101 Double_t GetCovariance( int i ) const ;
102 Double_t GetCovariance( int i, int j ) const ;
104 //* Accessors with calculations, value returned w/o error flag
106 Double_t GetP () const; //* momentum
107 Double_t GetPt () const; //* transverse momentum
108 Double_t GetEta () const; //* pseudorapidity
109 Double_t GetPhi () const; //* phi
110 Double_t GetMomentum () const; //* momentum (same as GetP() )
111 Double_t GetMass () const; //* mass
112 Double_t GetDecayLength () const; //* decay length
113 Double_t GetLifeTime () const; //* life time
114 Double_t GetR () const; //* distance to the origin
116 //* Accessors to estimated errors
118 Double_t GetErrX () const ; //* x of current position
119 Double_t GetErrY () const ; //* y of current position
120 Double_t GetErrZ () const ; //* z of current position
121 Double_t GetErrPx () const ; //* x-compoment of 3-momentum
122 Double_t GetErrPy () const ; //* y-compoment of 3-momentum
123 Double_t GetErrPz () const ; //* z-compoment of 3-momentum
124 Double_t GetErrE () const ; //* energy
125 Double_t GetErrS () const ; //* decay length / momentum
126 Double_t GetErrP () const ; //* momentum
127 Double_t GetErrPt () const ; //* transverse momentum
128 Double_t GetErrEta () const ; //* pseudorapidity
129 Double_t GetErrPhi () const ; //* phi
130 Double_t GetErrMomentum () const ; //* momentum
131 Double_t GetErrMass () const ; //* mass
132 Double_t GetErrDecayLength () const ; //* decay length
133 Double_t GetErrLifeTime () const ; //* life time
134 Double_t GetErrR () const ; //* distance to the origin
136 //* Accessors with calculations( &value, &estimated sigma )
137 //* error flag returned (0 means no error during calculations)
139 int GetP ( Double_t &P, Double_t &SigmaP ) const ; //* momentum
140 int GetPt ( Double_t &Pt, Double_t &SigmaPt ) const ; //* transverse momentum
141 int GetEta ( Double_t &Eta, Double_t &SigmaEta ) const ; //* pseudorapidity
142 int GetPhi ( Double_t &Phi, Double_t &SigmaPhi ) const ; //* phi
143 int GetMomentum ( Double_t &P, Double_t &SigmaP ) const ; //* momentum
144 int GetMass ( Double_t &M, Double_t &SigmaM ) const ; //* mass
145 int GetDecayLength ( Double_t &L, Double_t &SigmaL ) const ; //* decay length
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 Double_t GetDistanceFromVertexXY( const Double_t vtx[] ) const ;
271 Double_t GetDistanceFromVertexXY( const AliKFParticle &Vtx ) const ;
272 Double_t GetDistanceFromVertexXY( const AliVVertex &Vtx ) const ;
273 Double_t GetDistanceFromParticleXY( const AliKFParticle &p ) const ;
275 //* Calculate sqrt(Chi2/ndf) deviation from another object in XY plane
276 //* ( v = [xyz]-vertex, Cv=[Cxx,Cxy,Cyy,Cxz,Cyz,Czz]-covariance matrix )
278 Double_t GetDeviationFromVertexXY( const Double_t v[], const Double_t Cv[]=0 ) const ;
279 Double_t GetDeviationFromVertexXY( const AliKFParticle &Vtx ) const ;
280 Double_t GetDeviationFromVertexXY( const AliVVertex &Vtx ) const ;
281 Double_t GetDeviationFromParticleXY( const AliKFParticle &p ) const ;
283 //* Calculate opennig angle between two particles
285 Double_t GetAngle ( const AliKFParticle &p ) const ;
286 Double_t GetAngleXY( const AliKFParticle &p ) const ;
287 Double_t GetAngleRZ( const AliKFParticle &p ) const ;
289 //* Subtract the particle from the vertex
291 void SubtractFromVertex( AliKFParticle &v ) const ;
292 void SubtractFromVertex( AliESDVertex &v ) const ;
300 //* Method to access ALICE field
302 static Double_t GetFieldAlice();
304 //* Other methods required by the abstract AliKFParticleBase class
306 void GetFieldValue( const Double_t xyz[], Double_t B[] ) const ;
307 void GetDStoParticle( const AliKFParticleBase &p, Double_t &DS, Double_t &DSp )const ;
308 void Transport( Double_t dS, Double_t P[], Double_t C[] ) const ;
309 static void GetExternalTrackParam( const AliKFParticleBase &p, Double_t &X, Double_t &Alpha, Double_t P[5] ) ;
311 //void GetDStoParticleALICE( const AliKFParticleBase &p, Double_t &DS, Double_t &DS1 ) const;
316 static Double_t fgBz; //* Bz compoment of the magnetic field
318 ClassDef( AliKFParticle, 1 );
324 //---------------------------------------------------------------------
326 // Inline implementation of the AliKFParticle methods
328 //---------------------------------------------------------------------
331 inline void AliKFParticle::SetField( Double_t Bz )
337 inline AliKFParticle::AliKFParticle( const AliKFParticle &d1,
338 const AliKFParticle &d2 )
340 AliKFParticle mother;
346 inline AliKFParticle::AliKFParticle( const AliKFParticle &d1,
347 const AliKFParticle &d2,
348 const AliKFParticle &d3 )
350 AliKFParticle mother;
357 inline AliKFParticle::AliKFParticle( const AliKFParticle &d1,
358 const AliKFParticle &d2,
359 const AliKFParticle &d3,
360 const AliKFParticle &d4 )
362 AliKFParticle mother;
371 inline void AliKFParticle::Initialize()
373 AliKFParticleBase::Initialize();
376 inline void AliKFParticle::SetVtxGuess( Double_t x, Double_t y, Double_t z )
378 AliKFParticleBase::SetVtxGuess(x,y,z);
381 inline Double_t AliKFParticle::GetX () const
383 return AliKFParticleBase::GetX();
386 inline Double_t AliKFParticle::GetY () const
388 return AliKFParticleBase::GetY();
391 inline Double_t AliKFParticle::GetZ () const
393 return AliKFParticleBase::GetZ();
396 inline Double_t AliKFParticle::GetPx () const
398 return AliKFParticleBase::GetPx();
401 inline Double_t AliKFParticle::GetPy () const
403 return AliKFParticleBase::GetPy();
406 inline Double_t AliKFParticle::GetPz () const
408 return AliKFParticleBase::GetPz();
411 inline Double_t AliKFParticle::GetE () const
413 return AliKFParticleBase::GetE();
416 inline Double_t AliKFParticle::GetS () const
418 return AliKFParticleBase::GetS();
421 inline Int_t AliKFParticle::GetQ () const
423 return AliKFParticleBase::GetQ();
426 inline Double_t AliKFParticle::GetChi2 () const
428 return AliKFParticleBase::GetChi2();
431 inline Int_t AliKFParticle::GetNDF () const
433 return AliKFParticleBase::GetNDF();
436 inline Double_t AliKFParticle::GetParameter ( int i ) const
438 return AliKFParticleBase::GetParameter(i);
441 inline Double_t AliKFParticle::GetCovariance( int i ) const
443 return AliKFParticleBase::GetCovariance(i);
446 inline Double_t AliKFParticle::GetCovariance( int i, int j ) const
448 return AliKFParticleBase::GetCovariance(i,j);
452 inline Double_t AliKFParticle::GetP () const
455 if( AliKFParticleBase::GetMomentum( par, err ) ) return 0;
459 inline Double_t AliKFParticle::GetPt () const
462 if( AliKFParticleBase::GetPt( par, err ) ) return 0;
466 inline Double_t AliKFParticle::GetEta () const
469 if( AliKFParticleBase::GetEta( par, err ) ) return 0;
473 inline Double_t AliKFParticle::GetPhi () const
476 if( AliKFParticleBase::GetPhi( par, err ) ) return 0;
480 inline Double_t AliKFParticle::GetMomentum () const
483 if( AliKFParticleBase::GetMomentum( par, err ) ) return 0;
487 inline Double_t AliKFParticle::GetMass () const
490 if( AliKFParticleBase::GetMass( par, err ) ) return 0;
494 inline Double_t AliKFParticle::GetDecayLength () const
497 if( AliKFParticleBase::GetDecayLength( par, err ) ) return 0;
501 inline Double_t AliKFParticle::GetLifeTime () const
504 if( AliKFParticleBase::GetLifeTime( par, err ) ) return 0;
508 inline Double_t AliKFParticle::GetR () const
511 if( AliKFParticleBase::GetR( par, err ) ) return 0;
515 inline Double_t AliKFParticle::GetErrX () const
517 return TMath::Sqrt(TMath::Abs( GetCovariance(0,0) ));
520 inline Double_t AliKFParticle::GetErrY () const
522 return TMath::Sqrt(TMath::Abs( GetCovariance(1,1) ));
525 inline Double_t AliKFParticle::GetErrZ () const
527 return TMath::Sqrt(TMath::Abs( GetCovariance(2,2) ));
530 inline Double_t AliKFParticle::GetErrPx () const
532 return TMath::Sqrt(TMath::Abs( GetCovariance(3,3) ));
535 inline Double_t AliKFParticle::GetErrPy () const
537 return TMath::Sqrt(TMath::Abs( GetCovariance(4,4) ));
540 inline Double_t AliKFParticle::GetErrPz () const
542 return TMath::Sqrt(TMath::Abs( GetCovariance(5,5) ));
545 inline Double_t AliKFParticle::GetErrE () const
547 return TMath::Sqrt(TMath::Abs( GetCovariance(6,6) ));
550 inline Double_t AliKFParticle::GetErrS () const
552 return TMath::Sqrt(TMath::Abs( GetCovariance(7,7) ));
555 inline Double_t AliKFParticle::GetErrP () const
558 if( AliKFParticleBase::GetMomentum( par, err ) ) return 1.e10;
562 inline Double_t AliKFParticle::GetErrPt () const
565 if( AliKFParticleBase::GetPt( par, err ) ) return 1.e10;
569 inline Double_t AliKFParticle::GetErrEta () const
572 if( AliKFParticleBase::GetEta( par, err ) ) return 1.e10;
576 inline Double_t AliKFParticle::GetErrPhi () const
579 if( AliKFParticleBase::GetPhi( par, err ) ) return 1.e10;
583 inline Double_t AliKFParticle::GetErrMomentum () const
586 if( AliKFParticleBase::GetMomentum( par, err ) ) return 1.e10;
590 inline Double_t AliKFParticle::GetErrMass () const
593 if( AliKFParticleBase::GetMass( par, err ) ) return 1.e10;
597 inline Double_t AliKFParticle::GetErrDecayLength () const
600 if( AliKFParticleBase::GetDecayLength( par, err ) ) return 1.e10;
604 inline Double_t AliKFParticle::GetErrLifeTime () const
607 if( AliKFParticleBase::GetLifeTime( par, err ) ) return 1.e10;
611 inline Double_t AliKFParticle::GetErrR () const
614 if( AliKFParticleBase::GetR( par, err ) ) return 1.e10;
619 inline int AliKFParticle::GetP( Double_t &P, Double_t &SigmaP ) const
621 return AliKFParticleBase::GetMomentum( P, SigmaP );
624 inline int AliKFParticle::GetPt( Double_t &Pt, Double_t &SigmaPt ) const
626 return AliKFParticleBase::GetPt( Pt, SigmaPt );
629 inline int AliKFParticle::GetEta( Double_t &Eta, Double_t &SigmaEta ) const
631 return AliKFParticleBase::GetEta( Eta, SigmaEta );
634 inline int AliKFParticle::GetPhi( Double_t &Phi, Double_t &SigmaPhi ) const
636 return AliKFParticleBase::GetPhi( Phi, SigmaPhi );
639 inline int AliKFParticle::GetMomentum( Double_t &P, Double_t &SigmaP ) const
641 return AliKFParticleBase::GetMomentum( P, SigmaP );
644 inline int AliKFParticle::GetMass( Double_t &M, Double_t &SigmaM ) const
646 return AliKFParticleBase::GetMass( M, SigmaM );
649 inline int AliKFParticle::GetDecayLength( Double_t &L, Double_t &SigmaL ) const
651 return AliKFParticleBase::GetDecayLength( L, SigmaL );
654 inline int AliKFParticle::GetLifeTime( Double_t &T, Double_t &SigmaT ) const
656 return AliKFParticleBase::GetLifeTime( T, SigmaT );
659 inline int AliKFParticle::GetR( Double_t &R, Double_t &SigmaR ) const
661 return AliKFParticleBase::GetR( R, SigmaR );
664 inline Double_t & AliKFParticle::X()
666 return AliKFParticleBase::X();
669 inline Double_t & AliKFParticle::Y()
671 return AliKFParticleBase::Y();
674 inline Double_t & AliKFParticle::Z()
676 return AliKFParticleBase::Z();
679 inline Double_t & AliKFParticle::Px()
681 return AliKFParticleBase::Px();
684 inline Double_t & AliKFParticle::Py()
686 return AliKFParticleBase::Py();
689 inline Double_t & AliKFParticle::Pz()
691 return AliKFParticleBase::Pz();
694 inline Double_t & AliKFParticle::E()
696 return AliKFParticleBase::E();
699 inline Double_t & AliKFParticle::S()
701 return AliKFParticleBase::S();
704 inline Int_t & AliKFParticle::Q()
706 return AliKFParticleBase::Q();
709 inline Double_t & AliKFParticle::Chi2()
711 return AliKFParticleBase::Chi2();
714 inline Int_t & AliKFParticle::NDF()
716 return AliKFParticleBase::NDF();
719 inline Double_t & AliKFParticle::Parameter ( int i )
721 return AliKFParticleBase::Parameter(i);
724 inline Double_t & AliKFParticle::Covariance( int i )
726 return AliKFParticleBase::Covariance(i);
729 inline Double_t & AliKFParticle::Covariance( int i, int j )
731 return AliKFParticleBase::Covariance(i,j);
734 inline Double_t * AliKFParticle::Parameters ()
739 inline Double_t * AliKFParticle::CovarianceMatrix()
745 inline void AliKFParticle::operator +=( const AliKFParticle &Daughter )
747 AliKFParticleBase::operator +=( Daughter );
751 inline void AliKFParticle::AddDaughter( const AliKFParticle &Daughter )
753 AliKFParticleBase::AddDaughter( Daughter );
756 inline void AliKFParticle::SetProductionVertex( const AliKFParticle &Vtx )
758 AliKFParticleBase::SetProductionVertex( Vtx );
761 inline void AliKFParticle::SetMassConstraint( Double_t Mass, Double_t SigmaMass )
763 AliKFParticleBase::SetMassConstraint( Mass, SigmaMass );
766 inline void AliKFParticle::SetNoDecayLength()
768 AliKFParticleBase::SetNoDecayLength();
771 inline void AliKFParticle::Construct( const AliKFParticle *vDaughters[], int NDaughters,
772 const AliKFParticle *ProdVtx, Double_t Mass, Bool_t IsConstrained )
774 AliKFParticleBase::Construct( ( const AliKFParticleBase**)vDaughters, NDaughters,
775 ( const AliKFParticleBase*)ProdVtx, Mass, IsConstrained );
778 inline void AliKFParticle::TransportToDecayVertex()
780 AliKFParticleBase::TransportToDecayVertex();
783 inline void AliKFParticle::TransportToProductionVertex()
785 AliKFParticleBase::TransportToProductionVertex();
788 inline void AliKFParticle::TransportToPoint( const Double_t xyz[] )
790 TransportToDS( GetDStoPoint(xyz) );
793 inline void AliKFParticle::TransportToVertex( const AliVVertex &v )
795 TransportToPoint( AliKFParticle(v).fP );
798 inline void AliKFParticle::TransportToParticle( const AliKFParticle &p )
801 GetDStoParticle( p, dS, dSp );
805 inline void AliKFParticle::TransportToDS( Double_t dS )
807 AliKFParticleBase::TransportToDS( dS );
810 inline Double_t AliKFParticle::GetDStoPoint( const Double_t xyz[] ) const
812 return AliKFParticleBase::GetDStoPointBz( GetFieldAlice(), xyz );
816 inline void AliKFParticle::GetDStoParticle( const AliKFParticle &p,
817 Double_t &DS, Double_t &DSp ) const
819 GetDStoParticleXY( p, DS, DSp );
823 inline Double_t AliKFParticle::GetDistanceFromVertex( const Double_t vtx[] ) const
825 return AliKFParticleBase::GetDistanceFromVertex( vtx );
828 inline Double_t AliKFParticle::GetDeviationFromVertex( const Double_t v[],
829 const Double_t Cv[] ) const
831 return AliKFParticleBase::GetDeviationFromVertex( v, Cv);
834 inline Double_t AliKFParticle::GetDistanceFromVertex( const AliKFParticle &Vtx ) const
836 return AliKFParticleBase::GetDistanceFromVertex( Vtx );
839 inline Double_t AliKFParticle::GetDeviationFromVertex( const AliKFParticle &Vtx ) const
841 return AliKFParticleBase::GetDeviationFromVertex( Vtx );
844 inline Double_t AliKFParticle::GetDistanceFromVertex( const AliVVertex &Vtx ) const
846 return GetDistanceFromVertex( AliKFParticle(Vtx) );
849 inline Double_t AliKFParticle::GetDeviationFromVertex( const AliVVertex &Vtx ) const
851 return GetDeviationFromVertex( AliKFParticle(Vtx) );
854 inline Double_t AliKFParticle::GetDistanceFromParticle( const AliKFParticle &p ) const
856 return AliKFParticleBase::GetDistanceFromParticle( p );
859 inline Double_t AliKFParticle::GetDeviationFromParticle( const AliKFParticle &p ) const
861 return AliKFParticleBase::GetDeviationFromParticle( p );
864 inline void AliKFParticle::SubtractFromVertex( AliKFParticle &v ) const
866 AliKFParticleBase::SubtractFromVertex( v.fP, v.fC, v.fChi2, v.fNDF);
869 inline void AliKFParticle::SubtractFromVertex( AliESDVertex &v ) const
871 AliKFParticle vTmp(v);
872 SubtractFromVertex( vTmp );
873 v = AliESDVertex( vTmp.fP, vTmp.fC, vTmp.fChi2, (vTmp.fNDF +3)/2, v.GetName() );
876 inline void AliKFParticle::CopyToESDVertex( AliESDVertex &v ) const
878 AliKFParticle vTmp=*this;
879 v = AliESDVertex( vTmp.fP, vTmp.fC, vTmp.fChi2, (vTmp.fNDF +3)/2 );
882 inline Double_t AliKFParticle::GetFieldAlice()
887 inline void AliKFParticle::GetFieldValue( const Double_t * /*xyz*/, Double_t B[] ) const
890 B[2] = GetFieldAlice();
893 inline void AliKFParticle::GetDStoParticle( const AliKFParticleBase &p,
894 Double_t &DS, Double_t &DSp )const
896 GetDStoParticleXY( p, DS, DSp );
899 inline void AliKFParticle::GetDStoParticleXY( const AliKFParticleBase &p,
900 Double_t &DS, Double_t &DSp ) const
902 AliKFParticleBase::GetDStoParticleBz( GetFieldAlice(), p, DS, DSp ) ;
903 //GetDStoParticleALICE( p, DS, DSp ) ;
906 inline void AliKFParticle::Transport( Double_t dS, Double_t P[], Double_t C[] ) const
908 AliKFParticleBase::TransportBz( GetFieldAlice(), dS, P, C );