MarlinKinfit  0.4.0
ParticleFitObject.h
Go to the documentation of this file.
1 
9 #ifndef __PARTICLEFITOBJECT_H
10 #define __PARTICLEFITOBJECT_H
11 
12 #include "BaseFitObject.h"
13 #include "FourVector.h"
14 
15 
16 // Class ParticleFitObject
18 
64  public:
67 
68 
71  );
72 
75  );
76 
78  virtual ~ParticleFitObject();
79 
81  virtual ParticleFitObject& assign (const BaseFitObject& source
82  );
83 
85  virtual bool setMass (double mass_);
87  virtual double getMass () const;
88 
90  virtual std::ostream& print4Vector (std::ostream& os
91  ) const;
92 
93  virtual FourVector getFourMomentum() const;
94 
96  virtual double getE() const;
98  virtual double getPx() const;
100  virtual double getPy() const;
102  virtual double getPz() const;
103 
105  virtual double getP() const;
107  virtual double getP2() const;
109  virtual double getPt() const;
111  virtual double getPt2() const;
112 
114  virtual double getDPx (int ilocal
115  ) const = 0;
117  virtual double getDPy (int ilocal
118  ) const = 0;
120  virtual double getDPz (int ilocal
121  ) const = 0;
123  virtual double getDE (int ilocal
124  ) const = 0;
125 
126  virtual void getDerivatives (double der[], int idim) const;
127 
129  virtual void addToGlobalChi2DerMatrixNum (double *M,
130  int idim,
131  double eps
132  );
133 
135  virtual void addToGlobalChi2DerVectorNum (double *y,
136  int idim,
137  double eps
138  );
139 
141  virtual std::ostream& print (std::ostream& os
142  ) const;
143 
144  void test1stDerivatives ();
145  void test2ndDerivatives ();
146 
148  double num1stDerivative (int ilocal,
149  double eps
150  );
151 
153  double num2ndDerivative (int ilocal1,
154  double eps1,
155  int ilocal2,
156  double eps2
157  );
158 
159  virtual double getChi2 () const;
160 
161  protected:
163  double mass;
164 
165  mutable FourVector fourMomentum;
166 
167  // this is to flag phi angle, for example
168  double paramCycl[BaseDefs::MAXPAR];
169 
170 };
171 
172 #endif // __PARTICLEFITOBJECT_H
173 
virtual ~ParticleFitObject()
Virtual destructor.
Definition: ParticleFitObject.cc:73
virtual double getDPy(int ilocal) const =0
Return d p_y / d par_ilocal (derivative of py w.r.t. local parameter ilocal)
double mass
mass of particle
Definition: ParticleFitObject.h:163
virtual double getChi2() const
Get chi squared from measured and fitted parameters.
Definition: ParticleFitObject.cc:285
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()
Default constructor.
Definition: ParticleFitObject.cc:66
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:251
virtual std::ostream & print4Vector(std::ostream &os) const
print the four-momentum (E, px, py, pz)
Definition: ParticleFitObject.cc:119
virtual double getE() const
Return E.
Definition: ParticleFitObject.cc:129
virtual ParticleFitObject & assign(const BaseFitObject &source)
Assign from anther object, if of same type.
Definition: ParticleFitObject.cc:99
virtual double getP2() const
Return p (momentum) squared.
Definition: ParticleFitObject.cc:144
virtual double getMass() const
Get mass of particle.
Definition: ParticleFitObject.cc:115
virtual void addToGlobalChi2DerMatrixNum(double *M, int idim, double eps)
Add numerically determined derivatives of chi squared to global covariance matrix.
Definition: ParticleFitObject.cc:174
virtual double getP() const
Return p (momentum)
Definition: ParticleFitObject.cc:141
virtual std::ostream & print(std::ostream &os) const
print object to ostream
Definition: ParticleFitObject.cc:155
Abstract base class for particle objects of kinematic fits.
Definition: ParticleFitObject.h:63
virtual double getPz() const
Return pz.
Definition: ParticleFitObject.cc:138
virtual double getDPz(int ilocal) const =0
Return d p_z / d par_ilocal (derivative of pz w.r.t. local parameter ilocal)
virtual double getPy() const
Return py.
Definition: ParticleFitObject.cc:135
double num1stDerivative(int ilocal, double eps)
Evaluates numerically the 1st derivative of chi2 w.r.t. a parameter.
Definition: ParticleFitObject.cc:240
virtual double getPt() const
Return pt (transverse momentum)
Definition: ParticleFitObject.cc:147
Declares class BaseFitObject.
virtual void addToGlobalChi2DerVectorNum(double *y, int idim, double eps)
Add numerically determined derivatives of chi squared to global derivative vector.
Definition: ParticleFitObject.cc:166
ParticleFitObject & operator=(const ParticleFitObject &rhs)
Assignment.
Definition: ParticleFitObject.cc:83
Abstract base class for particle objects of kinematic fits.
Definition: BaseFitObject.h:110
virtual double getPx() const
Return px.
Definition: ParticleFitObject.cc:132
virtual bool setMass(double mass_)
Set mass of particle; return=success.
Definition: ParticleFitObject.cc:90
virtual double getPt2() const
Return pt (transverse momentum) squared.
Definition: ParticleFitObject.cc:150
Declares class FourVector.