![]() |
Object-oriented Scientific Computing Library: Version 0.910
|
Non-linear least-squares fitting class (GSL) More...
#include <gsl_fit.h>
The GSL-based fitting class using a Levenberg-Marquardt type algorithm. The algorithm stops when
where is the last step and
is the current position. If test_gradient is true, then additionally fit() requires that
where is the
-th component of the gradient of the function
where
Default template arguments
func_t
- fit_funct<>vec_t
- ovector_basealloc_vec_t
- ovectorbool_vec_t
- bool *
Properly generalize other vector types than ovector_base
Allow the user to specify the derivatives
Fix so that the user can specify automatic scaling of the fitting parameters, where the initial guess are used for scaling so that the fitting parameters are near unity.
Data Structures | |
struct | func_par |
A structure for passing to the functions func(), dfunc(), and fdfunc() More... | |
Public Member Functions | |
virtual int | fit (size_t ndat, vec_t &xdat, vec_t &ydat, vec_t &yerr, size_t npar, vec_t &par, mat_t &covar, double &chi2, func_t &fitfun) |
Fit the data specified in (xdat,ydat) to the function fitfun with the parameters in par . | |
virtual const char * | type () |
Return string denoting type ("gsl_fit") | |
Data Fields | |
int | max_iter |
(default 500) | |
double | epsabs |
(default 1.0e-4) | |
double | epsrel |
(default 1.0e-4) | |
bool | test_gradient |
If true, test the gradient also (default false) | |
bool | use_scaled |
Use the scaled routine if true (default true) | |
Protected Member Functions | |
virtual int | print_iter (size_t nv, gsl_vector *x, int iter, double val, double lim) |
Print the progress in the current iteration. | |
Static Protected Member Functions | |
static int | func (const gsl_vector *x, void *pa, gsl_vector *f) |
Evaluate the function. | |
static int | dfunc (const gsl_vector *x, void *pa, gsl_matrix *jac) |
Evaluate the jacobian. | |
static int | fdfunc (const gsl_vector *x, void *pa, gsl_vector *f, gsl_matrix *jac) |
Evaluate the function and the jacobian. |
virtual int gsl_fit< func_t, vec_t, mat_t, bool_vec_t >::fit | ( | size_t | ndat, |
vec_t & | xdat, | ||
vec_t & | ydat, | ||
vec_t & | yerr, | ||
size_t | npar, | ||
vec_t & | par, | ||
mat_t & | covar, | ||
double & | chi2, | ||
func_t & | fitfun | ||
) | [inline, virtual] |
The covariance matrix for the parameters is returned in covar
and the value of is returned in
chi2
.
Implements fit_base< func_t, vec_t, mat_t >.
Documentation generated with Doxygen. Provided under the GNU Free Documentation License (see License Information).