]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliKFParticle.h
06fe2bdd58418f9e0f0b609fe3829e0b636f970a
[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 "AliESDVertex.h"
23
24 class AliVTrack;
25
26 class AliKFParticle :public AliKFParticleBase
27 {
28   
29  public:
30
31   //*
32   //*  INITIALIZATION
33   //*
34
35   //* Set magnetic field for all particles
36   
37   static void SetField( Double_t Bz );
38
39   //* Constructor (empty)
40
41   AliKFParticle():AliKFParticleBase(){ ; }
42
43   //* Destructor (empty)
44
45   ~AliKFParticle(){ ; }
46
47   //* Construction of mother particle by its 2-3-4 daughters
48
49   AliKFParticle( const AliKFParticle &d1, const AliKFParticle &d2 );
50
51   AliKFParticle( const AliKFParticle &d1, const AliKFParticle &d2, 
52                  const AliKFParticle &d3 );
53
54   AliKFParticle( const AliKFParticle &d1, const AliKFParticle &d2, 
55                  const AliKFParticle &d3, const AliKFParticle &d4 );
56  
57  //* Initialisation from "cartesian" coordinates ( X Y Z Px Py Pz )
58  //* Parameters, covariance matrix, charge and PID hypothesis should be provided 
59
60   void Create( const Double_t Param[], const Double_t Cov[], Int_t Charge, Int_t PID );
61
62  //* Initialisation from ALICE track, PID hypothesis shoould be provided 
63
64   AliKFParticle( const AliVTrack &track, Int_t PID );
65
66   //* Initialisation from VVertex 
67
68   AliKFParticle( const AliVVertex &vertex );
69
70   //* Copy position part to ESD vertex 
71
72   void CopyToESDVertex( AliESDVertex &Vtx ) const ;
73
74   //* Initialise covariance matrix and set current parameters to 0.0 
75
76   void Initialize();
77
78   //* Set decay vertex parameters for linearisation 
79
80   void SetVtxGuess( Double_t x, Double_t y, Double_t z );
81
82   //*
83   //*  ACCESSORS
84   //*
85
86   //* Simple accessors 
87
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
99   
100   Double_t GetParameter ( int i ) const ;
101   Double_t GetCovariance( int i ) const ;
102   Double_t GetCovariance( int i, int j ) const ;
103
104   //* Accessors with calculations, value returned w/o error flag
105   
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
115
116   //* Accessors to estimated errors
117
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
135
136   //* Accessors with calculations( &value, &estimated sigma )
137   //* error flag returned (0 means no error during calculations) 
138
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
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   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 ;
274
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 )
277
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 ;
282
283   //* Calculate opennig angle between two particles
284
285   Double_t GetAngle  ( const AliKFParticle &p ) const ;
286   Double_t GetAngleXY( const AliKFParticle &p ) const ;
287   Double_t GetAngleRZ( const AliKFParticle &p ) const ;
288
289   //* Subtract the particle from the vertex  
290
291   void SubtractFromVertex( AliKFParticle &v ) const ;
292   void SubtractFromVertex( AliESDVertex &v ) const ;
293
294  protected: 
295   
296   //*
297   //*  INTERNAL STUFF
298   //* 
299
300   //* Method to access ALICE field 
301  
302   static Double_t GetFieldAlice();
303   
304   //* Other methods required by the abstract AliKFParticleBase class 
305   
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]  ) ;
310
311   //void GetDStoParticleALICE( const AliKFParticleBase &p, Double_t &DS, Double_t &DS1 ) const;
312
313
314  private:
315
316   static Double_t fgBz;  //* Bz compoment of the magnetic field
317
318   ClassDef( AliKFParticle, 1 );
319
320 };
321
322
323
324 //---------------------------------------------------------------------
325 //
326 //     Inline implementation of the AliKFParticle methods
327 //
328 //---------------------------------------------------------------------
329
330
331 inline void AliKFParticle::SetField( Double_t Bz )
332
333   fgBz = -Bz;//!!!
334 }
335
336
337 inline AliKFParticle::AliKFParticle( const AliKFParticle &d1, 
338                                      const AliKFParticle &d2 )
339 {
340   AliKFParticle mother;
341   mother+= d1;
342   mother+= d2;
343   *this = mother;
344 }
345
346 inline AliKFParticle::AliKFParticle( const AliKFParticle &d1, 
347                                      const AliKFParticle &d2, 
348                                      const AliKFParticle &d3 )
349 {
350   AliKFParticle mother;
351   mother+= d1;
352   mother+= d2;
353   mother+= d3;
354   *this = mother;
355 }
356
357 inline AliKFParticle::AliKFParticle( const AliKFParticle &d1, 
358                                      const AliKFParticle &d2, 
359                                      const AliKFParticle &d3, 
360                                      const AliKFParticle &d4 )
361 {
362   AliKFParticle mother;
363   mother+= d1;
364   mother+= d2;
365   mother+= d3;
366   mother+= d4;
367   *this = mother;
368 }
369
370
371 inline void AliKFParticle::Initialize()
372
373   AliKFParticleBase::Initialize(); 
374 }
375
376 inline void AliKFParticle::SetVtxGuess( Double_t x, Double_t y, Double_t z )
377 {
378   AliKFParticleBase::SetVtxGuess(x,y,z);
379 }  
380
381 inline Double_t AliKFParticle::GetX    () const 
382
383   return AliKFParticleBase::GetX();    
384 }
385
386 inline Double_t AliKFParticle::GetY    () const 
387
388   return AliKFParticleBase::GetY();    
389 }
390
391 inline Double_t AliKFParticle::GetZ    () const 
392
393   return AliKFParticleBase::GetZ();    
394 }
395
396 inline Double_t AliKFParticle::GetPx   () const 
397
398   return AliKFParticleBase::GetPx();   
399 }
400
401 inline Double_t AliKFParticle::GetPy   () const 
402
403   return AliKFParticleBase::GetPy();   
404 }
405
406 inline Double_t AliKFParticle::GetPz   () const 
407
408   return AliKFParticleBase::GetPz();   
409 }
410
411 inline Double_t AliKFParticle::GetE    () const 
412
413   return AliKFParticleBase::GetE();    
414 }
415
416 inline Double_t AliKFParticle::GetS    () const 
417
418   return AliKFParticleBase::GetS();    
419 }
420
421 inline Int_t    AliKFParticle::GetQ    () const 
422
423   return AliKFParticleBase::GetQ();    
424 }
425
426 inline Double_t AliKFParticle::GetChi2 () const 
427
428   return AliKFParticleBase::GetChi2(); 
429 }
430
431 inline Int_t    AliKFParticle::GetNDF  () const 
432
433   return AliKFParticleBase::GetNDF();  
434 }
435
436 inline Double_t AliKFParticle::GetParameter ( int i ) const 
437
438   return AliKFParticleBase::GetParameter(i);  
439 }
440
441 inline Double_t AliKFParticle::GetCovariance( int i ) const 
442
443   return AliKFParticleBase::GetCovariance(i); 
444 }
445
446 inline Double_t AliKFParticle::GetCovariance( int i, int j ) const 
447
448   return AliKFParticleBase::GetCovariance(i,j);
449 }
450
451
452 inline Double_t AliKFParticle::GetP    () const
453 {
454   Double_t par, err;
455   if( AliKFParticleBase::GetMomentum( par, err ) ) return 0;
456   else return par;
457 }
458
459 inline Double_t AliKFParticle::GetPt   () const
460 {
461   Double_t par, err;
462   if( AliKFParticleBase::GetPt( par, err ) ) return 0;
463   else return par;
464 }
465
466 inline Double_t AliKFParticle::GetEta   () const
467 {
468   Double_t par, err;
469   if( AliKFParticleBase::GetEta( par, err ) ) return 0;
470   else return par;
471 }
472
473 inline Double_t AliKFParticle::GetPhi   () const
474 {
475   Double_t par, err;
476   if( AliKFParticleBase::GetPhi( par, err ) ) return 0;
477   else return par;
478 }
479
480 inline Double_t AliKFParticle::GetMomentum    () const
481 {
482   Double_t par, err;
483   if( AliKFParticleBase::GetMomentum( par, err ) ) return 0;
484   else return par;
485 }
486
487 inline Double_t AliKFParticle::GetMass        () const
488 {
489   Double_t par, err;
490   if( AliKFParticleBase::GetMass( par, err ) ) return 0;
491   else return par;
492 }
493
494 inline Double_t AliKFParticle::GetDecayLength () const
495 {
496   Double_t par, err;
497   if( AliKFParticleBase::GetDecayLength( par, err ) ) return 0;
498   else return par;
499 }
500
501 inline Double_t AliKFParticle::GetLifeTime    () const
502 {
503   Double_t par, err;
504   if( AliKFParticleBase::GetLifeTime( par, err ) ) return 0;
505   else return par;
506 }
507
508 inline Double_t AliKFParticle::GetR   () const
509 {
510   Double_t par, err;
511   if( AliKFParticleBase::GetR( par, err ) ) return 0;
512   else return par;
513 }
514
515 inline Double_t AliKFParticle::GetErrX           () const 
516 {
517   return TMath::Sqrt(TMath::Abs( GetCovariance(0,0) ));
518 }
519
520 inline Double_t AliKFParticle::GetErrY           () const 
521 {
522   return TMath::Sqrt(TMath::Abs( GetCovariance(1,1) ));
523 }
524
525 inline Double_t AliKFParticle::GetErrZ           () const 
526 {
527   return TMath::Sqrt(TMath::Abs( GetCovariance(2,2) ));
528 }
529
530 inline Double_t AliKFParticle::GetErrPx          () const 
531 {
532   return TMath::Sqrt(TMath::Abs( GetCovariance(3,3) ));
533 }
534
535 inline Double_t AliKFParticle::GetErrPy          () const 
536 {
537   return TMath::Sqrt(TMath::Abs( GetCovariance(4,4) ));
538 }
539
540 inline Double_t AliKFParticle::GetErrPz          () const 
541 {
542   return TMath::Sqrt(TMath::Abs( GetCovariance(5,5) ));
543 }
544
545 inline Double_t AliKFParticle::GetErrE           () const 
546 {
547   return TMath::Sqrt(TMath::Abs( GetCovariance(6,6) ));
548 }
549
550 inline Double_t AliKFParticle::GetErrS           () const 
551 {
552   return TMath::Sqrt(TMath::Abs( GetCovariance(7,7) ));
553 }
554
555 inline Double_t AliKFParticle::GetErrP    () const
556 {
557   Double_t par, err;
558   if( AliKFParticleBase::GetMomentum( par, err ) ) return 1.e10;
559   else return err;
560 }
561
562 inline Double_t AliKFParticle::GetErrPt    () const
563 {
564   Double_t par, err;
565   if( AliKFParticleBase::GetPt( par, err ) ) return 1.e10;
566   else return err;
567 }
568
569 inline Double_t AliKFParticle::GetErrEta    () const
570 {
571   Double_t par, err;
572   if( AliKFParticleBase::GetEta( par, err ) ) return 1.e10;
573   else return err;
574 }
575
576 inline Double_t AliKFParticle::GetErrPhi    () const
577 {
578   Double_t par, err;
579   if( AliKFParticleBase::GetPhi( par, err ) ) return 1.e10;
580   else return err;
581 }
582
583 inline Double_t AliKFParticle::GetErrMomentum    () const
584 {
585   Double_t par, err;
586   if( AliKFParticleBase::GetMomentum( par, err ) ) return 1.e10;
587   else return err;
588 }
589
590 inline Double_t AliKFParticle::GetErrMass        () const
591 {
592   Double_t par, err;
593   if( AliKFParticleBase::GetMass( par, err ) ) return 1.e10;
594   else return err;
595 }
596
597 inline Double_t AliKFParticle::GetErrDecayLength () const
598 {
599   Double_t par, err;
600   if( AliKFParticleBase::GetDecayLength( par, err ) ) return 1.e10;
601   else return err;
602 }
603
604 inline Double_t AliKFParticle::GetErrLifeTime    () const
605 {
606   Double_t par, err;
607   if( AliKFParticleBase::GetLifeTime( par, err ) ) return 1.e10;
608   else return err;
609 }
610
611 inline Double_t AliKFParticle::GetErrR    () const
612 {
613   Double_t par, err;
614   if( AliKFParticleBase::GetR( par, err ) ) return 1.e10;
615   else return err;
616 }
617
618
619 inline int AliKFParticle::GetP( Double_t &P, Double_t &SigmaP ) const 
620 {
621   return AliKFParticleBase::GetMomentum( P, SigmaP );
622 }
623
624 inline int AliKFParticle::GetPt( Double_t &Pt, Double_t &SigmaPt ) const 
625 {
626   return AliKFParticleBase::GetPt( Pt, SigmaPt );
627 }
628
629 inline int AliKFParticle::GetEta( Double_t &Eta, Double_t &SigmaEta ) const 
630 {
631   return AliKFParticleBase::GetEta( Eta, SigmaEta );
632 }
633
634 inline int AliKFParticle::GetPhi( Double_t &Phi, Double_t &SigmaPhi ) const 
635 {
636   return AliKFParticleBase::GetPhi( Phi, SigmaPhi );
637 }
638
639 inline int AliKFParticle::GetMomentum( Double_t &P, Double_t &SigmaP ) const 
640 {
641   return AliKFParticleBase::GetMomentum( P, SigmaP );
642 }
643
644 inline int AliKFParticle::GetMass( Double_t &M, Double_t &SigmaM ) const 
645 {
646   return AliKFParticleBase::GetMass( M, SigmaM );
647 }
648
649 inline int AliKFParticle::GetDecayLength( Double_t &L, Double_t &SigmaL ) const 
650 {
651   return AliKFParticleBase::GetDecayLength( L, SigmaL );
652 }
653
654 inline int AliKFParticle::GetLifeTime( Double_t &T, Double_t &SigmaT ) const 
655 {
656   return AliKFParticleBase::GetLifeTime( T, SigmaT );
657 }
658
659 inline int AliKFParticle::GetR( Double_t &R, Double_t &SigmaR ) const 
660 {
661   return AliKFParticleBase::GetR( R, SigmaR );
662 }
663
664 inline Double_t & AliKFParticle::X() 
665
666   return AliKFParticleBase::X();    
667 }
668
669 inline Double_t & AliKFParticle::Y()
670
671   return AliKFParticleBase::Y();    
672 }
673
674 inline Double_t & AliKFParticle::Z() 
675
676   return AliKFParticleBase::Z();    
677 }
678
679 inline Double_t & AliKFParticle::Px() 
680
681   return AliKFParticleBase::Px();   
682 }
683
684 inline Double_t & AliKFParticle::Py() 
685
686   return AliKFParticleBase::Py();   
687 }
688
689 inline Double_t & AliKFParticle::Pz() 
690
691   return AliKFParticleBase::Pz();   
692 }
693
694 inline Double_t & AliKFParticle::E() 
695
696   return AliKFParticleBase::E();    
697 }
698
699 inline Double_t & AliKFParticle::S() 
700
701   return AliKFParticleBase::S();    
702 }
703
704 inline Int_t    & AliKFParticle::Q() 
705
706   return AliKFParticleBase::Q();    
707 }
708
709 inline Double_t & AliKFParticle::Chi2() 
710
711   return AliKFParticleBase::Chi2(); 
712 }
713
714 inline Int_t    & AliKFParticle::NDF() 
715
716   return AliKFParticleBase::NDF();  
717 }
718
719 inline Double_t & AliKFParticle::Parameter ( int i )        
720
721   return AliKFParticleBase::Parameter(i);
722 }
723
724 inline Double_t & AliKFParticle::Covariance( int i )        
725
726   return AliKFParticleBase::Covariance(i);
727 }
728
729 inline Double_t & AliKFParticle::Covariance( int i, int j ) 
730
731   return AliKFParticleBase::Covariance(i,j); 
732 }
733
734 inline Double_t * AliKFParticle::Parameters ()
735 {
736   return fP;
737 }
738
739 inline Double_t * AliKFParticle::CovarianceMatrix()
740 {
741   return fC;
742 }
743
744
745 inline void AliKFParticle::operator +=( const AliKFParticle &Daughter )
746 {
747   AliKFParticleBase::operator +=( Daughter );
748 }
749   
750
751 inline void AliKFParticle::AddDaughter( const AliKFParticle &Daughter )
752 {
753   AliKFParticleBase::AddDaughter( Daughter );
754 }
755
756 inline void AliKFParticle::SetProductionVertex( const AliKFParticle &Vtx )
757 {
758   AliKFParticleBase::SetProductionVertex( Vtx );
759 }
760
761 inline void AliKFParticle::SetMassConstraint( Double_t Mass, Double_t SigmaMass )
762 {
763   AliKFParticleBase::SetMassConstraint( Mass, SigmaMass );
764 }
765     
766 inline void AliKFParticle::SetNoDecayLength()
767 {
768   AliKFParticleBase::SetNoDecayLength();
769 }
770
771 inline void AliKFParticle::Construct( const AliKFParticle *vDaughters[], int NDaughters, 
772                                const AliKFParticle *ProdVtx,   Double_t Mass, Bool_t IsConstrained  )
773 {    
774   AliKFParticleBase::Construct( ( const AliKFParticleBase**)vDaughters, NDaughters, 
775                          ( const AliKFParticleBase*)ProdVtx, Mass, IsConstrained );
776 }
777
778 inline void AliKFParticle::TransportToDecayVertex()
779
780   AliKFParticleBase::TransportToDecayVertex(); 
781 }
782
783 inline void AliKFParticle::TransportToProductionVertex()
784 {
785   AliKFParticleBase::TransportToProductionVertex();
786 }
787
788 inline void AliKFParticle::TransportToPoint( const Double_t xyz[] )
789
790   TransportToDS( GetDStoPoint(xyz) );
791 }
792
793 inline void AliKFParticle::TransportToVertex( const AliVVertex &v )
794 {       
795   TransportToPoint( AliKFParticle(v).fP );
796 }
797
798 inline void AliKFParticle::TransportToParticle( const AliKFParticle &p )
799
800   Double_t dS, dSp;
801   GetDStoParticle( p, dS, dSp );
802   TransportToDS( dS );
803 }
804
805 inline void AliKFParticle::TransportToDS( Double_t dS )
806 {
807   AliKFParticleBase::TransportToDS( dS );
808
809
810 inline Double_t AliKFParticle::GetDStoPoint( const Double_t xyz[] ) const 
811 {
812   return AliKFParticleBase::GetDStoPointBz( GetFieldAlice(), xyz );
813 }
814
815   
816 inline void AliKFParticle::GetDStoParticle( const AliKFParticle &p, 
817                                             Double_t &DS, Double_t &DSp ) const 
818 {
819   GetDStoParticleXY( p, DS, DSp );
820 }
821
822
823 inline Double_t AliKFParticle::GetDistanceFromVertex( const Double_t vtx[] ) const
824 {
825   return AliKFParticleBase::GetDistanceFromVertex( vtx );
826 }
827
828 inline Double_t AliKFParticle::GetDeviationFromVertex( const Double_t v[], 
829                                                        const Double_t Cv[] ) const
830 {
831   return AliKFParticleBase::GetDeviationFromVertex( v, Cv);
832 }
833
834 inline Double_t AliKFParticle::GetDistanceFromVertex( const AliKFParticle &Vtx ) const
835 {
836   return AliKFParticleBase::GetDistanceFromVertex( Vtx );
837 }
838
839 inline Double_t AliKFParticle::GetDeviationFromVertex( const AliKFParticle &Vtx ) const
840 {
841   return AliKFParticleBase::GetDeviationFromVertex( Vtx );
842 }
843
844 inline Double_t AliKFParticle::GetDistanceFromVertex( const AliVVertex &Vtx ) const
845 {
846   return GetDistanceFromVertex( AliKFParticle(Vtx) );
847 }
848
849 inline Double_t AliKFParticle::GetDeviationFromVertex( const AliVVertex &Vtx ) const
850 {
851   return GetDeviationFromVertex( AliKFParticle(Vtx) );
852 }
853   
854 inline Double_t AliKFParticle::GetDistanceFromParticle( const AliKFParticle &p ) const 
855 {
856   return AliKFParticleBase::GetDistanceFromParticle( p );
857 }
858
859 inline Double_t AliKFParticle::GetDeviationFromParticle( const AliKFParticle &p ) const 
860 {
861   return AliKFParticleBase::GetDeviationFromParticle( p );
862 }
863
864 inline void AliKFParticle::SubtractFromVertex( AliKFParticle &v ) const 
865 {
866   AliKFParticleBase::SubtractFromVertex( v.fP, v.fC, v.fChi2, v.fNDF);
867 }
868
869 inline void AliKFParticle::SubtractFromVertex( AliESDVertex &v ) const 
870 {
871   AliKFParticle vTmp(v);
872   SubtractFromVertex( vTmp );
873   v = AliESDVertex( vTmp.fP, vTmp.fC, vTmp.fChi2, (vTmp.fNDF +3)/2, v.GetName() );
874 }
875  
876 inline void AliKFParticle::CopyToESDVertex( AliESDVertex &v ) const 
877 {
878   AliKFParticle vTmp=*this;
879   v = AliESDVertex( vTmp.fP, vTmp.fC, vTmp.fChi2, (vTmp.fNDF +3)/2 );
880 }
881
882 inline Double_t AliKFParticle::GetFieldAlice()
883
884   return fgBz; 
885 }
886
887 inline void AliKFParticle::GetFieldValue( const Double_t * /*xyz*/, Double_t B[] ) const 
888 {    
889   B[0] = B[1] = 0;
890   B[2] = GetFieldAlice();
891 }
892
893 inline void AliKFParticle::GetDStoParticle( const AliKFParticleBase &p, 
894                                             Double_t &DS, Double_t &DSp )const
895 {
896   GetDStoParticleXY( p, DS, DSp );
897 }
898
899 inline void AliKFParticle::GetDStoParticleXY( const AliKFParticleBase &p, 
900                                        Double_t &DS, Double_t &DSp ) const
901
902   AliKFParticleBase::GetDStoParticleBz( GetFieldAlice(), p, DS, DSp ) ;
903   //GetDStoParticleALICE( p, DS, DSp ) ;
904 }
905
906 inline void AliKFParticle::Transport( Double_t dS, Double_t P[], Double_t C[] ) const 
907 {
908   AliKFParticleBase::TransportBz( GetFieldAlice(), dS, P, C );
909 }
910
911 #endif