Belle II Software  release-05-02-19
ParticleFitObject.h
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * See https://github.com/tferber/OrcaKinfit, forked from *
4  * https://github.com/iLCSoft/MarlinKinfit *
5  * *
6  * Further information about the fit engine and the user interface *
7  * provided in MarlinKinfit can be found at *
8  * https://www.desy.de/~blist/kinfit/doc/html/ *
9  * and in the LCNotes LC-TOOL-2009-001 and LC-TOOL-2009-004 available *
10  * from http://www-flc.desy.de/lcnotes/ *
11  * *
12  * Adopted by: Torben Ferber (torben.ferber@desy.de) (TF) *
13  * *
14  * This software is provided "as is" without any warranty. *
15  **************************************************************************/
16 
17 #ifndef __PARTICLEFITOBJECT_H
18 #define __PARTICLEFITOBJECT_H
19 
20 #include "analysis/OrcaKinFit/BaseFitObject.h"
21 #include "analysis/OrcaKinFit/FourVector.h"
22 
23 
24 // Class ParticleFitObject
26 
70 namespace Belle2 {
76  namespace OrcaKinFit {
77 
78 
79  class ParticleFitObject: public BaseFitObject {
80  public:
83 
84 
87  );
88 
91  );
92 
94  virtual ~ParticleFitObject();
95 
97  virtual ParticleFitObject& assign(const BaseFitObject& source
98  ) override;
99 
101  virtual bool setMass(double mass_);
103  virtual double getMass() const;
104 
106  virtual std::ostream& print4Vector(std::ostream& os
107  ) const;
108 
109  virtual FourVector getFourMomentum() const;
110 
112  virtual double getE() const;
114  virtual double getPx() const;
116  virtual double getPy() const;
118  virtual double getPz() const;
119 
121  virtual double getP() const;
123  virtual double getP2() const;
125  virtual double getPt() const;
127  virtual double getPt2() const;
128 
130  virtual double getDPx(int ilocal
131  ) const = 0;
133  virtual double getDPy(int ilocal
134  ) const = 0;
136  virtual double getDPz(int ilocal
137  ) const = 0;
139  virtual double getDE(int ilocal
140  ) const = 0;
141 
142  virtual void getDerivatives(double der[], int idim) const override;
143 
145  virtual void addToGlobalChi2DerMatrixNum(double* M,
146  int idim,
147  double eps
148  );
149 
151  virtual void addToGlobalChi2DerVectorNum(double* y,
152  int idim,
153  double eps
154  );
155 
157  virtual std::ostream& print(std::ostream& os
158  ) const override;
159 
160  void test1stDerivatives();
161  void test2ndDerivatives();
162 
164  double num1stDerivative(int ilocal,
165  double eps
166  );
167 
169  double num2ndDerivative(int ilocal1,
170  double eps1,
171  int ilocal2,
172  double eps2
173  );
174 
175  virtual double getChi2() const override;
176 
177  protected:
179  double mass;
180 
181  mutable FourVector fourMomentum;
182 
183  // this is to flag phi angle, for example
184  double paramCycl[BaseDefs::MAXPAR];
185 
186  };
187 
188  }// end OrcaKinFit namespace
190 } // end Belle2 namespace
191 
192 #endif // __PARTICLEFITOBJECT_H
193 
Belle2::OrcaKinFit::ParticleFitObject::print
virtual std::ostream & print(std::ostream &os) const override
print object to ostream
Definition: ParticleFitObject.cc:140
Belle2::OrcaKinFit::ParticleFitObject::getE
virtual double getE() const
Return E.
Definition: ParticleFitObject.cc:106
Belle2::OrcaKinFit::ParticleFitObject::getP
virtual double getP() const
Return p (momentum)
Definition: ParticleFitObject.cc:122
Belle2::OrcaKinFit::ParticleFitObject::mass
double mass
mass of particle
Definition: ParticleFitObject.h:193
Belle2::OrcaKinFit::ParticleFitObject::setMass
virtual bool setMass(double mass_)
Set mass of particle; return=success.
Definition: ParticleFitObject.cc:63
Belle2::OrcaKinFit::ParticleFitObject::print4Vector
virtual std::ostream & print4Vector(std::ostream &os) const
print the four-momentum (E, px, py, pz)
Definition: ParticleFitObject.cc:94
Belle2::OrcaKinFit::ParticleFitObject::getPx
virtual double getPx() const
Return px.
Definition: ParticleFitObject.cc:110
Belle2::OrcaKinFit::ParticleFitObject::getP2
virtual double getP2() const
Return p (momentum) squared.
Definition: ParticleFitObject.cc:126
Belle2::OrcaKinFit::ParticleFitObject::getMass
virtual double getMass() const
Get mass of particle.
Definition: ParticleFitObject.cc:89
Belle2::OrcaKinFit::ParticleFitObject::~ParticleFitObject
virtual ~ParticleFitObject()
Virtual destructor.
Belle2::OrcaKinFit::ParticleFitObject::addToGlobalChi2DerMatrixNum
virtual void addToGlobalChi2DerMatrixNum(double *M, int idim, double eps)
Add numerically determined derivatives of chi squared to global covariance matrix.
Definition: ParticleFitObject.cc:162
Belle2::OrcaKinFit::ParticleFitObject::num1stDerivative
double num1stDerivative(int ilocal, double eps)
Evaluates numerically the 1st derivative of chi2 w.r.t. a parameter.
Definition: ParticleFitObject.cc:230
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::OrcaKinFit::ParticleFitObject::ParticleFitObject
ParticleFitObject()
Default constructor.
Definition: ParticleFitObject.cc:39
Belle2::OrcaKinFit::ParticleFitObject::getPt
virtual double getPt() const
Return pt (transverse momentum)
Definition: ParticleFitObject.cc:130
Belle2::OrcaKinFit::ParticleFitObject::getPt2
virtual double getPt2() const
Return pt (transverse momentum) squared.
Definition: ParticleFitObject.cc:134
Belle2::OrcaKinFit::ParticleFitObject::getPy
virtual double getPy() const
Return py.
Definition: ParticleFitObject.cc:114
Belle2::OrcaKinFit::ParticleFitObject::getDE
virtual double getDE(int ilocal) const =0
Return d E / d par_ilocal (derivative of E w.r.t. local parameter ilocal)
Belle2::OrcaKinFit::ParticleFitObject::getDPx
virtual double getDPx(int ilocal) const =0
Return d p_x / d par_ilocal (derivative of px w.r.t. local parameter ilocal)
Belle2::OrcaKinFit::ParticleFitObject::operator=
ParticleFitObject & operator=(const ParticleFitObject &rhs)
Assignment.
Definition: ParticleFitObject.cc:55
Belle2::OrcaKinFit::ParticleFitObject::num2ndDerivative
double num2ndDerivative(int ilocal1, double eps1, int ilocal2, double eps2)
Evaluates numerically the 2nd derivative of chi2 w.r.t. 2 parameters.
Definition: ParticleFitObject.cc:242
Belle2::OrcaKinFit::ParticleFitObject::getDPz
virtual double getDPz(int ilocal) const =0
Return d p_z / d par_ilocal (derivative of pz w.r.t. local parameter ilocal)
Belle2::OrcaKinFit::ParticleFitObject::addToGlobalChi2DerVectorNum
virtual void addToGlobalChi2DerVectorNum(double *y, int idim, double eps)
Add numerically determined derivatives of chi squared to global derivative vector.
Definition: ParticleFitObject.cc:152
Belle2::OrcaKinFit::ParticleFitObject::getDPy
virtual double getDPy(int ilocal) const =0
Return d p_y / d par_ilocal (derivative of py w.r.t. local parameter ilocal)
Belle2::OrcaKinFit::ParticleFitObject::getPz
virtual double getPz() const
Return pz.
Definition: ParticleFitObject.cc:118
Belle2::OrcaKinFit::ParticleFitObject::assign
virtual ParticleFitObject & assign(const BaseFitObject &source) override
Assign from anther object, if of same type.
Definition: ParticleFitObject.cc:73