Analysis Software
Documentation for sPHENIX simulation software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
KFParticle.h
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file KFParticle.h
1 /*
2  * This file is part of KFParticle package
3  * Copyright (C) 2007-2019 FIAS Frankfurt Institute for Advanced Studies
4  * 2007-2019 Goethe University of Frankfurt
5  * 2007-2019 Ivan Kisel <I.Kisel@compeng.uni-frankfurt.de>
6  * 2007-2019 Maksym Zyzak
7  * 2007-2019 Sergey Gorbunov
8  *
9  * KFParticle is free software: you can redistribute it and/or modify
10  * it under the terms of the GNU General Public License as published by
11  * the Free Software Foundation, either version 3 of the License, or
12  * (at your option) any later version.
13  *
14  * KFParticle is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17  * GNU General Public License for more details.
18  *
19  * You should have received a copy of the GNU General Public License
20  * along with this program. If not, see <https://www.gnu.org/licenses/>.
21  */
22 
23 //#define NonhomogeneousField
24 // #define HomogeneousField
25 
26 #ifndef KFPARTICLE_H
27 #define KFPARTICLE_H
28 
29 #include "KFParticleBase.h"
30 #include <cmath>
31 
32 class KFPTrack;
33 class KFPVertex;
34 
53 {
54 
55  public:
56 
57  //*
58  //* INITIALIZATION
59  //*
60 
61  //* Set magnetic field for all particles
62 #ifdef HomogeneousField
63  static void SetField( float Bz );
64 #endif
65  //* Constructor (empty)
66 
68 
69  //* Destructor (empty)
70 
71  ~KFParticle(){ ; }
72 
73  //* Construction of mother particle by its 2-3-4 daughters
74 
75  KFParticle( const KFParticle &d1, const KFParticle &d2 );
76 
77  KFParticle( const KFParticle &d1, const KFParticle &d2,
78  const KFParticle &d3 );
79 
80  KFParticle( const KFParticle &d1, const KFParticle &d2,
81  const KFParticle &d3, const KFParticle &d4 );
82 
83  //* Initialisation from "cartesian" coordinates ( X Y Z Px Py Pz )
84  //* Parameters, covariance matrix, charge and PID hypothesis should be provided
85 
86  void Create( const float Param[], const float Cov[], Int_t Charge, float mass /*Int_t PID*/ );
87  void Create( const Double_t Param[], const Double_t Cov[], Int_t Charge, float mass /*Int_t PID*/ );
88 
89  //* Initialisation from ALICE track, PID hypothesis shoould be provided
90 
91  KFParticle( const KFPTrack &track, const int PID );
92 
93 
94  //* Initialisation from VVertex
95 
96  KFParticle( const KFPVertex &vertex );
97 
98  //* Initialise covariance matrix and set current parameters to 0.0
99 
100  void Initialize();
101 
102  //*
103  //* ACCESSORS
104  //*
105 
106  //* Simple accessors
107 
108  float GetX () const ;
109  float GetY () const ;
110  float GetZ () const ;
111  float GetPx () const ;
112  float GetPy () const ;
113  float GetPz () const ;
114  float GetE () const ;
115  float GetS () const ;
116  char GetQ () const ;
117  float GetChi2 () const ;
118  Int_t GetNDF () const ;
119 
120  Bool_t GetAtProductionVertex() const { return fAtProductionVertex; }
122 
123 #ifdef NonhomogeneousField
124  const float* GetFieldCoeff() const { return fieldRegion; }
125  void SetFieldCoeff(float c, int i) { fieldRegion[i] = c; }
126 #endif
127 
128  const float& X () const { return fP[0]; }
129  const float& Y () const { return fP[1]; }
130  const float& Z () const { return fP[2]; }
131  const float& Px () const { return fP[3]; }
132  const float& Py () const { return fP[4]; }
133  const float& Pz () const { return fP[5]; }
134  const float& E () const { return fP[6]; }
135  const float& S () const { return fP[7]; }
136  const char& Q () const { return fQ; }
137  const float& Chi2 () const { return fChi2; }
138  const Int_t& NDF () const { return fNDF; }
139 
140  float GetParameter ( int i ) const ;
141  float GetCovariance( int i ) const ;
142  float GetCovariance( int i, int j ) const ;
143 
144  //* Accessors with calculations, value returned w/o error flag
145 
146  float GetP () const;
147  float GetPt () const;
148  float GetEta () const;
149  float GetPhi () const;
150  float GetMomentum () const;
151  float GetMass () const;
152  float GetDecayLength () const;
153  float GetDecayLengthXY () const;
154  float GetLifeTime () const;
155  float GetR () const;
156 
157  //* Accessors to estimated errors
158 
159  float GetErrX () const ;
160  float GetErrY () const ;
161  float GetErrZ () const ;
162  float GetErrPx () const ;
163  float GetErrPy () const ;
164  float GetErrPz () const ;
165  float GetErrE () const ;
166  float GetErrS () const ;
167  float GetErrP () const ;
168  float GetErrPt () const ;
169  float GetErrEta () const ;
170  float GetErrPhi () const ;
171  float GetErrMomentum () const ;
172  float GetErrMass () const ;
173  float GetErrDecayLength () const ;
174  float GetErrDecayLengthXY () const ;
175  float GetErrLifeTime () const ;
176  float GetErrR () const ;
177 
178  //* Accessors with calculations( &value, &estimated sigma )
179  //* error flag returned (0 means no error during calculations)
180 
181  int GetP ( float &P, float &SigmaP ) const ; //* momentum
182  int GetPt ( float &Pt, float &SigmaPt ) const ; //* transverse momentum
183  int GetEta ( float &Eta, float &SigmaEta ) const ; //* pseudorapidity
184  int GetPhi ( float &Phi, float &SigmaPhi ) const ; //* phi
185  int GetMomentum ( float &P, float &SigmaP ) const ; //* momentum
186  int GetMass ( float &M, float &SigmaM ) const ; //* mass
187  int GetDecayLength ( float &L, float &SigmaL ) const ; //* decay length
188  int GetDecayLengthXY ( float &L, float &SigmaL ) const ; //* decay length in XY
189  int GetLifeTime ( float &T, float &SigmaT ) const ; //* life time
190  int GetR ( float &R, float &SigmaR ) const ; //* R
191  float GetRapidity() const { return 0.5*log((fP[6] + fP[5])/(fP[6] - fP[5])); }
192  float GetTheta() const { return atan2(GetPt(),fP[5]); }
193 
194 
195  //*
196  //* MODIFIERS
197  //*
198 
199  float & X () ;
200  float & Y () ;
201  float & Z () ;
202  float & Px () ;
203  float & Py () ;
204  float & Pz () ;
205  float & E () ;
206  float & S () ;
207  char & Q () ;
208  float & Chi2 () ;
209  Int_t & NDF () ;
210 
211  float & Parameter ( int i ) ;
212  float & Covariance( int i ) ;
213  float & Covariance( int i, int j ) ;
214  float * Parameters () ;
215  float * CovarianceMatrix() ;
216 
217  //*
218  //* CONSTRUCTION OF THE PARTICLE BY ITS DAUGHTERS AND MOTHER
219  //* USING THE KALMAN FILTER METHOD
220  //*
221 
222 
223  //* Add daughter to the particle
224 
225  void AddDaughter( const KFParticle &Daughter );
226 
227  //* Add daughter via += operator: ex.{ D0; D0+=Pion; D0+= Kaon; }
228 
229  void operator +=( const KFParticle &Daughter );
230 
231  //* Everything in one go
232 
233  void Construct( const KFParticle *vDaughters[], int nDaughters,
234  const KFParticle *ProdVtx=0, float Mass=-1 );
235 
236  //*
237  //* TRANSPORT
238  //*
239  //* ( main transportation parameter is S = SignedPath/Momentum )
240  //* ( parameters of decay & production vertices are stored locally )
241  //*
242 
243  //* Transport the particle close to xyz[] point
244 
245  void TransportToPoint( const float xyz[] );
246 
247  //* Transport the particle close to VVertex
248 #ifdef HomogeneousField
249  void TransportToVertex( const KFPVertex &v );
250 #endif
251  //* Transport the particle close to another particle p
252  void TransportToParticle( const KFParticle &p );
253 
254  //* Get dS to a certain space point
255  float GetDStoPoint( const float xyz[3], float dsdr[6] ) const ;
256 
257  //* Get dS to other particle p (dSp for particle p also returned)
258  void GetDStoParticle( const KFParticleBase &p, float dS[2], float dsdr[4][6] ) const ;
259 
260 
261  //*
262  //* OTHER UTILITIES
263  //*
264 
265  //* Calculate distance from another object [cm] in XY-plane
266 
267  Bool_t GetDistanceFromVertexXY( const float vtx[], float &val, float &err ) const ;
268  Bool_t GetDistanceFromVertexXY( const float vtx[], const float Cv[], float &val, float &err ) const ;
269  Bool_t GetDistanceFromVertexXY( const KFParticle &Vtx, float &val, float &err ) const ;
270 #ifdef HomogeneousField
271  Bool_t GetDistanceFromVertexXY( const KFPVertex &Vtx, float &val, float &err ) const ;
272 #endif
273 
274  float GetDistanceFromVertexXY( const float vtx[] ) const ;
275  float GetDistanceFromVertexXY( const KFParticle &Vtx ) const ;
276 #ifdef HomogeneousField
277  float GetDistanceFromVertexXY( const KFPVertex &Vtx ) const ;
278 #endif
279  float GetDistanceFromParticleXY( const KFParticle &p ) const ;
280 
281  //* Calculate sqrt(Chi2/ndf) deviation from another object in XY plane
282  //* ( v = [xyz]-vertex, Cv=[Cxx,Cxy,Cyy,Cxz,Cyz,Czz]-covariance matrix )
283 
284  float GetDeviationFromVertexXY( const float v[], const float Cv[]=0 ) const ;
285  float GetDeviationFromVertexXY( const KFParticle &Vtx ) const ;
286 #ifdef HomogeneousField
287  float GetDeviationFromVertexXY( const KFPVertex &Vtx ) const ;
288 #endif
289  float GetDeviationFromParticleXY( const KFParticle &p ) const ;
290 
291  //* Get parameters at an arbitrary reconstructed point taking into account its errors
292  void GetParametersAtPoint(const float* point, const float* pointCov, float* m, float* mV);
293 
294  //* Calculate opennig angle between two particles
295 
296  float GetAngle ( const KFParticle &p ) const ;
297  float GetAngleXY( const KFParticle &p ) const ;
298  float GetAngleRZ( const KFParticle &p ) const ;
299 
300  float GetPseudoProperDecayTime( const KFParticle &primVertex, const float& mass, float* timeErr2 = 0 ) const;
301 
302  void GetFieldValue( const float xyz[], float B[] ) const ;
303 
304  void Transport( float dS, const float* dsdr, float P[], float C[], float* dsdr1=0, float* F=0, float* F1=0 ) const ;
305 
306  protected:
307 
308  //*
309  //* INTERNAL STUFF
310  //*
311 
312  //* Method to access ALICE field
313 #ifdef HomogeneousField
314  static float GetFieldAlice();
315 #endif
316 
317  private:
318 #ifdef HomogeneousField
319  static float fgBz;
320 #endif
321 #ifdef NonhomogeneousField
322 
325  float fieldRegion[10];
326 #endif
327 
328 #ifndef KFParticleStandalone
329  ClassDef( KFParticle, 3 )
330 #endif
331 };
332 
333 
334 
335 //---------------------------------------------------------------------
336 //
337 // Inline implementation of the KFParticle methods
338 //
339 //---------------------------------------------------------------------
340 
341 #ifdef HomogeneousField
342 inline void KFParticle::SetField( float Bz )
343 {
347  fgBz = Bz;
348 }
349 #endif
350 
352  const KFParticle &d2,
353  const KFParticle &d3 )
354 {
360  KFParticle mother;
361  mother+= d1;
362  mother+= d2;
363  mother+= d3;
364  *this = mother;
365 }
366 
368  const KFParticle &d2,
369  const KFParticle &d3,
370  const KFParticle &d4 )
371 {
378  KFParticle mother;
379  mother+= d1;
380  mother+= d2;
381  mother+= d3;
382  mother+= d4;
383  *this = mother;
384 }
385 
386 
388 {
391 }
392 
393 inline float KFParticle::GetX () const
394 {
395  return KFParticleBase::GetX();
396 }
397 
398 inline float KFParticle::GetY () const
399 {
400  return KFParticleBase::GetY();
401 }
402 
403 inline float KFParticle::GetZ () const
404 {
405  return KFParticleBase::GetZ();
406 }
407 
408 inline float KFParticle::GetPx () const
409 {
410  return KFParticleBase::GetPx();
411 }
412 
413 inline float KFParticle::GetPy () const
414 {
415  return KFParticleBase::GetPy();
416 }
417 
418 inline float KFParticle::GetPz () const
419 {
420  return KFParticleBase::GetPz();
421 }
422 
423 inline float KFParticle::GetE () const
424 {
425  return KFParticleBase::GetE();
426 }
427 
428 inline float KFParticle::GetS () const
429 {
430  return KFParticleBase::GetS();
431 }
432 
433 inline char KFParticle::GetQ () const
434 {
435  return KFParticleBase::GetQ();
436 }
437 
438 inline float KFParticle::GetChi2 () const
439 {
440  return KFParticleBase::GetChi2();
441 }
442 
443 inline Int_t KFParticle::GetNDF () const
444 {
445  return KFParticleBase::GetNDF();
446 }
447 
448 inline float KFParticle::GetParameter ( int i ) const
449 {
450  return KFParticleBase::GetParameter(i);
451 }
452 
453 inline float KFParticle::GetCovariance( int i ) const
454 {
456 }
457 
458 inline float KFParticle::GetCovariance( int i, int j ) const
459 {
460  return KFParticleBase::GetCovariance(i,j);
461 }
462 
463 
464 inline float KFParticle::GetP () const
465 {
466  float par, err;
467  if( KFParticleBase::GetMomentum( par, err ) ) return 0;
468  else return par;
469 }
470 
471 inline float KFParticle::GetPt () const
472 {
473  float par, err;
474  if( KFParticleBase::GetPt( par, err ) ) return 0;
475  else return par;
476 }
477 
478 inline float KFParticle::GetEta () const
479 {
480  float par, err;
481  if( KFParticleBase::GetEta( par, err ) ) return 0;
482  else return par;
483 }
484 
485 inline float KFParticle::GetPhi () const
486 {
487  float par, err;
488  if( KFParticleBase::GetPhi( par, err ) ) return 0;
489  else return par;
490 }
491 
492 inline float KFParticle::GetMomentum () const
493 {
494  float par, err;
495  if( KFParticleBase::GetMomentum( par, err ) ) return 0;
496  else return par;
497 }
498 
499 inline float KFParticle::GetMass () const
500 {
501  float par, err;
502  if( KFParticleBase::GetMass( par, err ) ) return 0;
503  else return par;
504 }
505 
506 inline float KFParticle::GetDecayLength () const
507 {
508  float par, err;
509  if( KFParticleBase::GetDecayLength( par, err ) ) return 0;
510  else return par;
511 }
512 
513 inline float KFParticle::GetDecayLengthXY () const
514 {
515  float par, err;
516  if( KFParticleBase::GetDecayLengthXY( par, err ) ) return 0;
517  else return par;
518 }
519 
520 inline float KFParticle::GetLifeTime () const
521 {
522  float par, err;
523  if( KFParticleBase::GetLifeTime( par, err ) ) return 0;
524  else return par;
525 }
526 
527 inline float KFParticle::GetR () const
528 {
529  float par, err;
530  if( KFParticleBase::GetR( par, err ) ) return 0;
531  else return par;
532 }
533 
534 inline float KFParticle::GetErrX () const
535 {
536  return sqrt(fabs( GetCovariance(0,0) ));
537 }
538 
539 inline float KFParticle::GetErrY () const
540 {
541  return sqrt(fabs( GetCovariance(1,1) ));
542 }
543 
544 inline float KFParticle::GetErrZ () const
545 {
546  return sqrt(fabs( GetCovariance(2,2) ));
547 }
548 
549 inline float KFParticle::GetErrPx () const
550 {
551  return sqrt(fabs( GetCovariance(3,3) ));
552 }
553 
554 inline float KFParticle::GetErrPy () const
555 {
556  return sqrt(fabs( GetCovariance(4,4) ));
557 }
558 
559 inline float KFParticle::GetErrPz () const
560 {
561  return sqrt(fabs( GetCovariance(5,5) ));
562 }
563 
564 inline float KFParticle::GetErrE () const
565 {
566  return sqrt(fabs( GetCovariance(6,6) ));
567 }
568 
569 inline float KFParticle::GetErrS () const
570 {
571  return sqrt(fabs( GetCovariance(7,7) ));
572 }
573 
574 inline float KFParticle::GetErrP () const
575 {
576  float par, err;
577  if( KFParticleBase::GetMomentum( par, err ) ) return 1.e10;
578  else return err;
579 }
580 
581 inline float KFParticle::GetErrPt () const
582 {
583  float par, err;
584  if( KFParticleBase::GetPt( par, err ) ) return 1.e10;
585  else return err;
586 }
587 
588 inline float KFParticle::GetErrEta () const
589 {
590  float par, err;
591  if( KFParticleBase::GetEta( par, err ) ) return 1.e10;
592  else return err;
593 }
594 
595 inline float KFParticle::GetErrPhi () const
596 {
597  float par, err;
598  if( KFParticleBase::GetPhi( par, err ) ) return 1.e10;
599  else return err;
600 }
601 
602 inline float KFParticle::GetErrMomentum () const
603 {
604  float par, err;
605  if( KFParticleBase::GetMomentum( par, err ) ) return 1.e10;
606  else return err;
607 }
608 
609 inline float KFParticle::GetErrMass () const
610 {
611  float par, err;
612  if( KFParticleBase::GetMass( par, err ) ) return 1.e10;
613  else return err;
614 }
615 
616 inline float KFParticle::GetErrDecayLength () const
617 {
618  float par, err;
619  if( KFParticleBase::GetDecayLength( par, err ) ) return 1.e10;
620  else return err;
621 }
622 
623 inline float KFParticle::GetErrDecayLengthXY () const
624 {
625  float par, err;
626  if( KFParticleBase::GetDecayLengthXY( par, err ) ) return 1.e10;
627  else return err;
628 }
629 
630 inline float KFParticle::GetErrLifeTime () const
631 {
632  float par, err;
633  if( KFParticleBase::GetLifeTime( par, err ) ) return 1.e10;
634  else return err;
635 }
636 
637 inline float KFParticle::GetErrR () const
638 {
639  float par, err;
640  if( KFParticleBase::GetR( par, err ) ) return 1.e10;
641  else return err;
642 }
643 
644 
645 inline int KFParticle::GetP( float &P, float &SigmaP ) const
646 {
651  return KFParticleBase::GetMomentum( P, SigmaP );
652 }
653 
654 inline int KFParticle::GetPt( float &Pt, float &SigmaPt ) const
655 {
660  return KFParticleBase::GetPt( Pt, SigmaPt );
661 }
662 
663 inline int KFParticle::GetEta( float &Eta, float &SigmaEta ) const
664 {
669  return KFParticleBase::GetEta( Eta, SigmaEta );
670 }
671 
672 inline int KFParticle::GetPhi( float &Phi, float &SigmaPhi ) const
673 {
678  return KFParticleBase::GetPhi( Phi, SigmaPhi );
679 }
680 
681 inline int KFParticle::GetMomentum( float &P, float &SigmaP ) const
682 {
687  return KFParticleBase::GetMomentum( P, SigmaP );
688 }
689 
690 inline int KFParticle::GetMass( float &M, float &SigmaM ) const
691 {
696  return KFParticleBase::GetMass( M, SigmaM );
697 }
698 
699 inline int KFParticle::GetDecayLength( float &L, float &SigmaL ) const
700 {
706  return KFParticleBase::GetDecayLength( L, SigmaL );
707 }
708 
709 inline int KFParticle::GetDecayLengthXY( float &L, float &SigmaL ) const
710 {
717  return KFParticleBase::GetDecayLengthXY( L, SigmaL );
718 }
719 
720 inline int KFParticle::GetLifeTime( float &T, float &SigmaT ) const
721 {
728  return KFParticleBase::GetLifeTime( T, SigmaT );
729 }
730 
731 inline int KFParticle::GetR( float &R, float &SigmaR ) const
732 {
737  return KFParticleBase::GetR( R, SigmaR );
738 }
739 
740 inline float & KFParticle::X()
741 {
742  return KFParticleBase::X();
743 }
744 
745 inline float & KFParticle::Y()
746 {
747  return KFParticleBase::Y();
748 }
749 
750 inline float & KFParticle::Z()
751 {
752  return KFParticleBase::Z();
753 }
754 
755 inline float & KFParticle::Px()
756 {
757  return KFParticleBase::Px();
758 }
759 
760 inline float & KFParticle::Py()
761 {
762  return KFParticleBase::Py();
763 }
764 
765 inline float & KFParticle::Pz()
766 {
767  return KFParticleBase::Pz();
768 }
769 
770 inline float & KFParticle::E()
771 {
772  return KFParticleBase::E();
773 }
774 
775 inline float & KFParticle::S()
776 {
777  return KFParticleBase::S();
778 }
779 
780 inline char & KFParticle::Q()
781 {
782  return KFParticleBase::Q();
783 }
784 
785 inline float & KFParticle::Chi2()
786 {
787  return KFParticleBase::Chi2();
788 }
789 
790 inline Int_t & KFParticle::NDF()
791 {
792  return KFParticleBase::NDF();
793 }
794 
795 inline float & KFParticle::Parameter ( int i )
796 {
797  return KFParticleBase::Parameter(i);
798 }
799 
800 inline float & KFParticle::Covariance( int i )
801 {
802  return KFParticleBase::Covariance(i);
803 }
804 
805 inline float & KFParticle::Covariance( int i, int j )
806 {
807  return KFParticleBase::Covariance(i,j);
808 }
809 
810 inline float * KFParticle::Parameters ()
811 {
812  return fP;
813 }
814 
816 {
817  return fC;
818 }
819 
820 
821 inline void KFParticle::operator +=( const KFParticle &Daughter )
822 {
826 #ifdef NonhomogeneousField
827  for(int i=0; i<10; i++)
828  SetFieldCoeff(Daughter.GetFieldCoeff()[i], i);
829 #endif
830  KFParticleBase::operator +=( Daughter );
831 }
832 
833 
834 inline void KFParticle::AddDaughter( const KFParticle &Daughter )
835 {
848 #ifdef NonhomogeneousField
849  for(int i=0; i<10; i++)
850  SetFieldCoeff(Daughter.GetFieldCoeff()[i], i);
851 #endif
852  KFParticleBase::AddDaughter( Daughter );
853 }
854 
855 inline void KFParticle::Construct( const KFParticle *vDaughters[], int nDaughters,
856  const KFParticle *ProdVtx, float Mass )
857 {
868 #ifdef NonhomogeneousField
869  for(int i=0; i<10; i++)
870  SetFieldCoeff(vDaughters[0]->GetFieldCoeff()[i], i);
871 #endif
872  KFParticleBase::Construct( ( const KFParticleBase**)vDaughters, nDaughters,
873  ( const KFParticleBase*)ProdVtx, Mass );
874 }
875 
876 inline void KFParticle::TransportToPoint( const float xyz[] )
877 {
881  float dsdr[6] = {0.f};
882  float dS = GetDStoPoint(xyz, dsdr);
883  TransportToDS( dS, dsdr );
884 }
885 #ifdef HomogeneousField
886 inline void KFParticle::TransportToVertex( const KFPVertex &v )
887 {
892 }
893 #endif
895 {
899  float dsdr[4][6];
900  float dS[2];
901  GetDStoParticle( p, dS, dsdr );
902  TransportToDS( dS[0], dsdr[0] );
903 }
904 
905 inline float KFParticle::GetDStoPoint( const float xyz[], float* dsdr ) const
906 {
918 #ifdef HomogeneousField
919  return KFParticleBase::GetDStoPointBz( GetFieldAlice(), xyz, dsdr );
920 #endif
921 #ifdef NonhomogeneousField
922  return KFParticleBase::GetDStoPointCBM( xyz, dsdr );
923 #endif
924 }
925 
926 #ifdef HomogeneousField
927 inline float KFParticle::GetFieldAlice()
928 {
930  return fgBz;
931 }
932 #endif
933 
934 #ifdef HomogeneousField
935 inline void KFParticle::GetFieldValue( const float * /*xyz*/, float B[] ) const
936 {
943  B[0] = B[1] = 0;
944  B[2] = GetFieldAlice();
945 }
946 #endif
947 
948 #ifdef NonhomogeneousField
949 inline void KFParticle::GetFieldValue( const float xyz[], float B[] ) const
950 {
957  const float dz = (xyz[2]-fieldRegion[9]);
958  const float dz2 = dz*dz;
959 
960  B[0] = fieldRegion[0] + fieldRegion[1]*dz + fieldRegion[2]*dz2;
961  B[1] = fieldRegion[3] + fieldRegion[4]*dz + fieldRegion[5]*dz2;
962  B[2] = fieldRegion[6] + fieldRegion[7]*dz + fieldRegion[8]*dz2;
963 }
964 #endif
965 
966 inline void KFParticle::GetDStoParticle( const KFParticleBase &p, float dS[2], float dsdr[4][6] ) const
967 {
984 #ifdef HomogeneousField
985  KFParticleBase::GetDStoParticleBz( GetFieldAlice(), p, dS, dsdr ) ;
986 #endif
987 #ifdef NonhomogeneousField
988  KFParticleBase::GetDStoParticleCBM( p, dS, dsdr ) ;
989 #endif
990 }
991 
992 inline void KFParticle::Transport( float dS, const float* dsdr, float P[], float C[], float* dsdr1, float* F, float* F1 ) const
993 {
1016 #ifdef HomogeneousField
1017  KFParticleBase::TransportBz( GetFieldAlice(), dS, dsdr, P, C, dsdr1, F, F1 );
1018 #endif
1019 #ifdef NonhomogeneousField
1020  KFParticleBase::TransportCBM( dS, dsdr, P, C, dsdr1, F, F1 );
1021 #endif
1022 }
1023 
1024 #endif