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" 65 #if GSL_MAJOR_VERSION > 1
70 if (
fType == 0)
fType = gsl_multifit_fdfsolver_lmsder;
78 if (
fVec != 0) gsl_vector_free(
fVec);
79 if (
fTmp != 0) gsl_vector_free(
fTmp);
80 if (
fCov != 0) gsl_matrix_free(
fCov);
81 #if GSL_MAJOR_VERSION > 1 82 if (fJac != 0) gsl_matrix_free(fJac);
98 if (
this == &rhs)
return *
this;
108 fSolver = gsl_multifit_fdfsolver_alloc(
fType, npoints, npar);
109 if (
fVec != 0) gsl_vector_free(
fVec);
110 fVec = gsl_vector_alloc( npar );
111 if (
fTmp != 0) gsl_vector_free(
fTmp);
112 fTmp = gsl_vector_alloc( npar );
113 if (
fCov != 0) gsl_matrix_free(
fCov);
114 fCov = gsl_matrix_alloc( npar, npar );
115 #if GSL_MAJOR_VERSION > 1 116 if (fJac != 0) gsl_matrix_free(fJac);
117 fJac = gsl_matrix_alloc( npoints, npar );
123 int Set(
const std::vector<Func> & funcVec,
const double *
x) {
126 unsigned int npts = funcVec.size();
127 if (npts == 0)
return -1;
129 unsigned int npar = funcVec[0].NDim();
140 std::copy(x,x+npar,
fVec->data);
142 gsl_vector_set_zero(
fTmp);
144 gsl_matrix_set_zero(
fCov);
145 #if GSL_MAJOR_VERSION > 1 147 gsl_matrix_set_zero(fJac);
153 if (
fSolver == 0)
return "undefined";
154 return std::string(gsl_multifit_fdfsolver_name(
fSolver) );
159 return gsl_multifit_fdfsolver_iterate(
fSolver);
163 const double *
X()
const {
165 gsl_vector *
x = gsl_multifit_fdfsolver_position(
fSolver);
172 #if GSL_MAJOR_VERSION > 1 183 static double kEpsrel = 0.0001;
184 #if GSL_MAJOR_VERSION > 1 185 gsl_multifit_fdfsolver_jac (
fSolver, fJac);
186 int ret = gsl_multifit_covar(fJac, kEpsrel,
fCov);
188 int ret = gsl_multifit_covar(
fSolver->J, kEpsrel,
fCov);
198 return gsl_multifit_test_gradient(
fVec, absTol);
205 return gsl_multifit_test_delta(
fSolver->dx,
fSolver->x, absTol, relTol);
213 if (g == 0)
return edm;
215 if (c == 0)
return edm;
216 if (
fTmp == 0)
return edm;
217 int status = gsl_blas_dgemv(CblasNoTrans, 1.0,
fCov,
fVec, 0.,
fTmp);
218 if (status == 0) status |= gsl_blas_ddot(
fVec,
fTmp, &edm);
219 if (status != 0)
return -1;
234 #if GSL_MAJOR_VERSION > 1 235 mutable gsl_matrix * fJac;
237 const gsl_multifit_fdfsolver_type *
fType;
const gsl_multifit_fdfsolver_type * fType
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...
int TestDelta(double absTol, double relTol) const
test using abs and relative tolerance |dx| < absTol + relTol*|x| for every component ...
int Set(const std::vector< Func > &funcVec, const double *x)
set the solver parameters
gsl_multifit_fdfsolver * fSolver
GSLMultiFit(const GSLMultiFit &)
Copy constructor.
const double * CovarMatrix() const
return covariance matrix of the parameters
const double * Gradient() const
gradient value at the minimum
GSLMultiFit & operator=(const GSLMultiFit &rhs)
Assignment operator.
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.
const double * X() const
parameter values at the minimum
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)
int TestGradient(double absTol) const
test gradient (ask from solver gradient vector)