27 #ifndef ROOT_Math_GSLMultiFit
28 #define ROOT_Math_GSLMultiFit
30 #include "gsl/gsl_vector.h"
31 #include "gsl/gsl_matrix.h"
32 #include "gsl/gsl_multifit_nlin.h"
33 #include "gsl/gsl_blas.h"
34 #include "gsl/gsl_version.h"
66 if (
fType == 0)
fType = gsl_multifit_fdfsolver_lmsder;
74 if (
fVec != 0) gsl_vector_free(
fVec);
75 if (
fCov != 0) gsl_matrix_free(
fCov);
90 if (
this == &rhs)
return *
this;
100 fSolver = gsl_multifit_fdfsolver_alloc(
fType, npoints, npar);
105 int Set(
const std::vector<Func> & funcVec,
const double *
x) {
108 unsigned int npts = funcVec.size();
109 if (npts == 0)
return -1;
111 unsigned int npar = funcVec[0].NDim();
120 if (
fVec != 0) gsl_vector_free(
fVec);
121 fVec = gsl_vector_alloc( npar );
122 std::copy(x,x+npar,
fVec->data);
128 if (
fSolver == 0)
return "undefined";
129 return std::string(gsl_multifit_fdfsolver_name(
fSolver) );
134 return gsl_multifit_fdfsolver_iterate(
fSolver);
138 const double *
X()
const {
140 gsl_vector *
x = gsl_multifit_fdfsolver_position(
fSolver);
147 #if GSL_MAJOR_VERSION > 1
158 if (
fCov != 0) gsl_matrix_free(
fCov);
159 unsigned int npar =
fSolver->fdf->p;
160 fCov = gsl_matrix_alloc( npar, npar );
161 static double kEpsrel = 0.0001;
162 #if GSL_MAJOR_VERSION > 1
163 gsl_matrix* J = gsl_matrix_alloc(npar,npar);
164 gsl_multifit_fdfsolver_jac (
fSolver, J);
165 int ret = gsl_multifit_covar(J, kEpsrel,
fCov);
168 int ret = gsl_multifit_covar(
fSolver->J, kEpsrel,
fCov);
178 return gsl_multifit_test_gradient(
fVec, absTol);
185 return gsl_multifit_test_delta(
fSolver->dx,
fSolver->x, absTol, relTol);
193 if (g == 0)
return edm;
195 if (c == 0)
return edm;
196 gsl_vector * tmp = gsl_vector_alloc(
fSolver->fdf->p );
197 int status = gsl_blas_dgemv(CblasNoTrans, 1.0,
fCov,
fVec, 0.,tmp);
198 if (status == 0) status |= gsl_blas_ddot(
fVec, tmp, &edm);
199 gsl_vector_free(tmp);
200 if (status != 0)
return -1;
214 const gsl_multifit_fdfsolver_type *
fType;
const gsl_multifit_fdfsolver_type * fType
int TestDelta(double absTol, double relTol) const
test using abs and relative tolerance |dx| < absTol + relTol*|x| for every component ...
Namespace for new ROOT classes and functions.
GSLMultiFit, internal class for implementing GSL non linear least square GSL fitting.
GSLMultiFit(const gsl_multifit_fdfsolver_type *type=0)
Default constructor No need to specify the type so far since only one solver exists so far...
const double * Gradient() const
gradient value at the minimum
int Set(const std::vector< Func > &funcVec, const double *x)
set the solver parameters
gsl_multifit_fdfsolver * fSolver
const double * CovarMatrix() const
return covariance matrix of the parameters
GSLMultiFit(const GSLMultiFit &)
Copy constructor.
const double * X() const
parameter values at the minimum
GSLMultiFit & operator=(const GSLMultiFit &rhs)
Assignment operator.
int TestGradient(double absTol) const
test gradient (ask from solver gradient vector)
GSLMultiFitFunctionWrapper fFunc
void SetFunction(const FuncVector &f, unsigned int nres, unsigned int npar)
Fill gsl function structure from a C++ function iterator and size and number of residuals.
Namespace for new Math classes and functions.
gsl_multifit_function_fdf * GetFunc()
wrapper to a multi-dim function withtout derivatives for multi-dimensional minimization algorithm ...
void CreateSolver(unsigned int npoints, unsigned int npar)
create the minimizer from the type and size of number of fitting points and number of parameters ...
~GSLMultiFit()
Destructor (no operations)