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 GSL_MAJOR_VERSION > 1
71 if (
fType ==
nullptr)
fType = gsl_multifit_fdfsolver_lmsder;
79 if (
fVec !=
nullptr) gsl_vector_free(
fVec);
80 if (
fTmp !=
nullptr) gsl_vector_free(
fTmp);
81 if (
fCov !=
nullptr) gsl_matrix_free(
fCov);
82#if GSL_MAJOR_VERSION > 1
83 if (fJac !=
nullptr) gsl_matrix_free(fJac);
96 fSolver = gsl_multifit_fdfsolver_alloc(
fType, npoints, npar);
97 if (
fVec !=
nullptr) gsl_vector_free(
fVec);
98 fVec = gsl_vector_alloc( npar );
99 if (
fTmp !=
nullptr) gsl_vector_free(
fTmp);
100 fTmp = gsl_vector_alloc( npar );
101 if (
fCov !=
nullptr) gsl_matrix_free(
fCov);
102 fCov = gsl_matrix_alloc( npar, npar );
103#if GSL_MAJOR_VERSION > 1
104 if (fJac !=
nullptr) gsl_matrix_free(fJac);
105 fJac = gsl_matrix_alloc( npoints, npar );
111 int Set(
const std::vector<Func> & funcVec,
const double *
x) {
114 unsigned int npts = funcVec.size();
115 if (npts == 0)
return -1;
117 unsigned int npar = funcVec[0].NDim();
127 assert(
fVec !=
nullptr);
128 std::copy(
x,
x+npar,
fVec->data);
129 assert(
fTmp !=
nullptr);
130 gsl_vector_set_zero(
fTmp);
131 assert(
fCov !=
nullptr);
132 gsl_matrix_set_zero(
fCov);
133#if GSL_MAJOR_VERSION > 1
134 assert(fJac !=
nullptr);
135 gsl_matrix_set_zero(fJac);
141 if (
fSolver ==
nullptr)
return "undefined";
142 return std::string(gsl_multifit_fdfsolver_name(
fSolver) );
146 if (
fSolver ==
nullptr)
return -1;
147 return gsl_multifit_fdfsolver_iterate(
fSolver);
151 const double *
X()
const {
152 if (
fSolver ==
nullptr)
return nullptr;
153 gsl_vector *
x = gsl_multifit_fdfsolver_position(
fSolver);
159 if (
fSolver ==
nullptr)
return nullptr;
160#if GSL_MAJOR_VERSION > 1
170 if (
fSolver ==
nullptr)
return nullptr;
171 static double kEpsrel = 0.0001;
172#if GSL_MAJOR_VERSION > 1
173 gsl_multifit_fdfsolver_jac (
fSolver, fJac);
174 int ret = gsl_multifit_covar(fJac, kEpsrel,
fCov);
176 int ret = gsl_multifit_covar(
fSolver->J, kEpsrel,
fCov);
184 if (
fSolver ==
nullptr)
return -1;
186 return gsl_multifit_test_gradient(
fVec, absTol);
192 if (
fSolver ==
nullptr)
return -1;
193 return gsl_multifit_test_delta(
fSolver->dx,
fSolver->x, absTol, relTol);
201 if (
g ==
nullptr)
return edm;
203 if (
c ==
nullptr)
return edm;
204 if (
fTmp ==
nullptr)
return edm;
205 int status = gsl_blas_dgemv(CblasNoTrans, 1.0,
fCov,
fVec, 0.,
fTmp);
206 if (status == 0) status |= gsl_blas_ddot(
fVec,
fTmp, &edm);
207 if (status != 0)
return -1;
222#if GSL_MAJOR_VERSION > 1
223 mutable gsl_matrix * fJac;
225 const gsl_multifit_fdfsolver_type *
fType;
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t type
wrapper to a multi-dim function withtout derivatives for multi-dimensional minimization algorithm
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.
gsl_multifit_function_fdf * GetFunc()
GSLMultiFit, internal class for implementing GSL non linear least square GSL fitting.
GSLMultiFit(GSLMultiFit &&)=delete
int TestGradient(double absTol) const
test gradient (ask from solver gradient vector)
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(const GSLMultiFit &)=delete
~GSLMultiFit()
Destructor (no operations)
int TestDelta(double absTol, double relTol) const
test using abs and relative tolerance |dx| < absTol + relTol*|x| for every component
GSLMultiFitFunctionWrapper fFunc
gsl_multifit_fdfsolver * fSolver
const double * Gradient() const
gradient value at the minimum
GSLMultiFit(const gsl_multifit_fdfsolver_type *type=nullptr)
Default constructor No need to specify the type so far since only one solver exists so far.
GSLMultiFit & operator=(const GSLMultiFit &rhs)=delete
const gsl_multifit_fdfsolver_type * fType
const double * CovarMatrix() const
return covariance matrix of the parameters
const double * X() const
parameter values at the minimum
int Set(const std::vector< Func > &funcVec, const double *x)
set the solver parameters
Namespace for new Math classes and functions.
tbb::task_arena is an alias of tbb::interface7::task_arena, which doesn't allow to forward declare tb...