Analysis Software
Documentation for sPHENIX simulation software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
KFParticleSIMD Class Reference

The main vectorised class of KF Particle pacakge, describes particle objects. More...

#include <KFParticle/blob/master/KFParticle/KFParticleSIMD.h>

+ Inheritance diagram for KFParticleSIMD:
+ Collaboration diagram for KFParticleSIMD:

Public Member Functions

void * operator new (size_t size)
 new operator for allocation of the SIMD-alligned dynamic memory allocation
 
void * operator new[] (size_t size)
 new operator for allocation of the SIMD-alligned dynamic memory allocation
 
void * operator new (size_t size, void *ptr)
 new operator for allocation of the SIMD-alligned dynamic memory allocation
 
void * operator new[] (size_t size, void *ptr)
 new operator for allocation of the SIMD-alligned dynamic memory allocation
 
void operator delete (void *ptr, size_t)
 delete operator for the SIMD-alligned dynamic memory release
 
void operator delete[] (void *ptr, size_t)
 delete operator for the SIMD-alligned dynamic memory release
 
 KFParticleSIMD ()
 
 ~KFParticleSIMD ()
 
 KFParticleSIMD (const KFParticleSIMD &d1, const KFParticleSIMD &d2)
 
 KFParticleSIMD (const KFParticleSIMD &d1, const KFParticleSIMD &d2, const KFParticleSIMD &d3)
 
 KFParticleSIMD (const KFParticleSIMD &d1, const KFParticleSIMD &d2, const KFParticleSIMD &d3, const KFParticleSIMD &d4)
 
void Create (const float_v Param[], const float_v Cov[], int_v Charge, float_v mass)
 
void SetOneEntry (int iEntry, KFParticleSIMD &part, int iEntryPart)
 
 KFParticleSIMD (const KFPTrack *track, Int_t PID)
 
 KFParticleSIMD (KFPTrack *Track[], int NTracks, const Int_t *pdg=0)
 
 KFParticleSIMD (KFPTrackVector &track, uint_v &index, const int_v &pdg)
 
void Create (KFPTrack *Track[], int NTracks, const Int_t *pdg=0)
 
void Create (KFPTrackVector &track, uint_v &index, const int_v &pdg)
 
void Load (KFPTrackVector &track, int index, const int_v &pdg)
 
void Rotate ()
 
 KFParticleSIMD (KFPTrack &Track, const Int_t *pdg=0)
 
 KFParticleSIMD (KFPTrackVector &track, int n, const Int_t *pdg=0)
 
 KFParticleSIMD (KFPEmcCluster &track, uint_v &index, const KFParticleSIMD &vertexGuess)
 
 KFParticleSIMD (KFPEmcCluster &track, int index, const KFParticleSIMD &vertexGuess)
 
void Create (KFPEmcCluster &track, uint_v &index, const KFParticleSIMD &vertexGuess)
 
void Load (KFPEmcCluster &track, int index, const KFParticleSIMD &vertexGuess)
 
 KFParticleSIMD (const KFPVertex &vertex)
 
 KFParticleSIMD (KFParticle *part[], const int nPart=0)
 
 KFParticleSIMD (KFParticle &part)
 
float_v GetX () const
 Retruns X coordinate of the particle, fP[0].
 
float_v GetY () const
 Retruns Y coordinate of the particle, fP[1].
 
float_v GetZ () const
 Retruns Z coordinate of the particle, fP[2].
 
float_v GetPx () const
 Retruns X component of the momentum, fP[3].
 
float_v GetPy () const
 Retruns Y component of the momentum, fP[4].
 
float_v GetPz () const
 Retruns Z component of the momentum, fP[5].
 
float_v GetE () const
 Returns energy of the particle, fP[6].
 
float_v GetS () const
 Returns dS=l/p, l - decay length, fP[7], defined if production vertex is set.
 
int_v GetQ () const
 Returns charge of the particle.
 
float_v GetChi2 () const
 Returns Chi2 of the fit.
 
int_v GetNDF () const
 Returns number of decrease of freedom.
 
Bool_t GetAtProductionVertex () const
 Returns a flag which shows if the particle is located at the production point.
 
const float_v & X () const
 Retruns X coordinate of the particle, fP[0].
 
const float_v & Y () const
 Retruns Y coordinate of the particle, fP[1].
 
const float_v & Z () const
 Retruns Z coordinate of the particle, fP[2].
 
const float_v & Px () const
 Retruns X component of the momentum, fP[3].
 
const float_v & Py () const
 Retruns Y component of the momentum, fP[4].
 
const float_v & Pz () const
 Retruns Z component of the momentum, fP[5].
 
const float_v & E () const
 Returns energy of the particle, fP[6].
 
const float_v & S () const
 Returns dS=l/p, l - decay length, fP[7], defined if production vertex is set.
 
const int_v & Q () const
 Returns charge of the particle.
 
const float_v & Chi2 () const
 Returns Chi2 of the fit.
 
const int_v & NDF () const
 Returns number of decrease of freedom.
 
float_v GetParameter (int i) const
 Returns P[i] parameter.
 
float_v GetCovariance (int i) const
 Returns C[i] element of the covariance matrix in the lower triangular form.
 
float_v GetCovariance (int i, int j) const
 Returns C[i,j] element of the covariance matrix.
 
float_v GetP () const
 Returns momentum.
 
float_v GetPt () const
 Returns transverse momentum.
 
float_v GetEta () const
 Returns pseudorapidity.
 
float_v GetPhi () const
 Returns the azimuthal angle phi.
 
float_v GetMomentum () const
 Returns momentum.
 
float_v GetMass () const
 Returns mass.
 
float_v GetDecayLength () const
 Returns decay length.
 
float_v GetDecayLengthXY () const
 Returns decay length in XY.
 
float_v GetLifeTime () const
 Returns life time ctau [cm].
 
float_v GetR () const
 Returns distance to the origin of the coordinate system {0,0,0}.
 
float_v GetRapidity () const
 Returns rapidity of the particle.
 
float_v GetErrX () const
 Returns the error of X of current position.
 
float_v GetErrY () const
 Returns the error of Y of current position.
 
float_v GetErrZ () const
 Returns the error of Z of current position.
 
float_v GetErrPx () const
 Returns the error of X-compoment of the particle momentum.
 
float_v GetErrPy () const
 Returns the error of Y-compoment of the particle momentum.
 
float_v GetErrPz () const
 Returns the error of Z-compoment of the particle momentum.
 
float_v GetErrE () const
 Returns the error of energy.
 
float_v GetErrS () const
 Returns the error of decay length / momentum.
 
float_v GetErrP () const
 Returns the error of momentum.
 
float_v GetErrPt () const
 Returns the error of transverse momentum.
 
float_v GetErrEta () const
 Returns the error of pseudorapidity.
 
float_v GetErrPhi () const
 Returns the error of the azimuthal angle phi.
 
float_v GetErrMomentum () const
 Returns the error of momentum.
 
float_v GetErrMass () const
 Returns the error of mass.
 
float_v GetErrDecayLength () const
 Returns the error of decay length.
 
float_v GetErrDecayLengthXY () const
 Returns the error of decay length in XY.
 
float_v GetErrLifeTime () const
 Returns the error of life time.
 
float_v GetErrR () const
 Returns the error of distance to the origin of the coordinate system {0,0,0}.
 
float_m GetP (float_v &P, float_v &SigmaP) const
 
float_m GetPt (float_v &Pt, float_v &SigmaPt) const
 
float_m GetEta (float_v &Eta, float_v &SigmaEta) const
 
float_m GetPhi (float_v &Phi, float_v &SigmaPhi) const
 
float_m GetMomentum (float_v &P, float_v &SigmaP) const
 
float_m GetMass (float_v &M, float_v &SigmaM) const
 
float_m GetDecayLength (float_v &L, float_v &SigmaL) const
 
float_m GetDecayLengthXY (float_v &L, float_v &SigmaL) const
 
float_m GetLifeTime (float_v &T, float_v &SigmaT) const
 
float_m GetR (float_v &R, float_v &SigmaR) const
 
float_v & X ()
 Modifier of X coordinate of the particle, fP[0].
 
float_v & Y ()
 Modifier of Y coordinate of the particle, fP[1].
 
float_v & Z ()
 Modifier of Z coordinate of the particle, fP[2].
 
float_v & Px ()
 Modifier of X component of the momentum, fP[3].
 
float_v & Py ()
 Modifier of Y component of the momentum, fP[4].
 
float_v & Pz ()
 Modifier of Z component of the momentum, fP[5].
 
float_v & E ()
 Modifier of energy of the particle, fP[6].
 
float_v & S ()
 Modifier of dS=l/p, l - decay length, fP[7], defined if production vertex is set.
 
int_v & Q ()
 Modifier of charge of the particle.
 
float_v & Chi2 ()
 Modifier of Chi2 of the fit.
 
int_v & NDF ()
 Modifier of number of decrease of freedom.
 
float_v & Parameter (int i)
 Modifier of P[i] parameter.
 
float_v & Covariance (int i)
 Modifier of C[i] element of the covariance matrix in the lower triangular form.
 
float_v & Covariance (int i, int j)
 Modifier of C[i,j] element of the covariance matrix.
 
float_v * Parameters ()
 Returns pointer to the parameters fP.
 
float_v * CovarianceMatrix ()
 Returns pointer to the covariance matrix fC.
 
void GetKFParticle (KFParticle &Part, int iPart=0)
 
void GetKFParticle (KFParticle *Part, int nPart=0)
 
void AddDaughter (const KFParticleSIMD &Daughter)
 
void operator+= (const KFParticleSIMD &Daughter)
 
void Construct (const KFParticleSIMD *vDaughters[], int nDaughters, const KFParticleSIMD *ProdVtx=0, Float_t Mass=-1)
 
void TransportToPoint (const float_v xyz[])
 
void TransportToParticle (const KFParticleSIMD &p)
 
float_v GetDStoPoint (const float_v xyz[3], float_v dsdr[6]) const
 
void GetDStoParticle (const KFParticleBaseSIMD &p, float_v dS[2], float_v dsdr[4][6]) const
 
void GetDStoParticleFast (const KFParticleBaseSIMD &p, float_v dS[2]) const
 
float_m GetDistanceFromVertexXY (const float_v vtx[], float_v &val, float_v &err) const
 
float_m GetDistanceFromVertexXY (const float_v vtx[], const float_v Cv[], float_v &val, float_v &err) const
 
float_m GetDistanceFromVertexXY (const KFParticleSIMD &Vtx, float_v &val, float_v &err) const
 
float_v GetDistanceFromVertexXY (const float_v vtx[]) const
 
float_v GetDistanceFromVertexXY (const KFParticleSIMD &Vtx) const
 
float_v GetDistanceFromParticleXY (const KFParticleSIMD &p) const
 
float_v GetDeviationFromVertexXY (const float_v v[], const float_v Cv[]=0) const
 
float_v GetDeviationFromVertexXY (const KFParticleSIMD &Vtx) const
 
float_v GetDeviationFromParticleXY (const KFParticleSIMD &p) const
 
float_v GetAngle (const KFParticleSIMD &p) const
 
float_v GetAngleXY (const KFParticleSIMD &p) const
 
float_v GetAngleRZ (const KFParticleSIMD &p) const
 
float_v GetPseudoProperDecayTime (const KFParticleSIMD &primVertex, const float_v &mass, float_v *timeErr2=0) const
 
void GetFieldValue (const float_v xyz[], float_v B[]) const
 
void Transport (float_v dS, const float_v *dsdr, float_v P[], float_v C[], float_v *dsdr1=0, float_v *F=0, float_v *F1=0) const
 
void TransportFast (float_v dS, float_v P[]) const
 
- Public Member Functions inherited from KFParticleBaseSIMD
void * operator new (size_t size)
 new operator for allocation of the SIMD-alligned dynamic memory allocation
 
void * operator new[] (size_t size)
 new operator for allocation of the SIMD-alligned dynamic memory allocation
 
void * operator new (size_t size, void *ptr)
 new operator for allocation of the SIMD-alligned dynamic memory allocation
 
void * operator new[] (size_t size, void *ptr)
 new operator for allocation of the SIMD-alligned dynamic memory allocation
 
void operator delete (void *ptr, size_t)
 delete operator for the SIMD-alligned dynamic memory release
 
void operator delete[] (void *ptr, size_t)
 delete operator for the SIMD-alligned dynamic memory release
 
float_v GetDStoPointLine (const float_v xyz[3], float_v dsdr[6]) const
 
float_v GetDStoPointBz (float_v Bz, const float_v xyz[3], float_v dsdr[6], const float_v *param=0) const
 
float_v GetDStoPointBy (float_v By, const float_v xyz[3], float_v dsdr[6]) const
 
float_v GetDStoPointCBM (const float_v xyz[3], float_v dsdr[6]) const
 
void GetDStoParticleLine (const KFParticleBaseSIMD &p, float_v dS[2], float_v dsdr[4][6]) const
 
void GetDStoParticleLine (const KFParticleBaseSIMD &p, float_v dS[2]) const
 
void GetDStoParticleBz (float_v Bz, const KFParticleBaseSIMD &p, float_v dS[2], float_v dsdr[4][6], const float_v *param1=0, const float_v *param2=0) const
 
void GetDStoParticleBz (float_v Bz, const KFParticleBaseSIMD &p, float_v dS[2], const float_v *param1=0, const float_v *param2=0) const
 
void GetDStoParticleBy (float_v B, const KFParticleBaseSIMD &p, float_v dS[2], float_v dsdr[4][6]) const
 
void GetDStoParticleBy (float_v B, const KFParticleBaseSIMD &p, float_v dS[2]) const
 
void GetDStoParticleB (float_v B[3], const KFParticleBaseSIMD &p, float_v dS[2], float_v dsdr[4][6]) const
 
void GetDStoParticleB (float_v B[3], const KFParticleBaseSIMD &p, float_v dS[2]) const
 
void GetDStoParticleCBM (const KFParticleBaseSIMD &p, float_v dS[2], float_v dsdr[4][6]) const
 
void GetDStoParticleCBM (const KFParticleBaseSIMD &p, float_v dS[2]) const
 
virtual void Transport (float_v dS, const float_v dsdr[6], float_v P[], float_v C[], float_v *dsdr1=0, float_v *F=0, float_v *F1=0) const =0
 
 KFParticleBaseSIMD ()
 
virtual ~KFParticleBaseSIMD ()
 The default destructor.
 
void Initialize (const float_v Param[], const float_v Cov[], int_v Charge, float_v Mass)
 
void Initialize ()
 
void SetConstructMethod (Int_t m)
 Defines the construction method for the current particle (see description of fConstructMethod).
 
void SetMassHypo (float_v m)
 Sets the mass hypothesis to the particle, is used when fConstructMethod = 2.
 
const float_v & GetMassHypo () const
 Returns the mass hypothesis.
 
const float_v & GetSumDaughterMass () const
 Returns the sum of masses of the daughters.
 
float_v GetX () const
 Returns the sum of masses of the daughters.
 
float_v GetY () const
 Returns the sum of masses of the daughters.
 
float_v GetZ () const
 Returns the sum of masses of the daughters.
 
float_v GetPx () const
 Returns the sum of masses of the daughters.
 
float_v GetPy () const
 Returns the sum of masses of the daughters.
 
float_v GetPz () const
 Returns the sum of masses of the daughters.
 
float_v GetE () const
 Returns the sum of masses of the daughters.
 
float_v GetS () const
 Returns the sum of masses of the daughters.
 
int_v GetQ () const
 Returns the sum of masses of the daughters.
 
float_v GetChi2 () const
 Returns the sum of masses of the daughters.
 
int_v GetNDF () const
 Returns the sum of masses of the daughters.
 
const float_v & X () const
 Retruns X coordinate of the particle, fP[0].
 
const float_v & Y () const
 Retruns Y coordinate of the particle, fP[1].
 
const float_v & Z () const
 Retruns Z coordinate of the particle, fP[2].
 
const float_v & Px () const
 Retruns X component of the momentum, fP[3].
 
const float_v & Py () const
 Retruns Y component of the momentum, fP[4].
 
const float_v & Pz () const
 Retruns Z component of the momentum, fP[5].
 
const float_v & E () const
 Returns energy of the particle, fP[6].
 
const float_v & S () const
 Returns dS=l/p, l - decay length, fP[7], defined if production vertex is set.
 
const int_v & Q () const
 Returns charge of the particle.
 
const float_v & Chi2 () const
 Returns Chi2 of the fit.
 
const int_v & NDF () const
 Returns number of decrease of freedom.
 
float_v GetParameter (Int_t i) const
 Returns P[i] parameter.
 
float_v GetCovariance (Int_t i) const
 Returns C[i] element of the covariance matrix in the lower triangular form.
 
float_v GetCovariance (Int_t i, Int_t j) const
 Returns C[i,j] element of the covariance matrix.
 
float_m GetMomentum (float_v &p, float_v &error) const
 
float_m GetPt (float_v &pt, float_v &error) const
 
float_m GetEta (float_v &eta, float_v &error) const
 
float_m GetPhi (float_v &phi, float_v &error) const
 
float_m GetMass (float_v &m, float_v &error) const
 
float_m GetDecayLength (float_v &l, float_v &error) const
 
float_m GetDecayLengthXY (float_v &l, float_v &error) const
 
float_m GetLifeTime (float_v &ctau, float_v &error) const
 
float_m GetR (float_v &r, float_v &error) const
 
float_v & X ()
 Modifier of X coordinate of the particle, fP[0].
 
float_v & Y ()
 Modifier of Y coordinate of the particle, fP[1].
 
float_v & Z ()
 Modifier of Z coordinate of the particle, fP[2].
 
float_v & Px ()
 Modifier of X component of the momentum, fP[3].
 
float_v & Py ()
 Modifier of Y component of the momentum, fP[4].
 
float_v & Pz ()
 Modifier of Z component of the momentum, fP[5].
 
float_v & E ()
 Modifier of energy of the particle, fP[6].
 
float_v & S ()
 Modifier of dS=l/p, l - decay length, fP[7], defined if production vertex is set.
 
int_v & Q ()
 Modifier of charge of the particle.
 
float_v & Chi2 ()
 Modifier of Chi2 of the fit.
 
int_v & NDF ()
 Modifier of number of decrease of freedom.
 
float_v & Parameter (Int_t i)
 Modifier of P[i] parameter.
 
float_v & Covariance (Int_t i)
 Modifier of C[i] element of the covariance matrix in the lower triangular form.
 
float_v & Covariance (Int_t i, Int_t j)
 Modifier of C[i,j] element of the covariance matrix.
 
void operator+= (const KFParticleBaseSIMD &Daughter)
 
void AddDaughter (const KFParticleBaseSIMD &Daughter)
 
void SubtractDaughter (const KFParticleBaseSIMD &Daughter)
 
void AddDaughterWithEnergyFit (const KFParticleBaseSIMD &Daughter)
 
void AddDaughterWithEnergyFitMC (const KFParticleBaseSIMD &Daughter)
 
void SetProductionVertex (const KFParticleBaseSIMD &Vtx)
 
void SetNonlinearMassConstraint (float_v Mass)
 
void SetMassConstraint (float_v Mass, float_v SigmaMass=float_v(0.f))
 
void SetNoDecayLength ()
 
void Construct (const KFParticleBaseSIMD *vDaughters[], Int_t nDaughters, const KFParticleBaseSIMD *ProdVtx=0, Float_t Mass=-1)
 
void TransportToDecayVertex ()
 
void TransportToProductionVertex ()
 
void TransportToDS (float_v dS, const float_v *dsdr)
 
void TransportToDSLine (float_v dS, const float_v *dsdr)
 
void TransportBz (float_v Bz, float_v dS, const float_v *dsdr, float_v P[], float_v C[], float_v *dsdr1=0, float_v *F=0, float_v *F1=0) const
 
void TransportBz (float_v Bz, float_v dS, float_v P[]) const
 
void TransportCBM (float_v dS, const float_v *dsdr, float_v P[], float_v C[], float_v *dsdr1=0, float_v *F=0, float_v *F1=0) const
 
void TransportCBM (float_v dS, float_v P[]) const
 
float_v GetDistanceFromVertex (const float_v vtx[]) const
 
float_v GetDistanceFromVertex (const KFParticleBaseSIMD &Vtx) const
 
float_v GetDistanceFromParticle (const KFParticleBaseSIMD &p) const
 
float_v GetDeviationFromVertex (const float_v v[], const float_v Cv[]=0) const
 
float_v GetDeviationFromVertex (const KFParticleBaseSIMD &Vtx) const
 
float_v GetDeviationFromParticle (const KFParticleBaseSIMD &p) const
 
void SubtractFromVertex (KFParticleBaseSIMD &Vtx) const
 
void SubtractFromParticle (KFParticleBaseSIMD &Vtx) const
 
void RotateXY (float_v angle, float_v Vtx[3])
 
int_v Id () const
 Returns Id of the particle.
 
int NDaughters () const
 Returns number of daughter particles.
 
std::vector< int_v > & DaughterIds ()
 Returns the vector with the indices of daughter particles.
 
int_v GetDaughterId (int iD) const
 Returns the daughter Id with the index iD.
 
void SetId (int_v id)
 Sets the Id of the particle. After the construction of a particle should be set by user.
 
void SetNDaughters (int n)
 Reserves the size of the vector with daughter Ids to n.
 
void AddDaughterId (int_v id)
 Adds index of the daughter particle.
 
void CleanDaughtersId ()
 Cleans the vector with the indices of daughter particles.
 
void SetPDG (int pdg)
 Sets the PDG hypothesis common for all elements of the SIMD vector.
 
void SetPDG (int_v &pdg)
 Sets the PDG hypothesis individual for each entry of the SIMD vector.
 
const int_v & GetPDG () const
 Returns the PDG hypothesis.
 
const int_v & PDG () const
 Returns the PDG hypothesis.
 
void GetDistanceToVertexLine (const KFParticleBaseSIMD &Vertex, float_v &l, float_v &dl, float_m *isParticleFromVertex=0) const
 

Additional Inherited Members

- Static Public Member Functions inherited from KFParticleBaseSIMD
static void GetArmenterosPodolanski (KFParticleBaseSIMD &positive, KFParticleBaseSIMD &negative, float_v QtAlfa[2])
 
static void MultQSQt (const float_v Q[], const float_v S[], float_v SOut[], const int kN)
 
- Protected Member Functions inherited from KFParticleBaseSIMD
float_v & Cij (Int_t i, Int_t j)
 
void TransportLine (float_v S, const float_v *dsdr, float_v P[], float_v C[], float_v *dsdr1=0, float_v *F=0, float_v *F1=0) const
 
void TransportLine (float_v S, float_v P[]) const
 
void GetMeasurement (const KFParticleBaseSIMD &daughter, float_v m[], float_v V[], float_v D[3][3])
 
void SetMassConstraint (float_v *mP, float_v *mC, float_v mJ[7][7], float_v mass, float_m mask)
 
- Static Protected Member Functions inherited from KFParticleBaseSIMD
static Int_t IJ (Int_t i, Int_t j)
 
static void InvertCholetsky3 (float_v a[6])
 
- Protected Attributes inherited from KFParticleBaseSIMD
float_v fP [8]
 Particle parameters { X, Y, Z, Px, Py, Pz, E, S[=DecayLength/P]}.
 
float_v fC [36]
 Low-triangle covariance matrix of fP.
 
int_v fQ
 The charge of the particle in the units of the elementary charge.
 
int_v fNDF
 Number of degrees of freedom.
 
float_v fChi2
 Chi^2.
 
float_v fSFromDecay
 Distance from the decay vertex to the current position.
 
float_v SumDaughterMass
 Sum of the daughter particles masses. Needed to set the constraint on the minimum mass during particle construction.
 
float_v fMassHypo
 The mass hypothesis, used for the constraints during particle construction.
 
int_v fId
 Id of the particle.
 
Bool_t fAtProductionVertex
 Flag shows if particle is at the production point.
 
int_v fPDG
 The PDG hypothesis assigned to the particle.
 
Int_t fConstructMethod
 Determines the method for the particle construction.
0 - Energy considered as an independent veriable, fitted independently from momentum, without any constraints on mass
2 - Energy considered as an independent variable, fitted independently from momentum, with constraints on mass of daughter particle.
 
std::vector< int_v > fDaughterIds
 A vector with ids of the daughter particles:
1) if particle is created from a track - the index of the track, in this case the size of the vector is always equal to one;
2) if particle is constructed from other particles - indices of these particles in the same array.
 

Detailed Description

The main vectorised class of KF Particle pacakge, describes particle objects.

Author
M.Zyzak
Date
05.02.2019
Version
1.0

Vectorised version of the KF Particle class, describes particle objects. The particle is described with the state vector { X, Y, Z, Px, Py, Pz, E } and the corresponding covariance matrix. It contains functionality to create particle-object from track, to construct short-lived particles from other tracks or particles. The mathematics is based on the Kalman filter method. It also allows to subtract particles from the already constructed object, to transport partices, get parameters together with their erros, get distance to other particles and vertices, get deviations from them in terms of errors, etc.

Definition at line 59 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 59 of file KFParticleSIMD.h

Constructor & Destructor Documentation

KFParticleSIMD::KFParticleSIMD ( )
inline

Definition at line 89 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 89 of file KFParticleSIMD.h

KFParticleSIMD::~KFParticleSIMD ( )
inline

Definition at line 97 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 97 of file KFParticleSIMD.h

KFParticleSIMD::KFParticleSIMD ( const KFParticleSIMD d1,
const KFParticleSIMD d2 
)

Constructs a particle from two input daughter particles

Parameters
[in]d1- the first daughter particle
[in]d2- the second daughter particle

Definition at line 32 of file KFParticleSIMD.cxx.

View newest version in sPHENIX GitHub at line 32 of file KFParticleSIMD.cxx

KFParticleSIMD::KFParticleSIMD ( const KFParticleSIMD d1,
const KFParticleSIMD d2,
const KFParticleSIMD d3 
)
inline

Constructs a particle from three input daughter particles

Parameters
[in]d1- the first daughter particle
[in]d2- the second daughter particle
[in]d3- the third daughter particle

Definition at line 382 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 382 of file KFParticleSIMD.h

KFParticleSIMD::KFParticleSIMD ( const KFParticleSIMD d1,
const KFParticleSIMD d2,
const KFParticleSIMD d3,
const KFParticleSIMD d4 
)
inline

Constructs a particle from four input daughter particles

Parameters
[in]d1- the first daughter particle
[in]d2- the second daughter particle
[in]d3- the third daughter particle
[in]d4- the fourth daughter particle

Definition at line 402 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 402 of file KFParticleSIMD.h

KFParticleSIMD::KFParticleSIMD ( const KFPTrack track,
Int_t  PID 
)

Constructor of the particle from an array of tracks.

Parameters
[in]track- pointer to the array of n=float_vLen tracks
[in]PID- the PID hypothesis common for all elements of the SIMD vector

Definition at line 72 of file KFParticleSIMD.cxx.

View newest version in sPHENIX GitHub at line 72 of file KFParticleSIMD.cxx

References C, KFPTrack::Charge(), Create(), KFParticleBaseSIMD::fC, KFParticleBaseSIMD::fChi2, float_vLen, KFParticleBaseSIMD::fNDF, KFParticleBaseSIMD::fP, KFParticleBaseSIMD::fQ, KFPTrack::GetChi2(), KFPTrack::GetCovarianceXYZPxPyPz(), KFParticleDatabase::GetMass(), KFPTrack::GetNDF(), i, KFParticleDatabase::Instance(), mass, KFPTrack::PxPyPz(), physmon_track_finding_ttbar::r, and KFPTrack::XvYvZv().

+ Here is the call graph for this function:

KFParticleSIMD::KFParticleSIMD ( KFPTrack Track[],
int  NTracks,
const Int_t *  pdg = 0 
)

Constructor of the particle from an array tracks.

Parameters
[in]Track- an array of pointers to tracks in the KFPTrack format
[in]NTracks- number of tracks in the arry
[in]pdg- pointer to the pdg hypothesis common for all elements of the SIMD vector

Definition at line 167 of file KFParticleSIMD.cxx.

View newest version in sPHENIX GitHub at line 167 of file KFParticleSIMD.cxx

References Create().

+ Here is the call graph for this function:

KFParticleSIMD::KFParticleSIMD ( KFPTrackVector track,
uint_v &  index,
const int_v &  pdg 
)

Constructor of the particle from a set of tracks with random indices "index" stored in the KFPTrackVector format.

Parameters
[in]track- an array with tracks in the KFPTrackVector format
[in]index- indices of the tracks to be converted to the KFParticleSIMD object
[in]pdg- a SIMD vector with an individual pdg hypothesis for each element

Definition at line 217 of file KFParticleSIMD.cxx.

View newest version in sPHENIX GitHub at line 217 of file KFParticleSIMD.cxx

References Create().

+ Here is the call graph for this function:

KFParticleSIMD::KFParticleSIMD ( KFPTrack Track,
const Int_t *  pdg = 0 
)

Constructor of the particle from a single track. The same track is put to each element of the SIMD vector.

Parameters
[in]Track- track in the KFPTrack format
[in]PID- the PID hypothesis common for all elements of the SIMD vector

Definition at line 109 of file KFParticleSIMD.cxx.

View newest version in sPHENIX GitHub at line 109 of file KFParticleSIMD.cxx

References C, KFPTrack::Charge(), Create(), KFParticleBaseSIMD::fC, KFParticleBaseSIMD::fChi2, KFParticleBaseSIMD::fNDF, KFParticleBaseSIMD::fP, KFParticleBaseSIMD::fQ, KFPTrack::GetChi2(), KFPTrack::GetCovarianceXYZPxPyPz(), KFParticleDatabase::GetMass(), KFPTrack::GetNDF(), i, KFParticleDatabase::Instance(), mass, KFPTrack::PxPyPz(), physmon_track_finding_ttbar::r, and KFPTrack::XvYvZv().

+ Here is the call graph for this function:

KFParticleSIMD::KFParticleSIMD ( KFPTrackVector track,
int  n,
const Int_t *  pdg = 0 
)

Constructor of the particle from a single track with index "n" stored in the KFPTrackVector format. The same track is put to each element of the SIMD vector.

Parameters
[in]track- an array with tracks in the KFPTrackVector format
[in]n- index of the track to be used
[in]pdg- pointer to the pdg hypothesis

Definition at line 140 of file KFParticleSIMD.cxx.

View newest version in sPHENIX GitHub at line 140 of file KFParticleSIMD.cxx

References KFPTrackVector::Covariance(), Create(), KFParticleBaseSIMD::fC, KFParticleBaseSIMD::fP, KFParticleBaseSIMD::fQ, KFParticleDatabase::GetMass(), i, KFParticleDatabase::Instance(), mass, n, KFPTrackVector::Parameter(), and KFPTrackVector::Q().

+ Here is the call graph for this function:

KFParticleSIMD::KFParticleSIMD ( KFPEmcCluster track,
uint_v &  index,
const KFParticleSIMD vertexGuess 
)

Constructor of gamma particles from a set of clusters of the electromagnetic calorimeter (EMC) with random indices "index". The vertex hypothesis should be provided for the estimation of the momentum.

Parameters
[in]track- an array of EMC clusters
[in]index- indices of the clusters to be converted to the KFParticleSIMD object
[in]vertexGuess- vertex guess for estimation of the momentum of created gamma particles

Definition at line 298 of file KFParticleSIMD.cxx.

View newest version in sPHENIX GitHub at line 298 of file KFParticleSIMD.cxx

References Create().

+ Here is the call graph for this function:

KFParticleSIMD::KFParticleSIMD ( KFPEmcCluster track,
int  index,
const KFParticleSIMD vertexGuess 
)

Constructr gamma particles from a set of consequetive clusters of the electromagnetic calorimeter (EMC) starting from the index "index". The vertex hypothesis should be provided for the estimation of the momentum.

Parameters
[in]track- an array with tracks in the KFPTrackVector format
[in]index- index of the first EMC cluster
[in]vertexGuess- vertex guess for estimation of the momentum of created gamma particles

Definition at line 380 of file KFParticleSIMD.cxx.

View newest version in sPHENIX GitHub at line 380 of file KFParticleSIMD.cxx

References Load().

+ Here is the call graph for this function:

KFParticleSIMD::KFParticleSIMD ( const KFPVertex vertex)

Copies a vertex in KFPVertex into each element of the vectorised KFParticle

Parameters
[in]vertex- vertex to b converted

Definition at line 461 of file KFParticleSIMD.cxx.

View newest version in sPHENIX GitHub at line 461 of file KFParticleSIMD.cxx

References C, KFParticleBaseSIMD::fAtProductionVertex, KFParticleBaseSIMD::fC, KFParticleBaseSIMD::fChi2, KFParticleBaseSIMD::fNDF, KFParticleBaseSIMD::fP, KFParticleBaseSIMD::fQ, KFParticleBaseSIMD::fSFromDecay, KFPVertex::GetChi2(), KFPVertex::GetCovarianceMatrix(), KFPVertex::GetNContributors(), KFPVertex::GetXYZ(), i, and physmon_track_finding_ttbar::r.

+ Here is the call graph for this function:

KFParticleSIMD::KFParticleSIMD ( KFParticle part[],
const int  nPart = 0 
)

Constructs a vectoriesd particle from an array of scalar KFParticle objects.

Parameters
[in]parts- array of scalar KFParticle objects
[in]nPart- number of particles in the array

Definition at line 522 of file KFParticleSIMD.cxx.

View newest version in sPHENIX GitHub at line 522 of file KFParticleSIMD.cxx

References KFParticle::CovarianceMatrix(), KFParticleBase::DaughterIds(), KFParticleBaseSIMD::fAtProductionVertex, KFParticleBaseSIMD::fC, KFParticleBaseSIMD::fChi2, KFParticleBaseSIMD::fDaughterIds, KFParticleBaseSIMD::fId, float_vLen, KFParticleBaseSIMD::fNDF, KFParticleBaseSIMD::fP, KFParticleBaseSIMD::fPDG, KFParticleBaseSIMD::fQ, KFParticle::GetAtProductionVertex(), KFParticle::GetChi2(), KFParticle::GetNDF(), KFParticleBase::GetPDG(), KFParticle::GetQ(), i, KFParticleBase::Id(), KFParticleBase::NDaughters(), KFParticleBaseSIMD::NDaughters(), and KFParticle::Parameters().

+ Here is the call graph for this function:

KFParticleSIMD::KFParticleSIMD ( KFParticle part)

Constructs a vectoriesd particle from a single scalar KFParticle object. The same particle is copied to each element.

Parameters
[in]part- a scalar particle which should be copied to the current vectorised particle

Definition at line 572 of file KFParticleSIMD.cxx.

View newest version in sPHENIX GitHub at line 572 of file KFParticleSIMD.cxx

References KFParticle::CovarianceMatrix(), KFParticleBase::DaughterIds(), KFParticleBaseSIMD::fAtProductionVertex, KFParticleBaseSIMD::fC, KFParticleBaseSIMD::fChi2, KFParticleBaseSIMD::fDaughterIds, KFParticleBaseSIMD::fId, KFParticleBaseSIMD::fNDF, KFParticleBaseSIMD::fP, KFParticleBaseSIMD::fPDG, KFParticleBaseSIMD::fQ, KFParticle::GetAtProductionVertex(), KFParticle::GetChi2(), KFParticle::GetNDF(), KFParticleBase::GetPDG(), KFParticle::GetQ(), i, KFParticleBase::Id(), KFParticleBase::NDaughters(), KFParticle::Parameters(), and KFParticleBaseSIMD::SetNDaughters().

+ Here is the call graph for this function:

Member Function Documentation

void KFParticleSIMD::AddDaughter ( const KFParticleSIMD Daughter)
inline

Adds daughter to the current particle. Depending on the selected construction method uses:
1) Either simplifyed fast mathematics which consideres momentum and energy as independent variables and thus ignores constraint on the fixed mass (fConstructMethod = 0). In this case the mass of the daughter particle can be corrupted when the constructed vertex is added as the measurement and the mass of the output short-lived particle can become unphysical - smaller then the threshold. Implemented in the AddDaughterWithEnergyFit() function
2) Or slower but correct mathematics which requires that the masses of daughter particles stays fixed in the construction process (fConstructMethod = 2). Implemented in the AddDaughterWithEnergyFitMC() function.

Parameters
[in]Daughter- the daughter particle

Definition at line 885 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 885 of file KFParticleSIMD.h

References KFParticleBaseSIMD::AddDaughter().

+ Here is the call graph for this function:

const float_v& KFParticleSIMD::Chi2 ( ) const
inline

Returns Chi2 of the fit.

Definition at line 165 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 165 of file KFParticleSIMD.h

References KFParticleBaseSIMD::fChi2.

Referenced by KFParticleFinder::CombinePartPart(), KFParticleFinder::ConstructTrackV0Cand(), KFParticleFinder::ConstructV0(), KFParticleFinder::NeutralDaughterDecay(), and SetOneEntry().

+ Here is the caller graph for this function:

float_v & KFParticleSIMD::Chi2 ( )
inline

Modifier of Chi2 of the fit.

Definition at line 837 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 837 of file KFParticleSIMD.h

References KFParticleBaseSIMD::Chi2().

+ Here is the call graph for this function:

void KFParticleSIMD::Construct ( const KFParticleSIMD vDaughters[],
int  nDaughters,
const KFParticleSIMD ProdVtx = 0,
Float_t  Mass = -1 
)
inline

Constructs a short-lived particle from a set of daughter particles:
1) all parameters of the "this" objects are initialised;
2) daughters are added one after another;
3) if Parent pointer is not null, the production vertex is set to it;
4) if Mass hypothesis >=0 the mass constraint is set.

Parameters
[in]vDaughters- array of daughter particles
[in]nDaughters- number of daughter particles in the input array
[in]Parent- optional parrent particle
[in]Mass- optional mass hypothesis

Definition at line 905 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 905 of file KFParticleSIMD.h

References KFParticleBaseSIMD::Construct().

Referenced by KFParticleFinder::CombinePartPart(), KFParticleFinder::ConstructTrackV0Cand(), KFParticleFinder::ConstructV0(), and KFParticleFinder::NeutralDaughterDecay().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

float_v & KFParticleSIMD::Covariance ( int  i)
inline

Modifier of C[i] element of the covariance matrix in the lower triangular form.

Definition at line 852 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 852 of file KFParticleSIMD.h

References KFParticleBaseSIMD::Covariance().

+ Here is the call graph for this function:

float_v & KFParticleSIMD::Covariance ( int  i,
int  j 
)
inline

Modifier of C[i,j] element of the covariance matrix.

Definition at line 857 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 857 of file KFParticleSIMD.h

References KFParticleBaseSIMD::Covariance().

+ Here is the call graph for this function:

float_v * KFParticleSIMD::CovarianceMatrix ( )
inline

Returns pointer to the covariance matrix fC.

Definition at line 867 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 867 of file KFParticleSIMD.h

References KFParticleBaseSIMD::fC.

Referenced by GetKFParticle(), and SetOneEntry().

+ Here is the caller graph for this function:

void KFParticleSIMD::Create ( const float_v  Param[],
const float_v  Cov[],
int_v  Charge,
float_v  mass 
)

Constructor from a "cartesian" track, mass hypothesis should be provided

Parameters
[in]Param[6]= { X, Y, Z, Px, Py, Pz } - position and momentum
[in]Cov[21]- lower-triangular part of the covariance matrix:
          (  0  .  .  .  .  . )
          (  1  2  .  .  .  . )
Cov[21] = (  3  4  5  .  .  . )
          (  6  7  8  9  .  . )
          ( 10 11 12 13 14  . )
          ( 15 16 17 18 19 20 )
[in]Charge- charge of the particle in elementary charge units
[in]mass- the mass hypothesis

Definition at line 48 of file KFParticleSIMD.cxx.

View newest version in sPHENIX GitHub at line 48 of file KFParticleSIMD.cxx

References C, i, and KFParticleBaseSIMD::Initialize().

Referenced by Create(), KFParticleTopoReconstructor::GetChiToPrimVertex(), KFParticleSIMD(), and Load().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void KFParticleSIMD::Create ( KFPTrack Track[],
int  NTracks,
const Int_t *  pdg = 0 
)

Create a particle from from an array tracks.

Parameters
[in]Track- an array of pointers to tracks in the KFPTrack format
[in]NTracks- number of tracks in the arry
[in]pdg- pointer to the pdg hypothesis common for all elements of the SIMD vector

Definition at line 180 of file KFParticleSIMD.cxx.

View newest version in sPHENIX GitHub at line 180 of file KFParticleSIMD.cxx

References C, KFPTrack::Charge(), Create(), KFParticleBaseSIMD::fC, KFParticleBaseSIMD::fChi2, float_vLen, KFParticleBaseSIMD::fNDF, KFParticleBaseSIMD::fP, KFParticleBaseSIMD::fQ, KFPTrack::GetChi2(), KFPTrack::GetCovarianceXYZPxPyPz(), KFParticleDatabase::GetMass(), KFPTrack::GetNDF(), i, KFParticleDatabase::Instance(), mass, KFPTrack::PxPyPz(), physmon_track_finding_ttbar::r, and KFPTrack::XvYvZv().

+ Here is the call graph for this function:

void KFParticleSIMD::Create ( KFPTrackVector track,
uint_v &  index,
const int_v &  pdg 
)

Create a particle from a set of tracks with indices "index" stored in the KFPTrackVector format. The function should be used in case if indices are random. If they are aligned please use function Load() that will benefit of the aligned memory reading and result in a faster code.

Parameters
[in]track- an array with tracks in the KFPTrackVector format
[in]index- indices of the tracks to be converted to the KFParticleSIMD object
[in]pdg- a SIMD vector with an individual pdg hypothesis for each element

Definition at line 231 of file KFParticleSIMD.cxx.

View newest version in sPHENIX GitHub at line 231 of file KFParticleSIMD.cxx

References KFPTrackVector::Covariance(), Create(), KFParticleBaseSIMD::fC, KFParticleBaseSIMD::fP, KFParticleBaseSIMD::fQ, KFParticleDatabase::GetMass(), i, KFParticleDatabase::Instance(), mass, KFPTrackVector::Parameter(), and KFPTrackVector::Q().

+ Here is the call graph for this function:

void KFParticleSIMD::Create ( KFPEmcCluster track,
uint_v &  index,
const KFParticleSIMD vertexGuess 
)

Creates gamma particles from a set of clusters of the electromagnetic calorimeter (EMC) with random indices "index". The vertex hypothesis should be provided for the estimation of the momentum.

Parameters
[in]track- an array of EMC clusters
[in]index- indices of the clusters to be converted to the KFParticleSIMD object
[in]vertexGuess- vertex guess for estimation of the momentum of created gamma particles

Definition at line 313 of file KFParticleSIMD.cxx.

View newest version in sPHENIX GitHub at line 313 of file KFParticleSIMD.cxx

References KFPEmcCluster::Covariance(), dy, dz, f, KFParticleBaseSIMD::fC, KFParticleBaseSIMD::fChi2, KFParticleBaseSIMD::fMassHypo, KFParticleBaseSIMD::fNDF, KFParticleBaseSIMD::fP, KFParticleBaseSIMD::fQ, i, KFParticleBaseSIMD::IJ(), j, Acts::UnitConstants::J, k, KFPEmcCluster::Parameter(), and KFParticleBaseSIMD::SumDaughterMass.

+ Here is the call graph for this function:

const float_v& KFParticleSIMD::E ( ) const
inline

Returns energy of the particle, fP[6].

Definition at line 162 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 162 of file KFParticleSIMD.h

References KFParticleBaseSIMD::fP.

Referenced by KFParticleFinder::NeutralDaughterDecay().

+ Here is the caller graph for this function:

float_v & KFParticleSIMD::E ( )
inline

Modifier of energy of the particle, fP[6].

Definition at line 822 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 822 of file KFParticleSIMD.h

References KFParticleBaseSIMD::E().

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetAngle ( const KFParticleSIMD p) const

Returns the opening angle between the current and the second particle in 3D.

Parameters
[in]p- the second particle

Definition at line 863 of file KFParticleSIMD.cxx.

View newest version in sPHENIX GitHub at line 863 of file KFParticleSIMD.cxx

References KFPMath::a, f, GetDStoParticle(), mask, n, and Transport().

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetAngleRZ ( const KFParticleSIMD p) const

Returns the opening angle between the current and the second particle in the RZ plane, R = sqrt(X*X+Y*Y).

Parameters
[in]p- the second particle

Definition at line 916 of file KFParticleSIMD.cxx.

View newest version in sPHENIX GitHub at line 916 of file KFParticleSIMD.cxx

References KFPMath::a, f, GetDStoParticle(), mask, n, and Transport().

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetAngleXY ( const KFParticleSIMD p) const

Returns the opening angle between the current and the second particle in the XY plane.

Parameters
[in]p- the second particle

Definition at line 889 of file KFParticleSIMD.cxx.

View newest version in sPHENIX GitHub at line 889 of file KFParticleSIMD.cxx

References KFPMath::a, f, GetDStoParticle(), mask, n, and Transport().

+ Here is the call graph for this function:

Bool_t KFParticleSIMD::GetAtProductionVertex ( ) const
inline

Returns a flag which shows if the particle is located at the production point.

Definition at line 154 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 154 of file KFParticleSIMD.h

References KFParticleBaseSIMD::fAtProductionVertex.

float_v KFParticleSIMD::GetChi2 ( ) const
inline

Returns Chi2 of the fit.

Definition at line 470 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 470 of file KFParticleSIMD.h

References KFParticleBaseSIMD::GetChi2().

Referenced by KFParticleFinder::CombinePartPart(), KFParticleFinder::ConstructTrackV0Cand(), KFParticleFinder::ConstructV0(), GetKFParticle(), and KFParticleFinder::SaveV0PrimSecCand().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

float_v KFParticleSIMD::GetCovariance ( int  i) const
inline

Returns C[i] element of the covariance matrix in the lower triangular form.

Definition at line 485 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 485 of file KFParticleSIMD.h

References KFParticleBaseSIMD::GetCovariance().

Referenced by GetDistanceFromVertexXY(), GetErrE(), GetErrPx(), GetErrPy(), GetErrPz(), GetErrS(), GetErrX(), GetErrY(), GetErrZ(), GetPseudoProperDecayTime(), and KFParticleTopoReconstructor::TransportPVTracksToPrimVertex().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

float_v KFParticleSIMD::GetCovariance ( int  i,
int  j 
) const
inline

Returns C[i,j] element of the covariance matrix.

Definition at line 490 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 490 of file KFParticleSIMD.h

References KFParticleBaseSIMD::GetCovariance().

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetDecayLength ( ) const
inline

Returns decay length.

Definition at line 538 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 538 of file KFParticleSIMD.h

References check_license::err().

Referenced by KFParticleFinder::ConstructV0(), and GetErrDecayLength().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

float_m KFParticleSIMD::GetDecayLength ( float_v &  L,
float_v &  SigmaL 
) const
inline

Calculates the decay length of the particle in the laboratory system and its error. If they are well defined returns 0, otherwise 1. The production point should be set before calling this function.

Parameters
[out]L- the decay length
[out]SigmaL- its error

Definition at line 751 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 751 of file KFParticleSIMD.h

References KFParticleBaseSIMD::GetDecayLength().

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetDecayLengthXY ( ) const
inline

Returns decay length in XY.

Definition at line 545 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 545 of file KFParticleSIMD.h

References check_license::err().

Referenced by GetErrDecayLengthXY().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

float_m KFParticleSIMD::GetDecayLengthXY ( float_v &  L,
float_v &  SigmaL 
) const
inline

Calculates the projection in the XY plane of the decay length of the particle in the laboratory system and its error. If they are well defined returns 0, otherwise 1. The production point should be set before calling this function.

Parameters
[out]L- the decay length
[out]SigmaL- its error

Definition at line 761 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 761 of file KFParticleSIMD.h

References KFParticleBaseSIMD::GetDecayLengthXY().

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetDeviationFromParticleXY ( const KFParticleSIMD p) const

Returns sqrt(Chi2/ndf) deviation from other particle in the XY plane.

Parameters
[in]p- the second particle

Definition at line 741 of file KFParticleSIMD.cxx.

View newest version in sPHENIX GitHub at line 741 of file KFParticleSIMD.cxx

References F1, F2, F3, F4, KFParticleBaseSIMD::fC, GetDStoParticle(), KFParticleBaseSIMD::MultQSQt(), and Transport().

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetDeviationFromVertexXY ( const float_v  v[],
const float_v  Cv[] = 0 
) const

Returns sqrt(Chi2/ndf) deviation from the vertex in the XY plane.

Parameters
[in]vtx[2]- { X, Y } coordinates of the vertex
[in]Cv[3]- lower-triangular part of the covariance matrix of the vertex

Definition at line 782 of file KFParticleSIMD.cxx.

View newest version in sPHENIX GitHub at line 782 of file KFParticleSIMD.cxx

References F, F1, GetDStoPoint(), i, KFParticleBaseSIMD::IJ(), KFParticleBaseSIMD::InvertCholetsky3(), j, k, and Transport().

Referenced by GetDeviationFromVertexXY().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

float_v KFParticleSIMD::GetDeviationFromVertexXY ( const KFParticleSIMD Vtx) const

Returns sqrt(Chi2/ndf) deviation from the vertex in the KFParticle format in the XY plane.

Parameters
[in]Vtx- the vertex in the KFParticle format

Definition at line 842 of file KFParticleSIMD.cxx.

View newest version in sPHENIX GitHub at line 842 of file KFParticleSIMD.cxx

References KFParticleBaseSIMD::fC, KFParticleBaseSIMD::fP, and GetDeviationFromVertexXY().

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetDistanceFromParticleXY ( const KFParticleSIMD p) const

Returns the DCA distance between the current and the second particles in the XY plane.

Parameters
[in]p- the second particle

Definition at line 724 of file KFParticleSIMD.cxx.

View newest version in sPHENIX GitHub at line 724 of file KFParticleSIMD.cxx

References dy, GetDStoParticle(), and Transport().

+ Here is the call graph for this function:

float_m KFParticleSIMD::GetDistanceFromVertexXY ( const float_v  vtx[],
float_v &  val,
float_v &  err 
) const

Calculates the DCA distance from a vertex together with the error in the XY plane. Returns "true" if calculation is failed, "false" if both value and the error are well defined.

Parameters
[in]vtx[2]- { X, Y } coordinates of the vertex
[out]val- the distance in the XY plane to the vertex
[out]err- the error of the calculated distance, takes into account errors of the particle only

Definition at line 655 of file KFParticleSIMD.cxx.

View newest version in sPHENIX GitHub at line 655 of file KFParticleSIMD.cxx

Referenced by GetDistanceFromVertexXY().

+ Here is the caller graph for this function:

float_m KFParticleSIMD::GetDistanceFromVertexXY ( const float_v  vtx[],
const float_v  Cv[],
float_v &  val,
float_v &  err 
) const

Calculates the DCA distance from a vertex together with the error in the XY plane. Returns "true" if calculation is failed, "false" if both value and the error are well defined.

Parameters
[in]vtx[2]- { X, Y } coordinates of the vertex
[in]Cv[3]- lower-triangular part of the covariance matrix of the vertex
[out]val- the distance in the XY plane to the vertex
[out]err- the error of the calculated distance, takes into account errors of the particle and vertex

Definition at line 603 of file KFParticleSIMD.cxx.

View newest version in sPHENIX GitHub at line 603 of file KFParticleSIMD.cxx

References dy, Acts::UnitConstants::e, f, GetCovariance(), GetDStoPoint(), h1, mask, physmon_track_finding_ttbar::pt, and Transport().

+ Here is the call graph for this function:

float_m KFParticleSIMD::GetDistanceFromVertexXY ( const KFParticleSIMD Vtx,
float_v &  val,
float_v &  err 
) const

Calculates the DCA distance from a vertex in the KFParticle format together with the error in the XY plane. Returns "true" if calculation is failed, "false" if both value and the error are well defined.

Parameters
[in]Vtx- the vertex in the KFParticle format
[out]val- the distance in the XY plane to the vertex
[out]err- the error of the calculated distance, takes into account errors of the particle and vertex

Definition at line 667 of file KFParticleSIMD.cxx.

View newest version in sPHENIX GitHub at line 667 of file KFParticleSIMD.cxx

References KFParticleBaseSIMD::fC, KFParticleBaseSIMD::fP, and GetDistanceFromVertexXY().

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetDistanceFromVertexXY ( const float_v  vtx[]) const

Returns the DCA distance from a vertex in the XY plane.

Parameters
[in]vtx[2]- { X, Y } coordinates of the vertex

Definition at line 693 of file KFParticleSIMD.cxx.

View newest version in sPHENIX GitHub at line 693 of file KFParticleSIMD.cxx

References check_license::err(), and GetDistanceFromVertexXY().

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetDistanceFromVertexXY ( const KFParticleSIMD Vtx) const

Returns the DCA distance from a vertex in the KFParticle format in the XY plane.

Parameters
[in]Vtx- the vertex in the KFParticle format

Definition at line 704 of file KFParticleSIMD.cxx.

View newest version in sPHENIX GitHub at line 704 of file KFParticleSIMD.cxx

References KFParticleBaseSIMD::fP, and GetDistanceFromVertexXY().

+ Here is the call graph for this function:

void KFParticleSIMD::GetDStoParticle ( const KFParticleBaseSIMD p,
float_v  dS[2],
float_v  dsdr[4][6] 
) const
inlinevirtual

Virtual method to get extrapolation parameter dS=l/p and ds/dr partial derivatives to another particle. Is defined in KFParticleSIMD.

Calculates dS = l/p parameters for two particles, where
1) l - signed distance to the DCA point with the other particle;
2) p - momentum of the particleÑŽ
dS[0] is the transport parameter for the current particle, dS[1] - for the particle "p". Also calculates partial derivatives dsdr of the parameters dS[0] and dS[1] over the state vectors of the particles:
1) dsdr[0][6] = d(dS[0])/d(param1);
2) dsdr[1][6] = d(dS[0])/d(param2);
3) dsdr[2][6] = d(dS[1])/d(param1);
4) dsdr[3][6] = d(dS[1])/d(param2);
where param1 are parameters of the current particle fP and param2 are parameters of the second particle p.fP. If "HomogeneousField" is defined KFParticleBaseSIMD::GetDStoParticleBz() is called, if "NonhomogeneousField" is defined - KFParticleBaseSIMD::GetDStoParticleCBM()

Parameters
[in]p- second particle
[out]dS[2]- transport parameters dS for the current particle (dS[0]) and the second particle "p" (dS[1])
[out]dsdr[4][6]- partial derivatives of the parameters dS[0] and dS[1] over the state vectors of the both particles

Implements KFParticleBaseSIMD.

Definition at line 1013 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 1013 of file KFParticleSIMD.h

References KFParticleBaseSIMD::GetDStoParticleBz(), and KFParticleBaseSIMD::GetDStoParticleCBM().

Referenced by KFParticleFinder::ConstructTrackV0Cand(), KFParticleFinder::ConstructV0(), GetAngle(), GetAngleRZ(), GetAngleXY(), GetDeviationFromParticleXY(), GetDistanceFromParticleXY(), and TransportToParticle().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void KFParticleSIMD::GetDStoParticleFast ( const KFParticleBaseSIMD p,
float_v  dS[2] 
) const
inlinevirtual

Virtual method to get extrapolation parameter dS=l/p to another particle. Is defined in KFParticleSIMD.

Calculates dS = l/p parameters for two particles, where
1) l - signed distance to the DCA point with the other particle;
2) p - momentum of the particleÑŽ
dS[0] is the transport parameter for the current particle, dS[1] - for the particle "p". If "HomogeneousField" is defined KFParticleBaseSIMD::GetDStoParticleBz() is called, if "NonhomogeneousField" is defined - KFParticleBaseSIMD::GetDStoParticleCBM()

Parameters
[in]p- second particle
[out]dS[2]- transport parameters dS for the current particle (dS[0]) and the second particle "p" (dS[1])

Implements KFParticleBaseSIMD.

Definition at line 1039 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 1039 of file KFParticleSIMD.h

References KFParticleBaseSIMD::GetDStoParticleBz(), and KFParticleBaseSIMD::GetDStoParticleCBM().

Referenced by KFParticleFinder::Find2DaughterDecay().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

float_v KFParticleSIMD::GetDStoPoint ( const float_v  xyz[3],
float_v  dsdr[6] 
) const
inlinevirtual

Virtual method to get extrapolation parameter dS=l/p to . Is defined in KFParticleSIMD.

Returns dS = l/p parameter, where
1) l - signed distance to the DCA point with the input xyz point;
2) p - momentum of the particle;
Also calculates partial derivatives dsdr of the parameter dS over the state vector of the current particle. If "HomogeneousField" is defined KFParticleBaseSIMD::GetDStoPointBz() is called, if "NonhomogeneousField" is defined - KFParticleBaseSIMD::GetDStoPointCBM()

Parameters
[in]xyz[3]- point, to which particle should be transported
[out]dsdr[6]= ds/dr partial derivatives of the parameter dS over the state vector of the current particle
[in]param- optional parameter, is used in case if the parameters of the particle are rotated to other coordinate system (see GetDStoPointBy() function), otherwise fP are used

Implements KFParticleBaseSIMD.

Definition at line 955 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 955 of file KFParticleSIMD.h

References KFParticleBaseSIMD::GetDStoPointBz(), and KFParticleBaseSIMD::GetDStoPointCBM().

Referenced by GetDeviationFromVertexXY(), GetDistanceFromVertexXY(), and TransportToPoint().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

float_v KFParticleSIMD::GetE ( ) const
inline

Returns energy of the particle, fP[6].

Definition at line 455 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 455 of file KFParticleSIMD.h

References KFParticleBaseSIMD::GetE().

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetErrDecayLength ( ) const
inline

Returns the error of decay length.

Definition at line 660 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 660 of file KFParticleSIMD.h

References check_license::err(), GetDecayLength(), and mask.

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetErrDecayLengthXY ( ) const
inline

Returns the error of decay length in XY.

Definition at line 669 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 669 of file KFParticleSIMD.h

References check_license::err(), GetDecayLengthXY(), and mask.

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetErrE ( ) const
inline

Returns the error of energy.

Definition at line 596 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 596 of file KFParticleSIMD.h

References GetCovariance().

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetErrEta ( ) const
inline

Returns the error of pseudorapidity.

Definition at line 624 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 624 of file KFParticleSIMD.h

References check_license::err(), GetEta(), and mask.

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetErrLifeTime ( ) const
inline

Returns the error of life time.

Definition at line 678 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 678 of file KFParticleSIMD.h

References check_license::err(), GetLifeTime(), and mask.

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetErrMass ( ) const
inline

Returns the error of mass.

Definition at line 651 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 651 of file KFParticleSIMD.h

References check_license::err(), GetMass(), and mask.

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetErrMomentum ( ) const
inline

Returns the error of momentum.

Definition at line 642 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 642 of file KFParticleSIMD.h

References check_license::err(), GetMomentum(), and mask.

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetErrP ( ) const
inline

Returns the error of momentum.

Definition at line 606 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 606 of file KFParticleSIMD.h

References check_license::err(), GetMomentum(), and mask.

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetErrPhi ( ) const
inline

Returns the error of the azimuthal angle phi.

Definition at line 633 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 633 of file KFParticleSIMD.h

References check_license::err(), GetPhi(), and mask.

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetErrPt ( ) const
inline

Returns the error of transverse momentum.

Definition at line 615 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 615 of file KFParticleSIMD.h

References check_license::err(), GetPt(), and mask.

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetErrPx ( ) const
inline

Returns the error of X-compoment of the particle momentum.

Definition at line 581 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 581 of file KFParticleSIMD.h

References GetCovariance().

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetErrPy ( ) const
inline

Returns the error of Y-compoment of the particle momentum.

Definition at line 586 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 586 of file KFParticleSIMD.h

References GetCovariance().

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetErrPz ( ) const
inline

Returns the error of Z-compoment of the particle momentum.

Definition at line 591 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 591 of file KFParticleSIMD.h

References GetCovariance().

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetErrR ( ) const
inline

Returns the error of distance to the origin of the coordinate system {0,0,0}.

Definition at line 687 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 687 of file KFParticleSIMD.h

References check_license::err(), GetR(), and mask.

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetErrS ( ) const
inline

Returns the error of decay length / momentum.

Definition at line 601 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 601 of file KFParticleSIMD.h

References GetCovariance().

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetErrX ( ) const
inline

Returns the error of X of current position.

Definition at line 566 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 566 of file KFParticleSIMD.h

References GetCovariance().

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetErrY ( ) const
inline

Returns the error of Y of current position.

Definition at line 571 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 571 of file KFParticleSIMD.h

References GetCovariance().

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetErrZ ( ) const
inline

Returns the error of Z of current position.

Definition at line 576 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 576 of file KFParticleSIMD.h

References GetCovariance().

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetEta ( ) const
inline

Returns pseudorapidity.

Definition at line 510 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 510 of file KFParticleSIMD.h

References check_license::err(), and KFParticleBaseSIMD::GetEta().

Referenced by GetErrEta().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

float_m KFParticleSIMD::GetEta ( float_v &  Eta,
float_v &  SigmaEta 
) const
inline

Calculates particle pseudorapidity and its error. If they are well defined returns 0, otherwise 1.

Parameters
[out]Eta- pseudorapidity of the particle
[out]SigmaEta- its error

Definition at line 715 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 715 of file KFParticleSIMD.h

References KFParticleBaseSIMD::GetEta().

+ Here is the call graph for this function:

void KFParticleSIMD::GetFieldValue ( const float_v  xyz[],
float_v  B[] 
) const
virtual

Virtual method to access the magnetic field

Implements KFParticleBaseSIMD.

void KFParticleSIMD::GetKFParticle ( KFParticle Part,
int  iPart = 0 
)

Copies an entry "iPart" of the current vectorised particle to the scalar KFParticle object.

Parameters
[out]Part- an output scalar particle, where element "iPart" will be copied
[in]iPart- index of the element to be copied to the scalar particle

Definition at line 1007 of file KFParticleSIMD.cxx.

View newest version in sPHENIX GitHub at line 1007 of file KFParticleSIMD.cxx

References KFParticleBase::AddDaughterId(), KFParticle::Chi2(), KFParticleBase::CleanDaughtersId(), KFParticle::CovarianceMatrix(), CovarianceMatrix(), KFParticleBaseSIMD::DaughterIds(), KFParticleBaseSIMD::fAtProductionVertex, GetChi2(), GetNDF(), KFParticleBaseSIMD::GetPDG(), GetQ(), i, KFParticleBaseSIMD::Id(), KFParticle::NDF(), KFParticle::Parameters(), Parameters(), KFParticle::Q(), KFParticle::SetAtProductionVertex(), KFParticleBase::SetId(), and KFParticleBase::SetPDG().

Referenced by KFParticleFinder::CombinePartPart(), KFParticleFinder::ConstructTrackV0Cand(), KFParticleFinder::ConstructV0(), KFParticleFinder::ExtrapolateToPV(), KFParticleFinder::FindParticles(), GetKFParticle(), KFParticleFinder::NeutralDaughterDecay(), and KFParticleFinder::SaveV0PrimSecCand().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void KFParticleSIMD::GetKFParticle ( KFParticle Part,
int  nPart = 0 
)

Copies "nPart" elements of the current vectorised particle to the array of scalar KFParticle objects.

Parameters
[out]Part- an output array of scalar particles
[in]nPart- number of elements to be copied to the array of scalar objects

Definition at line 1037 of file KFParticleSIMD.cxx.

View newest version in sPHENIX GitHub at line 1037 of file KFParticleSIMD.cxx

References GetKFParticle(), and i.

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetLifeTime ( ) const
inline

Returns life time ctau [cm].

Definition at line 552 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 552 of file KFParticleSIMD.h

References check_license::err().

Referenced by GetErrLifeTime().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

float_m KFParticleSIMD::GetLifeTime ( float_v &  T,
float_v &  SigmaT 
) const
inline

Calculates the lifetime times speed of life (ctau) [cm] of the particle in the center of mass frame and its error. If they are well defined returns 0, otherwise 1. The production point should be set before calling this function.

Parameters
[out]T- lifetime of the particle [cm]
[out]SigmaT- its error

Definition at line 772 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 772 of file KFParticleSIMD.h

References KFParticleBaseSIMD::GetLifeTime().

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetMass ( ) const
inline

Returns mass.

Definition at line 531 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 531 of file KFParticleSIMD.h

References check_license::err().

Referenced by KFParticleFinder::ConstructV0(), and GetErrMass().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

float_m KFParticleSIMD::GetMass ( float_v &  M,
float_v &  SigmaM 
) const
inline

Calculates the mass of the particle and its error. If they are well defined returns 0, otherwise 1.

Parameters
[out]M- mass of the particle
[out]SigmaM- its error

Definition at line 742 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 742 of file KFParticleSIMD.h

References KFParticleBaseSIMD::GetMass().

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetMomentum ( ) const
inline

Returns momentum.

Definition at line 524 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 524 of file KFParticleSIMD.h

References check_license::err().

Referenced by GetErrMomentum(), and GetErrP().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

float_m KFParticleSIMD::GetMomentum ( float_v &  P,
float_v &  SigmaP 
) const
inline

Calculates particle momentum and its error. If they are well defined returns 0, otherwise 1.

Parameters
[out]P- momentum of the particle
[out]SigmaP- its error

Definition at line 733 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 733 of file KFParticleSIMD.h

References KFParticleBaseSIMD::GetMomentum().

+ Here is the call graph for this function:

int_v KFParticleSIMD::GetNDF ( ) const
inline

Returns number of decrease of freedom.

Definition at line 475 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 475 of file KFParticleSIMD.h

References KFParticleBaseSIMD::GetNDF().

Referenced by KFParticleFinder::CombinePartPart(), GetKFParticle(), and KFParticleFinder::SaveV0PrimSecCand().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

float_v KFParticleSIMD::GetP ( ) const
inline

Returns momentum.

Definition at line 496 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 496 of file KFParticleSIMD.h

References check_license::err(), and KFParticleBaseSIMD::GetMomentum().

+ Here is the call graph for this function:

float_m KFParticleSIMD::GetP ( float_v &  P,
float_v &  SigmaP 
) const
inline

Calculates particle momentum and its error. If they are well defined returns 0, otherwise 1.

Parameters
[out]P- momentum of the particle
[out]SigmaP- its error

Definition at line 697 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 697 of file KFParticleSIMD.h

References KFParticleBaseSIMD::GetMomentum().

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetParameter ( int  i) const
inline

Returns P[i] parameter.

Definition at line 480 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 480 of file KFParticleSIMD.h

References KFParticleBaseSIMD::GetParameter().

Referenced by KFParticleTopoReconstructor::TransportPVTracksToPrimVertex().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

float_v KFParticleSIMD::GetPhi ( ) const
inline

Returns the azimuthal angle phi.

Definition at line 517 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 517 of file KFParticleSIMD.h

References check_license::err().

Referenced by GetErrPhi().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

float_m KFParticleSIMD::GetPhi ( float_v &  Phi,
float_v &  SigmaPhi 
) const
inline

Calculates particle polar angle at the current point and its error. If they are well defined returns 0, otherwise 1.

Parameters
[out]Phi- polar angle of the particle
[out]SigmaPhi- its error

Definition at line 724 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 724 of file KFParticleSIMD.h

References KFParticleBaseSIMD::GetPhi().

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetPseudoProperDecayTime ( const KFParticleSIMD primVertex,
const float_v &  mass,
float_v *  timeErr2 = 0 
) const

Returns the Pseudo Proper Time of the decay = (r*pt) / |pt| * M/|pt|

Parameters
[in]pV- the creation point of the particle
[in]mass- the mass of the particle
[out]timeErr2- error of the returned value, if null pointer is provided - is not calculated

Definition at line 944 of file KFParticleSIMD.cxx.

View newest version in sPHENIX GitHub at line 944 of file KFParticleSIMD.cxx

References dy, std::tr1::f1, std::tr1::f2, std::tr1::f3, std::tr1::f4, std::tr1::f5, GetCovariance(), Px(), Py(), X(), and Y().

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetPt ( ) const
inline

Returns transverse momentum.

Definition at line 503 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 503 of file KFParticleSIMD.h

References check_license::err(), and KFParticleBaseSIMD::GetPt().

Referenced by GetErrPt(), and KFParticleFinder::NeutralDaughterDecay().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

float_m KFParticleSIMD::GetPt ( float_v &  Pt,
float_v &  SigmaPt 
) const
inline

Calculates particle transverse momentum and its error. If they are well defined returns 0, otherwise 1.

Parameters
[out]Pt- transverse momentum of the particle
[out]SigmaPt- its error

Definition at line 706 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 706 of file KFParticleSIMD.h

References KFParticleBaseSIMD::GetPt().

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetPx ( ) const
inline

Retruns X component of the momentum, fP[3].

Definition at line 440 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 440 of file KFParticleSIMD.h

References KFParticleBaseSIMD::GetPx().

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetPy ( ) const
inline

Retruns Y component of the momentum, fP[4].

Definition at line 445 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 445 of file KFParticleSIMD.h

References KFParticleBaseSIMD::GetPy().

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetPz ( ) const
inline

Retruns Z component of the momentum, fP[5].

Definition at line 450 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 450 of file KFParticleSIMD.h

References KFParticleBaseSIMD::GetPz().

+ Here is the call graph for this function:

int_v KFParticleSIMD::GetQ ( ) const
inline

Returns charge of the particle.

Definition at line 465 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 465 of file KFParticleSIMD.h

References KFParticleBaseSIMD::GetQ().

Referenced by GetKFParticle().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

float_v KFParticleSIMD::GetR ( ) const
inline

Returns distance to the origin of the coordinate system {0,0,0}.

Definition at line 559 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 559 of file KFParticleSIMD.h

References check_license::err().

Referenced by GetErrR().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

float_m KFParticleSIMD::GetR ( float_v &  R,
float_v &  SigmaR 
) const
inline

Calculates the distance to the point {0,0,0} and its error. If they are well defined returns 0, otherwise 1.

Parameters
[out]R- polar angle of the particle
[out]SigmaR- its error

Definition at line 783 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 783 of file KFParticleSIMD.h

References KFParticleBaseSIMD::GetR().

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetRapidity ( ) const
inline

Returns rapidity of the particle.

Definition at line 184 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 184 of file KFParticleSIMD.h

References f, and KFParticleBaseSIMD::fP.

Referenced by KFParticleFinder::NeutralDaughterDecay().

+ Here is the caller graph for this function:

float_v KFParticleSIMD::GetS ( ) const
inline

Returns dS=l/p, l - decay length, fP[7], defined if production vertex is set.

Definition at line 460 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 460 of file KFParticleSIMD.h

References KFParticleBaseSIMD::GetS().

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetX ( ) const
inline

Retruns X coordinate of the particle, fP[0].

Definition at line 425 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 425 of file KFParticleSIMD.h

References KFParticleBaseSIMD::GetX().

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetY ( ) const
inline

Retruns Y coordinate of the particle, fP[1].

Definition at line 430 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 430 of file KFParticleSIMD.h

References KFParticleBaseSIMD::GetY().

+ Here is the call graph for this function:

float_v KFParticleSIMD::GetZ ( ) const
inline

Retruns Z coordinate of the particle, fP[2].

Definition at line 435 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 435 of file KFParticleSIMD.h

References KFParticleBaseSIMD::GetZ().

+ Here is the call graph for this function:

void KFParticleSIMD::Load ( KFPTrackVector track,
int  index,
const int_v &  pdg 
)

Create a particle from a set of consequetive tracks stored in the KFPTrackVector format starting from the index "index".

Parameters
[in]track- an array with tracks in the KFPTrackVector format
[in]index- index of the first track
[in]pdg- a SIMD vector with an individual pdg hypothesis for each element

Definition at line 257 of file KFParticleSIMD.cxx.

View newest version in sPHENIX GitHub at line 257 of file KFParticleSIMD.cxx

References KFPTrackVector::Covariance(), Create(), KFParticleBaseSIMD::fC, KFParticleBaseSIMD::fP, KFParticleBaseSIMD::fQ, KFParticleDatabase::GetMass(), i, index, KFParticleDatabase::Instance(), mass, KFPTrackVector::Parameter(), and KFPTrackVector::Q().

Referenced by KFParticleFinder::Find2DaughterDecay(), KFParticleFinder::FindParticles(), KFParticleFinder::FindTrackV0Decay(), KFParticleSIMD(), KFParticleFinder::NeutralDaughterDecay(), and KFParticleTopoReconstructor::TransportPVTracksToPrimVertex().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void KFParticleSIMD::Load ( KFPEmcCluster track,
int  index,
const KFParticleSIMD vertexGuess 
)

Create gamma particles from a set of consequetive clusters of the electromagnetic calorimeter (EMC) starting from the index "index". The vertex hypothesis should be provided for the estimation of the momentum.

Parameters
[in]track- an array with tracks in the KFPTrackVector format
[in]index- index of the first EMC cluster
[in]vertexGuess- vertex guess for estimation of the momentum of created gamma particles

Definition at line 395 of file KFParticleSIMD.cxx.

View newest version in sPHENIX GitHub at line 395 of file KFParticleSIMD.cxx

References KFPEmcCluster::Covariance(), dy, dz, f, KFParticleBaseSIMD::fC, KFParticleBaseSIMD::fChi2, KFParticleBaseSIMD::fMassHypo, KFParticleBaseSIMD::fNDF, KFParticleBaseSIMD::fP, KFParticleBaseSIMD::fQ, i, KFParticleBaseSIMD::IJ(), index, j, Acts::UnitConstants::J, k, KFPEmcCluster::Parameter(), and KFParticleBaseSIMD::SumDaughterMass.

+ Here is the call graph for this function:

const int_v& KFParticleSIMD::NDF ( ) const
inline

Returns number of decrease of freedom.

Definition at line 166 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 166 of file KFParticleSIMD.h

References KFParticleBaseSIMD::fNDF.

Referenced by KFParticleFinder::CombinePartPart(), KFParticleFinder::ConstructTrackV0Cand(), KFParticleFinder::ConstructV0(), KFParticleFinder::NeutralDaughterDecay(), and SetOneEntry().

+ Here is the caller graph for this function:

int_v & KFParticleSIMD::NDF ( )
inline

Modifier of number of decrease of freedom.

Definition at line 842 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 842 of file KFParticleSIMD.h

References KFParticleBaseSIMD::NDF().

+ Here is the call graph for this function:

void KFParticleSIMD::operator delete ( void *  ptr,
size_t   
)
inline

delete operator for the SIMD-alligned dynamic memory release

Definition at line 68 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 68 of file KFParticleSIMD.h

void KFParticleSIMD::operator delete[] ( void *  ptr,
size_t   
)
inline

delete operator for the SIMD-alligned dynamic memory release

Definition at line 69 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 69 of file KFParticleSIMD.h

void* KFParticleSIMD::operator new ( size_t  size)
inline

new operator for allocation of the SIMD-alligned dynamic memory allocation

Definition at line 64 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 64 of file KFParticleSIMD.h

References size.

void* KFParticleSIMD::operator new ( size_t  size,
void *  ptr 
)
inline

new operator for allocation of the SIMD-alligned dynamic memory allocation

Definition at line 66 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 66 of file KFParticleSIMD.h

References size.

void* KFParticleSIMD::operator new[] ( size_t  size)
inline

new operator for allocation of the SIMD-alligned dynamic memory allocation

Definition at line 65 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 65 of file KFParticleSIMD.h

References size.

void* KFParticleSIMD::operator new[] ( size_t  size,
void *  ptr 
)
inline

new operator for allocation of the SIMD-alligned dynamic memory allocation

Definition at line 67 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 67 of file KFParticleSIMD.h

References size.

void KFParticleSIMD::operator+= ( const KFParticleSIMD Daughter)
inline

Operator to add daughter to the current particle. Calls AddDaughter() function.

Parameters
[in]Daughter- the daughter particle

Definition at line 873 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 873 of file KFParticleSIMD.h

References KFParticleBaseSIMD::operator+=().

+ Here is the call graph for this function:

float_v & KFParticleSIMD::Parameter ( int  i)
inline

Modifier of P[i] parameter.

Definition at line 847 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 847 of file KFParticleSIMD.h

References KFParticleBaseSIMD::Parameter().

+ Here is the call graph for this function:

float_v * KFParticleSIMD::Parameters ( )
inline

Returns pointer to the parameters fP.

Definition at line 862 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 862 of file KFParticleSIMD.h

References KFParticleBaseSIMD::fP.

Referenced by KFParticleFinder::ConstructTrackV0Cand(), KFParticleFinder::ExtrapolateToPV(), GetKFParticle(), and SetOneEntry().

+ Here is the caller graph for this function:

const float_v& KFParticleSIMD::Px ( ) const
inline

Retruns X component of the momentum, fP[3].

Definition at line 159 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 159 of file KFParticleSIMD.h

References KFParticleBaseSIMD::fP.

Referenced by KFParticleFinder::ConstructV0(), KFParticleFinder::Find2DaughterDecay(), KFParticleFinder::FindTrackV0Decay(), and GetPseudoProperDecayTime().

+ Here is the caller graph for this function:

float_v & KFParticleSIMD::Px ( )
inline

Modifier of X component of the momentum, fP[3].

Definition at line 807 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 807 of file KFParticleSIMD.h

References KFParticleBaseSIMD::Px().

+ Here is the call graph for this function:

const float_v& KFParticleSIMD::Py ( ) const
inline

Retruns Y component of the momentum, fP[4].

Definition at line 160 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 160 of file KFParticleSIMD.h

References KFParticleBaseSIMD::fP.

Referenced by KFParticleFinder::ConstructV0(), KFParticleFinder::Find2DaughterDecay(), KFParticleFinder::FindTrackV0Decay(), and GetPseudoProperDecayTime().

+ Here is the caller graph for this function:

float_v & KFParticleSIMD::Py ( )
inline

Modifier of Y component of the momentum, fP[4].

Definition at line 812 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 812 of file KFParticleSIMD.h

References KFParticleBaseSIMD::Py().

+ Here is the call graph for this function:

const float_v& KFParticleSIMD::Pz ( ) const
inline

Retruns Z component of the momentum, fP[5].

Definition at line 161 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 161 of file KFParticleSIMD.h

References KFParticleBaseSIMD::fP.

float_v & KFParticleSIMD::Pz ( )
inline

Modifier of Z component of the momentum, fP[5].

Definition at line 817 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 817 of file KFParticleSIMD.h

References KFParticleBaseSIMD::Pz().

+ Here is the call graph for this function:

const int_v& KFParticleSIMD::Q ( ) const
inline

Returns charge of the particle.

Definition at line 164 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 164 of file KFParticleSIMD.h

References KFParticleBaseSIMD::fQ.

Referenced by SetOneEntry().

+ Here is the caller graph for this function:

int_v & KFParticleSIMD::Q ( )
inline

Modifier of charge of the particle.

Definition at line 832 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 832 of file KFParticleSIMD.h

References KFParticleBaseSIMD::Q().

+ Here is the call graph for this function:

void KFParticleSIMD::Rotate ( )

Rotates the entries of each SIMD vector of the data members.

Definition at line 282 of file KFParticleSIMD.cxx.

View newest version in sPHENIX GitHub at line 282 of file KFParticleSIMD.cxx

References KFParticleBaseSIMD::fC, KFParticleBaseSIMD::fId, KFParticleBaseSIMD::fP, KFParticleBaseSIMD::fQ, and i.

Referenced by KFParticleFinder::Find2DaughterDecay(), and KFParticleFinder::NeutralDaughterDecay().

+ Here is the caller graph for this function:

const float_v& KFParticleSIMD::S ( ) const
inline

Returns dS=l/p, l - decay length, fP[7], defined if production vertex is set.

Definition at line 163 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 163 of file KFParticleSIMD.h

References KFParticleBaseSIMD::fP.

float_v & KFParticleSIMD::S ( )
inline

Modifier of dS=l/p, l - decay length, fP[7], defined if production vertex is set.

Definition at line 827 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 827 of file KFParticleSIMD.h

References KFParticleBaseSIMD::S().

+ Here is the call graph for this function:

void KFParticleSIMD::SetOneEntry ( int  iEntry,
KFParticleSIMD part,
int  iEntryPart 
)

Copies one element of the KFParticleSIMD to one element of another KFParticleSIMD.

Parameters
[in]iEntry- index of the element of the current track, where the data will be copied
[in]part- particle, element of which should be copied to the current particle
[in]iEntryPart- index of the element of particle part, which should be copied to the current particle

Definition at line 486 of file KFParticleSIMD.cxx.

View newest version in sPHENIX GitHub at line 486 of file KFParticleSIMD.cxx

References Chi2(), CovarianceMatrix(), KFParticleBaseSIMD::fC, KFParticleBaseSIMD::fChi2, KFParticleBaseSIMD::fDaughterIds, KFParticleBaseSIMD::fId, KFParticleBaseSIMD::fNDF, KFParticleBaseSIMD::fP, KFParticleBaseSIMD::fPDG, KFParticleBaseSIMD::fQ, KFParticleBaseSIMD::GetPDG(), i, KFParticleBaseSIMD::Id(), KFParticleBaseSIMD::NDaughters(), NDF(), Parameters(), and Q().

Referenced by KFParticleFinder::ConstructV0().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void KFParticleSIMD::Transport ( float_v  dS,
const float_v *  dsdr,
float_v  P[],
float_v  C[],
float_v *  dsdr1 = 0,
float_v *  F = 0,
float_v *  F1 = 0 
) const
inline

Transports the parameters and their covariance matrix of the current particle on a length defined by the transport parameter dS = l/p, where l is the signed distance and p is the momentum of the current particle. If "HomogeneousField" is defined KFParticleBaseSIMD::TransportBz() is called, if "NonhomogeneousField" - KFParticleBaseSIMD::TransportCBM(). The obtained parameters and covariance matrix are stored to the arrays P and C respectively. P and C can be set to the parameters fP and covariance matrix fC of the current particle. In this case the particle parameters will be modified. Dependence of the transport parameter dS on the state vector of the current particle is taken into account in the covariance matrix using partial derivatives dsdr = d(dS)/d(fP). If a pointer to F is initialised the transport jacobian F = d(fP new)/d(fP old) is stored. Since dS can depend on the state vector r1 of other particle or vertex, the corelation matrix F1 = d(fP new)/d(r1) can be optionally calculated if a pointer F1 is provided. Parameters F and F1 should be either both initialised or both set to null pointer.

Parameters
[in]dS- transport parameter which defines the distance to which particle should be transported
[in]dsdr[6]= ds/dr - partial derivatives of the parameter dS over the state vector of the current particle
[out]P[8]- array, where transported parameters should be stored
[out]C[36]- array, where transported covariance matrix (8x8) should be stored in the lower triangular form
[in]dsdr1[6]= ds/dr - partial derivatives of the parameter dS over the state vector of another particle or vertex
[out]F[36]- optional parameter, transport jacobian, 6x6 matrix F = d(fP new)/d(fP old)
[out]F1[36]- optional parameter, corelation 6x6 matrix betweeen the current particle and particle or vertex with the state vector r1, to which the current particle is being transported, F1 = d(fP new)/d(r1)

Definition at line 1058 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 1058 of file KFParticleSIMD.h

References KFParticleBaseSIMD::TransportBz(), and KFParticleBaseSIMD::TransportCBM().

Referenced by GetAngle(), GetAngleRZ(), GetAngleXY(), GetDeviationFromParticleXY(), GetDeviationFromVertexXY(), GetDistanceFromParticleXY(), and GetDistanceFromVertexXY().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void KFParticleSIMD::TransportFast ( float_v  dS,
float_v  P[] 
) const
inlinevirtual

Virtual method to transport a particle parameters on a certain distance along the trajectory. Is defined in KFParticleSIMD.

Transports the parametersof the current particle on a length defined by the transport parameter dS = l/p, where l is the signed distance and p is the momentum of the current particle. If "HomogeneousField" is defined KFParticleBaseSIMD::TransportBz() is called, if "NonhomogeneousField" - KFParticleBaseSIMD::TransportCBM(). The obtained parameters are stored to the array P. P can be set to the parameters fP of the current particle. In this case the particle parameters will be modified.

Parameters
[in]dS- transport parameter which defines the distance to which particle should be transported
[out]P[8]- array, where transported parameters should be stored

Implements KFParticleBaseSIMD.

Definition at line 1090 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 1090 of file KFParticleSIMD.h

References KFParticleBaseSIMD::TransportBz(), and KFParticleBaseSIMD::TransportCBM().

Referenced by KFParticleFinder::Find2DaughterDecay().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void KFParticleSIMD::TransportToParticle ( const KFParticleSIMD p)
inline

Transports particle to the distance of closest approach to the particle p.

Parameters
[in]p- particle, to which the current particle should be transported.

Definition at line 944 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 944 of file KFParticleSIMD.h

References GetDStoParticle(), and KFParticleBaseSIMD::TransportToDS().

+ Here is the call graph for this function:

void KFParticleSIMD::TransportToPoint ( const float_v  xyz[])
inline

Transports particle to the distance of closest approach to the point xyz.

Parameters
[in]xyz[3]- point, where particle should be transported

Definition at line 925 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 925 of file KFParticleSIMD.h

References GetDStoPoint(), and KFParticleBaseSIMD::TransportToDS().

Referenced by KFParticleFinder::ConstructTrackV0Cand(), KFParticleFinder::ExtrapolateToPV(), KFParticleTopoReconstructor::GetChiToPrimVertex(), and KFParticleTopoReconstructor::TransportPVTracksToPrimVertex().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

const float_v& KFParticleSIMD::X ( ) const
inline

Retruns X coordinate of the particle, fP[0].

Definition at line 156 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 156 of file KFParticleSIMD.h

References KFParticleBaseSIMD::fP.

Referenced by KFParticleTopoReconstructor::GetChiToPrimVertex(), and GetPseudoProperDecayTime().

+ Here is the caller graph for this function:

float_v & KFParticleSIMD::X ( )
inline

Modifier of X coordinate of the particle, fP[0].

Definition at line 792 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 792 of file KFParticleSIMD.h

References KFParticleBaseSIMD::X().

+ Here is the call graph for this function:

const float_v& KFParticleSIMD::Y ( ) const
inline

Retruns Y coordinate of the particle, fP[1].

Definition at line 157 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 157 of file KFParticleSIMD.h

References KFParticleBaseSIMD::fP.

Referenced by KFParticleTopoReconstructor::GetChiToPrimVertex(), and GetPseudoProperDecayTime().

+ Here is the caller graph for this function:

float_v & KFParticleSIMD::Y ( )
inline

Modifier of Y coordinate of the particle, fP[1].

Definition at line 797 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 797 of file KFParticleSIMD.h

References KFParticleBaseSIMD::Y().

+ Here is the call graph for this function:

const float_v& KFParticleSIMD::Z ( ) const
inline

Retruns Z coordinate of the particle, fP[2].

Definition at line 158 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 158 of file KFParticleSIMD.h

References KFParticleBaseSIMD::fP.

Referenced by KFParticleTopoReconstructor::GetChiToPrimVertex(), and KFParticleFinder::NeutralDaughterDecay().

+ Here is the caller graph for this function:

float_v & KFParticleSIMD::Z ( )
inline

Modifier of Z coordinate of the particle, fP[2].

Definition at line 802 of file KFParticleSIMD.h.

View newest version in sPHENIX GitHub at line 802 of file KFParticleSIMD.h

References KFParticleBaseSIMD::Z().

+ Here is the call graph for this function:


The documentation for this class was generated from the following files: