17#ifndef __NEWFITTERGSL_H
18#define __NEWFITTERGSL_H
20#include "analysis/OrcaKinFit/BaseFitter.h"
22#include <gsl/gsl_vector.h>
23#include <gsl/gsl_matrix.h>
24#include <gsl/gsl_permutation.h>
25#include <gsl/gsl_eigen.h>
33 namespace OrcaKinFit {
59 virtual double fit()
override;
62 virtual int getError()
const override;
67 virtual double getChi2()
const override;
69 virtual int getDoF()
const override;
89 virtual void setDebug(
int debuglevel);
93 const gsl_matrix* MatM,
94 const gsl_vector* vecx,
102 const gsl_vector* vecxnew,
103 const gsl_matrix* MatM,
110 virtual gsl_matrix_view
calcZ(
int& rankA,
115 gsl_permutation* permW,
122 const gsl_vector* vecx,
127 gsl_permutation* permW,
133 const gsl_vector* vecx,
138 gsl_permutation* permW,
139 gsl_eigen_symm_workspace* eigenws,
146 void fillx(gsl_vector* vecx);
147 void fillperr(gsl_vector* vece);
150 void assembleM(gsl_matrix* MatM,
const gsl_vector* vecx,
bool errorpropagation =
false);
153 void assembleG(gsl_matrix* MatM,
const gsl_vector* vecx);
156 void scaleM(gsl_matrix* MatMscal,
const gsl_matrix* MatM,
const gsl_vector* vece);
159 void assembley(gsl_vector* vecy,
const gsl_vector* vecx);
162 void scaley(gsl_vector* vecyscal,
const gsl_vector* vecy,
const gsl_vector* vece);
165 int assembleChi2Der(gsl_vector* vecy);
168 void addConstraints(gsl_vector* vecy);
171 void assembleConstDer(gsl_matrix* MatM);
175 gsl_vector* vecdxscal,
177 const gsl_vector* vece,
179 gsl_matrix* MatMscal,
181 gsl_vector* vecyscal,
184 gsl_permutation* permW,
194 gsl_vector* vecdxhat,
195 const gsl_vector* vecdx,
196 const gsl_vector* vecdxscal,
197 const gsl_vector* vece,
198 const gsl_matrix* MatM,
199 const gsl_matrix* MatMscal,
214 const gsl_vector* vecx,
215 const gsl_vector* vecdx,
216 const gsl_vector* vece,
221 double calcMu(
const gsl_vector* vecx,
222 const gsl_vector* vece,
223 const gsl_vector* vecdx,
224 const gsl_vector* vecdxscal,
225 const gsl_vector* xnew,
226 const gsl_matrix* MatM,
227 const gsl_matrix* MatMscal,
233 const gsl_vector* vecx,
234 const gsl_vector* vece
238 const gsl_vector* vecx,
239 const gsl_vector* vece,
240 const gsl_vector* vecdx,
245 bool updateParams(gsl_vector* vecx);
250 int calcCovMatrix(gsl_matrix* MatW, gsl_permutation* permW, gsl_vector* vecx);
252 enum {NPARMAX = 50, NCONMAX = 10, NUNMMAX = 10};
264 static void ini_gsl_permutation(gsl_permutation*& p,
unsigned int size);
265 static void ini_gsl_vector(gsl_vector*& v,
int unsigned size);
266 static void ini_gsl_matrix(gsl_matrix*& m,
int unsigned size1,
unsigned int size2);
268 static void debug_print(
const gsl_matrix* m,
const char* name);
269 static void debug_print(
const gsl_vector* v,
const char* name);
272 static void add(gsl_vector* vecz,
const gsl_vector* vecx,
double a,
const gsl_vector* vecy);
275 static bool isfinite(
const gsl_vector* vec);
278 static bool isfinite(
const gsl_matrix* mat);
287 double calcpTLp(
const gsl_vector* vecdx,
288 const gsl_matrix* MatM,
295 const gsl_vector* vecyscal,
296 const gsl_matrix* MatMscal,
307 const gsl_vector* vecyscal,
308 const gsl_matrix* MatMscal,
317 const gsl_vector* vecyscal,
318 const gsl_matrix* MatMscal,
357 gsl_permutation* permW;
358 gsl_eigen_symm_workspace* eigenws;
360 unsigned int eigenwsdim;
371 double scalevals[NITMAX];
372 double fvals[NITMAX];
375 bool try2ndOrderCorr;
Abstract base class for fitting engines of kinematic fits.
A kinematic fitter using the Newton-Raphson method to solve the equations.
int ncon
total number of hard constraints
int npar
total number of parameters
double calcMu(const gsl_vector *vecx, const gsl_vector *vece, const gsl_vector *vecdx, const gsl_vector *vecdxscal, const gsl_vector *xnew, const gsl_matrix *MatM, const gsl_matrix *MatMscal, gsl_vector *vecw)
virtual double calcChi2()
Calculate the chi2.
NewFitterGSL()
Constructor.
static void MoorePenroseInverse(gsl_matrix *Ainv, gsl_matrix *A, gsl_matrix *W, gsl_vector *w, double eps=0)
Compute the Moore-Penrose pseudo-inverse A+ of A, using SVD.
virtual int getNsoft() const
Get the number of soft constraints of the last fit.
int solveSystemLU(gsl_vector *vecdxscal, double &detW, const gsl_vector *vecyscal, const gsl_matrix *MatMscal, gsl_matrix *MatW, gsl_vector *vecw, double eps)
solve system of equations Mscal*dxscal = yscal using LU decomposition
virtual gsl_matrix_view calcReducedHessian(int &rankH, gsl_matrix *MatW1, const gsl_vector *vecx, gsl_matrix *MatW2, gsl_matrix *MatW3, gsl_vector *vecw1, gsl_vector *vecw2, gsl_permutation *permW, double eps=0)
Calculate reduced Hessian; return value is a matrix view of MatW1 giving the reduced Hessian.
virtual int determineLambdas(gsl_vector *vecxnew, const gsl_matrix *MatM, const gsl_vector *vecx, gsl_matrix *MatW, gsl_vector *vecw, double eps=0)
Determine best lambda values.
virtual double fit() override
The fit method, returns the fit probability.
virtual bool initialize() override
Initialize the fitter.
int solveSystemSVD(gsl_vector *vecdxscal, const gsl_vector *vecyscal, const gsl_matrix *MatMscal, gsl_matrix *MatW, gsl_matrix *MatW2, gsl_vector *vecw, double eps)
solve system of equations Mscal*dxscal = yscal using SVD decomposition
virtual int getNcon() const
Get the number of hard constraints of the last fit.
virtual ~NewFitterGSL()
Virtual destructor.
int doLineSearch(double &alpha, gsl_vector *vecxnew, int imode, double phi0, double dphi0, double eta, double zeta, double mu, const gsl_vector *vecx, const gsl_vector *vecdx, const gsl_vector *vece, gsl_vector *vecw)
int calcNewtonDx(gsl_vector *vecdx, gsl_vector *vecdxscal, gsl_vector *vecx, const gsl_vector *vece, gsl_matrix *MatM, gsl_matrix *MatMscal, gsl_vector *vecy, gsl_vector *vecyscal, gsl_matrix *MatW, gsl_matrix *MatW2, gsl_permutation *permW, gsl_vector *vecw)
double fitprob
fit probability
int nsoft
total number of soft constraints
virtual int getError() const override
Get the error code of the last fit: 0=OK, 1=failed.
int nunm
total number of unmeasured parameters
int nit
Number of iterations.
virtual int getNunm() const
Get the number of unmeasured parameters of the last fit.
int solveSystem(gsl_vector *vecdxscal, double &detW, const gsl_vector *vecyscal, const gsl_matrix *MatMscal, gsl_matrix *MatW, gsl_matrix *MatW2, gsl_vector *vecw, double epsLU, double epsSV)
solve system of equations Mscal*dxscal = yscal
double meritFunction(double mu, const gsl_vector *vecx, const gsl_vector *vece)
virtual gsl_vector_view calcReducedHessianEigenvalues(int &rankH, gsl_matrix *MatW1, const gsl_vector *vecx, gsl_matrix *MatW2, gsl_matrix *MatW3, gsl_vector *vecw1, gsl_vector *vecw2, gsl_permutation *permW, gsl_eigen_symm_workspace *eigenws, double eps=0)
virtual int getIterations() const override
Get the number of iterations of the last fit.
virtual gsl_matrix_view calcZ(int &rankA, gsl_matrix *MatW1, gsl_matrix *MatW2, gsl_vector *vecw1, gsl_vector *vecw2, gsl_permutation *permW, double eps=0)
Calculate null space of constraints; return value is a matrix view of MatW1, with column vectors span...
virtual void calc2ndOrderCorr(gsl_vector *vecdxhat, const gsl_vector *vecxnew, const gsl_matrix *MatM, gsl_matrix *MatW, gsl_vector *vecw, double eps=0)
Calculate 2nd order correction step.
double meritFunctionDeriv(double mu, const gsl_vector *vecx, const gsl_vector *vece, const gsl_vector *vecdx, gsl_vector *vecw)
virtual double getProbability() const override
Get the fit probability of the last fit.
virtual void setDebug(int debuglevel)
Set the Debug Level.
virtual int getDoF() const override
Get the number of degrees of freedom of the last fit.
virtual double getChi2() const override
Get the chi**2 of the last fit.
virtual int getNpar() const
Get the number of all parameters of the last fit.
double calcpTLp(const gsl_vector *vecdx, const gsl_matrix *MatM, gsl_vector *vecw)
int calcLimitedDx(double &alpha, double &mu, gsl_vector *vecxnew, int imode, gsl_vector *vecx, gsl_vector *vecdxhat, const gsl_vector *vecdx, const gsl_vector *vecdxscal, const gsl_vector *vece, const gsl_matrix *MatM, const gsl_matrix *MatMscal, gsl_matrix *MatW, gsl_vector *vecw)
Abstract base class for different kinds of events.