Next: , Previous: Nonlinear Least-Squares Function Definition, Up: Nonlinear Least-Squares Fitting   [Index]


39.7 Iteration

The following functions drive the iteration of each algorithm. Each function performs one iteration of the trust region method and updates the state of the solver.

Function: int gsl_multifit_nlinear_iterate (gsl_multifit_nlinear_workspace * w)
Function: int gsl_multilarge_nlinear_iterate (gsl_multilarge_nlinear_workspace * w)

These functions perform a single iteration of the solver w. If the iteration encounters an unexpected problem then an error code will be returned. The solver workspace maintains a current estimate of the best-fit parameters at all times.

The solver workspace w contains the following entries, which can be used to track the progress of the solution:

gsl_vector * x

The current position, length p.

gsl_vector * f

The function residual vector at the current position f(x), length n.

gsl_matrix * J

The Jacobian matrix at the current position J(x), size n-by-p (only for gsl_multifit_nlinear interface).

gsl_vector * dx

The difference between the current position and the previous position, i.e. the last step \delta, taken as a vector, length p.

These quantities can be accessed with the following functions,

Function: gsl_vector * gsl_multifit_nlinear_position (const gsl_multifit_nlinear_workspace * w)
Function: gsl_vector * gsl_multilarge_nlinear_position (const gsl_multilarge_nlinear_workspace * w)

These functions return the current position x (i.e. best-fit parameters) of the solver w.

Function: gsl_vector * gsl_multifit_nlinear_residual (const gsl_multifit_nlinear_workspace * w)
Function: gsl_vector * gsl_multilarge_nlinear_residual (const gsl_multilarge_nlinear_workspace * w)

These functions return the current residual vector f(x) of the solver w. For weighted systems, the residual vector includes the weighting factor \sqrt{W}.

Function: gsl_matrix * gsl_multifit_nlinear_jac (const gsl_multifit_nlinear_workspace * w)

This function returns a pointer to the n-by-p Jacobian matrix for the current iteration of the solver w. This function is available only for the gsl_multifit_nlinear interface.

Function: size_t gsl_multifit_nlinear_niter (const gsl_multifit_nlinear_workspace * w)
Function: size_t gsl_multilarge_nlinear_niter (const gsl_multilarge_nlinear_workspace * w)

These functions return the number of iterations performed so far. The iteration counter is updated on each call to the _iterate functions above, and reset to 0 in the _init functions.

Function: int gsl_multifit_nlinear_rcond (double * rcond, const gsl_multifit_nlinear_workspace * w)
Function: int gsl_multilarge_nlinear_rcond (double * rcond, const gsl_multilarge_nlinear_workspace * w)

This function estimates the reciprocal condition number of the Jacobian matrix at the current position x and stores it in rcond. The computed value is only an estimate to give the user a guideline as to the conditioning of their particular problem. Its calculation is based on which factorization method is used (Cholesky, QR, or SVD).


Next: , Previous: Nonlinear Least-Squares Function Definition, Up: Nonlinear Least-Squares Fitting   [Index]