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
|
<!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: Fitting large linear systems example</title>
<meta name="description" content="GNU Scientific Library – Reference Manual: Fitting large linear systems example">
<meta name="keywords" content="GNU Scientific Library – Reference Manual: Fitting large linear systems 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="Fitting-Examples.html#Fitting-Examples" rel="up" title="Fitting Examples">
<link href="Fitting-References-and-Further-Reading.html#Fitting-References-and-Further-Reading" rel="next" title="Fitting References and Further Reading">
<link href="Fitting-robust-linear-regression-example.html#Fitting-robust-linear-regression-example" rel="previous" title="Fitting robust linear regression 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="Fitting-large-linear-systems-example"></a>
<div class="header">
<p>
Previous: <a href="Fitting-robust-linear-regression-example.html#Fitting-robust-linear-regression-example" accesskey="p" rel="previous">Fitting robust linear regression example</a>, Up: <a href="Fitting-Examples.html#Fitting-Examples" accesskey="u" rel="up">Fitting Examples</a> [<a href="Function-Index.html#Function-Index" title="Index" rel="index">Index</a>]</p>
</div>
<hr>
<a name="Large-Dense-Linear-Regression-Example"></a>
<h4 class="subsection">38.8.6 Large Dense Linear Regression Example</h4>
<p>The following program demostrates the large dense linear least squares
solvers. This example is adapted from Trefethen and Bau,
and fits the function <em>f(t) = \exp{(\sin^3{(10t)}})</em> on
the interval <em>[0,1]</em> with a degree 15 polynomial. The
program generates <em>n = 50000</em> equally spaced points
<em>t_i</em> on this interval, calculates the function value
and adds random noise to determine the observation value
<em>y_i</em>. The entries of the least squares matrix are
<em>X_{ij} = t_i^j</em>, representing a polynomial fit. The
matrix is highly ill-conditioned, with a condition number
of about <em>1.4 \cdot 10^{11}</em>. The program accumulates the
matrix into the least squares system in 5 blocks, each with
10000 rows. This way the full matrix <em>X</em> is never
stored in memory. We solve the system with both the
normal equations and TSQR methods. The results are shown
in the plot below. In the top left plot, we see the unregularized
normal equations solution has larger error than TSQR due to
the ill-conditioning of the matrix. In the bottom left plot,
we show the L-curve, which exhibits multiple corners.
In the top right panel, we plot a regularized solution using
<em>\lambda = 10^{-6}</em>. The TSQR and normal solutions now agree,
however they are unable to provide a good fit due to the damping.
This indicates that for some ill-conditioned
problems, regularizing the normal equations does not improve the
solution. This is further illustrated in the bottom right panel,
where we plot the L-curve calculated from the normal equations.
The curve agrees with the TSQR curve for larger damping parameters,
but for small <em>\lambda</em>, the normal equations approach cannot
provide accurate solution vectors leading to numerical
inaccuracies in the left portion of the curve.
</p>
<div class="example">
<pre class="verbatim">#include <gsl/gsl_math.h>
#include <gsl/gsl_vector.h>
#include <gsl/gsl_matrix.h>
#include <gsl/gsl_rng.h>
#include <gsl/gsl_randist.h>
#include <gsl/gsl_multifit.h>
#include <gsl/gsl_multilarge.h>
#include <gsl/gsl_blas.h>
/* function to be fitted */
double
func(const double t)
{
double x = sin(10.0 * t);
return exp(x*x*x);
}
/* construct a row of the least squares matrix */
int
build_row(const double t, gsl_vector *row)
{
const size_t p = row->size;
double Xj = 1.0;
size_t j;
for (j = 0; j < p; ++j)
{
gsl_vector_set(row, j, Xj);
Xj *= t;
}
return 0;
}
int
solve_system(const int print_data, const gsl_multilarge_linear_type * T,
const double lambda, const size_t n, const size_t p,
gsl_vector * c)
{
const size_t nblock = 5; /* number of blocks to accumulate */
const size_t nrows = n / nblock; /* number of rows per block */
gsl_multilarge_linear_workspace * w =
gsl_multilarge_linear_alloc(T, p);
gsl_matrix *X = gsl_matrix_alloc(nrows, p);
gsl_vector *y = gsl_vector_alloc(nrows);
gsl_rng *r = gsl_rng_alloc(gsl_rng_default);
const size_t nlcurve = 200;
gsl_vector *reg_param = gsl_vector_alloc(nlcurve);
gsl_vector *rho = gsl_vector_alloc(nlcurve);
gsl_vector *eta = gsl_vector_alloc(nlcurve);
size_t rowidx = 0;
double rnorm, snorm, rcond;
double t = 0.0;
double dt = 1.0 / (n - 1.0);
while (rowidx < n)
{
size_t nleft = n - rowidx; /* number of rows left to accumulate */
size_t nr = GSL_MIN(nrows, nleft); /* number of rows in this block */
gsl_matrix_view Xv = gsl_matrix_submatrix(X, 0, 0, nr, p);
gsl_vector_view yv = gsl_vector_subvector(y, 0, nr);
size_t i;
/* build (X,y) block with 'nr' rows */
for (i = 0; i < nr; ++i)
{
gsl_vector_view row = gsl_matrix_row(&Xv.matrix, i);
double fi = func(t);
double ei = gsl_ran_gaussian (r, 0.1 * fi); /* noise */
double yi = fi + ei;
/* construct this row of LS matrix */
build_row(t, &row.vector);
/* set right hand side value with added noise */
gsl_vector_set(&yv.vector, i, yi);
if (print_data && (i % 100 == 0))
printf("%f %f\n", t, yi);
t += dt;
}
/* accumulate (X,y) block into LS system */
gsl_multilarge_linear_accumulate(&Xv.matrix, &yv.vector, w);
rowidx += nr;
}
if (print_data)
printf("\n\n");
/* compute L-curve */
gsl_multilarge_linear_lcurve(reg_param, rho, eta, w);
/* solve large LS system and store solution in c */
gsl_multilarge_linear_solve(lambda, c, &rnorm, &snorm, w);
/* compute reciprocal condition number */
gsl_multilarge_linear_rcond(&rcond, w);
fprintf(stderr, "=== Method %s ===\n", gsl_multilarge_linear_name(w));
fprintf(stderr, "condition number = %e\n", 1.0 / rcond);
fprintf(stderr, "residual norm = %e\n", rnorm);
fprintf(stderr, "solution norm = %e\n", snorm);
/* output L-curve */
{
size_t i;
for (i = 0; i < nlcurve; ++i)
{
printf("%.12e %.12e %.12e\n",
gsl_vector_get(reg_param, i),
gsl_vector_get(rho, i),
gsl_vector_get(eta, i));
}
printf("\n\n");
}
gsl_matrix_free(X);
gsl_vector_free(y);
gsl_multilarge_linear_free(w);
gsl_rng_free(r);
gsl_vector_free(reg_param);
gsl_vector_free(rho);
gsl_vector_free(eta);
return 0;
}
int
main(int argc, char *argv[])
{
const size_t n = 50000; /* number of observations */
const size_t p = 16; /* polynomial order + 1 */
double lambda = 0.0; /* regularization parameter */
gsl_vector *c_tsqr = gsl_vector_alloc(p);
gsl_vector *c_normal = gsl_vector_alloc(p);
if (argc > 1)
lambda = atof(argv[1]);
/* solve system with TSQR method */
solve_system(1, gsl_multilarge_linear_tsqr, lambda, n, p, c_tsqr);
/* solve system with Normal equations method */
solve_system(0, gsl_multilarge_linear_normal, lambda, n, p, c_normal);
/* output solutions */
{
gsl_vector *v = gsl_vector_alloc(p);
double t;
for (t = 0.0; t <= 1.0; t += 0.01)
{
double f_exact = func(t);
double f_tsqr, f_normal;
build_row(t, v);
gsl_blas_ddot(v, c_tsqr, &f_tsqr);
gsl_blas_ddot(v, c_normal, &f_normal);
printf("%f %e %e %e\n", t, f_exact, f_tsqr, f_normal);
}
gsl_vector_free(v);
}
gsl_vector_free(c_tsqr);
gsl_vector_free(c_normal);
return 0;
}
</pre></div>
<hr>
<div class="header">
<p>
Previous: <a href="Fitting-robust-linear-regression-example.html#Fitting-robust-linear-regression-example" accesskey="p" rel="previous">Fitting robust linear regression example</a>, Up: <a href="Fitting-Examples.html#Fitting-Examples" accesskey="u" rel="up">Fitting Examples</a> [<a href="Function-Index.html#Function-Index" title="Index" rel="index">Index</a>]</p>
</div>
</body>
</html>
|