Next: Nonlinear Least-Squares Function Definition, Previous: Nonlinear Least-Squares Tunable Parameters, Up: Nonlinear Least-Squares Fitting [Index]
These functions return a pointer to a newly allocated instance of a
derivative solver of type T for n observations and p
parameters. The params input specifies a tunable set of
parameters which will affect important details in each iteration
of the trust region subproblem algorithm. It is recommended to start
with the suggested default parameters (see
gsl_multifit_nlinear_default_parameters
and
gsl_multilarge_nlinear_default_parameters
) and then tune
the parameters once the code is working correctly. See
Nonlinear Least-Squares Tunable Parameters
for descriptions of the various parameters.
For example, the following code creates an instance of a
Levenberg-Marquardt solver for 100 data points and 3 parameters,
using suggested defaults:
const gsl_multifit_nlinear_type * T = gsl_multifit_nlinear_lm; gsl_multifit_nlinear_parameters params = gsl_multifit_nlinear_default_parameters(); gsl_multifit_nlinear_workspace * w = gsl_multifit_nlinear_alloc (T, ¶ms, 100, 3);
The number of observations n must be greater than or equal to parameters p.
If there is insufficient memory to create the solver then the function
returns a null pointer and the error handler is invoked with an error
code of GSL_ENOMEM
.
These functions return a set of recommended default parameters for use in solving nonlinear least squares problems. The user can tune each parameter to improve the performance on their particular problem, see Nonlinear Least-Squares Tunable Parameters.
These functions initialize, or reinitialize, an existing workspace w to use the system fdf and the initial guess x. See Nonlinear Least-Squares Function Definition for a description of the fdf structure.
Optionally, a weight vector wts can be given to perform a weighted nonlinear regression. Here, the weighting matrix is W = diag(w_1,w_2,...,w_n).
These functions free all the memory associated with the workspace w.
These functions return a pointer to the name of the solver. For example,
printf ("w is a '%s' solver\n", gsl_multifit_nlinear_name (w));
would print something like w is a 'trust-region' solver
.
These functions return a pointer to the name of the trust region subproblem method. For example,
printf ("w is a '%s' solver\n", gsl_multifit_nlinear_trs_name (w));
would print something like w is a 'levenberg-marquardt' solver
.