manpagez: man pages & more
info gsl-ref
Home | html | info | man
[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

35.4 Providing a function to minimize

You must provide a parametric function of n variables for the minimizers to operate on. You may also need to provide a routine which calculates the gradient of the function and a third routine which calculates both the function value and the gradient together. In order to allow for general parameters the functions are defined by the following data types:

Data Type: gsl_multimin_function_fdf

This data type defines a general function of n variables with parameters and the corresponding gradient vector of derivatives,

double (* f) (const gsl_vector * x, void * params)

this function should return the result f(x,params) for argument x and parameters params. If the function cannot be computed, an error value of GSL_NAN should be returned.

void (* df) (const gsl_vector * x, void * params, gsl_vector * g)

this function should store the n-dimensional gradient g_i = d f(x,params) / d x_i in the vector g for argument x and parameters params, returning an appropriate error code if the function cannot be computed.

void (* fdf) (const gsl_vector * x, void * params, double * f, gsl_vector * g)

This function should set the values of the f and g as above, for arguments x and parameters params. This function provides an optimization of the separate functions for f(x) and g(x)—it is always faster to compute the function and its derivative at the same time.

size_t n

the dimension of the system, i.e. the number of components of the vectors x.

void * params

a pointer to the parameters of the function.

Data Type: gsl_multimin_function

This data type defines a general function of n variables with parameters,

double (* f) (const gsl_vector * x, void * params)

this function should return the result f(x,params) for argument x and parameters params. If the function cannot be computed, an error value of GSL_NAN should be returned.

size_t n

the dimension of the system, i.e. the number of components of the vectors x.

void * params

a pointer to the parameters of the function.

The following example function defines a simple two-dimensional paraboloid with five parameters,

 
/* Paraboloid centered on (p[0],p[1]), with  
   scale factors (p[2],p[3]) and minimum p[4] */

double
my_f (const gsl_vector *v, void *params)
{
  double x, y;
  double *p = (double *)params;
  
  x = gsl_vector_get(v, 0);
  y = gsl_vector_get(v, 1);
 
  return p[2] * (x - p[0]) * (x - p[0]) +
           p[3] * (y - p[1]) * (y - p[1]) + p[4]; 
}

/* The gradient of f, df = (df/dx, df/dy). */
void 
my_df (const gsl_vector *v, void *params, 
       gsl_vector *df)
{
  double x, y;
  double *p = (double *)params;
  
  x = gsl_vector_get(v, 0);
  y = gsl_vector_get(v, 1);
 
  gsl_vector_set(df, 0, 2.0 * p[2] * (x - p[0]));
  gsl_vector_set(df, 1, 2.0 * p[3] * (y - p[1]));
}

/* Compute both f and df together. */
void 
my_fdf (const gsl_vector *x, void *params, 
        double *f, gsl_vector *df) 
{
  *f = my_f(x, params); 
  my_df(x, params, df);
}

The function can be initialized using the following code,

 
gsl_multimin_function_fdf my_func;

/* Paraboloid center at (1,2), scale factors (10, 20), 
   minimum value 30 */
double p[5] = { 1.0, 2.0, 10.0, 20.0, 30.0 }; 

my_func.n = 2;  /* number of function components */
my_func.f = &my_f;
my_func.df = &my_df;
my_func.fdf = &my_fdf;
my_func.params = (void *)p;

[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]
© manpagez.com 2000-2024
Individual documents may contain additional copyright information.