]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliKFParticle.h
Coverity fix.
[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 );
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
346 inline AliKFParticle::AliKFParticle( const AliKFParticle &d1, 
347                                      const AliKFParticle &d2 )
348 {
349   AliKFParticle mother;
350   mother+= d1;
351   mother+= d2;
352   *this = mother;
353 }
354
355 inline AliKFParticle::AliKFParticle( const AliKFParticle &d1, 
356                                      const AliKFParticle &d2, 
357                                      const AliKFParticle &d3 )
358 {
359   AliKFParticle mother;
360   mother+= d1;
361   mother+= d2;
362   mother+= d3;
363   *this = mother;
364 }
365
366 inline AliKFParticle::AliKFParticle( const AliKFParticle &d1, 
367                                      const AliKFParticle &d2, 
368                                      const AliKFParticle &d3, 
369                                      const AliKFParticle &d4 )
370 {
371   AliKFParticle mother;
372   mother+= d1;
373   mother+= d2;
374   mother+= d3;
375   mother+= d4;
376   *this = mother;
377 }
378
379
380 inline void AliKFParticle::Initialize()
381
382   AliKFParticleBase::Initialize(); 
383 }
384
385 inline void AliKFParticle::SetVtxGuess( Double_t x, Double_t y, Double_t z )
386 {
387   AliKFParticleBase::SetVtxGuess(x,y,z);
388 }  
389
390 inline Double_t AliKFParticle::GetX    () const 
391
392   return AliKFParticleBase::GetX();    
393 }
394
395 inline Double_t AliKFParticle::GetY    () const 
396
397   return AliKFParticleBase::GetY();    
398 }
399
400 inline Double_t AliKFParticle::GetZ    () const 
401
402   return AliKFParticleBase::GetZ();    
403 }
404
405 inline Double_t AliKFParticle::GetPx   () const 
406
407   return AliKFParticleBase::GetPx();   
408 }
409
410 inline Double_t AliKFParticle::GetPy   () const 
411
412   return AliKFParticleBase::GetPy();   
413 }
414
415 inline Double_t AliKFParticle::GetPz   () const 
416
417   return AliKFParticleBase::GetPz();   
418 }
419
420 inline Double_t AliKFParticle::GetE    () const 
421
422   return AliKFParticleBase::GetE();    
423 }
424
425 inline Double_t AliKFParticle::GetS    () const 
426
427   return AliKFParticleBase::GetS();    
428 }
429
430 inline Int_t    AliKFParticle::GetQ    () const 
431
432   return AliKFParticleBase::GetQ();    
433 }
434
435 inline Double_t AliKFParticle::GetChi2 () const 
436
437   return AliKFParticleBase::GetChi2(); 
438 }
439
440 inline Int_t    AliKFParticle::GetNDF  () const 
441
442   return AliKFParticleBase::GetNDF();  
443 }
444
445 inline Double_t AliKFParticle::GetParameter ( int i ) const 
446
447   return AliKFParticleBase::GetParameter(i);  
448 }
449
450 inline Double_t AliKFParticle::GetCovariance( int i ) const 
451
452   return AliKFParticleBase::GetCovariance(i); 
453 }
454
455 inline Double_t AliKFParticle::GetCovariance( int i, int j ) const 
456
457   return AliKFParticleBase::GetCovariance(i,j);
458 }
459
460
461 inline Double_t AliKFParticle::GetP    () const
462 {
463   Double_t par, err;
464   if( AliKFParticleBase::GetMomentum( par, err ) ) return 0;
465   else return par;
466 }
467
468 inline Double_t AliKFParticle::GetPt   () const
469 {
470   Double_t par, err;
471   if( AliKFParticleBase::GetPt( par, err ) ) return 0;
472   else return par;
473 }
474
475 inline Double_t AliKFParticle::GetEta   () const
476 {
477   Double_t par, err;
478   if( AliKFParticleBase::GetEta( par, err ) ) return 0;
479   else return par;
480 }
481
482 inline Double_t AliKFParticle::GetPhi   () const
483 {
484   Double_t par, err;
485   if( AliKFParticleBase::GetPhi( par, err ) ) return 0;
486   else return par;
487 }
488
489 inline Double_t AliKFParticle::GetMomentum    () const
490 {
491   Double_t par, err;
492   if( AliKFParticleBase::GetMomentum( par, err ) ) return 0;
493   else return par;
494 }
495
496 inline Double_t AliKFParticle::GetMass        () const
497 {
498   Double_t par, err;
499   if( AliKFParticleBase::GetMass( par, err ) ) return 0;
500   else return par;
501 }
502
503 inline Double_t AliKFParticle::GetDecayLength () const
504 {
505   Double_t par, err;
506   if( AliKFParticleBase::GetDecayLength( par, err ) ) return 0;
507   else return par;
508 }
509
510 inline Double_t AliKFParticle::GetDecayLengthXY () const
511 {
512   Double_t par, err;
513   if( AliKFParticleBase::GetDecayLengthXY( par, err ) ) return 0;
514   else return par;
515 }
516
517 inline Double_t AliKFParticle::GetLifeTime    () const
518 {
519   Double_t par, err;
520   if( AliKFParticleBase::GetLifeTime( par, err ) ) return 0;
521   else return par;
522 }
523
524 inline Double_t AliKFParticle::GetR   () const
525 {
526   Double_t par, err;
527   if( AliKFParticleBase::GetR( par, err ) ) return 0;
528   else return par;
529 }
530
531 inline Double_t AliKFParticle::GetErrX           () const 
532 {
533   return TMath::Sqrt(TMath::Abs( GetCovariance(0,0) ));
534 }
535
536 inline Double_t AliKFParticle::GetErrY           () const 
537 {
538   return TMath::Sqrt(TMath::Abs( GetCovariance(1,1) ));
539 }
540
541 inline Double_t AliKFParticle::GetErrZ           () const 
542 {
543   return TMath::Sqrt(TMath::Abs( GetCovariance(2,2) ));
544 }
545
546 inline Double_t AliKFParticle::GetErrPx          () const 
547 {
548   return TMath::Sqrt(TMath::Abs( GetCovariance(3,3) ));
549 }
550
551 inline Double_t AliKFParticle::GetErrPy          () const 
552 {
553   return TMath::Sqrt(TMath::Abs( GetCovariance(4,4) ));
554 }
555
556 inline Double_t AliKFParticle::GetErrPz          () const 
557 {
558   return TMath::Sqrt(TMath::Abs( GetCovariance(5,5) ));
559 }
560
561 inline Double_t AliKFParticle::GetErrE           () const 
562 {
563   return TMath::Sqrt(TMath::Abs( GetCovariance(6,6) ));
564 }
565
566 inline Double_t AliKFParticle::GetErrS           () const 
567 {
568   return TMath::Sqrt(TMath::Abs( GetCovariance(7,7) ));
569 }
570
571 inline Double_t AliKFParticle::GetErrP    () const
572 {
573   Double_t par, err;
574   if( AliKFParticleBase::GetMomentum( par, err ) ) return 1.e10;
575   else return err;
576 }
577
578 inline Double_t AliKFParticle::GetErrPt    () const
579 {
580   Double_t par, err;
581   if( AliKFParticleBase::GetPt( par, err ) ) return 1.e10;
582   else return err;
583 }
584
585 inline Double_t AliKFParticle::GetErrEta    () const
586 {
587   Double_t par, err;
588   if( AliKFParticleBase::GetEta( par, err ) ) return 1.e10;
589   else return err;
590 }
591
592 inline Double_t AliKFParticle::GetErrPhi    () const
593 {
594   Double_t par, err;
595   if( AliKFParticleBase::GetPhi( par, err ) ) return 1.e10;
596   else return err;
597 }
598
599 inline Double_t AliKFParticle::GetErrMomentum    () const
600 {
601   Double_t par, err;
602   if( AliKFParticleBase::GetMomentum( par, err ) ) return 1.e10;
603   else return err;
604 }
605
606 inline Double_t AliKFParticle::GetErrMass        () const
607 {
608   Double_t par, err;
609   if( AliKFParticleBase::GetMass( par, err ) ) return 1.e10;
610   else return err;
611 }
612
613 inline Double_t AliKFParticle::GetErrDecayLength () const
614 {
615   Double_t par, err;
616   if( AliKFParticleBase::GetDecayLength( par, err ) ) return 1.e10;
617   else return err;
618 }
619
620 inline Double_t AliKFParticle::GetErrDecayLengthXY () const
621 {
622   Double_t par, err;
623   if( AliKFParticleBase::GetDecayLengthXY( par, err ) ) return 1.e10;
624   else return err;
625 }
626
627 inline Double_t AliKFParticle::GetErrLifeTime    () const
628 {
629   Double_t par, err;
630   if( AliKFParticleBase::GetLifeTime( par, err ) ) return 1.e10;
631   else return err;
632 }
633
634 inline Double_t AliKFParticle::GetErrR    () const
635 {
636   Double_t par, err;
637   if( AliKFParticleBase::GetR( par, err ) ) return 1.e10;
638   else return err;
639 }
640
641
642 inline int AliKFParticle::GetP( Double_t &P, Double_t &SigmaP ) const 
643 {
644   return AliKFParticleBase::GetMomentum( P, SigmaP );
645 }
646
647 inline int AliKFParticle::GetPt( Double_t &Pt, Double_t &SigmaPt ) const 
648 {
649   return AliKFParticleBase::GetPt( Pt, SigmaPt );
650 }
651
652 inline int AliKFParticle::GetEta( Double_t &Eta, Double_t &SigmaEta ) const 
653 {
654   return AliKFParticleBase::GetEta( Eta, SigmaEta );
655 }
656
657 inline int AliKFParticle::GetPhi( Double_t &Phi, Double_t &SigmaPhi ) const 
658 {
659   return AliKFParticleBase::GetPhi( Phi, SigmaPhi );
660 }
661
662 inline int AliKFParticle::GetMomentum( Double_t &P, Double_t &SigmaP ) const 
663 {
664   return AliKFParticleBase::GetMomentum( P, SigmaP );
665 }
666
667 inline int AliKFParticle::GetMass( Double_t &M, Double_t &SigmaM ) const 
668 {
669   return AliKFParticleBase::GetMass( M, SigmaM );
670 }
671
672 inline int AliKFParticle::GetDecayLength( Double_t &L, Double_t &SigmaL ) const 
673 {
674   return AliKFParticleBase::GetDecayLength( L, SigmaL );
675 }
676
677 inline int AliKFParticle::GetDecayLengthXY( Double_t &L, Double_t &SigmaL ) const 
678 {
679   return AliKFParticleBase::GetDecayLengthXY( L, SigmaL );
680 }
681
682 inline int AliKFParticle::GetLifeTime( Double_t &T, Double_t &SigmaT ) const 
683 {
684   return AliKFParticleBase::GetLifeTime( T, SigmaT );
685 }
686
687 inline int AliKFParticle::GetR( Double_t &R, Double_t &SigmaR ) const 
688 {
689   return AliKFParticleBase::GetR( R, SigmaR );
690 }
691
692 inline Double_t & AliKFParticle::X() 
693
694   return AliKFParticleBase::X();    
695 }
696
697 inline Double_t & AliKFParticle::Y()
698
699   return AliKFParticleBase::Y();    
700 }
701
702 inline Double_t & AliKFParticle::Z() 
703
704   return AliKFParticleBase::Z();    
705 }
706
707 inline Double_t & AliKFParticle::Px() 
708
709   return AliKFParticleBase::Px();   
710 }
711
712 inline Double_t & AliKFParticle::Py() 
713
714   return AliKFParticleBase::Py();   
715 }
716
717 inline Double_t & AliKFParticle::Pz() 
718
719   return AliKFParticleBase::Pz();   
720 }
721
722 inline Double_t & AliKFParticle::E() 
723
724   return AliKFParticleBase::E();    
725 }
726
727 inline Double_t & AliKFParticle::S() 
728
729   return AliKFParticleBase::S();    
730 }
731
732 inline Int_t    & AliKFParticle::Q() 
733
734   return AliKFParticleBase::Q();    
735 }
736
737 inline Double_t & AliKFParticle::Chi2() 
738
739   return AliKFParticleBase::Chi2(); 
740 }
741
742 inline Int_t    & AliKFParticle::NDF() 
743
744   return AliKFParticleBase::NDF();  
745 }
746
747 inline Double_t & AliKFParticle::Parameter ( int i )        
748
749   return AliKFParticleBase::Parameter(i);
750 }
751
752 inline Double_t & AliKFParticle::Covariance( int i )        
753
754   return AliKFParticleBase::Covariance(i);
755 }
756
757 inline Double_t & AliKFParticle::Covariance( int i, int j ) 
758
759   return AliKFParticleBase::Covariance(i,j); 
760 }
761
762 inline Double_t * AliKFParticle::Parameters ()
763 {
764   return fP;
765 }
766
767 inline Double_t * AliKFParticle::CovarianceMatrix()
768 {
769   return fC;
770 }
771
772
773 inline void AliKFParticle::operator +=( const AliKFParticle &Daughter )
774 {
775   AliKFParticleBase::operator +=( Daughter );
776 }
777   
778
779 inline void AliKFParticle::AddDaughter( const AliKFParticle &Daughter )
780 {
781   AliKFParticleBase::AddDaughter( Daughter );
782 }
783
784 inline void AliKFParticle::SetProductionVertex( const AliKFParticle &Vtx )
785 {
786   AliKFParticleBase::SetProductionVertex( Vtx );
787 }
788
789 inline void AliKFParticle::SetMassConstraint( Double_t Mass, Double_t SigmaMass )
790 {
791   AliKFParticleBase::SetMassConstraint( Mass, SigmaMass );
792 }
793     
794 inline void AliKFParticle::SetNoDecayLength()
795 {
796   AliKFParticleBase::SetNoDecayLength();
797 }
798
799 inline void AliKFParticle::Construct( const AliKFParticle *vDaughters[], int NDaughters, 
800                                const AliKFParticle *ProdVtx,   Double_t Mass, Bool_t IsConstrained  )
801 {    
802   AliKFParticleBase::Construct( ( const AliKFParticleBase**)vDaughters, NDaughters, 
803                          ( const AliKFParticleBase*)ProdVtx, Mass, IsConstrained );
804 }
805
806 inline void AliKFParticle::TransportToDecayVertex()
807
808   AliKFParticleBase::TransportToDecayVertex(); 
809 }
810
811 inline void AliKFParticle::TransportToProductionVertex()
812 {
813   AliKFParticleBase::TransportToProductionVertex();
814 }
815
816 inline void AliKFParticle::TransportToPoint( const Double_t xyz[] )
817
818   TransportToDS( GetDStoPoint(xyz) );
819 }
820
821 inline void AliKFParticle::TransportToVertex( const AliVVertex &v )
822 {       
823   TransportToPoint( AliKFParticle(v).fP );
824 }
825
826 inline void AliKFParticle::TransportToParticle( const AliKFParticle &p )
827
828   Double_t dS, dSp;
829   GetDStoParticle( p, dS, dSp );
830   TransportToDS( dS );
831 }
832
833 inline void AliKFParticle::TransportToDS( Double_t dS )
834 {
835   AliKFParticleBase::TransportToDS( dS );
836
837
838 inline Double_t AliKFParticle::GetDStoPoint( const Double_t xyz[] ) const 
839 {
840   return AliKFParticleBase::GetDStoPointBz( GetFieldAlice(), xyz );
841 }
842
843   
844 inline void AliKFParticle::GetDStoParticle( const AliKFParticle &p, 
845                                             Double_t &DS, Double_t &DSp ) const 
846 {
847   GetDStoParticleXY( p, DS, DSp );
848 }
849
850
851 inline Double_t AliKFParticle::GetDistanceFromVertex( const Double_t vtx[] ) const
852 {
853   return AliKFParticleBase::GetDistanceFromVertex( vtx );
854 }
855
856 inline Double_t AliKFParticle::GetDeviationFromVertex( const Double_t v[], 
857                                                        const Double_t Cv[] ) const
858 {
859   return AliKFParticleBase::GetDeviationFromVertex( v, Cv);
860 }
861
862 inline Double_t AliKFParticle::GetDistanceFromVertex( const AliKFParticle &Vtx ) const
863 {
864   return AliKFParticleBase::GetDistanceFromVertex( Vtx );
865 }
866
867 inline Double_t AliKFParticle::GetDeviationFromVertex( const AliKFParticle &Vtx ) const
868 {
869   return AliKFParticleBase::GetDeviationFromVertex( Vtx );
870 }
871
872 inline Double_t AliKFParticle::GetDistanceFromVertex( const AliVVertex &Vtx ) const
873 {
874   return GetDistanceFromVertex( AliKFParticle(Vtx) );
875 }
876
877 inline Double_t AliKFParticle::GetDeviationFromVertex( const AliVVertex &Vtx ) const
878 {
879   return GetDeviationFromVertex( AliKFParticle(Vtx) );
880 }
881  
882 inline Double_t AliKFParticle::GetDistanceFromParticle( const AliKFParticle &p ) const 
883 {
884   return AliKFParticleBase::GetDistanceFromParticle( p );
885 }
886
887 inline Double_t AliKFParticle::GetDeviationFromParticle( const AliKFParticle &p ) const 
888 {
889   return AliKFParticleBase::GetDeviationFromParticle( p );
890 }
891
892 inline void AliKFParticle::SubtractFromVertex( AliKFParticle &v ) const 
893 {
894   AliKFParticleBase::SubtractFromVertex( v );
895 }
896
897 inline Double_t AliKFParticle::GetFieldAlice()
898
899   return fgBz; 
900 }
901
902 inline void AliKFParticle::GetFieldValue( const Double_t * /*xyz*/, Double_t B[] ) const 
903 {    
904   B[0] = B[1] = 0;
905   B[2] = GetFieldAlice();
906 }
907
908 inline void AliKFParticle::GetDStoParticle( const AliKFParticleBase &p, 
909                                             Double_t &DS, Double_t &DSp )const
910 {
911   GetDStoParticleXY( p, DS, DSp );
912 }
913
914 inline void AliKFParticle::GetDStoParticleXY( const AliKFParticleBase &p, 
915                                        Double_t &DS, Double_t &DSp ) const
916
917   AliKFParticleBase::GetDStoParticleBz( GetFieldAlice(), p, DS, DSp ) ;
918   //GetDStoParticleALICE( p, DS, DSp ) ;
919 }
920
921 inline void AliKFParticle::Transport( Double_t dS, Double_t P[], Double_t C[] ) const 
922 {
923   AliKFParticleBase::TransportBz( GetFieldAlice(), dS, P, C );
924 }
925
926 inline void AliKFParticle::ConstructGamma( const AliKFParticle &daughter1,
927                                            const AliKFParticle &daughter2  )
928 {
929   AliKFParticleBase::ConstructGammaBz( daughter1, daughter2, GetFieldAlice() );
930 }
931
932 #endif