File: Nonlinear-Least_002dSquares-Tunable-Parameters.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 (387 lines) | stat: -rw-r--r-- 19,166 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
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
<!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 Tunable Parameters</title>

<meta name="description" content="GNU Scientific Library &ndash; Reference Manual: Nonlinear Least-Squares Tunable Parameters">
<meta name="keywords" content="GNU Scientific Library &ndash; Reference Manual: Nonlinear Least-Squares Tunable Parameters">
<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-Initialization.html#Nonlinear-Least_002dSquares-Initialization" rel="next" title="Nonlinear Least-Squares Initialization">
<link href="Nonlinear-Least_002dSquares-Weighted-Overview.html#Nonlinear-Least_002dSquares-Weighted-Overview" rel="previous" title="Nonlinear Least-Squares Weighted Overview">
<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-Tunable-Parameters"></a>
<div class="header">
<p>
Next: <a href="Nonlinear-Least_002dSquares-Initialization.html#Nonlinear-Least_002dSquares-Initialization" accesskey="n" rel="next">Nonlinear Least-Squares Initialization</a>, Previous: <a href="Nonlinear-Least_002dSquares-Weighted-Overview.html#Nonlinear-Least_002dSquares-Weighted-Overview" accesskey="p" rel="previous">Nonlinear Least-Squares Weighted Overview</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="Tunable-Parameters"></a>
<h3 class="section">39.4 Tunable Parameters</h3>

<p>The user can tune nearly all aspects of the iteration at allocation
time. For the <code>gsl_multifit_nlinear</code> interface, the user may
modify the <code>gsl_multifit_nlinear_parameters</code> structure, which is
defined as follows:
</p>
<div class="example">
<pre class="example">typedef struct
{
  const gsl_multifit_nlinear_trs *trs;        /* trust region subproblem method */
  const gsl_multifit_nlinear_scale *scale;    /* scaling method */
  const gsl_multifit_nlinear_solver *solver;  /* solver method */
  gsl_multifit_nlinear_fdtype fdtype;         /* finite difference method */
  double factor_up;                           /* factor for increasing trust radius */
  double factor_down;                         /* factor for decreasing trust radius */
  double avmax;                               /* max allowed |a|/|v| */
  double h_df;                                /* step size for finite difference Jacobian */
  double h_fvv;                               /* step size for finite difference fvv */
} gsl_multifit_nlinear_parameters;
</pre></div>

<p>For the <code>gsl_multilarge_nlinear</code> interface, the user may
modify the <code>gsl_multilarge_nlinear_parameters</code> structure, which is
defined as follows:
</p>
<div class="example">
<pre class="example">typedef struct
{
  const gsl_multilarge_nlinear_trs *trs;       /* trust region subproblem method */
  const gsl_multilarge_nlinear_scale *scale;   /* scaling method */
  const gsl_multilarge_nlinear_solver *solver; /* solver method */
  gsl_multilarge_nlinear_fdtype fdtype;        /* finite difference method */
  double factor_up;                            /* factor for increasing trust radius */
  double factor_down;                          /* factor for decreasing trust radius */
  double avmax;                                /* max allowed |a|/|v| */
  double h_df;                                 /* step size for finite difference Jacobian */
  double h_fvv;                                /* step size for finite difference fvv */
  size_t max_iter;                             /* maximum iterations for trs method */
  double tol;                                  /* tolerance for solving trs */
} gsl_multilarge_nlinear_parameters;
</pre></div>

<p>Each of these parameters is discussed in further detail below.
</p>
<dl>
<dt><a name="index-trs"></a>Parameter: <em>const gsl_multifit_nlinear_trs *</em> <strong>trs</strong></dt>
<dt><a name="index-trs-1"></a>Parameter: <em>const gsl_multilarge_nlinear_trs *</em> <strong>trs</strong></dt>
<dd>
<p>This parameter determines the method used to solve the trust region
subproblem, and may be selected from the following choices,
</p>
<dl>
<dt><a name="index-gsl_005fmultifit_005fnlinear_005ftrs_005flm"></a>Default: <strong>gsl_multifit_nlinear_trs_lm</strong></dt>
<dt><a name="index-gsl_005fmultilarge_005fnlinear_005ftrs_005flm"></a>Default: <strong>gsl_multilarge_nlinear_trs_lm</strong></dt>
<dd><p>This selects the Levenberg-Marquardt algorithm.
</p></dd></dl>

<dl>
<dt><a name="index-gsl_005fmultifit_005fnlinear_005ftrs_005flmaccel"></a>Option: <strong>gsl_multifit_nlinear_trs_lmaccel</strong></dt>
<dt><a name="index-gsl_005fmultilarge_005fnlinear_005ftrs_005flmaccel"></a>Option: <strong>gsl_multilarge_nlinear_trs_lmaccel</strong></dt>
<dd><p>This selects the Levenberg-Marquardt algorithm with geodesic
acceleration.
</p></dd></dl>

<dl>
<dt><a name="index-gsl_005fmultifit_005fnlinear_005ftrs_005fdogleg"></a>Option: <strong>gsl_multifit_nlinear_trs_dogleg</strong></dt>
<dt><a name="index-gsl_005fmultilarge_005fnlinear_005ftrs_005fdogleg"></a>Option: <strong>gsl_multilarge_nlinear_trs_dogleg</strong></dt>
<dd><p>This selects the dogleg algorithm.
</p></dd></dl>

<dl>
<dt><a name="index-gsl_005fmultifit_005fnlinear_005ftrs_005fddogleg"></a>Option: <strong>gsl_multifit_nlinear_trs_ddogleg</strong></dt>
<dt><a name="index-gsl_005fmultilarge_005fnlinear_005ftrs_005fddogleg"></a>Option: <strong>gsl_multilarge_nlinear_trs_ddogleg</strong></dt>
<dd><p>This selects the double dogleg algorithm.
</p></dd></dl>

<dl>
<dt><a name="index-gsl_005fmultifit_005fnlinear_005ftrs_005fsubspace2D"></a>Option: <strong>gsl_multifit_nlinear_trs_subspace2D</strong></dt>
<dt><a name="index-gsl_005fmultilarge_005fnlinear_005ftrs_005fsubspace2D"></a>Option: <strong>gsl_multilarge_nlinear_trs_subspace2D</strong></dt>
<dd><p>This selects the 2D subspace algorithm.
</p></dd></dl>

<dl>
<dt><a name="index-gsl_005fmultilarge_005fnlinear_005ftrs_005fcgst"></a>Option: <strong>gsl_multilarge_nlinear_trs_cgst</strong></dt>
<dd><p>This selects the Steihaug-Toint conjugate gradient algorithm. This
method is available only for large systems.
</p></dd></dl>

</dd></dl>

<dl>
<dt><a name="index-scale"></a>Parameter: <em>const gsl_multifit_nlinear_scale *</em> <strong>scale</strong></dt>
<dt><a name="index-scale-1"></a>Parameter: <em>const gsl_multilarge_nlinear_scale *</em> <strong>scale</strong></dt>
<dd>
<p>This parameter determines the diagonal scaling matrix <em>D</em> and
may be selected from the following choices,
</p>
<dl>
<dt><a name="index-gsl_005fmultifit_005fnlinear_005fscale_005fmore"></a>Default: <strong>gsl_multifit_nlinear_scale_more</strong></dt>
<dt><a name="index-gsl_005fmultilarge_005fnlinear_005fscale_005fmore"></a>Default: <strong>gsl_multilarge_nlinear_scale_more</strong></dt>
<dd><p>This damping strategy was suggested by Mor&eacute;, and
corresponds to <em>D^T D = </em> max(diag(<em>J^T J</em>)),
in other words the maximum elements of
diag(<em>J^T J</em>) encountered thus far in the iteration.
This choice of <em>D</em> makes the problem scale-invariant,
so that if the model parameters <em>x_i</em> are each scaled
by an arbitrary constant, <em>\tilde{x}_i = a_i x_i</em>, then
the sequence of iterates produced by the algorithm would
be unchanged. This method can work very well in cases
where the model parameters have widely different scales
(ie: if some parameters are measured in nanometers, while others
are measured in degrees Kelvin). This strategy has been proven
effective on a large class of problems and so it is the library
default, but it may not be the best choice for all problems.
</p></dd></dl>

<dl>
<dt><a name="index-gsl_005fmultifit_005fnlinear_005fscale_005flevenberg"></a>Option: <strong>gsl_multifit_nlinear_scale_levenberg</strong></dt>
<dt><a name="index-gsl_005fmultilarge_005fnlinear_005fscale_005flevenberg"></a>Option: <strong>gsl_multilarge_nlinear_scale_levenberg</strong></dt>
<dd><p>This damping strategy was originally suggested by Levenberg, and
corresponds to <em>D^T D = I</em>. This method has also proven
effective on a large class of problems, but is not scale-invariant.
However, some authors (e.g. Transtrum and Sethna 2012) argue
that this choice is better for problems which are susceptible
to parameter evaporation (ie: parameters go to infinity)
</p></dd></dl>

<dl>
<dt><a name="index-gsl_005fmultifit_005fnlinear_005fscale_005fmarquardt"></a>Option: <strong>gsl_multifit_nlinear_scale_marquardt</strong></dt>
<dt><a name="index-gsl_005fmultilarge_005fnlinear_005fscale_005fmarquardt"></a>Option: <strong>gsl_multilarge_nlinear_scale_marquardt</strong></dt>
<dd><p>This damping strategy was suggested by Marquardt, and
corresponds to <em>D^T D = </em> diag(<em>J^T J</em>). This
method is scale-invariant, but it is generally considered
inferior to both the Levenberg and Mor&eacute; strategies, though
may work well on certain classes of problems.
</p></dd></dl>

</dd></dl>

<dl>
<dt><a name="index-solver"></a>Parameter: <em>const gsl_multifit_nlinear_solver *</em> <strong>solver</strong></dt>
<dt><a name="index-solver-1"></a>Parameter: <em>const gsl_multilarge_nlinear_solver *</em> <strong>solver</strong></dt>
<dd>
<p>Solving the trust region subproblem on each iteration almost always
requires the solution of the following linear least squares system
</p>
<div class="example">
<pre class="example">[J; sqrt(mu) D] \delta = - [f; 0]
</pre></div>


<p>The <var>solver</var> parameter determines how the system is
solved and can be selected from the following choices:
</p>
<dl>
<dt><a name="index-gsl_005fmultifit_005fnlinear_005fsolver_005fqr"></a>Default: <strong>gsl_multifit_nlinear_solver_qr</strong></dt>
<dd><p>This method solves the system using a rank revealing QR
decomposition of the Jacobian <em>J</em>. This method will
produce reliable solutions in cases where the Jacobian
is rank deficient or near-singular but does require about
twice as many operations as the Cholesky method discussed
below.
</p></dd></dl>

<dl>
<dt><a name="index-gsl_005fmultifit_005fnlinear_005fsolver_005fcholesky"></a>Option: <strong>gsl_multifit_nlinear_solver_cholesky</strong></dt>
<dt><a name="index-gsl_005fmultilarge_005fnlinear_005fsolver_005fcholesky"></a>Default: <strong>gsl_multilarge_nlinear_solver_cholesky</strong></dt>
<dd><p>This method solves the alternate normal equations problem
</p>
<div class="example">
<pre class="example">( J^T J + \mu D^T D ) \delta = -J^T f
</pre></div>


<p>by using a Cholesky decomposition of the matrix
<em>J^T J + \mu D^T D</em>. This method is faster than the
QR approach, however it is susceptible to numerical instabilities
if the Jacobian matrix is rank deficient or near-singular. In
these cases, an attempt is made to reduce the condition number
of the matrix using Jacobi preconditioning, but for highly
ill-conditioned problems the QR approach is better. If it is
known that the Jacobian matrix is well conditioned, this method
is accurate and will perform faster than the QR approach.
</p></dd></dl>

<dl>
<dt><a name="index-gsl_005fmultifit_005fnlinear_005fsolver_005fsvd"></a>Option: <strong>gsl_multifit_nlinear_solver_svd</strong></dt>
<dd><p>This method solves the system using a singular value
decomposition of the Jacobian <em>J</em>. This method will
produce the most reliable solutions for ill-conditioned Jacobians
but is also the slowest solver method.
</p></dd></dl>

</dd></dl>

<dl>
<dt><a name="index-fdtype"></a>Parameter: <em>gsl_multifit_nlinear_fdtype</em> <strong>fdtype</strong></dt>
<dd>
<p>This parameter specifies whether to use forward or centered
differences when approximating the Jacobian. This is only
used when an analytic Jacobian is not provided to the solver.
This parameter may be set to one of the following choices.
</p>
<dl>
<dt><a name="index-GSL_005fMULTIFIT_005fNLINEAR_005fFWDIFF"></a>Default: <strong>GSL_MULTIFIT_NLINEAR_FWDIFF</strong></dt>
<dd><p>This specifies a forward finite difference to approximate
the Jacobian matrix. The Jacobian matrix will be calculated as
</p>
<div class="example">
<pre class="example">J_ij = 1 / \Delta_j ( f_i(x + \Delta_j e_j) - f_i(x) )
</pre></div>

<p>where <em>\Delta_j = h |x_j|</em> and <em>e_j</em> is the standard
<em>j</em>th Cartesian unit basis vector so that
<em>x + \Delta_j e_j</em> represents a small (forward) perturbation of
the <em>j</em>th parameter by an amount <em>\Delta_j</em>. The perturbation
<em>\Delta_j</em> is proportional to the current value <em>|x_j|</em> which
helps to calculate an accurate Jacobian when the various parameters have
different scale sizes. The value of <em>h</em> is specified by the <code>h_df</code>
parameter. The accuracy of this method is <em>O(h)</em>, and evaluating this
matrix requires an additional <em>p</em> function evaluations.
</p></dd></dl>

<dl>
<dt><a name="index-GSL_005fMULTIFIT_005fNLINEAR_005fCTRDIFF"></a>Option: <strong>GSL_MULTIFIT_NLINEAR_CTRDIFF</strong></dt>
<dd><p>This specifies a centered finite difference to approximate
the Jacobian matrix. The Jacobian matrix will be calculated as
</p>
<div class="example">
<pre class="example">J_ij = 1 / \Delta_j ( f_i(x + 1/2 \Delta_j e_j) - f_i(x - 1/2 \Delta_j e_j) )
</pre></div>

<p>See above for a description of <em>\Delta_j</em>. The accuracy of this
method is <em>O(h^2)</em>, but evaluating this
matrix requires an additional <em>2p</em> function evaluations.
</p></dd></dl>

</dd></dl>

<dl>
<dt><a name="index-factor_005fup"></a>Parameter: <em>double</em> <strong>factor_up</strong></dt>
<dd>
<p>When a step is accepted, the trust region radius will be increased
by this factor. The default value is <em>3</em>.
</p>
</dd></dl>

<dl>
<dt><a name="index-factor_005fdown"></a>Parameter: <em>double</em> <strong>factor_down</strong></dt>
<dd>
<p>When a step is rejected, the trust region radius will be decreased
by this factor. The default value is <em>2</em>.
</p>
</dd></dl>

<dl>
<dt><a name="index-avmax"></a>Parameter: <em>double</em> <strong>avmax</strong></dt>
<dd>
<p>When using geodesic acceleration to solve a nonlinear least squares problem,
an important parameter to monitor is the ratio of the acceleration term
to the velocity term,
</p>
<div class="example">
<pre class="example">|a| / |v|
</pre></div>


<p>If this ratio is small, it means the acceleration correction
is contributing very little to the step. This could be because
the problem is not &ldquo;nonlinear&rdquo; enough to benefit from
the acceleration. If the ratio is large (<em>&gt; 1</em>) it
means that the acceleration is larger than the velocity,
which shouldn&rsquo;t happen since the step represents a truncated
series and so the second order term <em>a</em> should be smaller than
the first order term <em>v</em> to guarantee convergence.
Therefore any steps with a ratio larger than the parameter
<var>avmax</var> are rejected. <var>avmax</var> is set to 0.75 by default.
For problems which experience difficulty converging, this threshold
could be lowered.
</p>
</dd></dl>

<dl>
<dt><a name="index-h_005fdf"></a>Parameter: <em>double</em> <strong>h_df</strong></dt>
<dd>
<p>This parameter specifies the step size for approximating the
Jacobian matrix with finite differences. It is set to
<em>\sqrt{\epsilon}</em> by default, where <em>\epsilon</em>
is <code>GSL_DBL_EPSILON</code>.
</p>
</dd></dl>

<dl>
<dt><a name="index-h_005ffvv"></a>Parameter: <em>double</em> <strong>h_fvv</strong></dt>
<dd>
<p>When using geodesic acceleration, the user must either supply
a function to calculate <em>f_{vv}(x)</em> or the library
can estimate this second directional derivative using a finite
difference method. When using finite differences, the library
must calculate <em>f(x + h v)</em> where <em>h</em> represents
a small step in the velocity direction. The parameter
<var>h_fvv</var> defines this step size and is set to 0.02 by
default.
</p>
</dd></dl>

<hr>
<div class="header">
<p>
Next: <a href="Nonlinear-Least_002dSquares-Initialization.html#Nonlinear-Least_002dSquares-Initialization" accesskey="n" rel="next">Nonlinear Least-Squares Initialization</a>, Previous: <a href="Nonlinear-Least_002dSquares-Weighted-Overview.html#Nonlinear-Least_002dSquares-Weighted-Overview" accesskey="p" rel="previous">Nonlinear Least-Squares Weighted Overview</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>