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
|
<!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 Comparison Example</title>
<meta name="description" content="GNU Scientific Library – Reference Manual: Nonlinear Least-Squares Comparison Example">
<meta name="keywords" content="GNU Scientific Library – Reference Manual: Nonlinear Least-Squares Comparison Example">
<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-Examples.html#Nonlinear-Least_002dSquares-Examples" rel="up" title="Nonlinear Least-Squares Examples">
<link href="Nonlinear-Least_002dSquares-Large-Example.html#Nonlinear-Least_002dSquares-Large-Example" rel="next" title="Nonlinear Least-Squares Large Example">
<link href="Nonlinear-Least_002dSquares-Geodesic-Acceleration-Example.html#Nonlinear-Least_002dSquares-Geodesic-Acceleration-Example" rel="previous" title="Nonlinear Least-Squares Geodesic Acceleration Example">
<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-Comparison-Example"></a>
<div class="header">
<p>
Next: <a href="Nonlinear-Least_002dSquares-Large-Example.html#Nonlinear-Least_002dSquares-Large-Example" accesskey="n" rel="next">Nonlinear Least-Squares Large Example</a>, Previous: <a href="Nonlinear-Least_002dSquares-Geodesic-Acceleration-Example.html#Nonlinear-Least_002dSquares-Geodesic-Acceleration-Example" accesskey="p" rel="previous">Nonlinear Least-Squares Geodesic Acceleration Example</a>, Up: <a href="Nonlinear-Least_002dSquares-Examples.html#Nonlinear-Least_002dSquares-Examples" accesskey="u" rel="up">Nonlinear Least-Squares Examples</a> [<a href="Function-Index.html#Function-Index" title="Index" rel="index">Index</a>]</p>
</div>
<hr>
<a name="Comparing-TRS-Methods-Example"></a>
<h4 class="subsection">39.12.3 Comparing TRS Methods Example</h4>
<p>The following program compares all available nonlinear least squares
trust-region subproblem (TRS) methods on the Branin function, a common
optimization test problem. The cost function is given by
</p>
<div class="example">
<pre class="example">\Phi(x) &= 1/2 (f_1^2 + f_2^2)
f_1 &= x_2 + a_1 x_1^2 + a_2 x_1 + a_3
f_2 &= sqrt(a_4) sqrt(1 + (1 - a_5) cos(x_1))
</pre></div>
<p>with <em>a_1 = -{5.1 \over 4 \pi^2}, a_2 = {5 \over \pi}, a_3 = -6, a_4 = 10, a_5 = {1 \over 8\pi}</em>.
There are three minima of this function in the range
<em>(x_1,x_2) \in [-5,15] \times [-5,15]</em>. The program
below uses the starting point <em>(x_1,x_2) = (6,14.5)</em>
and calculates the solution with all available nonlinear
least squares TRS methods. The program output is shown below.
</p>
<div class="smallformat">
<pre class="verbatim">Method NITER NFEV NJEV Initial Cost Final cost Final cond(J) Final x
levenberg-marquardt 20 27 21 1.9874e+02 3.9789e-01 6.1399e+07 (-3.14e+00, 1.23e+01)
levenberg-marquardt+accel 27 36 28 1.9874e+02 3.9789e-01 1.4465e+07 (3.14e+00, 2.27e+00)
dogleg 23 64 23 1.9874e+02 3.9789e-01 5.0692e+08 (3.14e+00, 2.28e+00)
double-dogleg 24 69 24 1.9874e+02 3.9789e-01 3.4879e+07 (3.14e+00, 2.27e+00)
2D-subspace 23 54 24 1.9874e+02 3.9789e-01 2.5142e+07 (3.14e+00, 2.27e+00)
</pre></div>
<p>The first row of output above corresponds to standard Levenberg-Marquardt, while
the second row includes geodesic acceleration. We see that the standard LM method
converges to the minimum at <em>(-\pi,12.275)</em> and also uses the least number
of iterations and Jacobian evaluations. All other methods converge to the minimum
<em>(\pi,2.275)</em> and perform similarly in terms of number of Jacobian evaluations.
We see that <em>J</em> is fairly ill-conditioned
at both minima, indicating that the QR (or SVD) solver is the best choice for this problem.
Since there are only two parameters in this optimization problem, we can easily
visualize the paths taken by each method, which are shown in the figure below.
The figure shows contours of the cost function <em>\Phi(x_1,x_2)</em> which exhibits
three global minima in the range <em>[-5,15] \times [-5,15]</em>. The paths taken
by each solver are shown as colored lines.
</p>
<p>The program is given below.
</p>
<div class="example">
<pre class="verbatim">#include <stdlib.h>
#include <stdio.h>
#include <gsl/gsl_vector.h>
#include <gsl/gsl_matrix.h>
#include <gsl/gsl_blas.h>
#include <gsl/gsl_multifit_nlinear.h>
/* parameters to model */
struct model_params
{
double a1;
double a2;
double a3;
double a4;
double a5;
};
/* Branin function */
int
func_f (const gsl_vector * x, void *params, gsl_vector * f)
{
struct model_params *par = (struct model_params *) params;
double x1 = gsl_vector_get(x, 0);
double x2 = gsl_vector_get(x, 1);
double f1 = x2 + par->a1 * x1 * x1 + par->a2 * x1 + par->a3;
double f2 = sqrt(par->a4) * sqrt(1.0 + (1.0 - par->a5) * cos(x1));
gsl_vector_set(f, 0, f1);
gsl_vector_set(f, 1, f2);
return GSL_SUCCESS;
}
int
func_df (const gsl_vector * x, void *params, gsl_matrix * J)
{
struct model_params *par = (struct model_params *) params;
double x1 = gsl_vector_get(x, 0);
double f2 = sqrt(par->a4) * sqrt(1.0 + (1.0 - par->a5) * cos(x1));
gsl_matrix_set(J, 0, 0, 2.0 * par->a1 * x1 + par->a2);
gsl_matrix_set(J, 0, 1, 1.0);
gsl_matrix_set(J, 1, 0, -0.5 * par->a4 / f2 * (1.0 - par->a5) * sin(x1));
gsl_matrix_set(J, 1, 1, 0.0);
return GSL_SUCCESS;
}
int
func_fvv (const gsl_vector * x, const gsl_vector * v,
void *params, gsl_vector * fvv)
{
struct model_params *par = (struct model_params *) params;
double x1 = gsl_vector_get(x, 0);
double v1 = gsl_vector_get(v, 0);
double c = cos(x1);
double s = sin(x1);
double f2 = sqrt(par->a4) * sqrt(1.0 + (1.0 - par->a5) * c);
double t = 0.5 * par->a4 * (1.0 - par->a5) / f2;
gsl_vector_set(fvv, 0, 2.0 * par->a1 * v1 * v1);
gsl_vector_set(fvv, 1, -t * (c + s*s/f2) * v1 * v1);
return GSL_SUCCESS;
}
void
callback(const size_t iter, void *params,
const gsl_multifit_nlinear_workspace *w)
{
gsl_vector * x = gsl_multifit_nlinear_position(w);
double x1 = gsl_vector_get(x, 0);
double x2 = gsl_vector_get(x, 1);
/* print out current location */
printf("%f %f\n", x1, x2);
}
void
solve_system(gsl_vector *x0, gsl_multifit_nlinear_fdf *fdf,
gsl_multifit_nlinear_parameters *params)
{
const gsl_multifit_nlinear_type *T = gsl_multifit_nlinear_trust;
const size_t max_iter = 200;
const double xtol = 1.0e-8;
const double gtol = 1.0e-8;
const double ftol = 1.0e-8;
const size_t n = fdf->n;
const size_t p = fdf->p;
gsl_multifit_nlinear_workspace *work =
gsl_multifit_nlinear_alloc(T, params, n, p);
gsl_vector * f = gsl_multifit_nlinear_residual(work);
gsl_vector * x = gsl_multifit_nlinear_position(work);
int info;
double chisq0, chisq, rcond;
printf("# %s/%s\n",
gsl_multifit_nlinear_name(work),
gsl_multifit_nlinear_trs_name(work));
/* initialize solver */
gsl_multifit_nlinear_init(x0, fdf, work);
/* store initial cost */
gsl_blas_ddot(f, f, &chisq0);
/* iterate until convergence */
gsl_multifit_nlinear_driver(max_iter, xtol, gtol, ftol,
callback, NULL, &info, work);
/* store final cost */
gsl_blas_ddot(f, f, &chisq);
/* store cond(J(x)) */
gsl_multifit_nlinear_rcond(&rcond, work);
/* print summary */
fprintf(stderr, "%-25s %-6zu %-5zu %-5zu %-13.4e %-12.4e %-13.4e (%.2e, %.2e)\n",
gsl_multifit_nlinear_trs_name(work),
gsl_multifit_nlinear_niter(work),
fdf->nevalf,
fdf->nevaldf,
chisq0,
chisq,
1.0 / rcond,
gsl_vector_get(x, 0),
gsl_vector_get(x, 1));
printf("\n\n");
gsl_multifit_nlinear_free(work);
}
int
main (void)
{
const size_t n = 2;
const size_t p = 2;
gsl_vector *f = gsl_vector_alloc(n);
gsl_vector *x = gsl_vector_alloc(p);
gsl_multifit_nlinear_fdf fdf;
gsl_multifit_nlinear_parameters fdf_params =
gsl_multifit_nlinear_default_parameters();
struct model_params params;
params.a1 = -5.1 / (4.0 * M_PI * M_PI);
params.a2 = 5.0 / M_PI;
params.a3 = -6.0;
params.a4 = 10.0;
params.a5 = 1.0 / (8.0 * M_PI);
/* print map of Phi(x1, x2) */
{
double x1, x2, chisq;
for (x1 = -5.0; x1 < 15.0; x1 += 0.1)
{
for (x2 = -5.0; x2 < 15.0; x2 += 0.1)
{
gsl_vector_set(x, 0, x1);
gsl_vector_set(x, 1, x2);
func_f(x, &params, f);
gsl_blas_ddot(f, f, &chisq);
printf("%f %f %f\n", x1, x2, chisq);
}
printf("\n");
}
printf("\n\n");
}
/* define function to be minimized */
fdf.f = func_f;
fdf.df = func_df;
fdf.fvv = func_fvv;
fdf.n = n;
fdf.p = p;
fdf.params = &params;
/* starting point */
gsl_vector_set(x, 0, 6.0);
gsl_vector_set(x, 1, 14.5);
fprintf(stderr, "%-25s %-6s %-5s %-5s %-13s %-12s %-13s %-15s\n",
"Method", "NITER", "NFEV", "NJEV", "Initial Cost",
"Final cost", "Final cond(J)", "Final x");
fdf_params.trs = gsl_multifit_nlinear_trs_lm;
solve_system(x, &fdf, &fdf_params);
fdf_params.trs = gsl_multifit_nlinear_trs_lmaccel;
solve_system(x, &fdf, &fdf_params);
fdf_params.trs = gsl_multifit_nlinear_trs_dogleg;
solve_system(x, &fdf, &fdf_params);
fdf_params.trs = gsl_multifit_nlinear_trs_ddogleg;
solve_system(x, &fdf, &fdf_params);
fdf_params.trs = gsl_multifit_nlinear_trs_subspace2D;
solve_system(x, &fdf, &fdf_params);
gsl_vector_free(f);
gsl_vector_free(x);
return 0;
}
</pre></div>
<hr>
<div class="header">
<p>
Next: <a href="Nonlinear-Least_002dSquares-Large-Example.html#Nonlinear-Least_002dSquares-Large-Example" accesskey="n" rel="next">Nonlinear Least-Squares Large Example</a>, Previous: <a href="Nonlinear-Least_002dSquares-Geodesic-Acceleration-Example.html#Nonlinear-Least_002dSquares-Geodesic-Acceleration-Example" accesskey="p" rel="previous">Nonlinear Least-Squares Geodesic Acceleration Example</a>, Up: <a href="Nonlinear-Least_002dSquares-Examples.html#Nonlinear-Least_002dSquares-Examples" accesskey="u" rel="up">Nonlinear Least-Squares Examples</a> [<a href="Function-Index.html#Function-Index" title="Index" rel="index">Index</a>]</p>
</div>
</body>
</html>
|