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

The base of KFParticleSIMD class. More...

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

+ Inheritance diagram for KFParticleBaseSIMD:
+ Collaboration diagram for KFParticleBaseSIMD:

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
 
virtual void GetFieldValue (const float_v xyz[], float_v B[]) const =0
 
virtual float_v GetDStoPoint (const float_v xyz[3], float_v dsdr[6]) const =0
 
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
 
virtual void GetDStoParticle (const KFParticleBaseSIMD &p, float_v dS[2], float_v dsdr[4][6]) const =0
 
virtual void GetDStoParticleFast (const KFParticleBaseSIMD &p, float_v dS[2]) const =0
 
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
 
virtual void TransportFast (float_v dS, float_v P[]) 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
 

Static Public Member Functions

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

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

static Int_t IJ (Int_t i, Int_t j)
 
static void InvertCholetsky3 (float_v a[6])
 

Protected Attributes

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 base of KFParticleSIMD class.

Author
M.Zyzak
Date
05.02.2019
Version
1.0

Contains the main mathematics of the vectorised KF Particle. Will be merged with the KFParticleSIMD class. The class stores all the data in a format of SIMD vectors, the mathematics is fully vectorised. THe mathematics is implemented in single precision. The functionality of the vectorised and scalar classes is the same.

Definition at line 45 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 45 of file KFParticleBaseSIMD.h

Constructor & Destructor Documentation

KFParticleBaseSIMD::KFParticleBaseSIMD ( )

The default constructor, initialises the parameters by:
1) all parameters are set to 0;
2) all elements of the covariance matrix are set to 0 except Cxx=Cyy=Czz=100;
3) Q = 0;
4) chi2 is set to 0;
5) NDF = -3, since 3 parameters should be fitted: X, Y, Z.

Definition at line 28 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 28 of file KFParticleBaseSIMD.cxx

References Initialize().

+ Here is the call graph for this function:

virtual KFParticleBaseSIMD::~KFParticleBaseSIMD ( )
inlinevirtual

The default destructor.

Definition at line 116 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 116 of file KFParticleBaseSIMD.h

Member Function Documentation

void KFParticleBaseSIMD::AddDaughter ( const KFParticleBaseSIMD Daughter)

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 514 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 514 of file KFParticleBaseSIMD.cxx

References AddDaughterId(), AddDaughterWithEnergyFit(), AddDaughterWithEnergyFitMC(), fC, fConstructMethod, fMassHypo, fNDF, fP, fQ, fSFromDecay, GetQ(), i, Id(), and SumDaughterMass.

Referenced by KFParticleSIMD::AddDaughter(), Construct(), and operator+=().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void KFParticleBaseSIMD::AddDaughterId ( int_v  id)
inline

Adds index of the daughter particle.

Definition at line 291 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 291 of file KFParticleBaseSIMD.h

References fDaughterIds.

Referenced by AddDaughter(), KFParticleFinder::NeutralDaughterDecay(), and SubtractDaughter().

+ Here is the caller graph for this function:

void KFParticleBaseSIMD::AddDaughterWithEnergyFit ( const KFParticleBaseSIMD Daughter)

Adds daughter to the current particle. Uses simplifyed fast mathematics which consideres momentum and energy as independent variables and thus ignores constraint on the fixed mass. 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.

Parameters
[in]Daughter- the daughter particle

Definition at line 551 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 551 of file KFParticleBaseSIMD.cxx

References A, fC, fChi2, fNDF, fP, fQ, fSFromDecay, GetMeasurement(), GetQ(), i, IJ(), InvertCholetsky3(), j, k, and Acts::UnitConstants::m.

Referenced by AddDaughter().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void KFParticleBaseSIMD::AddDaughterWithEnergyFitMC ( const KFParticleBaseSIMD Daughter)

Adds daughter to the current particle. Uses slower but correct mathematics which requires that the masses of daughter particles stays fixed in the construction process.

Parameters
[in]Daughter- the daughter particle

Definition at line 688 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 688 of file KFParticleBaseSIMD.cxx

References A, f, fC, fChi2, fMassHypo, fNDF, fP, fQ, fSFromDecay, GetMeasurement(), GetQ(), i, IJ(), InvertCholetsky3(), j, k, Acts::UnitConstants::m, SetMassConstraint(), and SumDaughterMass.

Referenced by AddDaughter().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

const float_v& KFParticleBaseSIMD::Chi2 ( ) const
inline

Returns Chi2 of the fit.

Definition at line 153 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 153 of file KFParticleBaseSIMD.h

References fChi2.

Referenced by KFParticleSIMD::Chi2().

+ Here is the caller graph for this function:

float_v& KFParticleBaseSIMD::Chi2 ( )
inline

Modifier of Chi2 of the fit.

Definition at line 186 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 186 of file KFParticleBaseSIMD.h

References fChi2.

float_v& KFParticleBaseSIMD::Cij ( Int_t  i,
Int_t  j 
)
inlineprotected

Return an element of the covariance matrix with {i,j} indices.

Definition at line 309 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 309 of file KFParticleBaseSIMD.h

References fC, and IJ().

Referenced by SetMassConstraint().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void KFParticleBaseSIMD::CleanDaughtersId ( )
inline

Cleans the vector with the indices of daughter particles.

Definition at line 292 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 292 of file KFParticleBaseSIMD.h

References fDaughterIds.

Referenced by Construct().

+ Here is the caller graph for this function:

void KFParticleBaseSIMD::Construct ( const KFParticleBaseSIMD vDaughters[],
Int_t  nDaughters,
const KFParticleBaseSIMD ProdVtx = 0,
Float_t  Mass = -1 
)

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 1441 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 1441 of file KFParticleBaseSIMD.cxx

References AddDaughter(), CleanDaughtersId(), fAtProductionVertex, fC, fChi2, fNDF, fQ, fSFromDecay, i, SetMassConstraint(), SetNDaughters(), SetProductionVertex(), and SumDaughterMass.

Referenced by KFParticleSIMD::Construct().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

float_v& KFParticleBaseSIMD::Covariance ( Int_t  i)
inline

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

Definition at line 190 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 190 of file KFParticleBaseSIMD.h

References fC, and i.

Referenced by KFParticleSIMD::Covariance(), and RotateXY().

+ Here is the caller graph for this function:

float_v& KFParticleBaseSIMD::Covariance ( Int_t  i,
Int_t  j 
)
inline

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

Definition at line 191 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 191 of file KFParticleBaseSIMD.h

References fC, and IJ().

+ Here is the call graph for this function:

std::vector<int_v>& KFParticleBaseSIMD::DaughterIds ( )
inline

Returns the vector with the indices of daughter particles.

Definition at line 286 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 286 of file KFParticleBaseSIMD.h

References fDaughterIds.

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

+ Here is the caller graph for this function:

const float_v& KFParticleBaseSIMD::E ( ) const
inline

Returns energy of the particle, fP[6].

Definition at line 150 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 150 of file KFParticleBaseSIMD.h

References fP.

Referenced by KFParticleSIMD::E(), and GetArmenterosPodolanski().

+ Here is the caller graph for this function:

float_v& KFParticleBaseSIMD::E ( )
inline

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

Definition at line 183 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 183 of file KFParticleBaseSIMD.h

References fP.

void KFParticleBaseSIMD::GetArmenterosPodolanski ( KFParticleBaseSIMD positive,
KFParticleBaseSIMD negative,
float_v  QtAlfa[2] 
)
static

Calculates parameters for the Armenteros-Podolanski plot for two particles. Example how to use:
KFParticle PosParticle(...)
KFParticle NegParticle(...)
Gamma.ConstructGamma(PosParticle, NegParticle);
float VertexGamma[3] = {Gamma.GetX(), Gamma.GetY(), Gamma.GetZ()};
PosParticle.TransportToPoint(VertexGamma);
NegParticle.TransportToPoint(VertexGamma);
float armenterosQtAlfa[2] = {0.};
KFParticle::GetArmenterosPodolanski(PosParticle, NegParticle, armenterosQtAlfa );

Parameters
[in]positive- first particle, positive or neutral
[in]negative- second particle, negative or neutral
[out]QtAlfa[2]- parameters for the Armenteros-Podolanski plot: QtAlfa[0] = qt - projection of the momenta of the particles on the transverse direction with respect to the total momentum, same for both particles; QtAlfa[1] = (Pl+ - Pl-)/(Pl+ + Pl-) - combination of the longitudinal components.

Definition at line 3926 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 3926 of file KFParticleBaseSIMD.cxx

References alpha, E(), Acts::UnitConstants::e, f, GetPx(), GetPy(), GetPz(), and mask.

+ Here is the call graph for this function:

float_v KFParticleBaseSIMD::GetChi2 ( ) const
inline

Returns the sum of masses of the daughters.

Definition at line 141 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 141 of file KFParticleBaseSIMD.h

References fChi2.

Referenced by KFParticleSIMD::GetChi2().

+ Here is the caller graph for this function:

float_v KFParticleBaseSIMD::GetCovariance ( Int_t  i) const
inline

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

Definition at line 157 of file KFParticleBaseSIMD.h.

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

References fC, and i.

Referenced by KFParticleSIMD::GetCovariance(), and RotateXY().

+ Here is the caller graph for this function:

float_v KFParticleBaseSIMD::GetCovariance ( Int_t  i,
Int_t  j 
) const
inline

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

Definition at line 158 of file KFParticleBaseSIMD.h.

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

References fC, and IJ().

+ Here is the call graph for this function:

int_v KFParticleBaseSIMD::GetDaughterId ( int  iD) const
inline

Returns the daughter Id with the index iD.

Definition at line 287 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 287 of file KFParticleBaseSIMD.h

References fDaughterIds.

float_m KFParticleBaseSIMD::GetDecayLength ( float_v &  l,
float_v &  error 
) const

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

Parameters
[out]l- the decay length
[out]error- its error

Definition at line 272 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 272 of file KFParticleBaseSIMD.cxx

References acts::error, f, fC, fP, mask, t, ambiguity_solver_full_chain::x, y, and physmon_track_finding_ttbar::z.

Referenced by KFParticleSIMD::GetDecayLength().

+ Here is the caller graph for this function:

float_m KFParticleBaseSIMD::GetDecayLengthXY ( float_v &  l,
float_v &  error 
) const

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 the corresponding element of the return mask is set to 0, otherwise 1. The production point should be set before calling this function.

Parameters
[out]l- the decay length
[out]error- its error

Definition at line 303 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 303 of file KFParticleBaseSIMD.cxx

References acts::error, f, fC, fP, mask, t, ambiguity_solver_full_chain::x, and y.

Referenced by KFParticleSIMD::GetDecayLengthXY().

+ Here is the caller graph for this function:

float_v KFParticleBaseSIMD::GetDeviationFromParticle ( const KFParticleBaseSIMD p) const

Returns Chi2 deviation of the current particle from another particle in 3D.

Parameters
[in]p- the second particle

Definition at line 3639 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 3639 of file KFParticleBaseSIMD.cxx

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

+ Here is the call graph for this function:

float_v KFParticleBaseSIMD::GetDeviationFromVertex ( const float_v  v[],
const float_v  Cv[] = 0 
) const

Returns Chi2 deviation of the current particle from the vertex v with the covariance matrix Cv in 3D.

Parameters
[in]v[3]- coordinates of the vertex {X, Y, Z}
[in]Cv[6]- covariance matrix of the vertex {Cxx, Cxy, Cyy, Cxz, Czy, Czz}

Definition at line 3579 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 3579 of file KFParticleBaseSIMD.cxx

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

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

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

float_v KFParticleBaseSIMD::GetDeviationFromVertex ( const KFParticleBaseSIMD Vtx) const

Returns Chi2 deviation of the current particle from the vertex in the KFParticle format in 3D.

Parameters
[in]Vtx- the vertex in KFPartcile format

Definition at line 3569 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 3569 of file KFParticleBaseSIMD.cxx

References fC, fP, and GetDeviationFromVertex().

+ Here is the call graph for this function:

float_v KFParticleBaseSIMD::GetDistanceFromParticle ( const KFParticleBaseSIMD p) const

Returns the DCA distance from another particle p.

Parameters
[in]p- the second particle

Definition at line 3552 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 3552 of file KFParticleBaseSIMD.cxx

References dy, dz, GetDStoParticleFast(), and TransportFast().

Referenced by KFParticleFinder::FindTrackV0Decay().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

float_v KFParticleBaseSIMD::GetDistanceFromVertex ( const float_v  vtx[]) const

Returns the DCA distance from vertex in 3D.

Parameters
[in]vtx[3]- the vertex coordinates {X, Y, Z}

Definition at line 3538 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 3538 of file KFParticleBaseSIMD.cxx

References GetDStoPoint(), and Transport().

Referenced by GetDistanceFromVertex().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

float_v KFParticleBaseSIMD::GetDistanceFromVertex ( const KFParticleBaseSIMD Vtx) const

Returns the DCA distance from vertex in the KFParticle format in 3D.

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

Definition at line 3529 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 3529 of file KFParticleBaseSIMD.cxx

References fP, and GetDistanceFromVertex().

+ Here is the call graph for this function:

void KFParticleBaseSIMD::GetDistanceToVertexLine ( const KFParticleBaseSIMD Vertex,
float_v &  l,
float_v &  dl,
float_m *  isParticleFromVertex = 0 
) const

Calculates the distance between the particle position and the vertex together with the error. Errors of both particle and vertex are taken into account. Also optionally checks if partcile is pointing flying from the vertex, not in the direction to the vertex if the pointer to the mask isParticleFromVertex is provided.

Parameters
[in]Vertex- vertex to which the distance should be calculated
[out]l- distance between the current position of the particle and a vertex
[out]dl- the error of the calculated distance
[out]isParticleFromVertex- mask which shows if particle is flying in the direction from the vertex

Definition at line 1525 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 1525 of file KFParticleBaseSIMD.cxx

References Acts::PhysicalConstants::c, dy, dz, Acts::UnitConstants::e, f, fC, and fP.

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

+ Here is the caller graph for this function:

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

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

Implemented in KFParticleSIMD.

Referenced by GetDeviationFromParticle(), and GetMeasurement().

+ Here is the caller graph for this function:

void KFParticleBaseSIMD::GetDStoParticleB ( float_v  B[3],
const KFParticleBaseSIMD p,
float_v  dS[2],
float_v  dsdr[4][6] 
) const

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;
under the assumption of the constant homogeneous field By. 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. The particle parameters are transformed to the coordinate system, where the magnetic field B is directed along the Z axis and the function GetDStoPointBz() is called. Derivatives dsdr are transformed back to the coordinate system of the particle.

Parameters
[in]B- magnetic field By
[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

Definition at line 2770 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 2770 of file KFParticleBaseSIMD.cxx

References Acts::IntegrationTest::Bz, Acts::UnitConstants::e, f, fP, and GetDStoParticleBz().

+ Here is the call graph for this function:

void KFParticleBaseSIMD::GetDStoParticleB ( float_v  B[3],
const KFParticleBaseSIMD p,
float_v  dS[2] 
) const

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;
under the assumption of the constant homogeneous field By. dS[0] is the transport parameter for the current particle, dS[1] - for the particle "p". The particle parameters are transformed to the coordinate system, where the magnetic field B is directed along the Z axis and the function GetDStoPointBz() is called.

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

Definition at line 2841 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 2841 of file KFParticleBaseSIMD.cxx

References Acts::IntegrationTest::Bz, Acts::UnitConstants::e, f, fP, and GetDStoParticleBz().

+ Here is the call graph for this function:

void KFParticleBaseSIMD::GetDStoParticleBy ( float_v  B,
const KFParticleBaseSIMD p,
float_v  dS[2],
float_v  dsdr[4][6] 
) const

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;
under the assumption of the constant homogeneous field By. 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. The particle parameters are transformed to the coordinate system, where the main component of the magnetic field By is directed along the Z axis: x->x, y->-z, z->y, and the function GetDStoPointBz() is called. Derivatives dsdr are transformed back to the coordinate system of the particle.

Parameters
[in]B- magnetic field By
[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

Definition at line 2706 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 2706 of file KFParticleBaseSIMD.cxx

References f, fP, and GetDStoParticleBz().

Referenced by GetDStoParticleCBM().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void KFParticleBaseSIMD::GetDStoParticleBy ( float_v  B,
const KFParticleBaseSIMD p,
float_v  dS[2] 
) const

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;
under the assumption of the constant homogeneous field By. dS[0] is the transport parameter for the current particle, dS[1] - for the particle "p". The particle parameters are transformed to the coordinate system, where the main component of the magnetic field By is directed along the Z axis: x->x, y->-z, z->y, and the function GetDStoPointBz() is called.

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

Definition at line 2750 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 2750 of file KFParticleBaseSIMD.cxx

References fP, and GetDStoParticleBz().

+ Here is the call graph for this function:

void KFParticleBaseSIMD::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

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;
under the assumption of the constant homogeneous field Bz. 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 (if the pointer is not provided it is initialised with fP) and param2 are parameters of the second particle "p" (if the pointer is not provided it is initialised with p.fP). Parameters param1 and param2 should be either provided both or both set to null pointers.

Parameters
[in]B- magnetic field Bz
[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
[in]param1- optional parameter, is used in case if the parameters of the current particles are rotated to other coordinate system (see GetDStoParticleBy() function), otherwise fP are used
[in]param2- optional parameter, is used in case if the parameters of the second particles are rotated to other coordinate system (see GetDStoParticleBy() function), otherwise p.fP are used

Definition at line 1770 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 1770 of file KFParticleBaseSIMD.cxx

References KFPMath::a, KFPMath::a2, KFPMath::b, c2, physmon_vertexing::delta, dy, dz, Acts::UnitConstants::e, f, fP, fQ, GetDStoParticleLine(), i, and Transport().

Referenced by KFParticleSIMD::GetDStoParticle(), GetDStoParticleB(), GetDStoParticleBy(), and KFParticleSIMD::GetDStoParticleFast().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void KFParticleBaseSIMD::GetDStoParticleBz ( float_v  Bz,
const KFParticleBaseSIMD p,
float_v  dS[2],
const float_v *  param1 = 0,
const float_v *  param2 = 0 
) const

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;
under the assumption of the constant homogeneous field Bz. dS[0] is the transport parameter for the current particle, dS[1] - for the particle "p". Parameters param1 and param2 should be either provided both or both set to null pointers.

Parameters
[in]B- magnetic field Bz
[in]p- second particle
[out]dS[2]- transport parameters dS for the current particle (dS[0]) and the second particle "p" (dS[1])
[in]param1- optional parameter, is used in case if the parameters of the current particles are rotated to other coordinate system (see GetDStoParticleBy() function), otherwise fP are used
[in]param2- optional parameter, is used in case if the parameters of the second particles are rotated to other coordinate system (see GetDStoParticleBy() function), otherwise p.fP are used

Definition at line 2452 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 2452 of file KFParticleBaseSIMD.cxx

References c2, dy, dz, Acts::UnitConstants::e, f, fP, fQ, and GetDStoParticleLine().

+ Here is the call graph for this function:

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

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;
in case of the CBM-like nonhomogeneous magnetic field. 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. For this the y-component of the magnetic field at the position of the current particle is obtained and the GetDStoParticleBy() is called. It is assumed that particles are already close to each other and that the difference in magnetic field approximation between two particles can be neglected. If the charge of both particles is zero or if the magnetic field is zero the function GetDStoParticleLine() is called.

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

Definition at line 2993 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 2993 of file KFParticleBaseSIMD.cxx

References Acts::UnitConstants::e, f, fP, fQ, GetDStoParticleBy(), GetDStoParticleLine(), and GetFieldValue().

Referenced by KFParticleSIMD::GetDStoParticle(), and KFParticleSIMD::GetDStoParticleFast().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void KFParticleBaseSIMD::GetDStoParticleCBM ( const KFParticleBaseSIMD p,
float_v  dS[2] 
) const

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;
in case of the CBM-like nonhomogeneous magnetic field. dS[0] is the transport parameter for the current particle, dS[1] - for the particle "p". For this the y-component of the magnetic field at the position of the current particle is obtained and the GetDStoParticleBy() is called. It is assumed that particles are already close to each other and that the difference in magnetic field approximation between two particles can be neglected. If the charge of both particles is zero or if the magnetic field is zero the function GetDStoParticleLine() is called.

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

Definition at line 3030 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 3030 of file KFParticleBaseSIMD.cxx

References Acts::UnitConstants::e, f, fP, fQ, GetDStoParticleBy(), GetDStoParticleLine(), and GetFieldValue().

+ Here is the call graph for this function:

virtual void KFParticleBaseSIMD::GetDStoParticleFast ( const KFParticleBaseSIMD p,
float_v  dS[2] 
) const
pure virtual

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

Implemented in KFParticleSIMD.

Referenced by GetDistanceFromParticle().

+ Here is the caller graph for this function:

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

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;
under the assumption of the straight line trajectory. Is used for particles with charge 0 or in case of zero magnetic field. 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.

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

Definition at line 2888 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 2888 of file KFParticleBaseSIMD.cxx

References KFPMath::a2, Acts::UnitConstants::e, f, fP, and i.

Referenced by GetDStoParticleBz(), and GetDStoParticleCBM().

+ Here is the caller graph for this function:

void KFParticleBaseSIMD::GetDStoParticleLine ( const KFParticleBaseSIMD p,
float_v  dS[2] 
) const

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;
under the assumption of the straight line trajectory. Is used for particles with charge 0 or in case of zero magnetic field. dS[0] is the transport parameter for the current particle, dS[1] - for the particle "p".

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

Definition at line 2968 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 2968 of file KFParticleBaseSIMD.cxx

References Acts::UnitConstants::e, f, and fP.

virtual float_v KFParticleBaseSIMD::GetDStoPoint ( const float_v  xyz[3],
float_v  dsdr[6] 
) const
pure virtual

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

Implemented in KFParticleSIMD.

Referenced by GetDeviationFromVertex(), GetDistanceFromVertex(), GetMeasurement(), and SetProductionVertex().

+ Here is the caller graph for this function:

float_v KFParticleBaseSIMD::GetDStoPointBy ( float_v  By,
const float_v  xyz[3],
float_v  dsdr[6] 
) const

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;
under the assumption of the constant homogeneous field By. Also calculates partial derivatives dsdr of the parameter dS over the state vector of the current particle. The particle parameters are transformed to the coordinate system, where the main component of the magnetic field By is directed along the Z axis: x->x, y->-z, z->y, and the function GetDStoPointBz() is called. Derivatives dsdr are transformed back to the coordinate system of the particle.

Parameters
[in]By- magnetic field By
[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

Definition at line 1714 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 1714 of file KFParticleBaseSIMD.cxx

References fP, and GetDStoPointBz().

Referenced by GetDStoPointCBM().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

float_v KFParticleBaseSIMD::GetDStoPointBz ( float_v  Bz,
const float_v  xyz[3],
float_v  dsdr[6],
const float_v *  param = 0 
) const

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;
under the assumption of the constant homogeneous field Bz. Also calculates partial derivatives dsdr of the parameter dS over the state vector of the current particle.

Parameters
[in]B- magnetic field Bz
[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

Definition at line 1593 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 1593 of file KFParticleBaseSIMD.cxx

References KFPMath::a, Acts::PhysicalConstants::c, dy, dz, Acts::UnitConstants::e, f, fP, fQ, mask, merge_hashes::p, physmon_simulation::s, ambiguity_solver_full_chain::x, y, and physmon_track_finding_ttbar::z.

Referenced by KFParticleSIMD::GetDStoPoint(), and GetDStoPointBy().

+ Here is the caller graph for this function:

float_v KFParticleBaseSIMD::GetDStoPointCBM ( const float_v  xyz[3],
float_v  dsdr[6] 
) const

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;
in case of the CBM-like nonhomogeneous magnetic field. Also calculates partial derivatives dsdr of the parameter dS over the state vector of the current particle. For this the y-component of the magnetic field at the current position of the particle is obtained and the GetDStoPointBy() is called.

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

Definition at line 1746 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 1746 of file KFParticleBaseSIMD.cxx

References fP, GetDStoPointBy(), and GetFieldValue().

Referenced by KFParticleSIMD::GetDStoPoint().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

float_v KFParticleBaseSIMD::GetDStoPointLine ( const float_v  xyz[3],
float_v  dsdr[6] 
) const

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;
assuming the straigth line trajectory. Is used for particles with charge 0 or in case of zero magnetic field. Also calculates partial derivatives dsdr of the parameter dS over the state vector of the current particle.

Parameters
[in]xyz[3]- point where particle should be transported
[out]dsdr[6]= ds/dr partial derivatives of the parameter dS over the state vector of the current particle

Definition at line 1568 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 1568 of file KFParticleBaseSIMD.cxx

References KFPMath::a, Acts::UnitConstants::e, f, and fP.

float_v KFParticleBaseSIMD::GetE ( ) const
inline

Returns the sum of masses of the daughters.

Definition at line 138 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 138 of file KFParticleBaseSIMD.h

References fP.

Referenced by KFParticleSIMD::GetE().

+ Here is the caller graph for this function:

float_m KFParticleBaseSIMD::GetEta ( float_v &  eta,
float_v &  error 
) const

Calculates particle pseudorapidity and its error. If they are well defined the corresponding element of the return mask is set to 0, otherwise 1.

Parameters
[out]eta- pseudorapidity of the particle
[out]error- its error

Definition at line 160 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 160 of file KFParticleBaseSIMD.cxx

References KFPMath::a, KFPMath::b, Acts::PhysicalConstants::c, acts::error, eta, f, fC, fP, testing::internal::Log(), mask, and merge_hashes::p.

Referenced by KFParticleSIMD::GetEta().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

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

Virtual method to access the magnetic field

Implemented in KFParticleSIMD.

Referenced by GetDStoParticleCBM(), GetDStoPointCBM(), and TransportCBM().

+ Here is the caller graph for this function:

float_m KFParticleBaseSIMD::GetLifeTime ( float_v &  ctau,
float_v &  error 
) const

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 the corresponding element of the return mask is set to 0, otherwise 1. The production point should be set before calling this function.

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

Definition at line 331 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 331 of file KFParticleBaseSIMD.cxx

References acts::error, fC, fP, GetMass(), Acts::UnitConstants::m, and mask.

Referenced by KFParticleSIMD::GetLifeTime().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

float_m KFParticleBaseSIMD::GetMass ( float_v &  m,
float_v &  error 
) const

Calculates the mass of the particle and its error. If they are well defined the corresponding element of the return mask is set to 0, otherwise 1.

Parameters
[out]m- mass of the particle
[out]error- its error

Definition at line 241 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 241 of file KFParticleBaseSIMD.cxx

References acts::error, f, fC, fP, Acts::UnitConstants::m, Acts::UnitConstants::m2, mask, and physmon_simulation::s.

Referenced by GetLifeTime(), and KFParticleSIMD::GetMass().

+ Here is the caller graph for this function:

const float_v& KFParticleBaseSIMD::GetMassHypo ( ) const
inline

Returns the mass hypothesis.

Definition at line 123 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 123 of file KFParticleBaseSIMD.h

References fMassHypo.

void KFParticleBaseSIMD::GetMeasurement ( const KFParticleBaseSIMD daughter,
float_v  m[],
float_v  V[],
float_v  D[3][3] 
)
protected

Obtains the measurements from the current particle and the daughter to be added for the Kalman filter mathematics. If these are two first daughters they are transported to the point of the closest approach, if the third or higher daughter is added it is transported to the DCA point of the already constructed vertex. The correlations are taken into account in the covariance matrices of both measurements, the correlation matrix of two measurements is also calculated. Parameters of the current particle are modified by this function, the daughter is not changed, its parameters are stored to the output arrays after modifications.

Parameters
[in]daughter- the daughter particle to be added, stays unchanged
[out]m[8]- the output parameters of the daughter particle at the DCA point
[out]V[36]- the output covariance matrix of the daughter parameters, takes into account the correlation
[out]D[3][3]- the correlation matrix between the current and daughter particles

Definition at line 364 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 364 of file KFParticleBaseSIMD.cxx

References C, F, F1, F2, F3, F4, fC, fNDF, fP, GetDStoParticle(), GetDStoPoint(), i, IJ(), j, k, MultQSQt(), and Transport().

Referenced by AddDaughterWithEnergyFit(), AddDaughterWithEnergyFitMC(), SubtractDaughter(), SubtractFromParticle(), and SubtractFromVertex().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

float_m KFParticleBaseSIMD::GetMomentum ( float_v &  p,
float_v &  error 
) const

Calculates particle momentum and its error. If they are well defined the corresponding element of the return mask is set to 0, otherwise 1.

Parameters
[out]p- momentum of the particle
[out]error- its error

Definition at line 114 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 114 of file KFParticleBaseSIMD.cxx

References acts::error, f, fC, fP, mask, ambiguity_solver_full_chain::x, y, and physmon_track_finding_ttbar::z.

Referenced by KFParticleSIMD::GetMomentum(), and KFParticleSIMD::GetP().

+ Here is the caller graph for this function:

int_v KFParticleBaseSIMD::GetNDF ( ) const
inline

Returns the sum of masses of the daughters.

Definition at line 142 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 142 of file KFParticleBaseSIMD.h

References fNDF.

Referenced by KFParticleSIMD::GetNDF().

+ Here is the caller graph for this function:

float_v KFParticleBaseSIMD::GetParameter ( Int_t  i) const
inline

Returns P[i] parameter.

Definition at line 156 of file KFParticleBaseSIMD.h.

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

References fP, and i.

Referenced by KFParticleSIMD::GetParameter().

+ Here is the caller graph for this function:

const int_v& KFParticleBaseSIMD::GetPDG ( ) const
inline

Returns the PDG hypothesis.

Definition at line 296 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 296 of file KFParticleBaseSIMD.h

References fPDG.

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

+ Here is the caller graph for this function:

float_m KFParticleBaseSIMD::GetPhi ( float_v &  phi,
float_v &  error 
) const

Calculates particle polar angle at the current point and its error. If they are well defined the corresponding element of the return mask is set to 0, otherwise 1.

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

Definition at line 198 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 198 of file KFParticleBaseSIMD.cxx

References Acts::UnitConstants::e, acts::error, f, fC, fP, and mask.

Referenced by KFParticleSIMD::GetPhi().

+ Here is the caller graph for this function:

float_m KFParticleBaseSIMD::GetPt ( float_v &  pt,
float_v &  error 
) const

Calculates particle transverse momentum and its error. If they are well defined the corresponding element of the return mask is set to 0, otherwise 1.

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

Definition at line 138 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 138 of file KFParticleBaseSIMD.cxx

References acts::error, f, fC, fP, and mask.

Referenced by KFParticleSIMD::GetPt().

+ Here is the caller graph for this function:

float_v KFParticleBaseSIMD::GetPx ( ) const
inline

Returns the sum of masses of the daughters.

Definition at line 135 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 135 of file KFParticleBaseSIMD.h

References fP.

Referenced by GetArmenterosPodolanski(), and KFParticleSIMD::GetPx().

+ Here is the caller graph for this function:

float_v KFParticleBaseSIMD::GetPy ( ) const
inline

Returns the sum of masses of the daughters.

Definition at line 136 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 136 of file KFParticleBaseSIMD.h

References fP.

Referenced by GetArmenterosPodolanski(), and KFParticleSIMD::GetPy().

+ Here is the caller graph for this function:

float_v KFParticleBaseSIMD::GetPz ( ) const
inline

Returns the sum of masses of the daughters.

Definition at line 137 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 137 of file KFParticleBaseSIMD.h

References fP.

Referenced by GetArmenterosPodolanski(), and KFParticleSIMD::GetPz().

+ Here is the caller graph for this function:

int_v KFParticleBaseSIMD::GetQ ( ) const
inline

Returns the sum of masses of the daughters.

Definition at line 140 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 140 of file KFParticleBaseSIMD.h

References fQ.

Referenced by AddDaughter(), AddDaughterWithEnergyFit(), AddDaughterWithEnergyFitMC(), KFParticleSIMD::GetQ(), SubtractDaughter(), and SubtractFromParticle().

+ Here is the caller graph for this function:

float_m KFParticleBaseSIMD::GetR ( float_v &  r,
float_v &  error 
) const

Calculates the distance to the point {0,0,0} and its error. If they are well defined the corresponding element of the return mask is set to 0, otherwise 1.

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

Definition at line 220 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 220 of file KFParticleBaseSIMD.cxx

References Acts::UnitConstants::e, acts::error, f, fC, fP, mask, physmon_track_finding_ttbar::r, ambiguity_solver_full_chain::x, and y.

Referenced by KFParticleSIMD::GetR().

+ Here is the caller graph for this function:

float_v KFParticleBaseSIMD::GetS ( ) const
inline

Returns the sum of masses of the daughters.

Definition at line 139 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 139 of file KFParticleBaseSIMD.h

References fP.

Referenced by KFParticleSIMD::GetS().

+ Here is the caller graph for this function:

const float_v& KFParticleBaseSIMD::GetSumDaughterMass ( ) const
inline

Returns the sum of masses of the daughters.

Definition at line 124 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 124 of file KFParticleBaseSIMD.h

References SumDaughterMass.

float_v KFParticleBaseSIMD::GetX ( ) const
inline

Returns the sum of masses of the daughters.

Definition at line 132 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 132 of file KFParticleBaseSIMD.h

References fP.

Referenced by KFParticleSIMD::GetX(), and RotateXY().

+ Here is the caller graph for this function:

float_v KFParticleBaseSIMD::GetY ( ) const
inline

Returns the sum of masses of the daughters.

Definition at line 133 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 133 of file KFParticleBaseSIMD.h

References fP.

Referenced by KFParticleSIMD::GetY(), and RotateXY().

+ Here is the caller graph for this function:

float_v KFParticleBaseSIMD::GetZ ( ) const
inline

Returns the sum of masses of the daughters.

Definition at line 134 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 134 of file KFParticleBaseSIMD.h

References fP.

Referenced by KFParticleSIMD::GetZ(), and RotateXY().

+ Here is the caller graph for this function:

int_v KFParticleBaseSIMD::Id ( ) const
inline

Returns Id of the particle.

Definition at line 284 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 284 of file KFParticleBaseSIMD.h

References fId.

Referenced by AddDaughter(), KFParticleSIMD::GetKFParticle(), KFParticleFinder::NeutralDaughterDecay(), KFParticleSIMD::SetOneEntry(), and SubtractDaughter().

+ Here is the caller graph for this function:

static Int_t KFParticleBaseSIMD::IJ ( Int_t  i,
Int_t  j 
)
inlinestaticprotected

Converts a pair of indices {i,j} of the covariance matrix to one index corresponding to the triangular form.

Definition at line 305 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 305 of file KFParticleBaseSIMD.h

References i.

Referenced by AddDaughterWithEnergyFit(), AddDaughterWithEnergyFitMC(), Cij(), Covariance(), KFParticleSIMD::Create(), GetCovariance(), GetDeviationFromVertex(), KFParticleSIMD::GetDeviationFromVertexXY(), GetMeasurement(), KFParticleSIMD::Load(), SetMassConstraint(), SetProductionVertex(), and SubtractDaughter().

+ Here is the caller graph for this function:

void KFParticleBaseSIMD::Initialize ( const float_v  Param[],
const float_v  Cov[],
int_v  Charge,
float_v  Mass 
)

Sets the parameters of the particle:

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 41 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 41 of file KFParticleBaseSIMD.cxx

References energy, f, fAtProductionVertex, fC, fChi2, fMassHypo, fNDF, fP, fQ, fSFromDecay, h1, h2, i, and SumDaughterMass.

void KFParticleBaseSIMD::Initialize ( )

Initialises the parameters by default:
1) all parameters are set to 0;
2) all elements of the covariance matrix are set to 0 except Cxx=Cyy=Czz=100;
3) Q = 0;
4) chi2 is set to 0;
5) NDF = -3, since 3 parameters should be fitted: X, Y, Z.

Definition at line 91 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 91 of file KFParticleBaseSIMD.cxx

References f, fAtProductionVertex, fC, fChi2, fMassHypo, fNDF, fP, fQ, fSFromDecay, i, and SumDaughterMass.

Referenced by KFParticleSIMD::Create(), and KFParticleBaseSIMD().

+ Here is the caller graph for this function:

void KFParticleBaseSIMD::InvertCholetsky3 ( float_v  a[6])
staticprotected

Inverts symmetric 3x3 matrix a using modified Choletsky decomposition. The result is stored to the same matrix a.

Parameters
[in,out]a- 3x3 symmetric matrix

Definition at line 4035 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 4035 of file KFParticleBaseSIMD.cxx

References Acts::UnitConstants::e, f, i, j, k, and physmon_ckf_tracking::u.

Referenced by AddDaughterWithEnergyFit(), AddDaughterWithEnergyFitMC(), GetDeviationFromVertex(), KFParticleSIMD::GetDeviationFromVertexXY(), SetProductionVertex(), SubtractDaughter(), SubtractFromParticle(), and SubtractFromVertex().

+ Here is the caller graph for this function:

void KFParticleBaseSIMD::MultQSQt ( const float_v  Q[],
const float_v  S[],
float_v  SOut[],
const int  kN 
)
static

Matrix multiplication SOut = Q*S*Q^T, where Q - square matrix, S - symmetric matrix.

Parameters
[in]Q- square matrix
[in]S- input symmetric matrix
[out]SOut- output symmetric matrix
[in]kN- dimensionality of the matrices

Definition at line 4093 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 4093 of file KFParticleBaseSIMD.cxx

References i, j, and k.

Referenced by GetDeviationFromParticle(), KFParticleSIMD::GetDeviationFromParticleXY(), GetMeasurement(), SetProductionVertex(), TransportBz(), TransportCBM(), and TransportLine().

+ Here is the caller graph for this function:

int KFParticleBaseSIMD::NDaughters ( ) const
inline

Returns number of daughter particles.

Definition at line 285 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 285 of file KFParticleBaseSIMD.h

References fDaughterIds.

Referenced by KFParticleSIMD::KFParticleSIMD(), and KFParticleSIMD::SetOneEntry().

+ Here is the caller graph for this function:

const int_v& KFParticleBaseSIMD::NDF ( ) const
inline

Returns number of decrease of freedom.

Definition at line 154 of file KFParticleBaseSIMD.h.

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

References fNDF.

Referenced by KFParticleSIMD::NDF().

+ Here is the caller graph for this function:

int_v& KFParticleBaseSIMD::NDF ( )
inline

Modifier of number of decrease of freedom.

Definition at line 187 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 187 of file KFParticleBaseSIMD.h

References fNDF.

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

delete operator for the SIMD-alligned dynamic memory release

Definition at line 60 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 60 of file KFParticleBaseSIMD.h

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

delete operator for the SIMD-alligned dynamic memory release

Definition at line 61 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 61 of file KFParticleBaseSIMD.h

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

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

Definition at line 56 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 56 of file KFParticleBaseSIMD.h

References size.

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

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

Definition at line 58 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 58 of file KFParticleBaseSIMD.h

References size.

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

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

Definition at line 57 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 57 of file KFParticleBaseSIMD.h

References size.

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

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

Definition at line 59 of file KFParticleBaseSIMD.h.

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

References size.

void KFParticleBaseSIMD::operator+= ( const KFParticleBaseSIMD Daughter)

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

Parameters
[in]Daughter- the daughter particle

Definition at line 355 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 355 of file KFParticleBaseSIMD.cxx

References AddDaughter().

Referenced by KFParticleSIMD::operator+=().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

float_v& KFParticleBaseSIMD::Parameter ( Int_t  i)
inline

Modifier of P[i] parameter.

Definition at line 189 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 189 of file KFParticleBaseSIMD.h

References fP, and i.

Referenced by KFParticleSIMD::Parameter().

+ Here is the caller graph for this function:

const int_v& KFParticleBaseSIMD::PDG ( ) const
inline

Returns the PDG hypothesis.

Definition at line 297 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 297 of file KFParticleBaseSIMD.h

References fPDG.

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

+ Here is the caller graph for this function:

const float_v& KFParticleBaseSIMD::Px ( ) const
inline

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

Definition at line 147 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 147 of file KFParticleBaseSIMD.h

References fP.

Referenced by KFParticleSIMD::Px().

+ Here is the caller graph for this function:

float_v& KFParticleBaseSIMD::Px ( )
inline

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

Definition at line 180 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 180 of file KFParticleBaseSIMD.h

References fP.

const float_v& KFParticleBaseSIMD::Py ( ) const
inline

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

Definition at line 148 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 148 of file KFParticleBaseSIMD.h

References fP.

Referenced by KFParticleSIMD::Py().

+ Here is the caller graph for this function:

float_v& KFParticleBaseSIMD::Py ( )
inline

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

Definition at line 181 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 181 of file KFParticleBaseSIMD.h

References fP.

const float_v& KFParticleBaseSIMD::Pz ( ) const
inline

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

Definition at line 149 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 149 of file KFParticleBaseSIMD.h

References fP.

Referenced by KFParticleSIMD::Pz().

+ Here is the caller graph for this function:

float_v& KFParticleBaseSIMD::Pz ( )
inline

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

Definition at line 182 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 182 of file KFParticleBaseSIMD.h

References fP.

const int_v& KFParticleBaseSIMD::Q ( ) const
inline

Returns charge of the particle.

Definition at line 152 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 152 of file KFParticleBaseSIMD.h

References fQ.

Referenced by KFParticleSIMD::Q().

+ Here is the caller graph for this function:

int_v& KFParticleBaseSIMD::Q ( )
inline

Modifier of charge of the particle.

Definition at line 185 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 185 of file KFParticleBaseSIMD.h

References fQ.

void KFParticleBaseSIMD::RotateXY ( float_v  angle,
float_v  Vtx[3] 
)

Rotates the KFParticle object around OZ axis, OZ axis is set by the vertex position.

Parameters
[in]angle- angle of rotation in XY plane in [rad]
[in]Vtx[3]- position of the vertex in [cm]

Definition at line 3967 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 3967 of file KFParticleBaseSIMD.cxx

References Acts::PhysicalConstants::c, Covariance(), fP, GetCovariance(), GetX(), GetY(), GetZ(), i, j, k, physmon_simulation::s, X(), Y(), and Z().

+ Here is the call graph for this function:

const float_v& KFParticleBaseSIMD::S ( ) const
inline

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

Definition at line 151 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 151 of file KFParticleBaseSIMD.h

References fP.

Referenced by KFParticleSIMD::S().

+ Here is the caller graph for this function:

float_v& KFParticleBaseSIMD::S ( )
inline

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

Definition at line 184 of file KFParticleBaseSIMD.h.

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

References fP.

void KFParticleBaseSIMD::SetConstructMethod ( Int_t  m)
inline

Defines the construction method for the current particle (see description of fConstructMethod).

Definition at line 121 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 121 of file KFParticleBaseSIMD.h

References fConstructMethod, and Acts::UnitConstants::m.

void KFParticleBaseSIMD::SetId ( int_v  id)
inline

Sets the Id of the particle. After the construction of a particle should be set by user.

Definition at line 289 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 289 of file KFParticleBaseSIMD.h

References fId, and train_ambiguity_solver::id.

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

+ Here is the caller graph for this function:

void KFParticleBaseSIMD::SetMassConstraint ( float_v  Mass,
float_v  SigmaMass = float_v(0.f) 
)

Sets linearised mass constraint on the current particle. The constraint can be set with an uncertainty.

Parameters
[in]Mass- the mass to be set on the state vector mP
[in]SigmaMass- uncertainty of the constraint

Definition at line 1360 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 1360 of file KFParticleBaseSIMD.cxx

References Cij(), fC, fChi2, fMassHypo, fNDF, fP, i, j, Acts::UnitConstants::m2, and SumDaughterMass.

Referenced by AddDaughterWithEnergyFitMC(), Construct(), and SetNonlinearMassConstraint().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void KFParticleBaseSIMD::SetMassConstraint ( float_v *  mP,
float_v *  mC,
float_v  mJ[7][7],
float_v  mass,
float_m  mask 
)
protected

Sets the exact nonlinear mass constraint on the state vector mP with the covariance matrix mC.

Parameters
[in,out]mP- the state vector to be modified
[in,out]mC- the corresponding covariance matrix
[in,out]mJ- the Jacobian between initial and modified parameters
[in]mass- the mass to be set on the state vector mP
[in]mask- mask defines entries of the SIMD vector, for which the constraint should be applied

Definition at line 1236 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 1236 of file KFParticleBaseSIMD.cxx

References KFPMath::a, KFPMath::b, Acts::PhysicalConstants::c, perf_headwind::df, Acts::UnitConstants::e, f, i, IJ(), j, k, mask, and mass.

+ Here is the call graph for this function:

void KFParticleBaseSIMD::SetMassHypo ( float_v  m)
inline

Sets the mass hypothesis to the particle, is used when fConstructMethod = 2.

Definition at line 122 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 122 of file KFParticleBaseSIMD.h

References fMassHypo, and Acts::UnitConstants::m.

void KFParticleBaseSIMD::SetNDaughters ( int  n)
inline

Reserves the size of the vector with daughter Ids to n.

Definition at line 290 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 290 of file KFParticleBaseSIMD.h

References fDaughterIds.

Referenced by Construct(), and KFParticleSIMD::KFParticleSIMD().

+ Here is the caller graph for this function:

void KFParticleBaseSIMD::SetNoDecayLength ( )

Sets constraint on the zero decay length. When the production point is set the measurement from this particle is created at the decay point.

Definition at line 1410 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 1410 of file KFParticleBaseSIMD.cxx

References Hydro::f, fC, fChi2, fNDF, fP, h, i, j, physmon_simulation::s, and TransportToDecayVertex().

+ Here is the call graph for this function:

void KFParticleBaseSIMD::SetNonlinearMassConstraint ( float_v  Mass)

Sets the exact nonlinear mass constraint on the current particle.

Parameters
[in]mass- the mass to be set on the particle

Definition at line 1336 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 1336 of file KFParticleBaseSIMD.cxx

References energy, f, fC, fChi2, fMassHypo, fNDF, fP, mass, SetMassConstraint(), and SumDaughterMass.

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

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void KFParticleBaseSIMD::SetPDG ( int  pdg)
inline

Sets the PDG hypothesis common for all elements of the SIMD vector.

Definition at line 294 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 294 of file KFParticleBaseSIMD.h

References fPDG, and pdg.

Referenced by KFParticleFinder::CombinePartPart(), KFParticleFinder::ConstructPrimaryBG(), KFParticleFinder::ConstructTrackV0Cand(), KFParticleFinder::Find2DaughterDecay(), and KFParticleFinder::FindTrackV0Decay().

+ Here is the caller graph for this function:

void KFParticleBaseSIMD::SetPDG ( int_v &  pdg)
inline

Sets the PDG hypothesis individual for each entry of the SIMD vector.

Definition at line 295 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 295 of file KFParticleBaseSIMD.h

References fPDG, and pdg.

void KFParticleBaseSIMD::SetProductionVertex ( const KFParticleBaseSIMD Vtx)

Adds a vertex as a point-like measurement to the current particle. The eights parameter of the state vector is filled with the decay length to the momentum ratio (s = l/p). The corresponding covariances are calculated as well. The parameters of the particle are stored at the position of the production vertex.

Parameters
[in]Vtx- the assumed producation vertex

Definition at line 1050 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 1050 of file KFParticleBaseSIMD.cxx

References A, f, F, F1, fAtProductionVertex, fC, fChi2, fNDF, fP, fSFromDecay, GetDStoPoint(), i, IJ(), InvertCholetsky3(), j, k, Acts::UnitConstants::m, MultQSQt(), Transport(), TransportToDecayVertex(), X(), Y(), and Z().

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

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void KFParticleBaseSIMD::SubtractDaughter ( const KFParticleBaseSIMD Daughter)

Subtracts a daughter particle from the mother particle. The mathematics is similar to AddDaughterWithEnergyFit() but momentum is subtracted.

Parameters
[in]Daughter- the daughter particle

Definition at line 919 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 919 of file KFParticleBaseSIMD.cxx

References A, AddDaughterId(), fC, fChi2, fNDF, fP, fQ, fSFromDecay, GetMeasurement(), GetQ(), i, Id(), IJ(), InvertCholetsky3(), j, k, and Acts::UnitConstants::m.

Referenced by KFParticleFinder::NeutralDaughterDecay().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void KFParticleBaseSIMD::SubtractFromParticle ( KFParticleBaseSIMD Vtx) const

Subtract the current particle from another particle Vtx using the Kalman filter mathematics. The function is depricated and is kept for compatibility reasons. Should be replaced with SubtractDaughter().

Parameters
[in]Vtx- particle from which the current particle should be subtracted

Definition at line 3741 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 3741 of file KFParticleBaseSIMD.cxx

References fC, fChi2, fNDF, fP, fQ, fSFromDecay, GetMeasurement(), GetQ(), i, InvertCholetsky3(), j, k, and Acts::UnitConstants::m.

+ Here is the call graph for this function:

void KFParticleBaseSIMD::SubtractFromVertex ( KFParticleBaseSIMD Vtx) const

Subtract the current particle from vertex Vtx using the Kalman filter mathematics.

Parameters
[in]Vtx- vertex from which particle should be subtracted

Definition at line 3680 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 3680 of file KFParticleBaseSIMD.cxx

References fC, fChi2, fNDF, fP, GetMeasurement(), i, InvertCholetsky3(), j, k, and Acts::UnitConstants::m.

+ Here is the call graph for this function:

virtual void KFParticleBaseSIMD::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
pure virtual

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

Referenced by GetDeviationFromParticle(), GetDeviationFromVertex(), GetDistanceFromVertex(), GetDStoParticleBz(), GetMeasurement(), SetProductionVertex(), and TransportToDS().

+ Here is the caller graph for this function:

void KFParticleBaseSIMD::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

Transports the parameters and their covariance matrix of the current particle assuming constant homogeneous magnetic field Bz on the length defined by the transport parameter dS = l/p, where l is the signed distance and p is the momentum of the current particle. 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]Bz- z-component of the constant homogeneous magnetic field Bz
[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 3395 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 3395 of file KFParticleBaseSIMD.cxx

References Acts::IntegrationTest::Bz, Acts::PhysicalConstants::c, f, fC, fP, fQ, i, j, MultQSQt(), and physmon_simulation::s.

Referenced by KFParticleSIMD::Transport(), and KFParticleSIMD::TransportFast().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void KFParticleBaseSIMD::TransportBz ( float_v  Bz,
float_v  dS,
float_v  P[] 
) const

Transports the parameters of the current particle assuming constant homogeneous magnetic field Bz on the length defined by the transport parameter dS = l/p, where l is the signed distance and p is the momentum of the current particle. The obtained parameters are stored to the array P. P be set to the parameters fP of the current particle. In this case the particle parameters will be modified.

Parameters
[in]Bz- z-component of the constant homogeneous magnetic field Bz
[in]dS- transport parameter which defines the distance to which particle should be transported
[out]P[8]- array, where transported parameters should be stored

Definition at line 3485 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 3485 of file KFParticleBaseSIMD.cxx

References Acts::IntegrationTest::Bz, Acts::PhysicalConstants::c, f, fP, fQ, and physmon_simulation::s.

void KFParticleBaseSIMD::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

Transports the parameters and their covariance matrix of the current particle assuming CBM-like nonhomogeneous magnetic field on the length defined by the transport parameter dS = l/p, where l is the signed distance and p is the momentum of the current particle. 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 3059 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 3059 of file KFParticleBaseSIMD.cxx

References Acts::PhysicalConstants::c, f, fC, fP, fQ, GetFieldValue(), i, j, Acts::UnitConstants::m, MultQSQt(), n, and TransportLine().

Referenced by KFParticleSIMD::Transport(), and KFParticleSIMD::TransportFast().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void KFParticleBaseSIMD::TransportCBM ( float_v  dS,
float_v  P[] 
) const

Transports the parameters of the current particle assuming CBM-like nonhomogeneous magnetic field on the length defined by the transport parameter dS = l/p, where l is the signed distance and p is the momentum of the current particle. The obtained parameters and covariance matrix are stored to the array P. P can be set to the parameters fP. 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

Definition at line 3263 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 3263 of file KFParticleBaseSIMD.cxx

References Acts::PhysicalConstants::c, f, fP, fQ, GetFieldValue(), i, j, Acts::UnitConstants::m, n, and TransportLine().

+ Here is the call graph for this function:

virtual void KFParticleBaseSIMD::TransportFast ( float_v  dS,
float_v  P[] 
) const
pure virtual

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

Implemented in KFParticleSIMD.

Referenced by GetDistanceFromParticle().

+ Here is the caller graph for this function:

void KFParticleBaseSIMD::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
protected

Transports the parameters and their covariance matrix of the current particle assuming the straight line trajectory on the length defined by the transport parameter dS = l/p, where l is the signed distance and p is the momentum of the current particle. 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 3835 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 3835 of file KFParticleBaseSIMD.cxx

References fC, fP, i, j, and MultQSQt().

Referenced by TransportCBM(), and TransportToDSLine().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void KFParticleBaseSIMD::TransportLine ( float_v  S,
float_v  P[] 
) const
protected

Transports the parameters of the current particle assuming the straight line trajectory on the length defined by the transport parameter dS = l/p, where l is the signed distance and p is the momentum of the current particle. 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

Definition at line 3905 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 3905 of file KFParticleBaseSIMD.cxx

References fP.

void KFParticleBaseSIMD::TransportToDecayVertex ( )

Transports the particle to its decay vertex

Definition at line 1481 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 1481 of file KFParticleBaseSIMD.cxx

References Acts::UnitConstants::e, f, fAtProductionVertex, fSFromDecay, and TransportToDS().

Referenced by SetNoDecayLength(), and SetProductionVertex().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void KFParticleBaseSIMD::TransportToDS ( float_v  dS,
const float_v *  dsdr 
)

Transport the particle on a certain distane. The distance is defined by the dS=l/p parameter, where
1) l - signed distance;
2) p - momentum of the particle.

Parameters
[in]dS= l/p - distance normalised to the momentum of the particle to be transported on
[in]dsdr[6]= ds/dr partial derivatives of the parameter dS over the state vector of the current particle

Definition at line 1498 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 1498 of file KFParticleBaseSIMD.cxx

References fC, fP, fSFromDecay, and Transport().

Referenced by KFParticleFinder::ConstructTrackV0Cand(), KFParticleFinder::ConstructV0(), TransportToDecayVertex(), KFParticleSIMD::TransportToParticle(), KFParticleSIMD::TransportToPoint(), and TransportToProductionVertex().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void KFParticleBaseSIMD::TransportToDSLine ( float_v  dS,
const float_v *  dsdr 
)

Transport the particle on a certain distane assuming the linear trajectory. The distance is defined by the dS=l/p parameter, where
1) l - signed distance;
2) p - momentum of the particle.

Parameters
[in]dS= l/p - distance normalised to the momentum of the particle to be transported on
[in]dsdr[6]= ds/dr partial derivatives of the parameter dS over the state vector of the current particle

Definition at line 1511 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 1511 of file KFParticleBaseSIMD.cxx

References fC, fP, fSFromDecay, and TransportLine().

+ Here is the call graph for this function:

void KFParticleBaseSIMD::TransportToProductionVertex ( )

Transports the particle to its production vertex

Definition at line 1489 of file KFParticleBaseSIMD.cxx.

View newest version in sPHENIX GitHub at line 1489 of file KFParticleBaseSIMD.cxx

References Acts::UnitConstants::e, f, fAtProductionVertex, fP, fSFromDecay, and TransportToDS().

+ Here is the call graph for this function:

const float_v& KFParticleBaseSIMD::X ( ) const
inline

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

Definition at line 144 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 144 of file KFParticleBaseSIMD.h

References fP.

Referenced by RotateXY(), SetProductionVertex(), and KFParticleSIMD::X().

+ Here is the caller graph for this function:

float_v& KFParticleBaseSIMD::X ( )
inline

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

Definition at line 177 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 177 of file KFParticleBaseSIMD.h

References fP.

const float_v& KFParticleBaseSIMD::Y ( ) const
inline

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

Definition at line 145 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 145 of file KFParticleBaseSIMD.h

References fP.

Referenced by RotateXY(), SetProductionVertex(), and KFParticleSIMD::Y().

+ Here is the caller graph for this function:

float_v& KFParticleBaseSIMD::Y ( )
inline

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

Definition at line 178 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 178 of file KFParticleBaseSIMD.h

References fP.

const float_v& KFParticleBaseSIMD::Z ( ) const
inline

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

Definition at line 146 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 146 of file KFParticleBaseSIMD.h

References fP.

Referenced by RotateXY(), SetProductionVertex(), and KFParticleSIMD::Z().

+ Here is the caller graph for this function:

float_v& KFParticleBaseSIMD::Z ( )
inline

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

Definition at line 179 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 179 of file KFParticleBaseSIMD.h

References fP.

Member Data Documentation

Bool_t KFParticleBaseSIMD::fAtProductionVertex
protected

Flag shows if particle is at the production point.

Definition at line 330 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 330 of file KFParticleBaseSIMD.h

Referenced by Construct(), KFParticleSIMD::GetAtProductionVertex(), KFParticleSIMD::GetKFParticle(), Initialize(), KFParticleSIMD::KFParticleSIMD(), SetProductionVertex(), TransportToDecayVertex(), and TransportToProductionVertex().

Int_t KFParticleBaseSIMD::fConstructMethod
protected

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.

Definition at line 336 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 336 of file KFParticleBaseSIMD.h

Referenced by AddDaughter(), and SetConstructMethod().

std::vector<int_v> KFParticleBaseSIMD::fDaughterIds
protected

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.

Definition at line 341 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 341 of file KFParticleBaseSIMD.h

Referenced by AddDaughterId(), CleanDaughtersId(), DaughterIds(), GetDaughterId(), KFParticleSIMD::KFParticleSIMD(), NDaughters(), SetNDaughters(), and KFParticleSIMD::SetOneEntry().

int_v KFParticleBaseSIMD::fId
protected

Id of the particle.

Definition at line 329 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 329 of file KFParticleBaseSIMD.h

Referenced by Id(), KFParticleSIMD::KFParticleSIMD(), KFParticleSIMD::Rotate(), SetId(), and KFParticleSIMD::SetOneEntry().

float_v KFParticleBaseSIMD::fMassHypo
protected

The mass hypothesis, used for the constraints during particle construction.

Definition at line 328 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 328 of file KFParticleBaseSIMD.h

Referenced by AddDaughter(), AddDaughterWithEnergyFitMC(), KFParticleSIMD::Create(), GetMassHypo(), Initialize(), KFParticleSIMD::Load(), SetMassConstraint(), SetMassHypo(), and SetNonlinearMassConstraint().

float_v KFParticleBaseSIMD::fP[8]
protected

Particle parameters { X, Y, Z, Px, Py, Pz, E, S[=DecayLength/P]}.

Definition at line 321 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 321 of file KFParticleBaseSIMD.h

Referenced by AddDaughter(), AddDaughterWithEnergyFit(), AddDaughterWithEnergyFitMC(), KFParticleSIMD::Create(), E(), KFParticleSIMD::E(), GetDecayLength(), GetDecayLengthXY(), GetDeviationFromVertex(), KFParticleSIMD::GetDeviationFromVertexXY(), GetDistanceFromVertex(), KFParticleSIMD::GetDistanceFromVertexXY(), GetDistanceToVertexLine(), GetDStoParticleB(), GetDStoParticleBy(), GetDStoParticleBz(), GetDStoParticleCBM(), GetDStoParticleLine(), GetDStoPointBy(), GetDStoPointBz(), GetDStoPointCBM(), GetDStoPointLine(), GetE(), GetEta(), GetLifeTime(), GetMass(), GetMeasurement(), GetMomentum(), GetParameter(), GetPhi(), GetPt(), GetPx(), GetPy(), GetPz(), GetR(), KFParticleSIMD::GetRapidity(), GetS(), GetX(), GetY(), GetZ(), Initialize(), KFParticleSIMD::KFParticleSIMD(), KFParticleSIMD::Load(), Parameter(), KFParticleSIMD::Parameters(), Px(), KFParticleSIMD::Px(), Py(), KFParticleSIMD::Py(), Pz(), KFParticleSIMD::Pz(), KFParticleSIMD::Rotate(), RotateXY(), S(), KFParticleSIMD::S(), SetMassConstraint(), SetNoDecayLength(), SetNonlinearMassConstraint(), KFParticleSIMD::SetOneEntry(), SetProductionVertex(), SubtractDaughter(), SubtractFromParticle(), SubtractFromVertex(), TransportBz(), TransportCBM(), TransportLine(), TransportToDS(), TransportToDSLine(), TransportToProductionVertex(), X(), KFParticleSIMD::X(), Y(), KFParticleSIMD::Y(), Z(), and KFParticleSIMD::Z().

int_v KFParticleBaseSIMD::fPDG
protected

The PDG hypothesis assigned to the particle.

Definition at line 331 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 331 of file KFParticleBaseSIMD.h

Referenced by GetPDG(), KFParticleSIMD::KFParticleSIMD(), PDG(), KFParticleSIMD::SetOneEntry(), and SetPDG().

float_v KFParticleBaseSIMD::fSFromDecay
protected
float_v KFParticleBaseSIMD::SumDaughterMass
protected

Sum of the daughter particles masses. Needed to set the constraint on the minimum mass during particle construction.

Definition at line 327 of file KFParticleBaseSIMD.h.

View newest version in sPHENIX GitHub at line 327 of file KFParticleBaseSIMD.h

Referenced by AddDaughter(), AddDaughterWithEnergyFitMC(), Construct(), KFParticleSIMD::Create(), GetSumDaughterMass(), Initialize(), KFParticleSIMD::Load(), SetMassConstraint(), and SetNonlinearMassConstraint().


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