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
|
<!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 the Function to be Minimized</title>
<meta name="description" content="GNU Scientific Library – Reference Manual: Providing the Function to be Minimized">
<meta name="keywords" content="GNU Scientific Library – Reference Manual: Providing the Function to be Minimized">
<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="Nonlinear-Least_002dSquares-Fitting.html#Nonlinear-Least_002dSquares-Fitting" rel="up" title="Nonlinear Least-Squares Fitting">
<link href="Finite-Difference-Jacobian.html#Finite-Difference-Jacobian" rel="next" title="Finite Difference Jacobian">
<link href="Initializing-the-Nonlinear-Least_002dSquares-Solver.html#Initializing-the-Nonlinear-Least_002dSquares-Solver" rel="previous" title="Initializing the Nonlinear Least-Squares Solver">
<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-the-Function-to-be-Minimized"></a>
<div class="header">
<p>
Next: <a href="Finite-Difference-Jacobian.html#Finite-Difference-Jacobian" accesskey="n" rel="next">Finite Difference Jacobian</a>, Previous: <a href="Initializing-the-Nonlinear-Least_002dSquares-Solver.html#Initializing-the-Nonlinear-Least_002dSquares-Solver" accesskey="p" rel="previous">Initializing the Nonlinear Least-Squares Solver</a>, Up: <a href="Nonlinear-Least_002dSquares-Fitting.html#Nonlinear-Least_002dSquares-Fitting" accesskey="u" rel="up">Nonlinear Least-Squares Fitting</a> [<a href="Function-Index.html#Function-Index" title="Index" rel="index">Index</a>]</p>
</div>
<hr>
<a name="Providing-the-Function-to-be-Minimized-1"></a>
<h3 class="section">38.3 Providing the Function to be Minimized</h3>
<p>You must provide <em>n</em> functions of <em>p</em> variables for the
minimization algorithms to operate on. In order to allow for
arbitrary parameters the functions are defined by the following data
types:
</p>
<dl>
<dt><a name="index-gsl_005fmultifit_005ffunction"></a>Data Type: <strong>gsl_multifit_function</strong></dt>
<dd><p>This data type defines a general system of functions with arbitrary parameters.
</p>
<dl compact="compact">
<dt><code>int (* f) (const gsl_vector * <var>x</var>, void * <var>params</var>, gsl_vector * <var>f</var>)</code></dt>
<dd><p>this function should store the vector result
<em>f(x,params)</em> in <var>f</var> for argument <var>x</var> and arbitrary parameters <var>params</var>,
returning an appropriate error code if the function cannot be computed.
</p>
</dd>
<dt><code>size_t n</code></dt>
<dd><p>the number of functions, i.e. the number of components of the
vector <var>f</var>.
</p>
</dd>
<dt><code>size_t p</code></dt>
<dd><p>the number of independent variables, i.e. the number of components of
the vector <var>x</var>.
</p>
</dd>
<dt><code>void * params</code></dt>
<dd><p>a pointer to the arbitrary parameters of the function.
</p></dd>
</dl>
</dd></dl>
<dl>
<dt><a name="index-gsl_005fmultifit_005ffunction_005ffdf"></a>Data Type: <strong>gsl_multifit_function_fdf</strong></dt>
<dd><p>This data type defines a general system of functions with arbitrary parameters and
the corresponding Jacobian matrix of derivatives,
</p>
<dl compact="compact">
<dt><code>int (* f) (const gsl_vector * <var>x</var>, void * <var>params</var>, gsl_vector * <var>f</var>)</code></dt>
<dd><p>this function should store the vector result
<em>f(x,params)</em> in <var>f</var> for argument <var>x</var> and arbitrary parameters <var>params</var>,
returning an appropriate error code if the function cannot be computed.
</p>
</dd>
<dt><code>int (* df) (const gsl_vector * <var>x</var>, void * <var>params</var>, gsl_matrix * <var>J</var>)</code></dt>
<dd><p>this function should store the <var>n</var>-by-<var>p</var> matrix result
<em>J_ij = d f_i(x,params) / d x_j</em> in <var>J</var> for argument <var>x</var>
and arbitrary parameters <var>params</var>, returning an appropriate error code if the
function cannot be computed. If an analytic Jacobian is unavailable, or too expensive
to compute, this function pointer may be set to NULL, in which
case the Jacobian will be internally computed using finite difference approximations
of the function <var>f</var>.
</p>
</dd>
<dt><code>int (* fdf) (const gsl_vector * <var>x</var>, void * <var>params</var>, gsl_vector * <var>f</var>, gsl_matrix * <var>J</var>)</code></dt>
<dd><p>This function should set the values of the <var>f</var> and <var>J</var> as above,
for arguments <var>x</var> and arbitrary parameters <var>params</var>. This function
provides an optimization of the separate functions for <em>f(x)</em> and
<em>J(x)</em>—it is always faster to compute the function and its
derivative at the same time. If an analytic Jacobian is unavailable, or too expensive
to compute, this function pointer may be set to NULL, in which
case the Jacobian will be internally computed using finite difference approximations
of the function <var>f</var>.
</p>
</dd>
<dt><code>size_t n</code></dt>
<dd><p>the number of functions, i.e. the number of components of the
vector <var>f</var>.
</p>
</dd>
<dt><code>size_t p</code></dt>
<dd><p>the number of independent variables, i.e. the number of components of
the vector <var>x</var>.
</p>
</dd>
<dt><code>void * params</code></dt>
<dd><p>a pointer to the arbitrary parameters of the function.
</p></dd>
</dl>
</dd></dl>
<p>Note that when fitting a non-linear model against experimental data,
the data is passed to the functions above using the
<var>params</var> argument and the trial best-fit parameters through the
<var>x</var> argument.
</p>
<hr>
<div class="header">
<p>
Next: <a href="Finite-Difference-Jacobian.html#Finite-Difference-Jacobian" accesskey="n" rel="next">Finite Difference Jacobian</a>, Previous: <a href="Initializing-the-Nonlinear-Least_002dSquares-Solver.html#Initializing-the-Nonlinear-Least_002dSquares-Solver" accesskey="p" rel="previous">Initializing the Nonlinear Least-Squares Solver</a>, Up: <a href="Nonlinear-Least_002dSquares-Fitting.html#Nonlinear-Least_002dSquares-Fitting" accesskey="u" rel="up">Nonlinear Least-Squares Fitting</a> [<a href="Function-Index.html#Function-Index" title="Index" rel="index">Index</a>]</p>
</div>
</body>
</html>
|