File: Nonlinear-Least_002dSquares-Function-Definition.html

package info (click to toggle)
gsl-ref-html 2.3-1
  • links: PTS
  • area: non-free
  • in suites: bullseye, buster, sid
  • size: 6,876 kB
  • ctags: 4,574
  • sloc: makefile: 35
file content (237 lines) | stat: -rw-r--r-- 12,501 bytes parent folder | download
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
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
<!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, 2014, 2015, 2016 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 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 the freedom to copy and modify this
GNU Manual." -->
<!-- Created by GNU Texinfo 5.1, http://www.gnu.org/software/texinfo/ -->
<head>
<title>GNU Scientific Library &ndash; Reference Manual: Nonlinear Least-Squares Function Definition</title>

<meta name="description" content="GNU Scientific Library &ndash; Reference Manual: Nonlinear Least-Squares Function Definition">
<meta name="keywords" content="GNU Scientific Library &ndash; Reference Manual: Nonlinear Least-Squares Function Definition">
<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="Nonlinear-Least_002dSquares-Iteration.html#Nonlinear-Least_002dSquares-Iteration" rel="next" title="Nonlinear Least-Squares Iteration">
<link href="Nonlinear-Least_002dSquares-Initialization.html#Nonlinear-Least_002dSquares-Initialization" rel="previous" title="Nonlinear Least-Squares Initialization">
<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="Nonlinear-Least_002dSquares-Function-Definition"></a>
<div class="header">
<p>
Next: <a href="Nonlinear-Least_002dSquares-Iteration.html#Nonlinear-Least_002dSquares-Iteration" accesskey="n" rel="next">Nonlinear Least-Squares Iteration</a>, Previous: <a href="Nonlinear-Least_002dSquares-Initialization.html#Nonlinear-Least_002dSquares-Initialization" accesskey="p" rel="previous">Nonlinear Least-Squares Initialization</a>, Up: <a href="Nonlinear-Least_002dSquares-Fitting.html#Nonlinear-Least_002dSquares-Fitting" accesskey="u" rel="up">Nonlinear Least-Squares Fitting</a> &nbsp; [<a href="Function-Index.html#Function-Index" title="Index" rel="index">Index</a>]</p>
</div>
<hr>
<a name="Providing-the-Function-to-be-Minimized"></a>
<h3 class="section">39.6 Providing the Function to be Minimized</h3>

<p>The user must provide <em>n</em> functions of <em>p</em> variables for the
minimization algorithm 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_005fnlinear_005ffdf"></a>Data Type: <strong>gsl_multifit_nlinear_fdf</strong></dt>
<dd><p>This data type defines a general system of functions with arbitrary parameters,
the corresponding Jacobian matrix of derivatives, and optionally the
second directional derivative of the functions for geodesic acceleration.
</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 <em>n</em> components of the vector
<em>f(x)</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) / 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
matrix 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 (* fvv) (const gsl_vector * <var>x</var>, const gsl_vector * <var>v</var>, void * <var>params</var>, gsl_vector * <var>fvv</var>)</code></dt>
<dd><p>When geodesic acceleration is enabled, this function should store the
<em>n</em> components of the vector
<em>f_{vv}(x) = \sum_{\alpha\beta} v_{\alpha} v_{\beta} {\partial \over \partial x_{\alpha}} {\partial \over \partial x_{\beta}} f(x)</em>,
representing second directional derivatives of the function to be minimized,
into the output <var>fvv</var>. The parameter vector is provided in <var>x</var> and
the velocity vector is provided in <var>v</var>, both of which have <em>p</em>
components. The arbitrary parameters are given in <var>params</var>. If
analytic expressions for <em>f_{vv}(x)</em> are unavailable or too difficult
to compute, this function pointer may be set to NULL, in which case
<em>f_{vv}(x)</em> will be computed internally using a finite difference
approximation.
</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>
<dt><code>size_t nevalf</code></dt>
<dd><p>This does not need to be set by the user. It counts the number of
function evaluations and is initialized by the <code>_init</code> function.
</p>
</dd>
<dt><code>size_t nevaldf</code></dt>
<dd><p>This does not need to be set by the user. It counts the number of
Jacobian evaluations and is initialized by the <code>_init</code> function.
</p>
</dd>
<dt><code>size_t nevalfvv</code></dt>
<dd><p>This does not need to be set by the user. It counts the number of
<em>f_{vv}(x)</em> evaluations and is initialized by the <code>_init</code> function.
</p></dd>
</dl>
</dd></dl>

<dl>
<dt><a name="index-gsl_005fmultilarge_005fnlinear_005ffdf"></a>Data Type: <strong>gsl_multilarge_nlinear_fdf</strong></dt>
<dd><p>This data type defines a general system of functions with arbitrary parameters,
a function to compute <em>J u</em> or <em>J^T u</em> for a given vector <em>u</em>,
the normal equations matrix <em>J^T J</em>,
and optionally the second directional derivative of the functions for geodesic acceleration.
</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 <em>n</em> components of the vector
<em>f(x)</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) (CBLAS_TRANSPOSE_t <var>TransJ</var>, const gsl_vector * <var>x</var>, const gsl_vector * <var>u</var>, void * <var>params</var>, gsl_vector * <var>v</var>, gsl_matrix * <var>JTJ</var>)</code></dt>
<dd><p>If <var>TransJ</var> is equal to <code>CblasNoTrans</code>, then this function should
compute the matrix-vector product <em>J u</em> and store the result in <var>v</var>.
If <var>TransJ</var> is equal to <code>CblasTrans</code>, then this function should
compute the matrix-vector product <em>J^T u</em> and store the result in <var>v</var>.
Additionally, the normal equations matrix <em>J^T J</em> should be stored in the
lower half of <var>JTJ</var>. The input matrix <var>JTJ</var> could be set to NULL,
for example by iterative methods which do not require this matrix, so the user
should check for this prior to constructing the matrix.
The input <var>params</var> contains the arbitrary parameters.
</p>
</dd>
<dt><code>int (* fvv) (const gsl_vector * <var>x</var>, const gsl_vector * <var>v</var>, void * <var>params</var>, gsl_vector * <var>fvv</var>)</code></dt>
<dd><p>When geodesic acceleration is enabled, this function should store the
<em>n</em> components of the vector
<em>f_{vv}(x) = \sum_{\alpha\beta} v_{\alpha} v_{\beta} {\partial \over \partial x_{\alpha}} {\partial \over \partial x_{\beta}} f(x)</em>,
representing second directional derivatives of the function to be minimized,
into the output <var>fvv</var>. The parameter vector is provided in <var>x</var> and
the velocity vector is provided in <var>v</var>, both of which have <em>p</em>
components. The arbitrary parameters are given in <var>params</var>. If
analytic expressions for <em>f_{vv}(x)</em> are unavailable or too difficult
to compute, this function pointer may be set to NULL, in which case
<em>f_{vv}(x)</em> will be computed internally using a finite difference
approximation.
</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>
<dt><code>size_t nevalf</code></dt>
<dd><p>This does not need to be set by the user. It counts the number of
function evaluations and is initialized by the <code>_init</code> function.
</p>
</dd>
<dt><code>size_t nevaldfu</code></dt>
<dd><p>This does not need to be set by the user. It counts the number of
Jacobian matrix-vector evaluations (<em>J u</em> or <em>J^T u</em>) and
is initialized by the <code>_init</code> function.
</p>
</dd>
<dt><code>size_t nevaldf2</code></dt>
<dd><p>This does not need to be set by the user. It counts the number of
<em>J^T J</em> evaluations and is initialized by the <code>_init</code> function.
</p>
</dd>
<dt><code>size_t nevalfvv</code></dt>
<dd><p>This does not need to be set by the user. It counts the number of
<em>f_{vv}(x)</em> evaluations and is initialized by the <code>_init</code> 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="Nonlinear-Least_002dSquares-Iteration.html#Nonlinear-Least_002dSquares-Iteration" accesskey="n" rel="next">Nonlinear Least-Squares Iteration</a>, Previous: <a href="Nonlinear-Least_002dSquares-Initialization.html#Nonlinear-Least_002dSquares-Initialization" accesskey="p" rel="previous">Nonlinear Least-Squares Initialization</a>, Up: <a href="Nonlinear-Least_002dSquares-Fitting.html#Nonlinear-Least_002dSquares-Fitting" accesskey="u" rel="up">Nonlinear Least-Squares Fitting</a> &nbsp; [<a href="Function-Index.html#Function-Index" title="Index" rel="index">Index</a>]</p>
</div>



</body>
</html>