18 #ifndef ALIKFPARTICLEBASE_H 19 #define ALIKFPARTICLEBASE_H 22 #include "RootTypesDef.h" 35 virtual void GetFieldValue(
const Double_t xyz[], Double_t B[])
const = 0;
44 virtual Double_t
GetDStoPoint(
const Double_t xyz[])
const = 0;
52 virtual void Transport(Double_t dS, Double_t P[], Double_t C[])
const = 0;
69 void Initialize(
const Double_t Param[],
const Double_t Cov[], Int_t Charge, Double_t
Mass);
77 void SetVtxGuess(Double_t x, Double_t y, Double_t z);
96 Double_t
GetX()
const {
return fP[0]; }
97 Double_t
GetY()
const {
return fP[1]; }
98 Double_t
GetZ()
const {
return fP[2]; }
108 const Double_t &
X()
const {
return fP[0]; }
109 const Double_t &
Y()
const {
return fP[1]; }
110 const Double_t &
Z()
const {
return fP[2]; }
111 const Double_t &
Px()
const {
return fP[3]; }
112 const Double_t &
Py()
const {
return fP[4]; }
113 const Double_t &
Pz()
const {
return fP[5]; }
114 const Double_t &
E()
const {
return fP[6]; }
115 const Double_t &
S()
const {
return fP[7]; }
116 const Int_t &
Q()
const {
return fQ; }
127 Int_t
GetMomentum(Double_t &P, Double_t &SigmaP)
const;
128 Int_t
GetPt(Double_t &Pt, Double_t &SigmaPt)
const;
129 Int_t
GetEta(Double_t &Eta, Double_t &SigmaEta)
const;
130 Int_t
GetPhi(Double_t &Phi, Double_t &SigmaPhi)
const;
131 Int_t
GetMass(Double_t &M, Double_t &SigmaM)
const;
134 Int_t
GetLifeTime(Double_t &T, Double_t &SigmaT)
const;
135 Int_t
GetR(Double_t &R, Double_t &SigmaR)
const;
141 Double_t &
X() {
return fP[0]; }
142 Double_t &
Y() {
return fP[1]; }
143 Double_t &
Z() {
return fP[2]; }
144 Double_t &
Px() {
return fP[3]; }
145 Double_t &
Py() {
return fP[4]; }
146 Double_t &
Pz() {
return fP[5]; }
147 Double_t &
E() {
return fP[6]; }
148 Double_t &
S() {
return fP[7]; }
221 void TransportBz(Double_t Bz, Double_t dS, Double_t P[], Double_t C[])
const;
222 void TransportCBM(Double_t dS, Double_t P[], Double_t C[])
const;
253 void RotateXY(Double_t angle, Double_t Vtx[3]);
266 static Int_t
IJ(Int_t
i, Int_t j) {
return (j <= i) ? i * (i + 1) / 2 + j : j * (j + 1) / 2 +
i; }
268 Double_t &
Cij(Int_t
i, Int_t j) {
return fC[
IJ(i, j)]; }
270 void Convert(
bool ToProduction);
275 void GetDSIter(
const KFParticleBase &p, Double_t
const &dS, Double_t x[3], Double_t dx[3], Double_t ddx[3])
const;
277 static Bool_t
InvertSym3(
const Double_t A[], Double_t Ainv[]);
279 static void MultQSQt(
const Double_t
Q[],
const Double_t S[], Double_t SOut[]);
281 static Double_t
GetSCorrection(
const Double_t Part[],
const Double_t XYZ[]);
283 void GetMeasurement(
const Double_t XYZ[], Double_t
m[], Double_t V[])
const;
286 void SetMassConstraint(Double_t *mP, Double_t *mC, Double_t mJ[7][7], Double_t mass);
Double_t GetDistanceFromVertex(const Double_t vtx[]) const
void SetConstructMethod(Int_t m)
void TransportToProductionVertex()
void TransportToDecayVertex()
static void GetArmenterosPodolanski(KFParticleBase &positive, KFParticleBase &negative, Double_t QtAlfa[2])
void GetDStoParticleLine(const KFParticleBase &p, Double_t &dS, Double_t &dS1) const
void SetNonlinearMassConstraint(Double_t Mass)
Int_t GetMomentum(Double_t &P, Double_t &SigmaP) const
void TransportCBM(Double_t dS, Double_t P[], Double_t C[]) const
Int_t GetPhi(Double_t &Phi, Double_t &SigmaPhi) const
Double_t GetDStoPointBz(Double_t Bz, const Double_t xyz[]) const
const Double_t & Pz() const
void GetDStoParticleCBM(const KFParticleBase &p, Double_t &dS, Double_t &dS1) const
Double_t GetDStoPointBy(Double_t By, const Double_t xyz[]) const
Int_t GetMass(Double_t &M, Double_t &SigmaM) const
const Double_t & Py() const
void SubtractFromVertex(KFParticleBase &Vtx) const
void SetVtxGuess(Double_t x, Double_t y, Double_t z)
Double_t GetParameter(Int_t i) const
void operator+=(const KFParticleBase &Daughter)
const Double_t & Chi2() const
void AddDaughterWithEnergyCalc(const KFParticleBase &Daughter)
static Int_t IJ(Int_t i, Int_t j)
const Double_t & Y() const
static Double_t GetSCorrection(const Double_t Part[], const Double_t XYZ[])
virtual void GetDStoParticle(const KFParticleBase &p, Double_t &DS, Double_t &DSp) const =0
Double_t & Cij(Int_t i, Int_t j)
Int_t GetEta(Double_t &Eta, Double_t &SigmaEta) const
Int_t GetPt(Double_t &Pt, Double_t &SigmaPt) const
void GetDStoParticleBy(Double_t B, const KFParticleBase &p, Double_t &dS, Double_t &dS1) const
void SetMassHypo(Double_t m)
void TransportLine(Double_t S, Double_t P[], Double_t C[]) const
void SetProductionVertex(const KFParticleBase &Vtx)
virtual void Transport(Double_t dS, Double_t P[], Double_t C[]) const =0
virtual void GetFieldValue(const Double_t xyz[], Double_t B[]) const =0
const Double_t & GetMassHypo() const
void GetDStoParticleBz(Double_t Bz, const KFParticleBase &p, Double_t &dS, Double_t &dS1) const
Double_t GetCovariance(Int_t i, Int_t j) const
void SetMassConstraint(Double_t Mass, Double_t SigmaMass=0)
Int_t GetDecayLengthXY(Double_t &L, Double_t &SigmaL) const
const Int_t & NDF() const
void RotateXY(Double_t angle, Double_t Vtx[3])
Double_t GetDeviationFromVertex(const Double_t v[], const Double_t Cv[]=0) const
virtual ~KFParticleBase()
const std::vector< int > & DaughterIds() const
static void MultQSQt(const Double_t Q[], const Double_t S[], Double_t SOut[])
Double_t GetCovariance(Int_t i) const
Int_t GetR(Double_t &R, Double_t &SigmaR) const
void AddDaughterWithEnergyFit(const KFParticleBase &Daughter)
const Double_t & Px() const
Double_t GetDStoPointCBM(const Double_t xyz[]) const
const Double_t & GetSumDaughterMass() const
void AddDaughterId(int id)
void Convert(bool ToProduction)
void AddDaughterWithEnergyFitMC(const KFParticleBase &Daughter)
void Construct(const KFParticleBase *vDaughters[], Int_t nDaughters, const KFParticleBase *ProdVtx=nullptr, Double_t Mass=-1, Bool_t IsConstrained=0)
void GetMeasurement(const Double_t XYZ[], Double_t m[], Double_t V[]) const
void TransportToDS(Double_t dS)
void GetDSIter(const KFParticleBase &p, Double_t const &dS, Double_t x[3], Double_t dx[3], Double_t ddx[3]) const
const Double_t & X() const
Double_t GetDeviationFromParticle(const KFParticleBase &p) const
static Bool_t InvertSym3(const Double_t A[], Double_t Ainv[])
Bool_t fAtProductionVertex
Double_t & Covariance(Int_t i, Int_t j)
Double_t GetDStoPointLine(const Double_t xyz[]) const
void AddDaughter(const KFParticleBase &Daughter)
Double_t & Parameter(Int_t i)
std::vector< int > fDaughtersIds
Double_t GetDistanceFromParticle(const KFParticleBase &p) const
const Double_t & S() const
Int_t GetDecayLength(Double_t &L, Double_t &SigmaL) const
const Double_t & E() const
const Double_t & Z() const
void ConstructGammaBz(const KFParticleBase &daughter1, const KFParticleBase &daughter2, double Bz)
virtual Double_t GetDStoPoint(const Double_t xyz[]) const =0
Double_t & Covariance(Int_t i)
void TransportBz(Double_t Bz, Double_t dS, Double_t P[], Double_t C[]) const
Int_t GetLifeTime(Double_t &T, Double_t &SigmaT) const