]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliKFParticle.h
Add -Wconversion to Mac warnings
[u/mrichter/AliRoot.git] / STEER / AliKFParticle.h
1 //---------------------------------------------------------------------------------
2 // The AliKFParticle 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 AliKFParticleBase
14 // 
15 //  -= Copyright &copy ALICE HLT Group =-
16 //_________________________________________________________________________________
17
18 #ifndef ALIKFPARTICLE_H
19 #define ALIKFPARTICLE_H
20
21 #include "AliKFParticleBase.h"
22 #include "TMath.h"
23
24 class AliVTrack;
25 class AliVVertex;
26
27 class AliKFParticle :public AliKFParticleBase
28 {
29   
30  public:
31
32   //*
33   //*  INITIALIZATION
34   //*
35
36   //* Set magnetic field for all particles
37   
38   static void SetField( Double_t Bz );
39
40   //* Constructor (empty)
41
42   AliKFParticle():AliKFParticleBase(){ ; }
43
44   //* Destructor (empty)
45
46   ~AliKFParticle(){ ; }
47
48   //* Construction of mother particle by its 2-3-4 daughters
49
50   AliKFParticle( const AliKFParticle &d1, const AliKFParticle &d2, Bool_t gamma = kFALSE );
51
52   AliKFParticle( const AliKFParticle &d1, const AliKFParticle &d2, 
53                  const AliKFParticle &d3 );
54
55   AliKFParticle( const AliKFParticle &d1, const AliKFParticle &d2, 
56                  const AliKFParticle &d3, const AliKFParticle &d4 );
57  
58  //* Initialisation from "cartesian" coordinates ( X Y Z Px Py Pz )
59  //* Parameters, covariance matrix, charge and PID hypothesis should be provided 
60
61   void Create( const Double_t Param[], const Double_t Cov[], Int_t Charge, Int_t PID );
62
63  //* Initialisation from ALICE track, PID hypothesis shoould be provided 
64
65   AliKFParticle( const AliVTrack &track, Int_t PID );
66
67   //* Initialisation from VVertex 
68
69   AliKFParticle( const AliVVertex &vertex );
70
71   //* Initialise covariance matrix and set current parameters to 0.0 
72
73   void Initialize();
74
75   //* Set decay vertex parameters for linearisation 
76
77   void SetVtxGuess( Double_t x, Double_t y, Double_t z );
78
79   //*
80   //*  ACCESSORS
81   //*
82
83   //* Simple accessors 
84
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
96   
97   Double_t GetParameter ( int i ) const ;
98   Double_t GetCovariance( int i ) const ;
99   Double_t GetCovariance( int i, int j ) const ;
100
101   //* Accessors with calculations, value returned w/o error flag
102   
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
113
114   //* Accessors to estimated errors
115
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
134
135   //* Accessors with calculations( &value, &estimated sigma )
136   //* error flag returned (0 means no error during calculations) 
137
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
148
149
150   //*
151   //*  MODIFIERS
152   //*
153   
154   Double_t & X    () ;
155   Double_t & Y    () ;
156   Double_t & Z    () ;
157   Double_t & Px   () ;
158   Double_t & Py   () ;
159   Double_t & Pz   () ;
160   Double_t & E    () ;
161   Double_t & S    () ;
162   Int_t    & Q    () ;
163   Double_t & Chi2 () ;
164   Int_t    & NDF  () ;
165
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() ;
171
172   //* 
173   //* CONSTRUCTION OF THE PARTICLE BY ITS DAUGHTERS AND MOTHER
174   //* USING THE KALMAN FILTER METHOD
175   //*
176
177
178   //* Add daughter to the particle 
179
180   void AddDaughter( const AliKFParticle &Daughter );
181
182   //* Add daughter via += operator: ex.{ D0; D0+=Pion; D0+= Kaon; }
183
184   void operator +=( const AliKFParticle &Daughter );  
185
186   //* Set production vertex 
187
188   void SetProductionVertex( const AliKFParticle &Vtx );
189
190   //* Set mass constraint 
191
192   void SetMassConstraint( Double_t Mass, Double_t SigmaMass = 0  );
193   
194   //* Set no decay length for resonances
195
196   void SetNoDecayLength();
197
198   //* Everything in one go  
199
200   void Construct( const AliKFParticle *vDaughters[], int NDaughters, 
201                   const AliKFParticle *ProdVtx=0,   Double_t Mass=-1, Bool_t IsConstrained=0  );
202
203   //*
204   //*                   TRANSPORT
205   //* 
206   //*  ( main transportation parameter is S = SignedPath/Momentum )
207   //*  ( parameters of decay & production vertices are stored locally )
208   //*
209
210   //* Transport the particle to its decay vertex 
211
212   void TransportToDecayVertex();
213
214   //* Transport the particle to its production vertex 
215
216   void TransportToProductionVertex();
217
218   //* Transport the particle close to xyz[] point 
219
220   void TransportToPoint( const Double_t xyz[] );
221
222   //* Transport the particle close to VVertex  
223
224   void TransportToVertex( const AliVVertex &v );
225
226   //* Transport the particle close to another particle p 
227
228   void TransportToParticle( const AliKFParticle &p );
229
230   //* Transport the particle on dS parameter (SignedPath/Momentum) 
231
232   void TransportToDS( Double_t dS );
233
234   //* Get dS to a certain space point 
235
236   Double_t GetDStoPoint( const Double_t xyz[] ) const ;
237   
238   //* Get dS to other particle p (dSp for particle p also returned) 
239
240   void GetDStoParticle( const AliKFParticle &p, 
241                         Double_t &DS, Double_t &DSp ) const ;
242   
243   //* Get dS to other particle p in XY-plane
244
245   void GetDStoParticleXY( const AliKFParticleBase &p, 
246                           Double_t &DS, Double_t &DSp ) const ;
247   
248   //* 
249   //* OTHER UTILITIES
250   //*
251
252
253   //* Calculate distance from another object [cm]
254
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 ;
259
260   //* Calculate sqrt(Chi2/ndf) deviation from another object
261   //* ( v = [xyz]-vertex, Cv=[Cxx,Cxy,Cyy,Cxz,Cyz,Czz]-covariance matrix )
262
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 ;
267  
268   //* Calculate distance from another object [cm] in XY-plane
269
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 ;
274
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 ;
279
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 )
282
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 ;
287
288   //* Calculate opennig angle between two particles
289
290   Double_t GetAngle  ( const AliKFParticle &p ) const ;
291   Double_t GetAngleXY( const AliKFParticle &p ) const ;
292   Double_t GetAngleRZ( const AliKFParticle &p ) const ;
293
294   //* Subtract the particle from the vertex  
295
296   void SubtractFromVertex( AliKFParticle &v ) const ;
297
298   //* Special method for creating gammas
299
300   void ConstructGamma( const AliKFParticle &daughter1,
301                        const AliKFParticle &daughter2  );
302
303  protected: 
304   
305   //*
306   //*  INTERNAL STUFF
307   //* 
308
309   //* Method to access ALICE field 
310  
311   static Double_t GetFieldAlice();
312   
313   //* Other methods required by the abstract AliKFParticleBase class 
314   
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]  ) ;
319
320   //void GetDStoParticleALICE( const AliKFParticleBase &p, Double_t &DS, Double_t &DS1 ) const;
321
322
323  private:
324
325   static Double_t fgBz;  //* Bz compoment of the magnetic field
326
327   ClassDef( AliKFParticle, 1 );
328
329 };
330
331
332
333 //---------------------------------------------------------------------
334 //
335 //     Inline implementation of the AliKFParticle methods
336 //
337 //---------------------------------------------------------------------
338
339
340 inline void AliKFParticle::SetField( Double_t Bz )
341
342   fgBz = Bz;
343 }
344
345 inline AliKFParticle::AliKFParticle( const AliKFParticle &d1, 
346                                      const AliKFParticle &d2, 
347                                      const AliKFParticle &d3 )
348 {
349   AliKFParticle mother;
350   mother+= d1;
351   mother+= d2;
352   mother+= d3;
353   *this = mother;
354 }
355
356 inline AliKFParticle::AliKFParticle( const AliKFParticle &d1, 
357                                      const AliKFParticle &d2, 
358                                      const AliKFParticle &d3, 
359                                      const AliKFParticle &d4 )
360 {
361   AliKFParticle mother;
362   mother+= d1;
363   mother+= d2;
364   mother+= d3;
365   mother+= d4;
366   *this = mother;
367 }
368
369
370 inline void AliKFParticle::Initialize()
371
372   AliKFParticleBase::Initialize(); 
373 }
374
375 inline void AliKFParticle::SetVtxGuess( Double_t x, Double_t y, Double_t z )
376 {
377   AliKFParticleBase::SetVtxGuess(x,y,z);
378 }  
379
380 inline Double_t AliKFParticle::GetX    () const 
381
382   return AliKFParticleBase::GetX();    
383 }
384
385 inline Double_t AliKFParticle::GetY    () const 
386
387   return AliKFParticleBase::GetY();    
388 }
389
390 inline Double_t AliKFParticle::GetZ    () const 
391
392   return AliKFParticleBase::GetZ();    
393 }
394
395 inline Double_t AliKFParticle::GetPx   () const 
396
397   return AliKFParticleBase::GetPx();   
398 }
399
400 inline Double_t AliKFParticle::GetPy   () const 
401
402   return AliKFParticleBase::GetPy();   
403 }
404
405 inline Double_t AliKFParticle::GetPz   () const 
406
407   return AliKFParticleBase::GetPz();   
408 }
409
410 inline Double_t AliKFParticle::GetE    () const 
411
412   return AliKFParticleBase::GetE();    
413 }
414
415 inline Double_t AliKFParticle::GetS    () const 
416
417   return AliKFParticleBase::GetS();    
418 }
419
420 inline Int_t    AliKFParticle::GetQ    () const 
421
422   return AliKFParticleBase::GetQ();    
423 }
424
425 inline Double_t AliKFParticle::GetChi2 () const 
426
427   return AliKFParticleBase::GetChi2(); 
428 }
429
430 inline Int_t    AliKFParticle::GetNDF  () const 
431
432   return AliKFParticleBase::GetNDF();  
433 }
434
435 inline Double_t AliKFParticle::GetParameter ( int i ) const 
436
437   return AliKFParticleBase::GetParameter(i);  
438 }
439
440 inline Double_t AliKFParticle::GetCovariance( int i ) const 
441
442   return AliKFParticleBase::GetCovariance(i); 
443 }
444
445 inline Double_t AliKFParticle::GetCovariance( int i, int j ) const 
446
447   return AliKFParticleBase::GetCovariance(i,j);
448 }
449
450
451 inline Double_t AliKFParticle::GetP    () const
452 {
453   Double_t par, err;
454   if( AliKFParticleBase::GetMomentum( par, err ) ) return 0;
455   else return par;
456 }
457
458 inline Double_t AliKFParticle::GetPt   () const
459 {
460   Double_t par, err;
461   if( AliKFParticleBase::GetPt( par, err ) ) return 0;
462   else return par;
463 }
464
465 inline Double_t AliKFParticle::GetEta   () const
466 {
467   Double_t par, err;
468   if( AliKFParticleBase::GetEta( par, err ) ) return 0;
469   else return par;
470 }
471
472 inline Double_t AliKFParticle::GetPhi   () const
473 {
474   Double_t par, err;
475   if( AliKFParticleBase::GetPhi( par, err ) ) return 0;
476   else return par;
477 }
478
479 inline Double_t AliKFParticle::GetMomentum    () const
480 {
481   Double_t par, err;
482   if( AliKFParticleBase::GetMomentum( par, err ) ) return 0;
483   else return par;
484 }
485
486 inline Double_t AliKFParticle::GetMass        () const
487 {
488   Double_t par, err;
489   if( AliKFParticleBase::GetMass( par, err ) ) return 0;
490   else return par;
491 }
492
493 inline Double_t AliKFParticle::GetDecayLength () const
494 {
495   Double_t par, err;
496   if( AliKFParticleBase::GetDecayLength( par, err ) ) return 0;
497   else return par;
498 }
499
500 inline Double_t AliKFParticle::GetDecayLengthXY () const
501 {
502   Double_t par, err;
503   if( AliKFParticleBase::GetDecayLengthXY( par, err ) ) return 0;
504   else return par;
505 }
506
507 inline Double_t AliKFParticle::GetLifeTime    () const
508 {
509   Double_t par, err;
510   if( AliKFParticleBase::GetLifeTime( par, err ) ) return 0;
511   else return par;
512 }
513
514 inline Double_t AliKFParticle::GetR   () const
515 {
516   Double_t par, err;
517   if( AliKFParticleBase::GetR( par, err ) ) return 0;
518   else return par;
519 }
520
521 inline Double_t AliKFParticle::GetErrX           () const 
522 {
523   return TMath::Sqrt(TMath::Abs( GetCovariance(0,0) ));
524 }
525
526 inline Double_t AliKFParticle::GetErrY           () const 
527 {
528   return TMath::Sqrt(TMath::Abs( GetCovariance(1,1) ));
529 }
530
531 inline Double_t AliKFParticle::GetErrZ           () const 
532 {
533   return TMath::Sqrt(TMath::Abs( GetCovariance(2,2) ));
534 }
535
536 inline Double_t AliKFParticle::GetErrPx          () const 
537 {
538   return TMath::Sqrt(TMath::Abs( GetCovariance(3,3) ));
539 }
540
541 inline Double_t AliKFParticle::GetErrPy          () const 
542 {
543   return TMath::Sqrt(TMath::Abs( GetCovariance(4,4) ));
544 }
545
546 inline Double_t AliKFParticle::GetErrPz          () const 
547 {
548   return TMath::Sqrt(TMath::Abs( GetCovariance(5,5) ));
549 }
550
551 inline Double_t AliKFParticle::GetErrE           () const 
552 {
553   return TMath::Sqrt(TMath::Abs( GetCovariance(6,6) ));
554 }
555
556 inline Double_t AliKFParticle::GetErrS           () const 
557 {
558   return TMath::Sqrt(TMath::Abs( GetCovariance(7,7) ));
559 }
560
561 inline Double_t AliKFParticle::GetErrP    () const
562 {
563   Double_t par, err;
564   if( AliKFParticleBase::GetMomentum( par, err ) ) return 1.e10;
565   else return err;
566 }
567
568 inline Double_t AliKFParticle::GetErrPt    () const
569 {
570   Double_t par, err;
571   if( AliKFParticleBase::GetPt( par, err ) ) return 1.e10;
572   else return err;
573 }
574
575 inline Double_t AliKFParticle::GetErrEta    () const
576 {
577   Double_t par, err;
578   if( AliKFParticleBase::GetEta( par, err ) ) return 1.e10;
579   else return err;
580 }
581
582 inline Double_t AliKFParticle::GetErrPhi    () const
583 {
584   Double_t par, err;
585   if( AliKFParticleBase::GetPhi( par, err ) ) return 1.e10;
586   else return err;
587 }
588
589 inline Double_t AliKFParticle::GetErrMomentum    () const
590 {
591   Double_t par, err;
592   if( AliKFParticleBase::GetMomentum( par, err ) ) return 1.e10;
593   else return err;
594 }
595
596 inline Double_t AliKFParticle::GetErrMass        () const
597 {
598   Double_t par, err;
599   if( AliKFParticleBase::GetMass( par, err ) ) return 1.e10;
600   else return err;
601 }
602
603 inline Double_t AliKFParticle::GetErrDecayLength () const
604 {
605   Double_t par, err;
606   if( AliKFParticleBase::GetDecayLength( par, err ) ) return 1.e10;
607   else return err;
608 }
609
610 inline Double_t AliKFParticle::GetErrDecayLengthXY () const
611 {
612   Double_t par, err;
613   if( AliKFParticleBase::GetDecayLengthXY( par, err ) ) return 1.e10;
614   else return err;
615 }
616
617 inline Double_t AliKFParticle::GetErrLifeTime    () const
618 {
619   Double_t par, err;
620   if( AliKFParticleBase::GetLifeTime( par, err ) ) return 1.e10;
621   else return err;
622 }
623
624 inline Double_t AliKFParticle::GetErrR    () const
625 {
626   Double_t par, err;
627   if( AliKFParticleBase::GetR( par, err ) ) return 1.e10;
628   else return err;
629 }
630
631
632 inline int AliKFParticle::GetP( Double_t &P, Double_t &SigmaP ) const 
633 {
634   return AliKFParticleBase::GetMomentum( P, SigmaP );
635 }
636
637 inline int AliKFParticle::GetPt( Double_t &Pt, Double_t &SigmaPt ) const 
638 {
639   return AliKFParticleBase::GetPt( Pt, SigmaPt );
640 }
641
642 inline int AliKFParticle::GetEta( Double_t &Eta, Double_t &SigmaEta ) const 
643 {
644   return AliKFParticleBase::GetEta( Eta, SigmaEta );
645 }
646
647 inline int AliKFParticle::GetPhi( Double_t &Phi, Double_t &SigmaPhi ) const 
648 {
649   return AliKFParticleBase::GetPhi( Phi, SigmaPhi );
650 }
651
652 inline int AliKFParticle::GetMomentum( Double_t &P, Double_t &SigmaP ) const 
653 {
654   return AliKFParticleBase::GetMomentum( P, SigmaP );
655 }
656
657 inline int AliKFParticle::GetMass( Double_t &M, Double_t &SigmaM ) const 
658 {
659   return AliKFParticleBase::GetMass( M, SigmaM );
660 }
661
662 inline int AliKFParticle::GetDecayLength( Double_t &L, Double_t &SigmaL ) const 
663 {
664   return AliKFParticleBase::GetDecayLength( L, SigmaL );
665 }
666
667 inline int AliKFParticle::GetDecayLengthXY( Double_t &L, Double_t &SigmaL ) const 
668 {
669   return AliKFParticleBase::GetDecayLengthXY( L, SigmaL );
670 }
671
672 inline int AliKFParticle::GetLifeTime( Double_t &T, Double_t &SigmaT ) const 
673 {
674   return AliKFParticleBase::GetLifeTime( T, SigmaT );
675 }
676
677 inline int AliKFParticle::GetR( Double_t &R, Double_t &SigmaR ) const 
678 {
679   return AliKFParticleBase::GetR( R, SigmaR );
680 }
681
682 inline Double_t & AliKFParticle::X() 
683
684   return AliKFParticleBase::X();    
685 }
686
687 inline Double_t & AliKFParticle::Y()
688
689   return AliKFParticleBase::Y();    
690 }
691
692 inline Double_t & AliKFParticle::Z() 
693
694   return AliKFParticleBase::Z();    
695 }
696
697 inline Double_t & AliKFParticle::Px() 
698
699   return AliKFParticleBase::Px();   
700 }
701
702 inline Double_t & AliKFParticle::Py() 
703
704   return AliKFParticleBase::Py();   
705 }
706
707 inline Double_t & AliKFParticle::Pz() 
708
709   return AliKFParticleBase::Pz();   
710 }
711
712 inline Double_t & AliKFParticle::E() 
713
714   return AliKFParticleBase::E();    
715 }
716
717 inline Double_t & AliKFParticle::S() 
718
719   return AliKFParticleBase::S();    
720 }
721
722 inline Int_t    & AliKFParticle::Q() 
723
724   return AliKFParticleBase::Q();    
725 }
726
727 inline Double_t & AliKFParticle::Chi2() 
728
729   return AliKFParticleBase::Chi2(); 
730 }
731
732 inline Int_t    & AliKFParticle::NDF() 
733
734   return AliKFParticleBase::NDF();  
735 }
736
737 inline Double_t & AliKFParticle::Parameter ( int i )        
738
739   return AliKFParticleBase::Parameter(i);
740 }
741
742 inline Double_t & AliKFParticle::Covariance( int i )        
743
744   return AliKFParticleBase::Covariance(i);
745 }
746
747 inline Double_t & AliKFParticle::Covariance( int i, int j ) 
748
749   return AliKFParticleBase::Covariance(i,j); 
750 }
751
752 inline Double_t * AliKFParticle::Parameters ()
753 {
754   return fP;
755 }
756
757 inline Double_t * AliKFParticle::CovarianceMatrix()
758 {
759   return fC;
760 }
761
762
763 inline void AliKFParticle::operator +=( const AliKFParticle &Daughter )
764 {
765   AliKFParticleBase::operator +=( Daughter );
766 }
767   
768
769 inline void AliKFParticle::AddDaughter( const AliKFParticle &Daughter )
770 {
771   AliKFParticleBase::AddDaughter( Daughter );
772 }
773
774 inline void AliKFParticle::SetProductionVertex( const AliKFParticle &Vtx )
775 {
776   AliKFParticleBase::SetProductionVertex( Vtx );
777 }
778
779 inline void AliKFParticle::SetMassConstraint( Double_t Mass, Double_t SigmaMass )
780 {
781   AliKFParticleBase::SetMassConstraint( Mass, SigmaMass );
782 }
783     
784 inline void AliKFParticle::SetNoDecayLength()
785 {
786   AliKFParticleBase::SetNoDecayLength();
787 }
788
789 inline void AliKFParticle::Construct( const AliKFParticle *vDaughters[], int NDaughters, 
790                                const AliKFParticle *ProdVtx,   Double_t Mass, Bool_t IsConstrained  )
791 {    
792   AliKFParticleBase::Construct( ( const AliKFParticleBase**)vDaughters, NDaughters, 
793                          ( const AliKFParticleBase*)ProdVtx, Mass, IsConstrained );
794 }
795
796 inline void AliKFParticle::TransportToDecayVertex()
797
798   AliKFParticleBase::TransportToDecayVertex(); 
799 }
800
801 inline void AliKFParticle::TransportToProductionVertex()
802 {
803   AliKFParticleBase::TransportToProductionVertex();
804 }
805
806 inline void AliKFParticle::TransportToPoint( const Double_t xyz[] )
807
808   TransportToDS( GetDStoPoint(xyz) );
809 }
810
811 inline void AliKFParticle::TransportToVertex( const AliVVertex &v )
812 {       
813   TransportToPoint( AliKFParticle(v).fP );
814 }
815
816 inline void AliKFParticle::TransportToParticle( const AliKFParticle &p )
817
818   Double_t dS, dSp;
819   GetDStoParticle( p, dS, dSp );
820   TransportToDS( dS );
821 }
822
823 inline void AliKFParticle::TransportToDS( Double_t dS )
824 {
825   AliKFParticleBase::TransportToDS( dS );
826
827
828 inline Double_t AliKFParticle::GetDStoPoint( const Double_t xyz[] ) const 
829 {
830   return AliKFParticleBase::GetDStoPointBz( GetFieldAlice(), xyz );
831 }
832
833   
834 inline void AliKFParticle::GetDStoParticle( const AliKFParticle &p, 
835                                             Double_t &DS, Double_t &DSp ) const 
836 {
837   GetDStoParticleXY( p, DS, DSp );
838 }
839
840
841 inline Double_t AliKFParticle::GetDistanceFromVertex( const Double_t vtx[] ) const
842 {
843   return AliKFParticleBase::GetDistanceFromVertex( vtx );
844 }
845
846 inline Double_t AliKFParticle::GetDeviationFromVertex( const Double_t v[], 
847                                                        const Double_t Cv[] ) const
848 {
849   return AliKFParticleBase::GetDeviationFromVertex( v, Cv);
850 }
851
852 inline Double_t AliKFParticle::GetDistanceFromVertex( const AliKFParticle &Vtx ) const
853 {
854   return AliKFParticleBase::GetDistanceFromVertex( Vtx );
855 }
856
857 inline Double_t AliKFParticle::GetDeviationFromVertex( const AliKFParticle &Vtx ) const
858 {
859   return AliKFParticleBase::GetDeviationFromVertex( Vtx );
860 }
861
862 inline Double_t AliKFParticle::GetDistanceFromVertex( const AliVVertex &Vtx ) const
863 {
864   return GetDistanceFromVertex( AliKFParticle(Vtx) );
865 }
866
867 inline Double_t AliKFParticle::GetDeviationFromVertex( const AliVVertex &Vtx ) const
868 {
869   return GetDeviationFromVertex( AliKFParticle(Vtx) );
870 }
871  
872 inline Double_t AliKFParticle::GetDistanceFromParticle( const AliKFParticle &p ) const 
873 {
874   return AliKFParticleBase::GetDistanceFromParticle( p );
875 }
876
877 inline Double_t AliKFParticle::GetDeviationFromParticle( const AliKFParticle &p ) const 
878 {
879   return AliKFParticleBase::GetDeviationFromParticle( p );
880 }
881
882 inline void AliKFParticle::SubtractFromVertex( AliKFParticle &v ) const 
883 {
884   AliKFParticleBase::SubtractFromVertex( v );
885 }
886
887 inline Double_t AliKFParticle::GetFieldAlice()
888
889   return fgBz; 
890 }
891
892 inline void AliKFParticle::GetFieldValue( const Double_t * /*xyz*/, Double_t B[] ) const 
893 {    
894   B[0] = B[1] = 0;
895   B[2] = GetFieldAlice();
896 }
897
898 inline void AliKFParticle::GetDStoParticle( const AliKFParticleBase &p, 
899                                             Double_t &DS, Double_t &DSp )const
900 {
901   GetDStoParticleXY( p, DS, DSp );
902 }
903
904 inline void AliKFParticle::GetDStoParticleXY( const AliKFParticleBase &p, 
905                                        Double_t &DS, Double_t &DSp ) const
906
907   AliKFParticleBase::GetDStoParticleBz( GetFieldAlice(), p, DS, DSp ) ;
908   //GetDStoParticleALICE( p, DS, DSp ) ;
909 }
910
911 inline void AliKFParticle::Transport( Double_t dS, Double_t P[], Double_t C[] ) const 
912 {
913   AliKFParticleBase::TransportBz( GetFieldAlice(), dS, P, C );
914 }
915
916 inline void AliKFParticle::ConstructGamma( const AliKFParticle &daughter1,
917                                            const AliKFParticle &daughter2  )
918 {
919   AliKFParticleBase::ConstructGammaBz( daughter1, daughter2, GetFieldAlice() );
920 }
921
922 #endif