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
|
<!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: Large Dense Linear Systems Routines</title>
<meta name="description" content="GNU Scientific Library – Reference Manual: Large Dense Linear Systems Routines">
<meta name="keywords" content="GNU Scientific Library – Reference Manual: Large Dense Linear Systems Routines">
<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="Large-Dense-Linear-Systems.html#Large-Dense-Linear-Systems" rel="up" title="Large Dense Linear Systems">
<link href="Troubleshooting.html#Troubleshooting" rel="next" title="Troubleshooting">
<link href="Large-Dense-Linear-Systems-Solution-Steps.html#Large-Dense-Linear-Systems-Solution-Steps" rel="previous" title="Large Dense Linear Systems Solution Steps">
<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="Large-Dense-Linear-Systems-Routines"></a>
<div class="header">
<p>
Previous: <a href="Large-Dense-Linear-Systems-Solution-Steps.html#Large-Dense-Linear-Systems-Solution-Steps" accesskey="p" rel="previous">Large Dense Linear Systems Solution Steps</a>, Up: <a href="Large-Dense-Linear-Systems.html#Large-Dense-Linear-Systems" accesskey="u" rel="up">Large Dense Linear Systems</a> [<a href="Function-Index.html#Function-Index" title="Index" rel="index">Index</a>]</p>
</div>
<hr>
<a name="Large-Dense-Linear-Least-Squares-Routines"></a>
<h4 class="subsection">38.6.4 Large Dense Linear Least Squares Routines</h4>
<a name="index-large-linear-least-squares_002c-routines"></a>
<dl>
<dt><a name="index-gsl_005fmultilarge_005flinear_005falloc"></a>Function: <em>gsl_multilarge_linear_workspace *</em> <strong>gsl_multilarge_linear_alloc</strong> <em>(const gsl_multilarge_linear_type * <var>T</var>, const size_t <var>p</var>)</em></dt>
<dd><p>This function allocates a workspace for solving large linear least squares
systems. The least squares matrix <em>X</em> has <var>p</var> columns,
but may have any number of rows. The parameter <var>T</var> specifies
the method to be used for solving the large least squares system
and may be selected from the following choices
</p>
<dl>
<dt><a name="index-gsl_005fmultilarge_005flinear_005fnormal"></a>Multilarge type: <strong>gsl_multilarge_linear_normal</strong></dt>
<dd><p>This specifies the normal equations approach for
solving the least squares system. This method is suitable
in cases where performance is critical and it is known that the
least squares matrix <em>X</em> is well conditioned. The size
of this workspace is <em>O(p^2)</em>.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005fmultilarge_005flinear_005ftsqr"></a>Multilarge type: <strong>gsl_multilarge_linear_tsqr</strong></dt>
<dd><p>This specifies the sequential Tall Skinny QR (TSQR) approach for
solving the least squares system. This method is a good
general purpose choice for large systems, but requires about
twice as many operations as the normal equations method for
<em>n >> p</em>. The size of this workspace is <em>O(p^2)</em>.
</p></dd></dl>
</dd></dl>
<dl>
<dt><a name="index-gsl_005fmultilarge_005flinear_005ffree"></a>Function: <em>void</em> <strong>gsl_multilarge_linear_free</strong> <em>(gsl_multilarge_linear_workspace * <var>w</var>)</em></dt>
<dd><p>This function frees the memory associated with the
workspace <var>w</var>.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005fmultilarge_005flinear_005fname"></a>Function: <em>const char *</em> <strong>gsl_multilarge_linear_name</strong> <em>(gsl_multilarge_linear_workspace * <var>w</var>)</em></dt>
<dd><p>This function returns a string pointer to the name
of the multilarge solver.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005fmultilarge_005flinear_005freset"></a>Function: <em>int</em> <strong>gsl_multilarge_linear_reset</strong> <em>(gsl_multilarge_linear_workspace * <var>w</var>)</em></dt>
<dd><p>This function resets the workspace <var>w</var> so
it can begin to accumulate a new least squares
system.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005fmultilarge_005flinear_005fstdform1"></a>Function: <em>int</em> <strong>gsl_multilarge_linear_stdform1</strong> <em>(const gsl_vector * <var>L</var>, const gsl_matrix * <var>X</var>, const gsl_vector * <var>y</var>, gsl_matrix * <var>Xs</var>, gsl_vector * <var>ys</var>, gsl_multilarge_linear_workspace * <var>work</var>)</em></dt>
<dt><a name="index-gsl_005fmultilarge_005flinear_005fwstdform1"></a>Function: <em>int</em> <strong>gsl_multilarge_linear_wstdform1</strong> <em>(const gsl_vector * <var>L</var>, const gsl_matrix * <var>X</var>, const gsl_vector * <var>w</var>, const gsl_vector * <var>y</var>, gsl_matrix * <var>Xs</var>, gsl_vector * <var>ys</var>, gsl_multilarge_linear_workspace * <var>work</var>)</em></dt>
<dd><p>These functions define a regularization matrix
<em>L =</em> diag<em>(l_0,l_1,...,l_{p-1})</em>.
The diagonal matrix element <em>l_i</em> is provided by the
<em>i</em>th element of the input vector <var>L</var>.
The block (<var>X</var>,<var>y</var>) is converted to standard form and
the parameters (<em>\tilde{X}</em>,<em>\tilde{y}</em>) are stored in <var>Xs</var>
and <var>ys</var> on output. <var>Xs</var> and <var>ys</var> have the same dimensions as
<var>X</var> and <var>y</var>. Optional data weights may be supplied in the
vector <var>w</var>. In order to apply this transformation,
<em>L^{-1}</em> must exist and so none of the <em>l_i</em>
may be zero. After the standard form system has been solved,
use <code>gsl_multilarge_linear_genform1</code> to recover the original solution vector.
It is allowed to have <var>X</var> = <var>Xs</var> and <var>y</var> = <var>ys</var> for an in-place transform.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005fmultilarge_005flinear_005fL_005fdecomp"></a>Function: <em>int</em> <strong>gsl_multilarge_linear_L_decomp</strong> <em>(gsl_matrix * <var>L</var>, gsl_vector * <var>tau</var>)</em></dt>
<dd><p>This function calculates the QR decomposition of the <em>m</em>-by-<em>p</em> regularization matrix
<var>L</var>. <var>L</var> must have <em>m \ge p</em>. On output,
the Householder scalars are stored in the vector <var>tau</var> of size <em>p</em>.
These outputs will be used by <code>gsl_multilarge_linear_wstdform2</code> to complete the
transformation to standard form.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005fmultilarge_005flinear_005fstdform2"></a>Function: <em>int</em> <strong>gsl_multilarge_linear_stdform2</strong> <em>(const gsl_matrix * <var>LQR</var>, const gsl_vector * <var>Ltau</var>, const gsl_matrix * <var>X</var>, const gsl_vector * <var>y</var>, gsl_matrix * <var>Xs</var>, gsl_vector * <var>ys</var>, gsl_multilarge_linear_workspace * <var>work</var>)</em></dt>
<dt><a name="index-gsl_005fmultilarge_005flinear_005fwstdform2"></a>Function: <em>int</em> <strong>gsl_multilarge_linear_wstdform2</strong> <em>(const gsl_matrix * <var>LQR</var>, const gsl_vector * <var>Ltau</var>, const gsl_matrix * <var>X</var>, const gsl_vector * <var>w</var>, const gsl_vector * <var>y</var>, gsl_matrix * <var>Xs</var>, gsl_vector * <var>ys</var>, gsl_multilarge_linear_workspace * <var>work</var>)</em></dt>
<dd><p>These functions convert a block of rows (<var>X</var>,<var>y</var>,<var>w</var>) to standard
form (<em>\tilde{X}</em>,<em>\tilde{y}</em>) which are stored in <var>Xs</var> and <var>ys</var>
respectively. <var>X</var>, <var>y</var>, and <var>w</var> must all have the same number of rows.
The <em>m</em>-by-<em>p</em> regularization matrix <var>L</var> is specified by the inputs
<var>LQR</var> and <var>Ltau</var>, which are outputs from <code>gsl_multilarge_linear_L_decomp</code>.
<var>Xs</var> and <var>ys</var> have the same dimensions as <var>X</var> and <var>y</var>. After the
standard form system has been solved, use <code>gsl_multilarge_linear_genform2</code> to
recover the original solution vector. Optional data weights may be supplied in the
vector <var>w</var>, where <em>W =</em> diag(w).
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005fmultilarge_005flinear_005faccumulate"></a>Function: <em>int</em> <strong>gsl_multilarge_linear_accumulate</strong> <em>(gsl_matrix * <var>X</var>, gsl_vector * <var>y</var>, gsl_multilarge_linear_workspace * <var>w</var>)</em></dt>
<dd><p>This function accumulates the standard form block (<em>X,y</em>) into the
current least squares system. <var>X</var> and <var>y</var> have the same number
of rows, which can be arbitrary. <var>X</var> must have <em>p</em> columns.
For the TSQR method, <var>X</var> and <var>y</var> are destroyed on output.
For the normal equations method, they are both unchanged.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005fmultilarge_005flinear_005fsolve"></a>Function: <em>int</em> <strong>gsl_multilarge_linear_solve</strong> <em>(const double <var>lambda</var>, gsl_vector * <var>c</var>, double * <var>rnorm</var>, double * <var>snorm</var>, gsl_multilarge_linear_workspace * <var>w</var>)</em></dt>
<dd><p>After all blocks (<em>X_i,y_i</em>) have been accumulated into
the large least squares system, this function will compute
the solution vector which is stored in <var>c</var> on output.
The regularization parameter <em>\lambda</em> is provided in
<var>lambda</var>. On output, <var>rnorm</var> contains the residual norm
<em>||y - X c||_W</em> and <var>snorm</var> contains the solution
norm <em>||L c||</em>.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005fmultilarge_005flinear_005fgenform1"></a>Function: <em>int</em> <strong>gsl_multilarge_linear_genform1</strong> <em>(const gsl_vector * <var>L</var>, const gsl_vector * <var>cs</var>, gsl_vector * <var>c</var>, gsl_multilarge_linear_workspace * <var>work</var>)</em></dt>
<dd><p>After a regularized system has been solved with
<em>L =</em> diag<em>(\l_0,\l_1,...,\l_{p-1})</em>,
this function backtransforms the standard form solution vector <var>cs</var>
to recover the solution vector of the original problem <var>c</var>. The
diagonal matrix elements <em>l_i</em> are provided in
the vector <var>L</var>. It is allowed to have <var>c</var> = <var>cs</var> for an
in-place transform.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005fmultilarge_005flinear_005fgenform2"></a>Function: <em>int</em> <strong>gsl_multilarge_linear_genform2</strong> <em>(const gsl_matrix * <var>LQR</var>, const gsl_vector * <var>Ltau</var>, const gsl_vector * <var>cs</var>, gsl_vector * <var>c</var>, gsl_multilarge_linear_workspace * <var>work</var>)</em></dt>
<dd><p>After a regularized system has been solved with a regularization matrix <em>L</em>,
specified by (<var>LQR</var>,<var>Ltau</var>), this function backtransforms the standard form solution <var>cs</var>
to recover the solution vector of the original problem, which is stored in <var>c</var>,
of length <em>p</em>.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005fmultilarge_005flinear_005flcurve"></a>Function: <em>int</em> <strong>gsl_multilarge_linear_lcurve</strong> <em>(gsl_vector * <var>reg_param</var>, gsl_vector * <var>rho</var>, gsl_vector * <var>eta</var>, gsl_multilarge_linear_workspace * <var>work</var>)</em></dt>
<dd><p>This function computes the L-curve for a large least squares system
after it has been fully accumulated into the workspace <var>work</var>.
The output vectors <var>reg_param</var>, <var>rho</var>, and <var>eta</var> must all
be the same size, and will contain the regularization parameters
<em>\lambda_i</em>, residual norms <em>||y - X c_i||</em>, and solution
norms <em>|| L c_i ||</em> which compose the L-curve, where <em>c_i</em>
is the regularized solution vector corresponding to <em>\lambda_i</em>.
The user may determine the number of points on the L-curve by
adjusting the size of these input arrays. For the TSQR method,
the regularization parameters <em>\lambda_i</em> are estimated from the
singular values of the triangular <em>R</em> factor. For the normal
equations method, they are estimated from the eigenvalues of the
<em>X^T X</em> matrix.
</p></dd></dl>
<dl>
<dt><a name="index-gsl_005fmultilarge_005flinear_005frcond"></a>Function: <em>int</em> <strong>gsl_multilarge_linear_rcond</strong> <em>(double * <var>rcond</var>, gsl_multilarge_linear_workspace * <var>work</var>)</em></dt>
<dd><p>This function computes the reciprocal condition number, stored in
<var>rcond</var>, of the least squares matrix after it has been accumulated
into the workspace <var>work</var>. For the TSQR algorithm, this is
accomplished by calculating the SVD of the <em>R</em> factor, which
has the same singular values as the matrix <em>X</em>. For the normal
equations method, this is done by computing the eigenvalues of
<em>X^T X</em>, which could be inaccurate for ill-conditioned matrices
<em>X</em>.
</p></dd></dl>
<hr>
<div class="header">
<p>
Previous: <a href="Large-Dense-Linear-Systems-Solution-Steps.html#Large-Dense-Linear-Systems-Solution-Steps" accesskey="p" rel="previous">Large Dense Linear Systems Solution Steps</a>, Up: <a href="Large-Dense-Linear-Systems.html#Large-Dense-Linear-Systems" accesskey="u" rel="up">Large Dense Linear Systems</a> [<a href="Function-Index.html#Function-Index" title="Index" rel="index">Index</a>]</p>
</div>
</body>
</html>
|