1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206
|
<!DOCTYPE html PUBLIC "-//W3C//DTD HTML 4.01 Transitional//EN" "http://www.w3.org/TR/html4/loose.dtd">
<html>
<!-- Copyright (C) 1996, 1997, 1998, 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008, 2009, 2010, 2011, 2012, 2013 The GSL Team.
Permission is granted to copy, distribute and/or modify this document
under the terms of the GNU Free Documentation License, Version 1.3 or
any later version published by the Free Software Foundation; with no
Invariant Sections and no cover texts. A copy of the license is
included in the section entitled "GNU Free Documentation License". -->
<!-- Created by GNU Texinfo 5.1, http://www.gnu.org/software/texinfo/ -->
<head>
<title>GNU Scientific Library – Reference Manual: Providing a function to minimize</title>
<meta name="description" content="GNU Scientific Library – Reference Manual: Providing a function to minimize">
<meta name="keywords" content="GNU Scientific Library – Reference Manual: Providing a function to minimize">
<meta name="resource-type" content="document">
<meta name="distribution" content="global">
<meta name="Generator" content="makeinfo">
<meta http-equiv="Content-Type" content="text/html; charset=utf-8">
<link href="index.html#Top" rel="start" title="Top">
<link href="Function-Index.html#Function-Index" rel="index" title="Function Index">
<link href="Multidimensional-Minimization.html#Multidimensional-Minimization" rel="up" title="Multidimensional Minimization">
<link href="Multimin-Iteration.html#Multimin-Iteration" rel="next" title="Multimin Iteration">
<link href="Initializing-the-Multidimensional-Minimizer.html#Initializing-the-Multidimensional-Minimizer" rel="previous" title="Initializing the Multidimensional Minimizer">
<style type="text/css">
<!--
a.summary-letter {text-decoration: none}
blockquote.smallquotation {font-size: smaller}
div.display {margin-left: 3.2em}
div.example {margin-left: 3.2em}
div.indentedblock {margin-left: 3.2em}
div.lisp {margin-left: 3.2em}
div.smalldisplay {margin-left: 3.2em}
div.smallexample {margin-left: 3.2em}
div.smallindentedblock {margin-left: 3.2em; font-size: smaller}
div.smalllisp {margin-left: 3.2em}
kbd {font-style:oblique}
pre.display {font-family: inherit}
pre.format {font-family: inherit}
pre.menu-comment {font-family: serif}
pre.menu-preformatted {font-family: serif}
pre.smalldisplay {font-family: inherit; font-size: smaller}
pre.smallexample {font-size: smaller}
pre.smallformat {font-family: inherit; font-size: smaller}
pre.smalllisp {font-size: smaller}
span.nocodebreak {white-space:nowrap}
span.nolinebreak {white-space:nowrap}
span.roman {font-family:serif; font-weight:normal}
span.sansserif {font-family:sans-serif; font-weight:normal}
ul.no-bullet {list-style: none}
-->
</style>
</head>
<body lang="en" bgcolor="#FFFFFF" text="#000000" link="#0000FF" vlink="#800080" alink="#FF0000">
<a name="Providing-a-function-to-minimize"></a>
<div class="header">
<p>
Next: <a href="Multimin-Iteration.html#Multimin-Iteration" accesskey="n" rel="next">Multimin Iteration</a>, Previous: <a href="Initializing-the-Multidimensional-Minimizer.html#Initializing-the-Multidimensional-Minimizer" accesskey="p" rel="previous">Initializing the Multidimensional Minimizer</a>, Up: <a href="Multidimensional-Minimization.html#Multidimensional-Minimization" accesskey="u" rel="up">Multidimensional Minimization</a> [<a href="Function-Index.html#Function-Index" title="Index" rel="index">Index</a>]</p>
</div>
<hr>
<a name="Providing-a-function-to-minimize-1"></a>
<h3 class="section">36.4 Providing a function to minimize</h3>
<p>You must provide a parametric function of <em>n</em> 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:
</p>
<dl>
<dt><a name="index-gsl_005fmultimin_005ffunction_005ffdf"></a>Data Type: <strong>gsl_multimin_function_fdf</strong></dt>
<dd><p>This data type defines a general function of <em>n</em> variables with
parameters and the corresponding gradient vector of derivatives,
</p>
<dl compact="compact">
<dt><code>double (* f) (const gsl_vector * <var>x</var>, void * <var>params</var>)</code></dt>
<dd><p>this function should return the result
<em>f(x,params)</em> for argument <var>x</var> and parameters <var>params</var>.
If the function cannot be computed, an error value of <code>GSL_NAN</code>
should be returned.
</p>
</dd>
<dt><code>void (* df) (const gsl_vector * <var>x</var>, void * <var>params</var>, gsl_vector * <var>g</var>)</code></dt>
<dd><p>this function should store the <var>n</var>-dimensional gradient
<em>g_i = d f(x,params) / d x_i</em> in the vector <var>g</var> for argument <var>x</var>
and parameters <var>params</var>, returning an appropriate error code if the
function cannot be computed.
</p>
</dd>
<dt><code>void (* fdf) (const gsl_vector * <var>x</var>, void * <var>params</var>, double * f, gsl_vector * <var>g</var>)</code></dt>
<dd><p>This function should set the values of the <var>f</var> and <var>g</var> as above,
for arguments <var>x</var> and parameters <var>params</var>. This function
provides an optimization of the separate functions for <em>f(x)</em> and
<em>g(x)</em>—it is always faster to compute the function and its
derivative at the same time.
</p>
</dd>
<dt><code>size_t n</code></dt>
<dd><p>the dimension of the system, i.e. the number of components of the
vectors <var>x</var>.
</p>
</dd>
<dt><code>void * params</code></dt>
<dd><p>a pointer to the parameters of the function.
</p></dd>
</dl>
</dd></dl>
<dl>
<dt><a name="index-gsl_005fmultimin_005ffunction"></a>Data Type: <strong>gsl_multimin_function</strong></dt>
<dd><p>This data type defines a general function of <em>n</em> variables with
parameters,
</p>
<dl compact="compact">
<dt><code>double (* f) (const gsl_vector * <var>x</var>, void * <var>params</var>)</code></dt>
<dd><p>this function should return the result
<em>f(x,params)</em> for argument <var>x</var> and parameters <var>params</var>.
If the function cannot be computed, an error value of <code>GSL_NAN</code>
should be returned.
</p>
</dd>
<dt><code>size_t n</code></dt>
<dd><p>the dimension of the system, i.e. the number of components of the
vectors <var>x</var>.
</p>
</dd>
<dt><code>void * params</code></dt>
<dd><p>a pointer to the parameters of the function.
</p></dd>
</dl>
</dd></dl>
<p>The following example function defines a simple two-dimensional
paraboloid with five parameters,
</p>
<div class="example">
<pre class="verbatim">/* 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);
}
</pre></div>
<p>The function can be initialized using the following code,
</p>
<div class="example">
<pre class="example">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;
</pre></div>
<hr>
<div class="header">
<p>
Next: <a href="Multimin-Iteration.html#Multimin-Iteration" accesskey="n" rel="next">Multimin Iteration</a>, Previous: <a href="Initializing-the-Multidimensional-Minimizer.html#Initializing-the-Multidimensional-Minimizer" accesskey="p" rel="previous">Initializing the Multidimensional Minimizer</a>, Up: <a href="Multidimensional-Minimization.html#Multidimensional-Minimization" accesskey="u" rel="up">Multidimensional Minimization</a> [<a href="Function-Index.html#Function-Index" title="Index" rel="index">Index</a>]</p>
</div>
</body>
</html>
|