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
|
<html lang="en">
<head>
<title>Providing the Function to be Minimized - GNU Scientific Library -- Reference Manual</title>
<meta http-equiv="Content-Type" content="text/html">
<meta name="description" content="GNU Scientific Library -- Reference Manual">
<meta name="generator" content="makeinfo 4.8">
<link title="Top" rel="start" href="index.html#Top">
<link rel="up" href="Nonlinear-Least_002dSquares-Fitting.html" title="Nonlinear Least-Squares Fitting">
<link rel="prev" href="Initializing-the-Nonlinear-Least_002dSquares-Solver.html" title="Initializing the Nonlinear Least-Squares Solver">
<link rel="next" href="Iteration-of-the-Minimization-Algorithm.html" title="Iteration of the Minimization Algorithm">
<link href="http://www.gnu.org/software/texinfo/" rel="generator-home" title="Texinfo Homepage">
<!--
Copyright (C) 1996, 1997, 1998, 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007 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.2 or
any later version published by the Free Software Foundation; with the
Invariant Sections being ``GNU General Public License'' and ``Free Software
Needs Free Documentation'', the Front-Cover text being ``A GNU Manual'',
and with the Back-Cover Text being (a) (see below). A copy of the
license is included in the section entitled ``GNU Free Documentation
License''.
(a) The Back-Cover Text is: ``You have freedom to copy and modify this
GNU Manual, like GNU software.''-->
<meta http-equiv="Content-Style-Type" content="text/css">
<style type="text/css"><!--
pre.display { font-family:inherit }
pre.format { font-family:inherit }
pre.smalldisplay { font-family:inherit; font-size:smaller }
pre.smallformat { font-family:inherit; font-size:smaller }
pre.smallexample { font-size:smaller }
pre.smalllisp { font-size:smaller }
span.sc { font-variant:small-caps }
span.roman { font-family:serif; font-weight:normal; }
span.sansserif { font-family:sans-serif; font-weight:normal; }
--></style>
</head>
<body>
<div class="node">
<p>
<a name="Providing-the-Function-to-be-Minimized"></a>
Next: <a rel="next" accesskey="n" href="Iteration-of-the-Minimization-Algorithm.html">Iteration of the Minimization Algorithm</a>,
Previous: <a rel="previous" accesskey="p" href="Initializing-the-Nonlinear-Least_002dSquares-Solver.html">Initializing the Nonlinear Least-Squares Solver</a>,
Up: <a rel="up" accesskey="u" href="Nonlinear-Least_002dSquares-Fitting.html">Nonlinear Least-Squares Fitting</a>
<hr>
</div>
<h3 class="section">37.3 Providing the Function to be Minimized</h3>
<p>You must provide n functions of p variables for the
minimization algorithms to operate on. In order to allow for
arbitrary parameters the functions are defined by the following data
types:
<div class="defun">
— Data Type: <b>gsl_multifit_function</b><var><a name="index-gsl_005fmultifit_005ffunction-2403"></a></var><br>
<blockquote><p>This data type defines a general system of functions with arbitrary parameters.
<dl>
<dt><code>int (* f) (const gsl_vector * </code><var>x</var><code>, void * </code><var>params</var><code>, gsl_vector * </code><var>f</var><code>)</code><dd>this function should store the vector result
<!-- {$f(x,\hbox{\it params})$} -->
f(x,params) 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.
<br><dt><code>size_t n</code><dd>the number of functions, i.e. the number of components of the
vector <var>f</var>.
<br><dt><code>size_t p</code><dd>the number of independent variables, i.e. the number of components of
the vector <var>x</var>.
<br><dt><code>void * params</code><dd>a pointer to the arbitrary parameters of the function.
</dl>
</p></blockquote></div>
<div class="defun">
— Data Type: <b>gsl_multifit_function_fdf</b><var><a name="index-gsl_005fmultifit_005ffunction_005ffdf-2404"></a></var><br>
<blockquote><p>This data type defines a general system of functions with arbitrary parameters and
the corresponding Jacobian matrix of derivatives,
<dl>
<dt><code>int (* f) (const gsl_vector * </code><var>x</var><code>, void * </code><var>params</var><code>, gsl_vector * </code><var>f</var><code>)</code><dd>this function should store the vector result
<!-- {$f(x,\hbox{\it params})$} -->
f(x,params) 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.
<br><dt><code>int (* df) (const gsl_vector * </code><var>x</var><code>, void * </code><var>params</var><code>, gsl_matrix * </code><var>J</var><code>)</code><dd>this function should store the <var>n</var>-by-<var>p</var> matrix result
<!-- {$J_{ij} = \partial f_i(x,\hbox{\it params}) / \partial x_j$} -->
J_ij = d f_i(x,params) / d x_j 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.
<br><dt><code>int (* fdf) (const gsl_vector * </code><var>x</var><code>, void * </code><var>params</var><code>, gsl_vector * </code><var>f</var><code>, gsl_matrix * </code><var>J</var><code>)</code><dd>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 f(x) and
J(x)—it is always faster to compute the function and its
derivative at the same time.
<br><dt><code>size_t n</code><dd>the number of functions, i.e. the number of components of the
vector <var>f</var>.
<br><dt><code>size_t p</code><dd>the number of independent variables, i.e. the number of components of
the vector <var>x</var>.
<br><dt><code>void * params</code><dd>a pointer to the arbitrary parameters of the function.
</dl>
</p></blockquote></div>
<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.
<hr>The GNU Scientific Library - a free numerical library licensed under the GNU GPL<br>Back to the <a href="/software/gsl/">GNU Scientific Library Homepage</a></body></html>
|