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"
60 GSLMultiFit (
const gsl_multifit_fdfsolver_type * type = 0) :
65 #if GSL_MAJOR_VERSION > 1
70 if (fType == 0) fType = gsl_multifit_fdfsolver_lmsder;
77 if (fSolver) gsl_multifit_fdfsolver_free(fSolver);
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);
92 GSLMultiFit(
const GSLMultiFit &) {}
97 GSLMultiFit & operator = (
const GSLMultiFit & rhs) {
98 if (
this == &rhs)
return *
this;
106 void CreateSolver(
unsigned int npoints,
unsigned int npar) {
107 if (fSolver) gsl_multifit_fdfsolver_free(fSolver);
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();
134 fFunc.SetFunction(funcVec, npts, npar);
136 CreateSolver( npts, npar );
137 assert(fSolver != 0);
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);
149 return gsl_multifit_fdfsolver_set(fSolver, fFunc.GetFunc(), fVec);
152 std::string Name()
const {
153 if (fSolver == 0)
return "undefined";
154 return std::string(gsl_multifit_fdfsolver_name(fSolver) );
158 if (fSolver == 0)
return -1;
159 return gsl_multifit_fdfsolver_iterate(fSolver);
163 const double * X()
const {
164 if (fSolver == 0)
return 0;
165 gsl_vector * x = gsl_multifit_fdfsolver_position(fSolver);
170 const double * Gradient()
const {
171 if (fSolver == 0)
return 0;
172 #if GSL_MAJOR_VERSION > 1
173 fType->gradient(fSolver->state, fVec);
175 gsl_multifit_gradient(fSolver->J, fSolver->f,fVec);
181 const double * CovarMatrix()
const {
182 if (fSolver == 0)
return 0;
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);
190 if (ret != GSL_SUCCESS)
return 0;
195 int TestGradient(
double absTol)
const {
196 if (fSolver == 0)
return -1;
198 return gsl_multifit_test_gradient( fVec, absTol);
203 int TestDelta(
double absTol,
double relTol)
const {
204 if (fSolver == 0)
return -1;
205 return gsl_multifit_test_delta(fSolver->dx, fSolver->x, absTol, relTol);
212 const double * g = Gradient();
213 if (g == 0)
return edm;
214 const double * c = CovarMatrix();
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;
228 GSLMultiFitFunctionWrapper fFunc;
229 gsl_multifit_fdfsolver * fSolver;
231 mutable gsl_vector * fVec;
232 mutable gsl_vector * fTmp;
233 mutable gsl_matrix * fCov;
234 #if GSL_MAJOR_VERSION > 1
235 mutable gsl_matrix * fJac;
237 const gsl_multifit_fdfsolver_type * fType;