Next: Providing the Function to be Minimized, Previous: Overview of Regularized Nonlinear Least-Squares Fitting, Up: Nonlinear Least-Squares Fitting [Index]
This function returns a pointer to a newly allocated instance of a solver of type T for n observations and p parameters. 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
.
This function returns a pointer to a newly allocated instance of a derivative solver of type T for n observations and p parameters. For example, the following code creates an instance of a Levenberg-Marquardt solver for 100 data points and 3 parameters,
const gsl_multifit_fdfsolver_type * T = gsl_multifit_fdfsolver_lmder; gsl_multifit_fdfsolver * s = gsl_multifit_fdfsolver_alloc (T, 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
.
This function returns a pointer to a newly allocated instance of a
derivative solver of type T for n observations and p
parameters. The solver will automatically form the augmented
system \tilde{f}(x) and \tilde{J} for ridge (Tikhonov)
regression.
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
.
This function initializes, or reinitializes, an existing solver s to use the function f and the initial guess x.
These functions initialize, or reinitialize, an existing solver s to use the function and derivative fdf and the initial guess x.
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). The wts vector is referenced throughout the iteration so it should not be freed by the caller until the iteration terminates.
This function initializes, or reinitializes, an existing ridge solver s to use the function and derivative fdf and the initial guess x. Here, the regularization matrix is set to L = \lambda I, with \lambda specified in lambda.
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). The wts vector is referenced throughout the iteration so it should not be freed by the caller until the iteration terminates.
This function initializes, or reinitializes, an existing ridge solver s to use the function and derivative fdf and the initial guess x. Here, the regularization matrix is set to L = diag(\lambda_1,\lambda_2,...,\lambda_p), where the \lambda_i are given in lambda.
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). The wts vector is referenced throughout the iteration so it should not be freed by the caller until the iteration terminates.
This function initializes, or reinitializes, an existing ridge solver s to use the function and derivative fdf and the initial guess x. Here, the regularization matrix is set to L, which must have p columns but may have any number of rows.
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). The wts vector is referenced throughout the iteration so it should not be freed by the caller until the iteration terminates.
These functions free all the memory associated with the solver s.
These functions return a pointer to the name of the solver. For example,
printf ("s is a '%s' solver\n", gsl_multifit_fdfsolver_name (s));
would print something like s is a 'lmder' solver
.