Node:Initializing the Nonlinear Least-Squares Solver, Next:, Previous:Overview of Nonlinear Least-Squares Fitting, Up:Nonlinear Least-Squares Fitting



Initializing the Solver

gsl_multifit_fsolver * gsl_multifit_fsolver_alloc (const gsl_multifit_fsolver_type * T, size_t n, size_t p) Function
This function returns a pointer to a newly allocated instance of a solver of type T for n observations and p parameters.

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.

gsl_multifit_fdfsolver * gsl_multifit_fdfsolver_alloc (const gsl_multifit_fdfsolver_type * T, size_t n, size_t p) Function
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);
          

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.

int gsl_multifit_fsolver_set (gsl_multifit_fsolver * s, gsl_multifit_function * f, gsl_vector * x) Function
This function initializes, or reinitializes, an existing solver s to use the function f and the initial guess x.

int gsl_multifit_fdfsolver_set (gsl_multifit_fdfsolver * s, gsl_function_fdf * fdf, gsl_vector * x) Function
This function initializes, or reinitializes, an existing solver s to use the function and derivative fdf and the initial guess x.

void gsl_multifit_fsolver_free (gsl_multifit_fsolver * s) Function
void gsl_multifit_fdfsolver_free (gsl_multifit_fdfsolver * s) Function
These functions free all the memory associated with the solver s.

const char * gsl_multifit_fsolver_name (const gsl_multifit_fdfsolver * s) Function
const char * gsl_multifit_fdfsolver_name (const gsl_multifit_fdfsolver * s) Function
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.