 |
Belle II Software
release-05-02-19
|
17 #ifndef __PARTICLEFITOBJECT_H
18 #define __PARTICLEFITOBJECT_H
20 #include "analysis/OrcaKinFit/BaseFitObject.h"
21 #include "analysis/OrcaKinFit/FourVector.h"
76 namespace OrcaKinFit {
79 class ParticleFitObject:
public BaseFitObject {
101 virtual bool setMass(
double mass_);
103 virtual double getMass()
const;
109 virtual FourVector getFourMomentum()
const;
112 virtual double getE()
const;
114 virtual double getPx()
const;
116 virtual double getPy()
const;
118 virtual double getPz()
const;
121 virtual double getP()
const;
123 virtual double getP2()
const;
125 virtual double getPt()
const;
127 virtual double getPt2()
const;
130 virtual double getDPx(
int ilocal
133 virtual double getDPy(
int ilocal
136 virtual double getDPz(
int ilocal
139 virtual double getDE(
int ilocal
142 virtual void getDerivatives(
double der[],
int idim)
const override;
157 virtual std::ostream&
print(std::ostream& os
160 void test1stDerivatives();
161 void test2ndDerivatives();
175 virtual double getChi2()
const override;
181 mutable FourVector fourMomentum;
184 double paramCycl[BaseDefs::MAXPAR];
192 #endif // __PARTICLEFITOBJECT_H
virtual std::ostream & print(std::ostream &os) const override
print object to ostream
virtual double getE() const
Return E.
virtual double getP() const
Return p (momentum)
double mass
mass of particle
virtual bool setMass(double mass_)
Set mass of particle; return=success.
virtual std::ostream & print4Vector(std::ostream &os) const
print the four-momentum (E, px, py, pz)
virtual double getPx() const
Return px.
virtual double getP2() const
Return p (momentum) squared.
virtual double getMass() const
Get mass of particle.
virtual ~ParticleFitObject()
Virtual destructor.
virtual void addToGlobalChi2DerMatrixNum(double *M, int idim, double eps)
Add numerically determined derivatives of chi squared to global covariance matrix.
double num1stDerivative(int ilocal, double eps)
Evaluates numerically the 1st derivative of chi2 w.r.t. a parameter.
Abstract base class for different kinds of events.
ParticleFitObject()
Default constructor.
virtual double getPt() const
Return pt (transverse momentum)
virtual double getPt2() const
Return pt (transverse momentum) squared.
virtual double getPy() const
Return py.
virtual double getDE(int ilocal) const =0
Return d E / d par_ilocal (derivative of E w.r.t. local parameter ilocal)
virtual double getDPx(int ilocal) const =0
Return d p_x / d par_ilocal (derivative of px w.r.t. local parameter ilocal)
ParticleFitObject & operator=(const ParticleFitObject &rhs)
Assignment.
double num2ndDerivative(int ilocal1, double eps1, int ilocal2, double eps2)
Evaluates numerically the 2nd derivative of chi2 w.r.t. 2 parameters.
virtual double getDPz(int ilocal) const =0
Return d p_z / d par_ilocal (derivative of pz w.r.t. local parameter ilocal)
virtual void addToGlobalChi2DerVectorNum(double *y, int idim, double eps)
Add numerically determined derivatives of chi squared to global derivative vector.
virtual double getDPy(int ilocal) const =0
Return d p_y / d par_ilocal (derivative of py w.r.t. local parameter ilocal)
virtual double getPz() const
Return pz.
virtual ParticleFitObject & assign(const BaseFitObject &source) override
Assign from anther object, if of same type.