FairRoot/PandaRoot
KFParticleSIMD.h
Go to the documentation of this file.
1 //---------------------------------------------------------------------------------
2 // The KFParticleSIMD 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 KFParticleBaseSIMD
14 //
15 // -= Copyright &copy ALICE HLT and CBM L1 Groups =-
16 //_________________________________________________________________________________
17 
18 //#define NonhomogeneousField
19 #define HomogeneousField
20 
21 #ifndef KFPARTICLESIMD_H
22 #define KFPARTICLESIMD_H
23 
24 #include "KFParticleBaseSIMD.h"
25 
26 #ifdef HomogeneousField
27 class KFPTrack;
28 class KFPVertex;
29 #endif
30 
31 #ifdef NonhomogeneousField
32 class CbmKFTrackInterface;
33 class CbmKFVertexInterface;
34 #include "L1Field.h"
35 #endif
36 
37 class KFParticle;
38 
40 {
41 
42  public:
43 
44  //*
45  //* INITIALIZATION
46  //*
47 
48  //* Set magnetic field for all particles
49 #ifdef HomogeneousField
50  static void SetField( fvec Bz );
51 #endif
52 #ifdef NonhomogeneousField
53  void SetField(const L1FieldRegion &field, bool isOneEntry=0, const int iVec=0)
54  {
55  if(!isOneEntry)
56  fField = field;
57  else
58  fField.SetOneEntry(field,iVec);
59  }
60 #endif
61  //* Constructor (empty)
62 
64 
65  //* Destructor (empty)
66 
68 
69  //* Construction of mother particle by its 2-3-4 daughters
70 
71  KFParticleSIMD( const KFParticleSIMD &d1, const KFParticleSIMD &d2, Bool_t gamma = false );
72 
73  KFParticleSIMD( const KFParticleSIMD &d1, const KFParticleSIMD &d2, const KFParticleSIMD &d3 );
74 
75  KFParticleSIMD( const KFParticleSIMD &d1, const KFParticleSIMD &d2, const KFParticleSIMD &d3, const KFParticleSIMD &d4 );
76 
77  //* Initialisation from "cartesian" coordinates ( X Y Z Px Py Pz )
78  //* Parameters, covariance matrix, charge and PID hypothesis should be provided
79 
80  void Create( const fvec Param[], const fvec Cov[], fvec Charge, fvec mass /*Int_t PID*/ );
81 
82 #ifdef HomogeneousField
83  //* Initialisation from ALICE track, PID hypothesis shoould be provided
84 
85  KFParticleSIMD( const KFPTrack &track, Int_t PID );
86  KFParticleSIMD( const KFPTrack *track, Int_t PID );
87  KFParticleSIMD( KFPTrack* Track[], int NTracks, Int_t *qHypo=0, const Int_t *pdg=0 );
88  void Create(KFPTrack* Track[], int NTracks, Int_t *qHypo=0, const Int_t *pdg=0);
89  KFParticleSIMD(KFPTrack &Track, Int_t *qHypo=0, const Int_t *pdg=0);
90  //* Initialisation from VVertex
91 
92  KFParticleSIMD( const KFPVertex &vertex );
93 #endif
94 
95 #ifdef NonhomogeneousField
96  KFParticleSIMD( CbmKFTrackInterface* Track[], int NTracks, Int_t *qHypo=0, const Int_t *pdg=0 );
97  KFParticleSIMD( CbmKFTrackInterface &Track, Int_t *qHypo=0, const Int_t *pdg=0);
98  KFParticleSIMD( CbmKFVertexInterface &vertex );
99  void Create(CbmKFTrackInterface* Track[], int NTracks, Int_t *qHypo=0, const Int_t *pdg=0);
100 #endif
101  KFParticleSIMD( KFParticle* part[], const int nPart = 0 );
102  KFParticleSIMD( KFParticle &part );
103 
104  //* Initialise covariance matrix and set current parameters to 0.0
105 
106  void Initialize();
107 
108  //* Set decay vertex parameters for linearisation
109 
110  void SetVtxGuess( fvec x, fvec y, fvec z );
111 
112  //*
113  //* ACCESSORS
114  //*
115 
116  //* Simple accessors
117 
118  fvec GetX () const ; //* x of current position
119  fvec GetY () const ; //* y of current position
120  fvec GetZ () const ; //* z of current position
121  fvec GetPx () const ; //* x-compoment of 3-momentum
122  fvec GetPy () const ; //* y-compoment of 3-momentum
123  fvec GetPz () const ; //* z-compoment of 3-momentum
124  fvec GetE () const ; //* energy
125  fvec GetS () const ; //* decay length / momentum
126  fvec GetQ () const ; //* charge
127  fvec GetChi2 () const ; //* chi^2
128  fvec GetNDF () const ; //* Number of Degrees of Freedom
129 
131 
132  const fvec& X () const { return fP[0]; }
133  const fvec& Y () const { return fP[1]; }
134  const fvec& Z () const { return fP[2]; }
135  const fvec& Px () const { return fP[3]; }
136  const fvec& Py () const { return fP[4]; }
137  const fvec& Pz () const { return fP[5]; }
138  const fvec& E () const { return fP[6]; }
139  const fvec& S () const { return fP[7]; }
140  const fvec& Q () const { return fQ; }
141  const fvec& Chi2 () const { return fChi2; }
142  const fvec& NDF () const { return fNDF; }
143 
144  fvec GetParameter ( int i ) const ;
145  fvec GetCovariance( int i ) const ;
146  fvec GetCovariance( int i, int j ) const ;
147 
148  //* Accessors with calculations, value returned w/o error flag
149 
150  fvec GetP () const; //* momentum
151  fvec GetPt () const; //* transverse momentum
152  fvec GetEta () const; //* pseudorapidity
153  fvec GetPhi () const; //* phi
154  fvec GetMomentum () const; //* momentum (same as GetP() )
155  fvec GetMass () const; //* mass
156  fvec GetDecayLength () const; //* decay length
157  fvec GetDecayLengthXY () const; //* decay length in XY
158  fvec GetLifeTime () const; //* life time
159  fvec GetR () const; //* distance to the origin
160 
161  //* Accessors to estimated errors
162 
163  fvec GetErrX () const ; //* x of current position
164  fvec GetErrY () const ; //* y of current position
165  fvec GetErrZ () const ; //* z of current position
166  fvec GetErrPx () const ; //* x-compoment of 3-momentum
167  fvec GetErrPy () const ; //* y-compoment of 3-momentum
168  fvec GetErrPz () const ; //* z-compoment of 3-momentum
169  fvec GetErrE () const ; //* energy
170  fvec GetErrS () const ; //* decay length / momentum
171  fvec GetErrP () const ; //* momentum
172  fvec GetErrPt () const ; //* transverse momentum
173  fvec GetErrEta () const ; //* pseudorapidity
174  fvec GetErrPhi () const ; //* phi
175  fvec GetErrMomentum () const ; //* momentum
176  fvec GetErrMass () const ; //* mass
177  fvec GetErrDecayLength () const ; //* decay length
178  fvec GetErrDecayLengthXY () const ; //* decay length in XY
179  fvec GetErrLifeTime () const ; //* life time
180  fvec GetErrR () const ; //* distance to the origin
181 
182  //* Accessors with calculations( &value, &estimated sigma )
183  //* error flag returned (0 means no error during calculations)
184 
185  fvec GetP ( fvec &P, fvec &SigmaP ) const ; //* momentum
186  fvec GetPt ( fvec &Pt, fvec &SigmaPt ) const ; //* transverse momentum
187  fvec GetEta ( fvec &Eta, fvec &SigmaEta ) const ; //* pseudorapidity
188  fvec GetPhi ( fvec &Phi, fvec &SigmaPhi ) const ; //* phi
189  fvec GetMomentum ( fvec &P, fvec &SigmaP ) const ; //* momentum
190  fvec GetMass ( fvec &M, fvec &SigmaM ) const ; //* mass
191  fvec GetDecayLength ( fvec &L, fvec &SigmaL ) const ; //* decay length
192  fvec GetDecayLengthXY ( fvec &L, fvec &SigmaL ) const ; //* decay length in XY
193  fvec GetLifeTime ( fvec &T, fvec &SigmaT ) const ; //* life time
194  fvec GetR ( fvec &R, fvec &SigmaR ) const ; //* R
195 
196 
197  //*
198  //* MODIFIERS
199  //*
200 
201  fvec & X () ;
202  fvec & Y () ;
203  fvec & Z () ;
204  fvec & Px () ;
205  fvec & Py () ;
206  fvec & Pz () ;
207  fvec & E () ;
208  fvec & S () ;
209  fvec & Q () ;
210  fvec & Chi2 () ;
211  fvec & NDF () ;
212 
213  fvec & Parameter ( int i ) ;
214  fvec & Covariance( int i ) ;
215  fvec & Covariance( int i, int j ) ;
216  fvec * Parameters () ;
217  fvec * CovarianceMatrix() ;
218 
219  void GetKFParticle( KFParticle &Part, int iPart = 0);
220  void GetKFParticle( KFParticle *Part, int nPart = 0);
221 
222  //*
223  //* CONSTRUCTION OF THE PARTICLE BY ITS DAUGHTERS AND MOTHER
224  //* USING THE KALMAN FILTER METHOD
225  //*
226 
227 
228  //* Add daughter to the particle
229 
230  void AddDaughter( const KFParticleSIMD &Daughter );
231 
232  //* Add daughter via += operator: ex.{ D0; D0+=Pion; D0+= Kaon; }
233 
234  void operator +=( const KFParticleSIMD &Daughter );
235 
236  //* Set production vertex
237 
238  void SetProductionVertex( const KFParticleSIMD &Vtx );
239 
240  //* Set mass constraint
241 
242  void SetMassConstraint( fvec Mass, fvec SigmaMass = 0 );
243 
244  //* Set no decay length for resonances
245 
246  void SetNoDecayLength();
247 
248  //* Everything in one go
249 
250  void Construct( const KFParticleSIMD *vDaughters[], int nDaughters,
251  const KFParticleSIMD *ProdVtx=0, Float_t Mass=-1, Bool_t IsConstrained=0,
252  Bool_t isAtVtxGuess = 0);
253 
254  //*
255  //* TRANSPORT
256  //*
257  //* ( main transportation parameter is S = SignedPath/Momentum )
258  //* ( parameters of decay & production vertices are stored locally )
259  //*
260 
261  //* Transport the particle to its decay vertex
262 
263  void TransportToDecayVertex();
264 
265  //* Transport the particle to its production vertex
266 
268 
269  //* Transport the particle close to xyz[] point
270 
271  void TransportToPoint( const fvec xyz[] );
272 
273  //* Transport the particle close to VVertex
274 #ifdef HomogeneousField
275  void TransportToVertex( const KFPVertex &v );
276 #endif
277  //* Transport the particle close to another particle p
278 
279  void TransportToParticle( const KFParticleSIMD &p );
280 
281  //* Transport the particle on dS parameter (SignedPath/Momentum)
282 
283  void TransportToDS( fvec dS );
284 
285  //* Get dS to a certain space point
286 
287  fvec GetDStoPoint( const fvec xyz[] ) const ;
288 
289  //* Get dS to other particle p (dSp for particle p also returned)
290 
291  void GetDStoParticle( const KFParticleSIMD &p,
292  fvec &DS, fvec &DSp ) const ;
293 
294  //* Get dS to other particle p in XY-plane
295 
296  void GetDStoParticleXY( const KFParticleBaseSIMD &p,
297  fvec &DS, fvec &DSp ) const ;
298 
299  //*
300  //* OTHER UTILITIES
301  //*
302 
303 
304  //* Calculate distance from another object [cm]
305 
306  fvec GetDistanceFromVertex( const fvec vtx[] ) const ;
307  fvec GetDistanceFromVertex( const KFParticleSIMD &Vtx ) const ;
308 #ifdef HomogeneousField
309  fvec GetDistanceFromVertex( const KFPVertex &Vtx ) const ;
310 #endif
311  fvec GetDistanceFromParticle( const KFParticleSIMD &p ) const ;
312 
313  //* Calculate sqrt(Chi2/ndf) deviation from another object
314  //* ( v = [xyz]-vertex, Cv=[Cxx,Cxy,Cyy,Cxz,Cyz,Czz]-covariance matrix )
315 
316  fvec GetDeviationFromVertex( const fvec v[], const fvec Cv[]=0 ) const ;
317  fvec GetDeviationFromVertex( const KFParticleSIMD &Vtx ) const ;
318 #ifdef HomogeneousField
319  fvec GetDeviationFromVertex( const KFPVertex &Vtx ) const ;
320 #endif
322 
323  //* Calculate distance from another object [cm] in XY-plane
324 
325  fvec GetDistanceFromVertexXY( const fvec vtx[], fvec &val, fvec &err ) const ;
326  fvec GetDistanceFromVertexXY( const fvec vtx[], const fvec Cv[], fvec &val, fvec &err ) const ;
327  fvec GetDistanceFromVertexXY( const KFParticleSIMD &Vtx, fvec &val, fvec &err ) const ;
328 #ifdef HomogeneousField
329  fvec GetDistanceFromVertexXY( const KFPVertex &Vtx, fvec &val, fvec &err ) const ;
330 #endif
331 
332  fvec GetDistanceFromVertexXY( const fvec vtx[] ) const ;
333  fvec GetDistanceFromVertexXY( const KFParticleSIMD &Vtx ) const ;
334 #ifdef HomogeneousField
335  fvec GetDistanceFromVertexXY( const KFPVertex &Vtx ) const ;
336 #endif
338 
339  //* Calculate sqrt(Chi2/ndf) deviation from another object in XY plane
340  //* ( v = [xyz]-vertex, Cv=[Cxx,Cxy,Cyy,Cxz,Cyz,Czz]-covariance matrix )
341 
342  fvec GetDeviationFromVertexXY( const fvec v[], const fvec Cv[]=0 ) const ;
343  fvec GetDeviationFromVertexXY( const KFParticleSIMD &Vtx ) const ;
344 #ifdef HomogeneousField
345  fvec GetDeviationFromVertexXY( const KFPVertex &Vtx ) const ;
346 #endif
348 
349  //* Calculate opennig angle between two particles
350 
351  fvec GetAngle ( const KFParticleSIMD &p ) const ;
352  fvec GetAngleXY( const KFParticleSIMD &p ) const ;
353  fvec GetAngleRZ( const KFParticleSIMD &p ) const ;
354 
355  //* Subtract the particle from the vertex
356 
357  void SubtractFromVertex( KFParticleSIMD &v ) const ;
358  void SubtractFromParticle( KFParticleSIMD &v ) const ;
359 
360  //* Special method for creating gammas
361 
362  void ConstructGamma( const KFParticleSIMD &daughter1,
363  const KFParticleSIMD &daughter2 );
364 
365  // * Pseudo Proper Time of decay = (r*pt) / |pt| * M/|pt|
366  // @primVertex - primary vertex
367  // @mass - mass of the mother particle (in the case of "Hb -> JPsi" it would be JPsi mass)
368  // @*timeErr2 - squared error of the decay time. If timeErr2 = 0 it isn't calculated
369  fvec GetPseudoProperDecayTime( const KFParticleSIMD &primVertex, const fvec& mass, fvec* timeErr2 = 0 ) const;
370 
371  void GetFieldValue( const fvec xyz[], fvec B[] ) const ;
372 
373  protected:
374 
375  //*
376  //* INTERNAL STUFF
377  //*
378 
379  //* Method to access ALICE field
380 #ifdef HomogeneousField
381  static fvec GetFieldAlice();
382 #endif
383  //* Other methods required by the abstract KFParticleBaseSIMD class
384 
385  void GetDStoParticle( const KFParticleBaseSIMD &p, fvec &DS, fvec &DSp )const ;
386  void Transport( fvec dS, fvec P[], fvec C[] ) const ;
387  static void GetExternalTrackParam( const KFParticleBaseSIMD &p, Double_t X[fvecLen], Double_t Alpha[fvecLen], Double_t P[5][fvecLen] ) ;
388 
389  //void GetDStoParticleALICE( const KFParticleBaseSIMD &p, fvec &DS, fvec &DS1 ) const;
390 
391 
392  private:
393 #ifdef HomogeneousField
394  static fvec fgBz; //* Bz compoment of the magnetic field
395 #endif
396 #ifdef NonhomogeneousField
398 #endif
399 };
400 
401 
402 
403 //---------------------------------------------------------------------
404 //
405 // Inline implementation of the KFParticleSIMD methods
406 //
407 //---------------------------------------------------------------------
408 
409 #ifdef HomogeneousField
411 {
412  fgBz = Bz;
413 }
414 #endif
415 
417  const KFParticleSIMD &d2,
418  const KFParticleSIMD &d3 )
419 {
420  KFParticleSIMD mother;
421  mother+= d1;
422  mother+= d2;
423  mother+= d3;
424  *this = mother;
425 }
426 
428  const KFParticleSIMD &d2,
429  const KFParticleSIMD &d3,
430  const KFParticleSIMD &d4 )
431 {
432  KFParticleSIMD mother;
433  mother+= d1;
434  mother+= d2;
435  mother+= d3;
436  mother+= d4;
437  *this = mother;
438 }
439 
440 
442 {
444 }
445 
447 {
449 }
450 
451 inline fvec KFParticleSIMD::GetX () const
452 {
453  return KFParticleBaseSIMD::GetX();
454 }
455 
456 inline fvec KFParticleSIMD::GetY () const
457 {
458  return KFParticleBaseSIMD::GetY();
459 }
460 
461 inline fvec KFParticleSIMD::GetZ () const
462 {
463  return KFParticleBaseSIMD::GetZ();
464 }
465 
466 inline fvec KFParticleSIMD::GetPx () const
467 {
468  return KFParticleBaseSIMD::GetPx();
469 }
470 
471 inline fvec KFParticleSIMD::GetPy () const
472 {
473  return KFParticleBaseSIMD::GetPy();
474 }
475 
476 inline fvec KFParticleSIMD::GetPz () const
477 {
478  return KFParticleBaseSIMD::GetPz();
479 }
480 
481 inline fvec KFParticleSIMD::GetE () const
482 {
483  return KFParticleBaseSIMD::GetE();
484 }
485 
486 inline fvec KFParticleSIMD::GetS () const
487 {
488  return KFParticleBaseSIMD::GetS();
489 }
490 
491 inline fvec KFParticleSIMD::GetQ () const
492 {
493  return KFParticleBaseSIMD::GetQ();
494 }
495 
497 {
498  return KFParticleBaseSIMD::GetChi2();
499 }
500 
501 inline fvec KFParticleSIMD::GetNDF () const
502 {
503  return KFParticleBaseSIMD::GetNDF();
504 }
505 
506 inline fvec KFParticleSIMD::GetParameter ( int i ) const
507 {
509 }
510 
512 {
514 }
515 
516 inline fvec KFParticleSIMD::GetCovariance( int i, int j ) const
517 {
519 }
520 
521 
522 inline fvec KFParticleSIMD::GetP () const
523 {
524  fvec par, err;
525  fvec mask = KFParticleBaseSIMD::GetMomentum( par, err );
526  return ((!mask) & par);
527 }
528 
530 {
531  fvec par, err;
532  fvec mask = KFParticleBaseSIMD::GetPt( par, err ) ;
533  return ((!mask) & par);
534 }
535 
537 {
538  fvec par, err;
539  fvec mask = KFParticleBaseSIMD::GetEta( par, err );
540  return ((!mask) & par);
541 }
542 
544 {
545  fvec par, err;
546  fvec mask = KFParticleSIMD::GetPhi( par, err );
547  return ((!mask) & par);
548 }
549 
551 {
552  fvec par, err;
553  fvec mask = KFParticleSIMD::GetMomentum( par, err );
554  return ((!mask) & par);
555 }
556 
558 {
559  fvec par, err;
560  fvec mask = KFParticleSIMD::GetMass( par, err );
561  return ((!mask) & par);
562 }
563 
565 {
566  fvec par, err;
567  fvec mask = KFParticleSIMD::GetDecayLength( par, err );
568  return ((!mask) & par);
569 }
570 
572 {
573  fvec par, err;
574  fvec mask = KFParticleSIMD::GetDecayLengthXY( par, err );
575  return ((!mask) & par);
576 }
577 
579 {
580  fvec par, err;
581  fvec mask = KFParticleSIMD::GetLifeTime( par, err );
582  return ((!mask) & par);
583 }
584 
585 inline fvec KFParticleSIMD::GetR () const
586 {
587  fvec par, err;
588  fvec mask = KFParticleSIMD::GetR( par, err );
589  return ((!mask) & par);
590 }
591 
593 {
594  return sqrt(fabs( GetCovariance(0,0) ));
595 }
596 
598 {
599  return sqrt(fabs( GetCovariance(1,1) ));
600 }
601 
603 {
604  return sqrt(fabs( GetCovariance(2,2) ));
605 }
606 
608 {
609  return sqrt(fabs( GetCovariance(3,3) ));
610 }
611 
613 {
614  return sqrt(fabs( GetCovariance(4,4) ));
615 }
616 
618 {
619  return sqrt(fabs( GetCovariance(5,5) ));
620 }
621 
623 {
624  return sqrt(fabs( GetCovariance(6,6) ));
625 }
626 
628 {
629  return sqrt(fabs( GetCovariance(7,7) ));
630 }
631 
633 {
634  fvec par, err;
635  fvec mask = KFParticleSIMD::GetMomentum( par, err );
636  fvec ret(1.e10);
637  ret = (mask & ret) + ((!mask) & err);
638  return ret;
639 }
640 
642 {
643  fvec par, err;
644  fvec mask = KFParticleSIMD::GetPt( par, err );
645  fvec ret(1.e10);
646  ret = (mask & ret) + ((!mask) & err);
647  return ret;
648 }
649 
651 {
652  fvec par, err;
653  fvec mask = KFParticleSIMD::GetEta( par, err );
654  fvec ret(1.e10);
655  ret = (mask & ret) + ((!mask) & err);
656  return ret;
657 }
658 
660 {
661  fvec par, err;
662  fvec mask = KFParticleSIMD::GetPhi( par, err );
663  fvec ret(1.e10);
664  ret = (mask & ret) + ((!mask) & err);
665  return ret;
666 }
667 
669 {
670  fvec par, err;
671  fvec mask = KFParticleSIMD::GetMomentum( par, err );
672  fvec ret(1.e10);
673  ret = (mask & ret) + ((!mask) & err);
674  return ret;
675 }
676 
678 {
679  fvec par, err;
680  fvec mask = KFParticleSIMD::GetMass( par, err );
681  fvec ret(1.e10);
682  ret = (mask & ret) + ((!mask) & err);
683  return ret;
684 }
685 
687 {
688  fvec par, err;
689  fvec mask = KFParticleSIMD::GetDecayLength( par, err );
690  fvec ret(1.e10);
691  ret = (mask & ret) + ((!mask) & err);
692  return ret;
693 }
694 
696 {
697  fvec par, err;
698  fvec mask = KFParticleSIMD::GetDecayLengthXY( par, err );
699  fvec ret(1.e10);
700  ret = (mask & ret) + ((!mask) & err);
701  return ret;
702 }
703 
705 {
706  fvec par, err;
707  fvec mask = KFParticleSIMD::GetLifeTime( par, err );
708  fvec ret(1.e10);
709  ret = (mask & ret) + ((!mask) & err);
710  return ret;
711 }
712 
714 {
715  fvec par, err;
716  fvec mask = KFParticleSIMD::GetR( par, err );
717  fvec ret(1.e10);
718  ret = (mask & ret) + ((!mask) & err);
719  return ret;
720 }
721 
722 
723 inline fvec KFParticleSIMD::GetP( fvec &P, fvec &SigmaP ) const
724 {
725  return KFParticleBaseSIMD::GetMomentum( P, SigmaP );
726 }
727 
728 inline fvec KFParticleSIMD::GetPt( fvec &Pt, fvec &SigmaPt ) const
729 {
730  return KFParticleBaseSIMD::GetPt( Pt, SigmaPt );
731 }
732 
733 inline fvec KFParticleSIMD::GetEta( fvec &Eta, fvec &SigmaEta ) const
734 {
735  return KFParticleBaseSIMD::GetEta( Eta, SigmaEta );
736 }
737 
738 inline fvec KFParticleSIMD::GetPhi( fvec &Phi, fvec &SigmaPhi ) const
739 {
740  return KFParticleBaseSIMD::GetPhi( Phi, SigmaPhi );
741 }
742 
743 inline fvec KFParticleSIMD::GetMomentum( fvec &P, fvec &SigmaP ) const
744 {
745  return KFParticleBaseSIMD::GetMomentum( P, SigmaP );
746 }
747 
748 inline fvec KFParticleSIMD::GetMass( fvec &M, fvec &SigmaM ) const
749 {
750  return KFParticleBaseSIMD::GetMass( M, SigmaM );
751 }
752 
753 inline fvec KFParticleSIMD::GetDecayLength( fvec &L, fvec &SigmaL ) const
754 {
755  return KFParticleBaseSIMD::GetDecayLength( L, SigmaL );
756 }
757 
758 inline fvec KFParticleSIMD::GetDecayLengthXY( fvec &L, fvec &SigmaL ) const
759 {
760  return KFParticleBaseSIMD::GetDecayLengthXY( L, SigmaL );
761 }
762 
763 inline fvec KFParticleSIMD::GetLifeTime( fvec &T, fvec &SigmaT ) const
764 {
765  return KFParticleBaseSIMD::GetLifeTime( T, SigmaT );
766 }
767 
768 inline fvec KFParticleSIMD::GetR( fvec &R, fvec &SigmaR ) const
769 {
770  return KFParticleBaseSIMD::GetR( R, SigmaR );
771 }
772 
774 {
775  return KFParticleBaseSIMD::X();
776 }
777 
779 {
780  return KFParticleBaseSIMD::Y();
781 }
782 
784 {
785  return KFParticleBaseSIMD::Z();
786 }
787 
789 {
790  return KFParticleBaseSIMD::Px();
791 }
792 
794 {
795  return KFParticleBaseSIMD::Py();
796 }
797 
799 {
800  return KFParticleBaseSIMD::Pz();
801 }
802 
804 {
805  return KFParticleBaseSIMD::E();
806 }
807 
809 {
810  return KFParticleBaseSIMD::S();
811 }
812 
814 {
815  return KFParticleBaseSIMD::Q();
816 }
817 
819 {
820  return KFParticleBaseSIMD::Chi2();
821 }
822 
824 {
825  return KFParticleBaseSIMD::NDF();
826 }
827 
829 {
831 }
832 
834 {
836 }
837 
838 inline fvec & KFParticleSIMD::Covariance( int i, int j )
839 {
840  return KFParticleBaseSIMD::Covariance(i,j);
841 }
842 
844 {
845  return fP;
846 }
847 
849 {
850  return fC;
851 }
852 
853 
854 inline void KFParticleSIMD::operator +=( const KFParticleSIMD &Daughter )
855 {
857 }
858 
859 
860 inline void KFParticleSIMD::AddDaughter( const KFParticleSIMD &Daughter )
861 {
863 }
864 
866 {
868 }
869 
871 {
872  KFParticleBaseSIMD::SetMassConstraint( Mass, SigmaMass );
873 }
874 
876 {
878 }
879 
880 inline void KFParticleSIMD::Construct( const KFParticleSIMD *vDaughters[], int nDaughters,
881  const KFParticleSIMD *ProdVtx, Float_t Mass, Bool_t IsConstrained,
882  Bool_t isAtVtxGuess )
883 {
884  KFParticleBaseSIMD::Construct( ( const KFParticleBaseSIMD**)vDaughters, nDaughters,
885  ( const KFParticleBaseSIMD*)ProdVtx, Mass, IsConstrained, isAtVtxGuess );
886 
887  #ifdef NonhomogeneousField
888  // calculate a field region for the constructed particle
889  L1FieldValue field[3];
890  fvec zField[3] = {0, fP[2]/2, fP[2]};
891 
892  for(int iPoint=0; iPoint<3; iPoint++)
893  {
894  for(int iD=0; iD<nDaughters; ++iD)
895  {
896  L1FieldValue b = const_cast<KFParticleSIMD *>(vDaughters[iD])->fField.Get(zField[iPoint]);
897  field[iPoint].x += b.x;
898  field[iPoint].y += b.y;
899  field[iPoint].z += b.z;
900  }
901  field[iPoint].x /= nDaughters;
902  field[iPoint].y /= nDaughters;
903  field[iPoint].z /= nDaughters;
904  }
905 
906  fField.Set( field[2], zField[2], field[1], zField[1], field[0], zField[0] );
907  #endif
908 }
909 
911 {
913 }
914 
916 {
918 }
919 
920 inline void KFParticleSIMD::TransportToPoint( const fvec xyz[] )
921 {
922  TransportToDS( GetDStoPoint(xyz) );
923 }
924 #ifdef HomogeneousField
926 {
928 }
929 #endif
931 {
932  fvec dS, dSp;
933  GetDStoParticle( p, dS, dSp );
934  TransportToDS( dS );
935 }
936 
938 {
940 }
941 
942 inline fvec KFParticleSIMD::GetDStoPoint( const fvec xyz[] ) const
943 {
944 #ifdef HomogeneousField
946 #endif
947 #ifdef NonhomogeneousField
949 #endif
950 }
951 
952 
954  fvec &DS, fvec &DSp ) const
955 {
956  GetDStoParticleXY( p, DS, DSp );
957 }
958 
959 
961 {
963 }
964 
966  const fvec Cv[] ) const
967 {
969 }
970 
972 {
974 }
975 
977 {
979 }
980 #ifdef HomogeneousField
982 {
983  return GetDistanceFromVertex( KFParticleSIMD(Vtx) );
984 }
985 
987 {
988  return GetDeviationFromVertex( KFParticleSIMD(Vtx) );
989 }
990 #endif
992 {
994 }
995 
997 {
999 }
1000 
1002 {
1004 }
1005 
1007 {
1009 }
1010 
1011 #ifdef HomogeneousField
1013 {
1014  return fgBz;
1015 }
1016 #endif
1017 
1018 #ifdef HomogeneousField
1019 inline void KFParticleSIMD::GetFieldValue( const fvec * /*xyz*/, fvec B[] ) const
1020 {
1021  B[0] = B[1] = 0;
1022  B[2] = GetFieldAlice();
1023 }
1024 #endif
1025 
1026 #ifdef NonhomogeneousField
1027 inline void KFParticleSIMD::GetFieldValue( const fvec xyz[], fvec B[] ) const
1028 {
1029  L1FieldValue mB = const_cast<L1FieldRegion&>(fField).Get(xyz[2]);
1030  B[0] = mB.x;
1031  B[1] = mB.y;
1032  B[2] = mB.z;
1033 }
1034 #endif
1035 
1037  fvec &DS, fvec &DSp )const
1038 {
1039  GetDStoParticleXY( p, DS, DSp );
1040 }
1041 
1043  fvec &DS, fvec &DSp ) const
1044 {
1045 #ifdef HomogeneousField
1047 #endif
1048 #ifdef NonhomogeneousField
1050 #endif
1051  //GetDStoParticleALICE( p, DS, DSp ) ;
1052 }
1053 
1054 inline void KFParticleSIMD::Transport( fvec dS, fvec P[], fvec C[] ) const
1055 {
1056 #ifdef HomogeneousField
1058 #endif
1059 #ifdef NonhomogeneousField
1061 #endif
1062 }
1063 
1064 inline void KFParticleSIMD::ConstructGamma( const KFParticleSIMD &daughter1,
1065  const KFParticleSIMD &daughter2 )
1066 {
1067 #ifdef HomogeneousField
1068  KFParticleBaseSIMD::ConstructGammaBz( daughter1, daughter2, GetFieldAlice() );
1069 #endif
1070 }
1071 
1072 #endif
fvec GetR() const
void GetDStoParticle(const KFParticleSIMD &p, fvec &DS, fvec &DSp) const
PndMultiField * fField
Definition: sim_emc_apd.C:97
const fvec & X() const
Double_t p
Definition: anasim.C:58
fvec GetAngleXY(const KFParticleSIMD &p) const
fvec GetErrS() const
fvec GetEta(fvec &Eta, fvec &SigmaEta) const
fvec GetEta() const
fvec GetDeviationFromParticle(const KFParticleSIMD &p) const
static void SetField(fvec Bz)
fvec GetDStoPointCBM(const fvec xyz[]) const
fvec GetChi2() const
fvec GetDeviationFromParticle(const KFParticleBaseSIMD &p) const
void Construct(const KFParticleSIMD *vDaughters[], int nDaughters, const KFParticleSIMD *ProdVtx=0, Float_t Mass=-1, Bool_t IsConstrained=0, Bool_t isAtVtxGuess=0)
const fvec & NDF() const
const fvec & Py() const
Int_t i
Definition: run_full.C:25
TTree * b
fvec GetDecayLength(fvec &L, fvec &SigmaL) const
fvec GetMomentum(fvec &P, fvec &SigmaP) const
fvec GetDeviationFromVertex(const fvec v[], const fvec Cv[]=0) const
const fvec & S() const
void SubtractFromParticle(KFParticleSIMD &v) const
const fvec & Q() const
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
fvec GetErrPy() const
Double_t val[nBoxes][nFEBox]
Definition: createCalib.C:11
fvec GetErrMomentum() const
const fvec & Y() const
void TransportToDS(fvec dS)
fvec GetDistanceFromParticle(const KFParticleSIMD &p) const
fvec GetDeviationFromVertexXY(const fvec v[], const fvec Cv[]=0) const
fvec GetPt() const
Double_t par[3]
fvec GetErrLifeTime() const
fvec GetPz() const
void ConstructGammaBz(const KFParticleBaseSIMD &daughter1, const KFParticleBaseSIMD &daughter2, fvec Bz)
const fvec & Y() const
void SetVtxGuess(fvec x, fvec y, fvec z)
fvec GetErrPt() const
fvec GetErrDecayLength() const
void TransportToProductionVertex()
void SetNoDecayLength()
fvec GetLifeTime(fvec &T, fvec &SigmaT) const
fvec GetErrPx() const
void operator+=(const KFParticleBaseSIMD &Daughter)
fvec GetParameter(int i) const
void GetDStoParticleCBM(const KFParticleBaseSIMD &p, fvec &dS, fvec &dS1) const
fvec GetR(fvec &R, fvec &SigmaR) const
fvec GetS() const
fvec GetQ() const
fvec & Covariance(int i)
fvec GetP() const
fvec GetErrPhi() const
__m128 v
Definition: P4_F32vec4.h:4
void TransportToParticle(const KFParticleSIMD &p)
fvec GetNDF() const
fvec GetCovariance(Int_t i) const
void GetDStoParticleXY(const KFParticleBaseSIMD &p, fvec &DS, fvec &DSp) const
fvec GetPseudoProperDecayTime(const KFParticleSIMD &primVertex, const fvec &mass, fvec *timeErr2=0) const
fvec GetDStoPoint(const fvec xyz[]) const
void TransportCBM(fvec dS, fvec P[], fvec C[]) const
fvec GetLifeTime() const
int Pic_FED Eff_lEE C()
fvec * CovarianceMatrix()
TTree * T
Definition: anaLmdReco.C:32
fvec GetCovariance(int i) const
const fvec & Chi2() const
fvec GetErrPz() const
fvec & Covariance(Int_t i)
const fvec & Z() const
static fvec GetFieldAlice()
fvec GetErrZ() const
const fvec & Pz() const
FairMCTracks * Track
Definition: drawEveTracks.C:8
Double_t
const fvec & Px() const
void GetFieldValue(const fvec xyz[], fvec B[]) const
fvec GetDistanceFromParticle(const KFParticleBaseSIMD &p) const
void GetDStoParticleBz(fvec Bz, const KFParticleBaseSIMD &p, fvec &dS, fvec &dS1) const
PndMCTrack * track
Definition: anaLmdCluster.C:89
fvec GetDeviationFromParticleXY(const KFParticleSIMD &p) const
const fvec & E() const
const fvec & X() const
const fvec & NDF() const
void TransportToPoint(const fvec xyz[])
const fvec & Z() const
fvec GetE() const
fvec GetDecayLengthXY() const
Double_t z
fvec GetMomentum() const
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
const fvec & Pz() const
fvec GetPx() const
const fvec & Chi2() const
fvec GetDistanceFromVertex(const fvec vtx[]) const
fvec GetPt(fvec &Pt, fvec &SigmaPt) const
const fvec & Q() const
const fvec & S() const
Bool_t GetAtProductionVertex() const
void ConstructGamma(const KFParticleSIMD &daughter1, const KFParticleSIMD &daughter2)
fvec GetPhi() const
fvec GetDeviationFromVertex(const fvec v[], const fvec Cv[]=0) const
void SubtractFromParticle(KFParticleBaseSIMD &Vtx) const
fvec GetMass() const
const int fvecLen
Definition: P4_F32vec4.h:220
fvec GetErrP() const
fvec GetErrE() const
const fvec & E() const
fvec GetZ() const
void SetMassConstraint(fvec Mass, fvec SigmaMass=0)
GeV c P
void AddDaughter(const KFParticleBaseSIMD &Daughter, Bool_t isAtVtxGuess=0)
fvec GetErrEta() const
fvec GetDStoPointBz(fvec Bz, const fvec xyz[]) const
Double_t x
void TransportToVertex(const KFPVertex &v)
file Get("FairBaseParSet")
fvec GetDecayLength() const
fvec GetErrDecayLengthXY() const
fvec GetAngle(const KFParticleSIMD &p) const
fvec GetX() const
void AddDaughter(const KFParticleSIMD &Daughter)
void operator+=(const KFParticleSIMD &Daughter)
void Transport(fvec dS, fvec P[], fvec C[]) const
fvec GetMass(fvec &M, fvec &SigmaM) const
static fvec fgBz
fvec GetDistanceFromVertex(const fvec vtx[]) const
fvec GetPhi(fvec &Phi, fvec &SigmaPhi) const
fvec GetPy() const
fvec GetErrX() const
void TransportToDecayVertex()
void GetKFParticle(KFParticle &Part, int iPart=0)
Double_t y
const fvec & Px() const
fvec & Parameter(int i)
void SetProductionVertex(const KFParticleBaseSIMD &Vtx)
fvec GetParameter(Int_t i) const
fvec GetErrY() const
fvec GetAngleRZ(const KFParticleSIMD &p) const
fvec * Parameters()
Double_t R
Definition: checkhelixhit.C:61
fvec GetDistanceFromVertexXY(const fvec vtx[], fvec &val, fvec &err) const
void SetVtxGuess(fvec x, fvec y, fvec z)
void TransportBz(fvec Bz, fvec dS, fvec P[], fvec C[]) const
fvec GetErrR() const
void Create(const fvec Param[], const fvec Cov[], fvec Charge, fvec mass)
fvec GetDistanceFromParticleXY(const KFParticleSIMD &p) const
void SubtractFromVertex(KFParticleBaseSIMD &Vtx) const
fvec & Parameter(Int_t i)
fvec GetDecayLengthXY(fvec &L, fvec &SigmaL) const
static void GetExternalTrackParam(const KFParticleBaseSIMD &p, Double_t X[fvecLen], Double_t Alpha[fvecLen], Double_t P[5][fvecLen])
const fvec & Py() const
fvec GetErrMass() const
void SetMassConstraint(fvec Mass, fvec SigmaMass=0)
void SubtractFromVertex(KFParticleSIMD &v) const
void Construct(const KFParticleBaseSIMD *vDaughters[], Int_t nDaughters, const KFParticleBaseSIMD *ProdVtx=0, Float_t Mass=-1, Bool_t IsConstrained=0, Bool_t isAtVtxGuess=0)
fvec GetY() const
void SetProductionVertex(const KFParticleSIMD &Vtx)