File: QR-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 (227 lines) | stat: -rw-r--r-- 13,302 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
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
<!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: QR Decomposition</title>

<meta name="description" content="GNU Scientific Library &ndash; Reference Manual: QR Decomposition">
<meta name="keywords" content="GNU Scientific Library &ndash; Reference Manual: QR 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="QR-Decomposition-with-Column-Pivoting.html#QR-Decomposition-with-Column-Pivoting" rel="next" title="QR Decomposition with Column Pivoting">
<link href="LU-Decomposition.html#LU-Decomposition" rel="previous" title="LU Decomposition">
<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="QR-Decomposition"></a>
<div class="header">
<p>
Next: <a href="QR-Decomposition-with-Column-Pivoting.html#QR-Decomposition-with-Column-Pivoting" accesskey="n" rel="next">QR Decomposition with Column Pivoting</a>, Previous: <a href="LU-Decomposition.html#LU-Decomposition" accesskey="p" rel="previous">LU Decomposition</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="QR-Decomposition-1"></a>
<h3 class="section">14.2 QR Decomposition</h3>
<a name="index-QR-decomposition"></a>

<p>A general rectangular <em>M</em>-by-<em>N</em> matrix <em>A</em> has a
<em>QR</em> decomposition into the product of an orthogonal
<em>M</em>-by-<em>M</em> square matrix <em>Q</em> (where <em>Q^T Q = I</em>) and
an <em>M</em>-by-<em>N</em> right-triangular matrix <em>R</em>,
</p>
<div class="example">
<pre class="example">A = Q R
</pre></div>

<p>This decomposition can be used to convert the linear system <em>A x =
b</em> into the triangular system <em>R x = Q^T b</em>, which can be solved by
back-substitution. Another use of the <em>QR</em> decomposition is to
compute an orthonormal basis for a set of vectors. The first <em>N</em>
columns of <em>Q</em> form an orthonormal basis for the range of <em>A</em>,
<em>ran(A)</em>, when <em>A</em> has full column rank.
</p>
<dl>
<dt><a name="index-gsl_005flinalg_005fQR_005fdecomp"></a>Function: <em>int</em> <strong>gsl_linalg_QR_decomp</strong> <em>(gsl_matrix * <var>A</var>, gsl_vector * <var>tau</var>)</em></dt>
<dd><p>This function factorizes the <em>M</em>-by-<em>N</em> matrix <var>A</var> into
the <em>QR</em> decomposition <em>A = Q R</em>.  On output the diagonal and
upper triangular part of the input matrix contain the matrix
<em>R</em>. The vector <var>tau</var> and the columns of the lower triangular
part of the matrix <var>A</var> contain the Householder coefficients and
Householder vectors which encode the orthogonal matrix <var>Q</var>.  The
vector <var>tau</var> must be of length <em>k=\min(M,N)</em>. The matrix
<em>Q</em> is related to these components by, <em>Q = Q_k ... Q_2 Q_1</em>
where <em>Q_i = I - \tau_i v_i v_i^T</em> and <em>v_i</em> is the
Householder vector <em>v_i =
(0,...,1,A(i+1,i),A(i+2,i),...,A(m,i))</em>. This is the same storage scheme
as used by <small>LAPACK</small>.
</p>
<p>The algorithm used to perform the decomposition is Householder QR (Golub
&amp; Van Loan, <cite>Matrix Computations</cite>, Algorithm 5.2.1).
</p></dd></dl>

<dl>
<dt><a name="index-gsl_005flinalg_005fQR_005fsolve"></a>Function: <em>int</em> <strong>gsl_linalg_QR_solve</strong> <em>(const gsl_matrix * <var>QR</var>, const gsl_vector * <var>tau</var>, const gsl_vector * <var>b</var>, gsl_vector * <var>x</var>)</em></dt>
<dd><p>This function solves the square system <em>A x = b</em> using the <em>QR</em>
decomposition of <em>A</em> held in (<var>QR</var>, <var>tau</var>) which must 
have been computed previously with <code>gsl_linalg_QR_decomp</code>. 
The least-squares solution for 
rectangular systems can be found using <code>gsl_linalg_QR_lssolve</code>.
</p></dd></dl>

<dl>
<dt><a name="index-gsl_005flinalg_005fQR_005fsvx"></a>Function: <em>int</em> <strong>gsl_linalg_QR_svx</strong> <em>(const gsl_matrix * <var>QR</var>, const gsl_vector * <var>tau</var>, gsl_vector * <var>x</var>)</em></dt>
<dd><p>This function solves the square system <em>A x = b</em> in-place using
the <em>QR</em> decomposition of <em>A</em> held in (<var>QR</var>,<var>tau</var>)
which must have been computed previously by
<code>gsl_linalg_QR_decomp</code>.  On input <var>x</var> should contain the
right-hand side <em>b</em>, which is replaced by the solution on output.
</p></dd></dl>

<dl>
<dt><a name="index-gsl_005flinalg_005fQR_005flssolve"></a>Function: <em>int</em> <strong>gsl_linalg_QR_lssolve</strong> <em>(const gsl_matrix * <var>QR</var>, const gsl_vector * <var>tau</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>||Ax - b||</em>.The routine requires as input 
the <em>QR</em> decomposition
of <em>A</em> into (<var>QR</var>, <var>tau</var>) given by
<code>gsl_linalg_QR_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_005fQR_005fQTvec"></a>Function: <em>int</em> <strong>gsl_linalg_QR_QTvec</strong> <em>(const gsl_matrix * <var>QR</var>, const gsl_vector * <var>tau</var>, gsl_vector * <var>v</var>)</em></dt>
<dd><p>This function applies the matrix <em>Q^T</em> encoded in the decomposition
(<var>QR</var>,<var>tau</var>) to the vector <var>v</var>, storing the result <em>Q^T
v</em> in <var>v</var>.  The matrix multiplication is carried out directly using
the encoding of the Householder vectors without needing to form the full
matrix <em>Q^T</em>.
</p></dd></dl>

<dl>
<dt><a name="index-gsl_005flinalg_005fQR_005fQvec"></a>Function: <em>int</em> <strong>gsl_linalg_QR_Qvec</strong> <em>(const gsl_matrix * <var>QR</var>, const gsl_vector * <var>tau</var>, gsl_vector * <var>v</var>)</em></dt>
<dd><p>This function applies the matrix <em>Q</em> encoded in the decomposition
(<var>QR</var>,<var>tau</var>) to the vector <var>v</var>, storing the result <em>Q
v</em> in <var>v</var>.  The matrix multiplication is carried out directly using
the encoding of the Householder vectors without needing to form the full
matrix <em>Q</em>.
</p></dd></dl>

<dl>
<dt><a name="index-gsl_005flinalg_005fQR_005fQTmat"></a>Function: <em>int</em> <strong>gsl_linalg_QR_QTmat</strong> <em>(const gsl_matrix * <var>QR</var>, const gsl_vector * <var>tau</var>, gsl_matrix * <var>A</var>)</em></dt>
<dd><p>This function applies the matrix <em>Q^T</em> encoded in the decomposition
(<var>QR</var>,<var>tau</var>) to the matrix <var>A</var>, storing the result <em>Q^T
A</em> in <var>A</var>.  The matrix multiplication is carried out directly using
the encoding of the Householder vectors without needing to form the full
matrix <em>Q^T</em>.
</p></dd></dl>

<dl>
<dt><a name="index-gsl_005flinalg_005fQR_005fRsolve"></a>Function: <em>int</em> <strong>gsl_linalg_QR_Rsolve</strong> <em>(const gsl_matrix * <var>QR</var>, const gsl_vector * <var>b</var>, gsl_vector * <var>x</var>)</em></dt>
<dd><p>This function solves the triangular system <em>R x = b</em> for
<var>x</var>. It may be useful if the product <em>b' = Q^T b</em> has already
been computed using <code>gsl_linalg_QR_QTvec</code>.
</p></dd></dl>

<dl>
<dt><a name="index-gsl_005flinalg_005fQR_005fRsvx"></a>Function: <em>int</em> <strong>gsl_linalg_QR_Rsvx</strong> <em>(const gsl_matrix * <var>QR</var>, gsl_vector * <var>x</var>)</em></dt>
<dd><p>This function solves the triangular system <em>R x = b</em> for <var>x</var>
in-place. On input <var>x</var> should contain the right-hand side <em>b</em>
and is replaced by the solution on output. This function may be useful if
the product <em>b' = Q^T b</em> has already been computed using
<code>gsl_linalg_QR_QTvec</code>.
</p></dd></dl>

<dl>
<dt><a name="index-gsl_005flinalg_005fQR_005funpack"></a>Function: <em>int</em> <strong>gsl_linalg_QR_unpack</strong> <em>(const gsl_matrix * <var>QR</var>, const gsl_vector * <var>tau</var>, gsl_matrix * <var>Q</var>, gsl_matrix * <var>R</var>)</em></dt>
<dd><p>This function unpacks the encoded <em>QR</em> decomposition
(<var>QR</var>,<var>tau</var>) into the matrices <var>Q</var> and <var>R</var>, where
<var>Q</var> is <em>M</em>-by-<em>M</em> and <var>R</var> is <em>M</em>-by-<em>N</em>. 
</p></dd></dl>

<dl>
<dt><a name="index-gsl_005flinalg_005fQR_005fQRsolve"></a>Function: <em>int</em> <strong>gsl_linalg_QR_QRsolve</strong> <em>(gsl_matrix * <var>Q</var>, gsl_matrix * <var>R</var>, const gsl_vector * <var>b</var>, gsl_vector * <var>x</var>)</em></dt>
<dd><p>This function solves the system <em>R x = Q^T b</em> for <var>x</var>. It can
be used when the <em>QR</em> decomposition of a matrix is available in
unpacked form as (<var>Q</var>, <var>R</var>).
</p></dd></dl>

<dl>
<dt><a name="index-gsl_005flinalg_005fQR_005fupdate"></a>Function: <em>int</em> <strong>gsl_linalg_QR_update</strong> <em>(gsl_matrix * <var>Q</var>, gsl_matrix * <var>R</var>, gsl_vector * <var>w</var>, const gsl_vector * <var>v</var>)</em></dt>
<dd><p>This function performs a rank-1 update <em>w v^T</em> of the <em>QR</em>
decomposition (<var>Q</var>, <var>R</var>). The update is given by <em>Q'R' = Q
(R + w v^T)</em> where the output matrices <em>Q'</em> and <em>R'</em> are also
orthogonal and right triangular. Note that <var>w</var> is destroyed by the
update.
</p></dd></dl>

<dl>
<dt><a name="index-gsl_005flinalg_005fR_005fsolve"></a>Function: <em>int</em> <strong>gsl_linalg_R_solve</strong> <em>(const gsl_matrix * <var>R</var>, const gsl_vector * <var>b</var>, gsl_vector * <var>x</var>)</em></dt>
<dd><p>This function solves the triangular system <em>R x = b</em> for the
<em>N</em>-by-<em>N</em> matrix <var>R</var>.
</p></dd></dl>

<dl>
<dt><a name="index-gsl_005flinalg_005fR_005fsvx"></a>Function: <em>int</em> <strong>gsl_linalg_R_svx</strong> <em>(const gsl_matrix * <var>R</var>, gsl_vector * <var>x</var>)</em></dt>
<dd><p>This function solves the triangular system <em>R x = b</em> in-place. On
input <var>x</var> should contain the right-hand side <em>b</em>, which is
replaced by the solution on output.
</p></dd></dl>

<hr>
<div class="header">
<p>
Next: <a href="QR-Decomposition-with-Column-Pivoting.html#QR-Decomposition-with-Column-Pivoting" accesskey="n" rel="next">QR Decomposition with Column Pivoting</a>, Previous: <a href="LU-Decomposition.html#LU-Decomposition" accesskey="p" rel="previous">LU Decomposition</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>