File: Complete-Orthogonal-Decomposition.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 (147 lines) | stat: -rw-r--r-- 9,440 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
<!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: Complete Orthogonal Decomposition</title>

<meta name="description" content="GNU Scientific Library &ndash; Reference Manual: Complete Orthogonal Decomposition">
<meta name="keywords" content="GNU Scientific Library &ndash; Reference Manual: Complete Orthogonal Decomposition">
<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="Linear-Algebra.html#Linear-Algebra" rel="up" title="Linear Algebra">
<link href="Singular-Value-Decomposition.html#Singular-Value-Decomposition" rel="next" title="Singular Value Decomposition">
<link href="QR-Decomposition-with-Column-Pivoting.html#QR-Decomposition-with-Column-Pivoting" rel="previous" title="QR Decomposition with Column Pivoting">
<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="Complete-Orthogonal-Decomposition"></a>
<div class="header">
<p>
Next: <a href="Singular-Value-Decomposition.html#Singular-Value-Decomposition" accesskey="n" rel="next">Singular Value Decomposition</a>, Previous: <a href="QR-Decomposition-with-Column-Pivoting.html#QR-Decomposition-with-Column-Pivoting" accesskey="p" rel="previous">QR Decomposition with Column Pivoting</a>, Up: <a href="Linear-Algebra.html#Linear-Algebra" accesskey="u" rel="up">Linear Algebra</a> &nbsp; [<a href="Function-Index.html#Function-Index" title="Index" rel="index">Index</a>]</p>
</div>
<hr>
<a name="Complete-Orthogonal-Decomposition-1"></a>
<h3 class="section">14.4 Complete Orthogonal Decomposition</h3>

<p>The complete orthogonal decomposition of a <em>M</em>-by-<em>N</em> matrix
<em>A</em> is a generalization of the QR decomposition with column pivoting, given by
</p>
<div class="example">
<pre class="example">A P = Q [ R11 0 ] Z
        [  0  0 ]
</pre></div>

<p>where <em>P</em> is a <em>N</em>-by-<em>N</em> permutation matrix,
<em>Q</em> is <em>M</em>-by-<em>M</em> orthogonal, <em>R_{11}</em> is
<em>r</em>-by-<em>r</em> upper triangular, with <em>r = {\rm rank}(A)</em>,
and <em>Z</em> is <em>N</em>-by-<em>N</em> orthogonal. If <em>A</em>
has full rank, then <em>R_{11} = R</em>, <em>Z = I</em> and this reduces to the
QR decomposition with column pivoting. The advantage of using
the complete orthogonal decomposition for rank deficient matrices
is the ability to compute the minimum norm solution to the linear
least squares problem <em>Ax = b</em>, which is given by
</p>
<div class="example">
<pre class="example">x = P Z^T [ R11^-1 c1 ]
          [    0      ]
</pre></div>

<p>and the vector <em>c_1</em> is the first <em>r</em> elements of <em>Q^T b</em>.
</p>
<dl>
<dt><a name="index-gsl_005flinalg_005fCOD_005fdecomp"></a>Function: <em>int</em> <strong>gsl_linalg_COD_decomp</strong> <em>(gsl_matrix * <var>A</var>, gsl_vector * <var>tau_Q</var>, gsl_vector * <var>tau_Z</var>, gsl_permutation * <var>p</var>, size_t * <var>rank</var>, gsl_vector * <var>work</var>)</em></dt>
<dt><a name="index-gsl_005flinalg_005fCOD_005fdecomp_005fe"></a>Function: <em>int</em> <strong>gsl_linalg_COD_decomp_e</strong> <em>(gsl_matrix * <var>A</var>, gsl_vector * <var>tau_Q</var>, gsl_vector * <var>tau_Z</var>, gsl_permutation * <var>p</var>, double <var>tol</var>, size_t * <var>rank</var>, gsl_vector * <var>work</var>)</em></dt>
<dd><p>These functions factor the <em>M</em>-by-<em>N</em> matrix <var>A</var> into the decomposition <em>A = Q R Z P^T</em>. The rank of <var>A</var>
is computed as the number of diagonal elements of <em>R</em> greater than the tolerance <var>tol</var> and output in <var>rank</var>.
If <var>tol</var> is not specified, a default value is used (see <code>gsl_linalg_QRPT_rank</code>). On output, the permutation
matrix <em>P</em> is stored in <var>p</var>. The matrix <em>R_{11}</em> is stored in the upper <var>rank</var>-by-<var>rank</var> block of <var>A</var>.
The matrices <em>Q</em> and <em>Z</em> are encoded in packed storage in <var>A</var> on output. The vectors <var>tau_Q</var> and <var>tau_Z</var>
contain the Householder scalars corresponding to the matrices <em>Q</em> and <em>Z</em> respectively and must be
of length <em>k = \min(M,N)</em>. The vector <var>work</var> is additional workspace of length <em>N</em>.
</p></dd></dl>

<dl>
<dt><a name="index-gsl_005flinalg_005fCOD_005flssolve"></a>Function: <em>int</em> <strong>gsl_linalg_COD_lssolve</strong> <em>(const gsl_matrix * <var>QRZ</var>, const gsl_vector * <var>tau_Q</var>, const gsl_vector * <var>tau_Z</var>, const gsl_permutation * <var>p</var>, const size_t <var>rank</var>, const gsl_vector * <var>b</var>, gsl_vector * <var>x</var>, gsl_vector * <var>residual</var>)</em></dt>
<dd><p>This function finds the least squares solution to the overdetermined
system <em>A x = b</em> where the matrix <var>A</var> has more rows than
columns.  The least squares solution minimizes the Euclidean norm of the
residual, <em>||b - A x||</em>.  The routine requires as input 
the <em>QRZ</em> decomposition of <em>A</em> into (<var>QRZ</var>, <var>tau_Q</var>, <var>tau_Z</var>, <var>p</var>, <var>rank</var>)
given by <code>gsl_linalg_COD_decomp</code>.  The solution is returned in <var>x</var>.  The
residual is computed as a by-product and stored in <var>residual</var>.
</p></dd></dl>

<dl>
<dt><a name="index-gsl_005flinalg_005fCOD_005funpack"></a>Function: <em>int</em> <strong>gsl_linalg_COD_unpack</strong> <em>(const gsl_matrix * <var>QRZ</var>, const gsl_vector * <var>tau_Q</var>, const gsl_vector * <var>tau_Z</var>, const size_t <var>rank</var>, gsl_matrix * <var>Q</var>, gsl_matrix * <var>R</var>, gsl_matrix * <var>Z</var>)</em></dt>
<dd><p>This function unpacks the encoded <em>QRZ</em> decomposition
(<var>QRZ</var>, <var>tau_Q</var>, <var>tau_Z</var>, <var>rank</var>) into the matrices
<var>Q</var>, <var>R</var>, and <var>Z</var>, where <var>Q</var> is <em>M</em>-by-<em>M</em>,
<var>R</var> is <em>M</em>-by-<em>N</em>, and <var>Z</var> is <em>N</em>-by-<em>N</em>.
</p></dd></dl>

<dl>
<dt><a name="index-gsl_005flinalg_005fCOD_005fmatZ"></a>Function: <em>int</em> <strong>gsl_linalg_COD_matZ</strong> <em>(const gsl_matrix * <var>QRZ</var>, const gsl_vector * <var>tau_Z</var>, const size_t <var>rank</var>, gsl_matrix * <var>A</var>, gsl_vector * <var>work</var>)</em></dt>
<dd><p>This function multiplies the input matrix <var>A</var> on the right by <var>Z</var>,
<em>A' = A Z</em> using the encoded <em>QRZ</em> decomposition
(<var>QRZ</var>, <var>tau_Z</var>, <var>rank</var>). <var>A</var> must have <em>N</em> columns but may
have any number of rows. Additional workspace of length <em>M</em> is provided
in <var>work</var>.
</p></dd></dl>

<hr>
<div class="header">
<p>
Next: <a href="Singular-Value-Decomposition.html#Singular-Value-Decomposition" accesskey="n" rel="next">Singular Value Decomposition</a>, Previous: <a href="QR-Decomposition-with-Column-Pivoting.html#QR-Decomposition-with-Column-Pivoting" accesskey="p" rel="previous">QR Decomposition with Column Pivoting</a>, Up: <a href="Linear-Algebra.html#Linear-Algebra" accesskey="u" rel="up">Linear Algebra</a> &nbsp; [<a href="Function-Index.html#Function-Index" title="Index" rel="index">Index</a>]</p>
</div>



</body>
</html>