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
|
<!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 – Reference Manual: Nonlinear Least-Squares Initialization</title>
<meta name="description" content="GNU Scientific Library – Reference Manual: Nonlinear Least-Squares Initialization">
<meta name="keywords" content="GNU Scientific Library – Reference Manual: Nonlinear Least-Squares Initialization">
<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-Function-Definition.html#Nonlinear-Least_002dSquares-Function-Definition" rel="next" title="Nonlinear Least-Squares Function Definition">
<link href="Nonlinear-Least_002dSquares-Tunable-Parameters.html#Nonlinear-Least_002dSquares-Tunable-Parameters" rel="previous" title="Nonlinear Least-Squares Tunable Parameters">
<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-Initialization"></a>
<div class="header">
<p>
Next: <a href="Nonlinear-Least_002dSquares-Function-Definition.html#Nonlinear-Least_002dSquares-Function-Definition" accesskey="n" rel="next">Nonlinear Least-Squares Function Definition</a>, Previous: <a href="Nonlinear-Least_002dSquares-Tunable-Parameters.html#Nonlinear-Least_002dSquares-Tunable-Parameters" accesskey="p" rel="previous">Nonlinear Least-Squares Tunable Parameters</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="Initializing-the-Solver-3"></a>
<h3 class="section">39.5 Initializing the Solver</h3>
<dl>
<dt><a name="index-gsl_005fmultifit_005fnlinear_005falloc"></a>Function: <em>gsl_multifit_nlinear_workspace *</em> <strong>gsl_multifit_nlinear_alloc</strong> <em>(const gsl_multifit_nlinear_type * <var>T</var>, const gsl_multifit_nlinear_parameters * <var>params</var>, const size_t <var>n</var>, const size_t <var>p</var>)</em></dt>
<dt><a name="index-gsl_005fmultilarge_005fnlinear_005falloc"></a>Function: <em>gsl_multilarge_nlinear_workspace *</em> <strong>gsl_multilarge_nlinear_alloc</strong> <em>(const gsl_multilarge_nlinear_type * <var>T</var>, const gsl_multilarge_nlinear_parameters * <var>params</var>, const size_t <var>n</var>, const size_t <var>p</var>)</em></dt>
<dd><a name="index-gsl_005fmultifit_005fnlinear_005falloc-1"></a>
<a name="index-gsl_005fmultifit_005fnlinear_005ftype"></a>
<p>These functions return a pointer to a newly allocated instance of a
derivative solver of type <var>T</var> for <var>n</var> observations and <var>p</var>
parameters. The <var>params</var> input specifies a tunable set of
parameters which will affect important details in each iteration
of the trust region subproblem algorithm. It is recommended to start
with the suggested default parameters (see
<code>gsl_multifit_nlinear_default_parameters</code> and
<code>gsl_multilarge_nlinear_default_parameters</code>) and then tune
the parameters once the code is working correctly. See
<a href="Nonlinear-Least_002dSquares-Tunable-Parameters.html#Nonlinear-Least_002dSquares-Tunable-Parameters">Nonlinear Least-Squares Tunable Parameters</a>
for descriptions of the various parameters.
For example, the following code creates an instance of a
Levenberg-Marquardt solver for 100 data points and 3 parameters,
using suggested defaults:
</p>
<div class="example">
<pre class="example">const gsl_multifit_nlinear_type * T
= gsl_multifit_nlinear_lm;
gsl_multifit_nlinear_parameters params
= gsl_multifit_nlinear_default_parameters();
gsl_multifit_nlinear_workspace * w
= gsl_multifit_nlinear_alloc (T, &params, 100, 3);
</pre></div>
<p>The number of observations <var>n</var> must be greater than or equal to
parameters <var>p</var>.
</p>
<p>If there is insufficient memory to create the solver then the function
returns a null pointer and the error handler is invoked with an error
code of <code>GSL_ENOMEM</code>.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005fmultifit_005fnlinear_005fdefault_005fparameters"></a>Function: <em>gsl_multifit_nlinear_parameters</em> <strong>gsl_multifit_nlinear_default_parameters</strong> <em>(void)</em></dt>
<dt><a name="index-gsl_005fmultilarge_005fnlinear_005fdefault_005fparameters"></a>Function: <em>gsl_multilarge_nlinear_parameters</em> <strong>gsl_multilarge_nlinear_default_parameters</strong> <em>(void)</em></dt>
<dd><p>These functions return a set of recommended default parameters
for use in solving nonlinear least squares problems. The user
can tune each parameter to improve the performance on their
particular problem, see
<a href="Nonlinear-Least_002dSquares-Tunable-Parameters.html#Nonlinear-Least_002dSquares-Tunable-Parameters">Nonlinear Least-Squares Tunable Parameters</a>.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005fmultifit_005fnlinear_005finit"></a>Function: <em>int</em> <strong>gsl_multifit_nlinear_init</strong> <em>(const gsl_vector * <var>x</var>, gsl_multifit_nlinear_fdf * <var>fdf</var>, gsl_multifit_nlinear_workspace * <var>w</var>)</em></dt>
<dt><a name="index-gsl_005fmultifit_005fnlinear_005fwinit"></a>Function: <em>int</em> <strong>gsl_multifit_nlinear_winit</strong> <em>(const gsl_vector * <var>x</var>, const gsl_vector * <var>wts</var>, gsl_multifit_nlinear_fdf * <var>fdf</var>, gsl_multifit_nlinear_workspace * <var>w</var>)</em></dt>
<dt><a name="index-gsl_005fmultilarge_005fnlinear_005finit"></a>Function: <em>int</em> <strong>gsl_multilarge_nlinear_init</strong> <em>(const gsl_vector * <var>x</var>, gsl_multilarge_nlinear_fdf * <var>fdf</var>, gsl_multilarge_nlinear_workspace * <var>w</var>)</em></dt>
<dt><a name="index-gsl_005fmultilarge_005fnlinear_005fwinit"></a>Function: <em>int</em> <strong>gsl_multilarge_nlinear_winit</strong> <em>(const gsl_vector * <var>x</var>, const gsl_vector * <var>wts</var>, gsl_multilarge_nlinear_fdf * <var>fdf</var>, gsl_multilarge_nlinear_workspace * <var>w</var>)</em></dt>
<dd><p>These functions initialize, or reinitialize, an existing workspace <var>w</var>
to use the system <var>fdf</var> and the initial guess
<var>x</var>. See <a href="Nonlinear-Least_002dSquares-Function-Definition.html#Nonlinear-Least_002dSquares-Function-Definition">Nonlinear Least-Squares Function Definition</a>
for a description of the <var>fdf</var> structure.
</p>
<p>Optionally, a weight vector <var>wts</var> can be given to perform
a weighted nonlinear regression. Here, the weighting matrix is
<em>W = diag(w_1,w_2,...,w_n)</em>.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005fmultifit_005fnlinear_005ffree"></a>Function: <em>void</em> <strong>gsl_multifit_nlinear_free</strong> <em>(gsl_multifit_nlinear_workspace * <var>w</var>)</em></dt>
<dt><a name="index-gsl_005fmultilarge_005fnlinear_005ffree"></a>Function: <em>void</em> <strong>gsl_multilarge_nlinear_free</strong> <em>(gsl_multilarge_nlinear_workspace * <var>w</var>)</em></dt>
<dd><p>These functions free all the memory associated with the workspace <var>w</var>.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005fmultifit_005fnlinear_005fname"></a>Function: <em>const char *</em> <strong>gsl_multifit_nlinear_name</strong> <em>(const gsl_multifit_nlinear_workspace * <var>w</var>)</em></dt>
<dt><a name="index-gsl_005fmultilarge_005fnlinear_005fname"></a>Function: <em>const char *</em> <strong>gsl_multilarge_nlinear_name</strong> <em>(const gsl_multilarge_nlinear_workspace * <var>w</var>)</em></dt>
<dd><p>These functions return a pointer to the name of the solver. For example,
</p>
<div class="example">
<pre class="example">printf ("w is a '%s' solver\n",
gsl_multifit_nlinear_name (w));
</pre></div>
<p>would print something like <code>w is a 'trust-region' solver</code>.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005fmultifit_005fnlinear_005ftrs_005fname"></a>Function: <em>const char *</em> <strong>gsl_multifit_nlinear_trs_name</strong> <em>(const gsl_multifit_nlinear_workspace * <var>w</var>)</em></dt>
<dt><a name="index-gsl_005fmultilarge_005fnlinear_005ftrs_005fname"></a>Function: <em>const char *</em> <strong>gsl_multilarge_nlinear_trs_name</strong> <em>(const gsl_multilarge_nlinear_workspace * <var>w</var>)</em></dt>
<dd><p>These functions return a pointer to the name of the trust region subproblem
method. For example,
</p>
<div class="example">
<pre class="example">printf ("w is a '%s' solver\n",
gsl_multifit_nlinear_trs_name (w));
</pre></div>
<p>would print something like <code>w is a 'levenberg-marquardt' solver</code>.
</p></dd></dl>
<hr>
<div class="header">
<p>
Next: <a href="Nonlinear-Least_002dSquares-Function-Definition.html#Nonlinear-Least_002dSquares-Function-Definition" accesskey="n" rel="next">Nonlinear Least-Squares Function Definition</a>, Previous: <a href="Nonlinear-Least_002dSquares-Tunable-Parameters.html#Nonlinear-Least_002dSquares-Tunable-Parameters" accesskey="p" rel="previous">Nonlinear Least-Squares Tunable Parameters</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>
|