![]() |
Object-oriented Scientific Computing Library: Version 0.910
|
Non-linear least-squares fitting class with generic minimizer. More...
#include <min_fit.h>
This minimizes a generic fitting function using any multi_min object, and then uses the GSL routines to calculate the uncertainties in the parameters and the covariance matrix.
This can be useful for fitting problems which might be better handled by more complex minimizers than those that are used in gsl_fit. For problems with many local minima near the global minimum, using a sim_anneal object with this class can sometimes produce better results than gsl_fit.
Default template arguments
func_t
- fit_funct<>vec_t
- ovector_basealloc_vec_t
- ovector Data Structures | |
struct | func_par |
A structure for passing information to the GSL functions for the min_fit class. 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 . | |
int | set_multi_min (multi_min< multi_funct< vec_t > > &mm) |
Set the multi_min object to use (default is gsl_mmin_nmsimplex) | |
virtual const char * | type () |
Return string denoting type ("min_fit") | |
Data Fields | |
gsl_mmin_simp< multi_funct < vec_t > > | def_multi_min |
The default minimizer. | |
Protected Member Functions | |
double | min_func (size_t np, const vec_t &xp, func_par *&fp) |
The function to minimize. | |
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. | |
Protected Attributes | |
multi_min< multi_funct< vec_t > > * | mmp |
The minimizer. | |
bool | min_set |
True if the minimizer has been set by the user. |
virtual int min_fit< func_t, vec_t, mat_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).