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