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

39.7 Iteration

The following functions drive the iteration of each algorithm. Each function performs one iteration to update the state of any solver of the corresponding type. The same functions work for all solvers so that different methods can be substituted at runtime without modifications to the code.

Function: int gsl_multifit_fsolver_iterate (gsl_multifit_fsolver * s)
Function: int gsl_multifit_fdfsolver_iterate (gsl_multifit_fdfsolver * s)
Function: int gsl_multifit_fdfridge_iterate (gsl_multifit_fdfridge * s)

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

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

gsl_vector * x

The current position.

gsl_vector * f

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

gsl_vector * dx

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

The best-fit information also can be accessed with the following auxiliary functions,

Function: gsl_vector * gsl_multifit_fsolver_position (const gsl_multifit_fsolver * s)
Function: gsl_vector * gsl_multifit_fdfsolver_position (const gsl_multifit_fdfsolver * s)
Function: gsl_vector * gsl_multifit_fdfridge_position (const gsl_multifit_fdfridge * s)

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

Function: gsl_vector * gsl_multifit_fdfsolver_residual (const gsl_multifit_fdfsolver * s)
Function: gsl_vector * gsl_multifit_fdfridge_residual (const gsl_multifit_fdfridge * s)

These functions return the current residual vector f of the solver s. For weighted cases, the residual vector includes the weighting factor \sqrt{W}. For ridge regression, the residual vector is the augmented vector \tilde{f}.

Function: size_t gsl_multifit_fdfsolver_niter (const gsl_multifit_fdfsolver * s)
Function: size_t gsl_multifit_fdfridge_niter (const gsl_multifit_fdfridge * s)

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 _set functions.

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