Mercurial > hg > octave-nkf
annotate src/DLD-FUNCTIONS/qr.cc @ 7553:56be6f31dd4e
implementation of QR factorization updating
author | Jaroslav Hajek <highegg@gmail.com> |
---|---|
date | Tue, 04 Mar 2008 21:47:11 -0500 |
parents | f5005d9510f4 |
children | 40574114c514 |
rev | line source |
---|---|
2928 | 1 /* |
2 | |
7017 | 3 Copyright (C) 1996, 1997, 1999, 2000, 2005, 2006, 2007 John W. Eaton |
2928 | 4 |
5 This file is part of Octave. | |
6 | |
7 Octave is free software; you can redistribute it and/or modify it | |
8 under the terms of the GNU General Public License as published by the | |
7016 | 9 Free Software Foundation; either version 3 of the License, or (at your |
10 option) any later version. | |
2928 | 11 |
12 Octave is distributed in the hope that it will be useful, but WITHOUT | |
13 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or | |
14 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License | |
15 for more details. | |
16 | |
17 You should have received a copy of the GNU General Public License | |
7016 | 18 along with Octave; see the file COPYING. If not, see |
19 <http://www.gnu.org/licenses/>. | |
2928 | 20 |
21 */ | |
22 | |
7553
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
23 // The qrupdate, qrinsert, and qrdelete functions were written by |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
24 // Jaroslav Hajek <highegg@gmail.com>, Copyright (C) 2008 VZLU |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
25 // Prague, a.s., Czech Republic. |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
26 |
2928 | 27 #ifdef HAVE_CONFIG_H |
28 #include <config.h> | |
29 #endif | |
30 | |
31 #include "CmplxQR.h" | |
32 #include "CmplxQRP.h" | |
33 #include "dbleQR.h" | |
34 #include "dbleQRP.h" | |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
35 #include "SparseQR.h" |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
36 #include "SparseCmplxQR.h" |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
37 |
2928 | 38 |
39 #include "defun-dld.h" | |
40 #include "error.h" | |
41 #include "gripes.h" | |
42 #include "oct-obj.h" | |
43 #include "utils.h" | |
44 | |
3372 | 45 // [Q, R] = qr (X): form Q unitary and R upper triangular such |
46 // that Q * R = X | |
47 // | |
48 // [Q, R] = qr (X, 0): form the economy decomposition such that if X is | |
49 // m by n then only the first n columns of Q are | |
50 // computed. | |
51 // | |
52 // [Q, R, P] = qr (X): form QRP factorization of X where | |
53 // P is a permutation matrix such that | |
54 // A * P = Q * R | |
55 // | |
56 // [Q, R, P] = qr (X, 0): form the economy decomposition with | |
57 // permutation vector P such that Q * R = X (:, P) | |
58 // | |
59 // qr (X) alone returns the output of the LAPACK routine dgeqrf, such | |
60 // that R = triu (qr (X)) | |
61 | |
2928 | 62 DEFUN_DLD (qr, args, nargout, |
3548 | 63 "-*- texinfo -*-\n\ |
3372 | 64 @deftypefn {Loadable Function} {[@var{q}, @var{r}, @var{p}] =} qr (@var{a})\n\ |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
65 @deftypefnx {Loadable Function} {[@var{q}, @var{r}, @var{p}] =} qr (@var{a}, '0')\n\ |
3372 | 66 @cindex QR factorization\n\ |
67 Compute the QR factorization of @var{a}, using standard @sc{Lapack}\n\ | |
68 subroutines. For example, given the matrix @code{a = [1, 2; 3, 4]},\n\ | |
2928 | 69 \n\ |
3372 | 70 @example\n\ |
71 [q, r] = qr (a)\n\ | |
72 @end example\n\ | |
73 \n\ | |
74 @noindent\n\ | |
75 returns\n\ | |
76 \n\ | |
77 @example\n\ | |
78 q =\n\ | |
79 \n\ | |
80 -0.31623 -0.94868\n\ | |
81 -0.94868 0.31623\n\ | |
82 \n\ | |
83 r =\n\ | |
84 \n\ | |
85 -3.16228 -4.42719\n\ | |
86 0.00000 -0.63246\n\ | |
87 @end example\n\ | |
88 \n\ | |
89 The @code{qr} factorization has applications in the solution of least\n\ | |
90 squares problems\n\ | |
91 @iftex\n\ | |
92 @tex\n\ | |
93 $$\n\ | |
94 \\min_x \\left\\Vert A x - b \\right\\Vert_2\n\ | |
95 $$\n\ | |
96 @end tex\n\ | |
97 @end iftex\n\ | |
98 @ifinfo\n\ | |
2928 | 99 \n\ |
3372 | 100 @example\n\ |
101 @code{min norm(A x - b)}\n\ | |
102 @end example\n\ | |
103 \n\ | |
104 @end ifinfo\n\ | |
105 for overdetermined systems of equations (i.e.,\n\ | |
106 @iftex\n\ | |
107 @tex\n\ | |
108 $A$\n\ | |
109 @end tex\n\ | |
110 @end iftex\n\ | |
111 @ifinfo\n\ | |
112 @code{a}\n\ | |
113 @end ifinfo\n\ | |
114 is a tall, thin matrix). The QR factorization is\n\ | |
115 @iftex\n\ | |
116 @tex\n\ | |
117 $QR = A$ where $Q$ is an orthogonal matrix and $R$ is upper triangular.\n\ | |
118 @end tex\n\ | |
119 @end iftex\n\ | |
120 @ifinfo\n\ | |
121 @code{q * r = a} where @code{q} is an orthogonal matrix and @code{r} is\n\ | |
122 upper triangular.\n\ | |
123 @end ifinfo\n\ | |
2928 | 124 \n\ |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
125 If given a second argument of '0', @code{qr} returns an economy-sized\n\ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
126 QR factorization, omitting zero rows of @var{R} and the corresponding\n\ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
127 columns of @var{Q}.\n\ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
128 \n\ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
129 If the matrix @var{a} is full, the permuted QR factorization\n\ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
130 @code{[@var{q}, @var{r}, @var{p}] = qr (@var{a})} forms the QR factorization\n\ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
131 such that the diagonal entries of @code{r} are decreasing in magnitude\n\ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
132 order. For example,given the matrix @code{a = [1, 2; 3, 4]},\n\ |
3372 | 133 \n\ |
134 @example\n\ | |
3600 | 135 [q, r, p] = qr(a)\n\ |
3372 | 136 @end example\n\ |
137 \n\ | |
138 @noindent\n\ | |
139 returns\n\ | |
140 \n\ | |
141 @example\n\ | |
142 q = \n\ | |
2928 | 143 \n\ |
3372 | 144 -0.44721 -0.89443\n\ |
145 -0.89443 0.44721\n\ | |
146 \n\ | |
147 r =\n\ | |
148 \n\ | |
149 -4.47214 -3.13050\n\ | |
150 0.00000 0.44721\n\ | |
151 \n\ | |
152 p =\n\ | |
153 \n\ | |
154 0 1\n\ | |
155 1 0\n\ | |
156 @end example\n\ | |
157 \n\ | |
158 The permuted @code{qr} factorization @code{[q, r, p] = qr (a)}\n\ | |
159 factorization allows the construction of an orthogonal basis of\n\ | |
160 @code{span (a)}.\n\ | |
7491 | 161 \n\ |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
162 If the matrix @var{a} is sparse, then compute the sparse QR factorization\n\ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
163 of @var{a}, using @sc{CSparse}. As the matrix @var{Q} is in general a full\n\ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
164 matrix, this function returns the @var{Q}-less factorization @var{r} of\n\ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
165 @var{a}, such that @code{@var{r} = chol (@var{a}' * @var{a})}.\n\ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
166 \n\ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
167 If the final argument is the scalar @code{0} and the number of rows is\n\ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
168 larger than the number of columns, then an economy factorization is\n\ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
169 returned. That is @var{r} will have only @code{size (@var{a},1)} rows.\n\ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
170 \n\ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
171 If an additional matrix @var{b} is supplied, then @code{qr} returns\n\ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
172 @var{c}, where @code{@var{c} = @var{q}' * @var{b}}. This allows the\n\ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
173 least squares approximation of @code{@var{a} \\ @var{b}} to be calculated\n\ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
174 as\n\ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
175 \n\ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
176 @example\n\ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
177 [@var{c},@var{r}] = spqr (@var{a},@var{b})\n\ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
178 @var{x} = @var{r} \\ @var{c}\n\ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
179 @end example\n\ |
3372 | 180 @end deftypefn") |
2928 | 181 { |
182 octave_value_list retval; | |
183 | |
184 int nargin = args.length (); | |
185 | |
7553
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
186 if (nargin < 1 || nargin > (args(0).is_sparse_type() ? 3 : 2)) |
2928 | 187 { |
5823 | 188 print_usage (); |
2928 | 189 return retval; |
190 } | |
191 | |
192 octave_value arg = args(0); | |
193 | |
194 int arg_is_empty = empty_arg ("qr", arg.rows (), arg.columns ()); | |
195 | |
196 if (arg_is_empty < 0) | |
197 return retval; | |
198 else if (arg_is_empty > 0) | |
199 return octave_value_list (3, Matrix ()); | |
200 | |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
201 if (arg.is_sparse_type ()) |
2928 | 202 { |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
203 bool economy = false; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
204 bool is_cmplx = false; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
205 int have_b = 0; |
3068 | 206 |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
207 if (arg.is_complex_type ()) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
208 is_cmplx = true; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
209 if (nargin > 1) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
210 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
211 have_b = 1; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
212 if (args(nargin-1).is_scalar_type ()) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
213 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
214 int val = args(nargin-1).int_value (); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
215 if (val == 0) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
216 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
217 economy = true; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
218 have_b = (nargin > 2 ? 2 : 0); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
219 } |
2928 | 220 } |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
221 if (have_b > 0 && args(have_b).is_complex_type ()) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
222 is_cmplx = true; |
2928 | 223 } |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
224 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
225 if (!error_state) |
2928 | 226 { |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
227 if (have_b && nargout < 2) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
228 error ("qr: incorrect number of output arguments"); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
229 else if (is_cmplx) |
2928 | 230 { |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
231 SparseComplexQR q (arg.sparse_complex_matrix_value ()); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
232 if (!error_state) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
233 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
234 if (have_b > 0) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
235 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
236 retval(1) = q.R (economy); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
237 retval(0) = q.C (args(have_b).complex_matrix_value ()); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
238 if (arg.rows() < arg.columns()) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
239 warning ("qr: non minimum norm solution for under-determined problem"); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
240 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
241 else if (nargout > 1) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
242 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
243 retval(1) = q.R (economy); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
244 retval(0) = q.Q (); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
245 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
246 else |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
247 retval(0) = q.R (economy); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
248 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
249 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
250 else |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
251 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
252 SparseQR q (arg.sparse_matrix_value ()); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
253 if (!error_state) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
254 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
255 if (have_b > 0) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
256 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
257 retval(1) = q.R (economy); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
258 retval(0) = q.C (args(have_b).matrix_value ()); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
259 if (args(0).rows() < args(0).columns()) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
260 warning ("qr: non minimum norm solution for under-determined problem"); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
261 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
262 else if (nargout > 1) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
263 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
264 retval(1) = q.R (economy); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
265 retval(0) = q.Q (); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
266 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
267 else |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
268 retval(0) = q.R (economy); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
269 } |
2928 | 270 } |
271 } | |
272 } | |
273 else | |
274 { | |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
275 QR::type type = (nargout == 0 || nargout == 1) ? QR::raw |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
276 : (nargin == 2 ? QR::economy : QR::std); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
277 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
278 if (arg.is_real_type ()) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
279 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
280 Matrix m = arg.matrix_value (); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
281 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
282 if (! error_state) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
283 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
284 switch (nargout) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
285 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
286 case 0: |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
287 case 1: |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
288 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
289 QR fact (m, type); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
290 retval(0) = fact.R (); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
291 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
292 break; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
293 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
294 case 2: |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
295 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
296 QR fact (m, type); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
297 retval(1) = fact.R (); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
298 retval(0) = fact.Q (); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
299 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
300 break; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
301 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
302 default: |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
303 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
304 QRP fact (m, type); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
305 retval(2) = fact.P (); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
306 retval(1) = fact.R (); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
307 retval(0) = fact.Q (); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
308 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
309 break; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
310 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
311 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
312 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
313 else if (arg.is_complex_type ()) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
314 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
315 ComplexMatrix m = arg.complex_matrix_value (); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
316 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
317 if (! error_state) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
318 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
319 switch (nargout) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
320 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
321 case 0: |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
322 case 1: |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
323 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
324 ComplexQR fact (m, type); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
325 retval(0) = fact.R (); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
326 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
327 break; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
328 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
329 case 2: |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
330 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
331 ComplexQR fact (m, type); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
332 retval(1) = fact.R (); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
333 retval(0) = fact.Q (); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
334 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
335 break; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
336 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
337 default: |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
338 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
339 ComplexQRP fact (m, type); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
340 retval(2) = fact.P (); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
341 retval(1) = fact.R (); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
342 retval(0) = fact.Q (); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
343 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
344 break; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
345 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
346 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
347 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
348 else |
7553
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
349 gripe_wrong_type_arg ("qr", arg); |
2928 | 350 } |
351 | |
352 return retval; | |
353 } | |
354 | |
355 /* | |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
356 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
357 The deactivated tests below can't be tested till rectangular back-subs is |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
358 implemented for sparse matrices. |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
359 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
360 %!testif HAVE_CXSPARSE |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
361 %! n = 20; d= 0.2; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
362 %! a = sprandn(n,n,d)+speye(n,n); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
363 %! r = qr(a); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
364 %! assert(r'*r,a'*a,1e-10) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
365 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
366 %!testif HAVE_CXSPARSE |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
367 %! n = 20; d= 0.2; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
368 %! a = sprandn(n,n,d)+speye(n,n); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
369 %! q = symamd(a); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
370 %! a = a(q,q); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
371 %! r = qr(a); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
372 %! assert(r'*r,a'*a,1e-10) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
373 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
374 %!testif HAVE_CXSPARSE |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
375 %! n = 20; d= 0.2; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
376 %! a = sprandn(n,n,d)+speye(n,n); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
377 %! [c,r] = qr(a,ones(n,1)); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
378 %! assert (r\c,full(a)\ones(n,1),10e-10) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
379 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
380 %!testif HAVE_CXSPARSE |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
381 %! n = 20; d= 0.2; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
382 %! a = sprandn(n,n,d)+speye(n,n); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
383 %! b = randn(n,2); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
384 %! [c,r] = qr(a,b); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
385 %! assert (r\c,full(a)\b,10e-10) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
386 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
387 %% Test under-determined systems!! |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
388 %!#testif HAVE_CXSPARSE |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
389 %! n = 20; d= 0.2; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
390 %! a = sprandn(n,n+1,d)+speye(n,n+1); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
391 %! b = randn(n,2); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
392 %! [c,r] = qr(a,b); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
393 %! assert (r\c,full(a)\b,10e-10) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
394 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
395 %!testif HAVE_CXSPARSE |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
396 %! n = 20; d= 0.2; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
397 %! a = 1i*sprandn(n,n,d)+speye(n,n); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
398 %! r = qr(a); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
399 %! assert(r'*r,a'*a,1e-10) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
400 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
401 %!testif HAVE_CXSPARSE |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
402 %! n = 20; d= 0.2; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
403 %! a = 1i*sprandn(n,n,d)+speye(n,n); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
404 %! q = symamd(a); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
405 %! a = a(q,q); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
406 %! r = qr(a); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
407 %! assert(r'*r,a'*a,1e-10) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
408 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
409 %!testif HAVE_CXSPARSE |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
410 %! n = 20; d= 0.2; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
411 %! a = 1i*sprandn(n,n,d)+speye(n,n); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
412 %! [c,r] = qr(a,ones(n,1)); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
413 %! assert (r\c,full(a)\ones(n,1),10e-10) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
414 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
415 %!testif HAVE_CXSPARSE |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
416 %! n = 20; d= 0.2; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
417 %! a = 1i*sprandn(n,n,d)+speye(n,n); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
418 %! b = randn(n,2); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
419 %! [c,r] = qr(a,b); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
420 %! assert (r\c,full(a)\b,10e-10) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
421 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
422 %% Test under-determined systems!! |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
423 %!#testif HAVE_CXSPARSE |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
424 %! n = 20; d= 0.2; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
425 %! a = 1i*sprandn(n,n+1,d)+speye(n,n+1); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
426 %! b = randn(n,2); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
427 %! [c,r] = qr(a,b); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
428 %! assert (r\c,full(a)\b,10e-10) |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
429 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
430 %!error qr(sprandn(10,10,0.2),ones(10,1)); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
431 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
432 */ |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
433 |
7553
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
434 DEFUN_DLD (qrupdate, args, , |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
435 "-*- texinfo -*-\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
436 @deftypefn{Loadable Function} {[@var{Q1}, @var{R1}]} = qrupdate (@var{Q}, @var{R}, @var{u}, @var{v})\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
437 Given a QR@tie{}factorization of a real or complex matrix\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
438 @w{@var{A} = @var{Q}*@var{R}}, @var{Q}@tie{}unitary and\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
439 @var{R}@tie{}upper trapezoidal, return the QR@tie{}factorization\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
440 of @w{@var{A} + @var{u}*@var{v}'}, where @var{u} and @var{v} are\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
441 column vectors (rank-1 update).\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
442 \n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
443 If the matrix @var{Q} is not square, the matrix @var{A} is updated by\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
444 Q*Q'*u*v' instead of u*v'.\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
445 @seealso{qr, qrinsert, qrdelete}\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
446 @end deftypefn") |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
447 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
448 octave_idx_type nargin = args.length (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
449 octave_value_list retval; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
450 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
451 octave_value argQ,argR,argu,argv; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
452 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
453 if (nargin == 4 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
454 && (argQ = args(0),argQ.is_matrix_type ()) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
455 && (argR = args(1),argR.is_matrix_type ()) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
456 && (argu = args(2),argu.is_matrix_type ()) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
457 && (argv = args(3),argv.is_matrix_type ())) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
458 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
459 octave_idx_type m = argQ.rows (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
460 octave_idx_type n = argR.columns (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
461 octave_idx_type k = argQ.columns (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
462 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
463 if (argR.rows () == k |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
464 && argu.rows () == m && argu.columns () == 1 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
465 && argv.rows () == n && argv.columns () == 1) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
466 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
467 if (argQ.is_real_matrix () |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
468 && argR.is_real_matrix () |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
469 && argu.is_real_matrix () |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
470 && argv.is_real_matrix ()) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
471 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
472 // all real case |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
473 Matrix Q = argQ.matrix_value (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
474 Matrix R = argR.matrix_value (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
475 Matrix u = argu.matrix_value (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
476 Matrix v = argv.matrix_value (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
477 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
478 QR fact (Q, R); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
479 fact.update (u, v); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
480 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
481 retval(1) = fact.R (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
482 retval(0) = fact.Q (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
483 } |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
484 else |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
485 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
486 // complex case |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
487 ComplexMatrix Q = argQ.complex_matrix_value (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
488 ComplexMatrix R = argR.complex_matrix_value (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
489 ComplexMatrix u = argu.complex_matrix_value (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
490 ComplexMatrix v = argv.complex_matrix_value (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
491 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
492 ComplexQR fact (Q, R); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
493 fact.update (u, v); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
494 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
495 retval(1) = fact.R (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
496 retval(0) = fact.Q (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
497 } |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
498 } |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
499 else |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
500 error ("qrupdate: dimensions mismatch"); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
501 } |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
502 else |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
503 print_usage (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
504 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
505 return retval; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
506 } |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
507 /* |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
508 %!test |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
509 %! A = [0.091364 0.613038 0.999083; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
510 %! 0.594638 0.425302 0.603537; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
511 %! 0.383594 0.291238 0.085574; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
512 %! 0.265712 0.268003 0.238409; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
513 %! 0.669966 0.743851 0.445057 ]; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
514 %! |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
515 %! u = [0.85082; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
516 %! 0.76426; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
517 %! 0.42883; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
518 %! 0.53010; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
519 %! 0.80683 ]; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
520 %! |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
521 %! v = [0.98810; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
522 %! 0.24295; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
523 %! 0.43167 ]; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
524 %! |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
525 %! [Q,R] = qr(A); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
526 %! [Q,R] = qrupdate(Q,R,u,v); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
527 %! assert(norm(vec(Q'*Q - eye(5)),Inf) < 1e1*eps) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
528 %! assert(norm(vec(triu(R)-R),Inf) == 0) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
529 %! assert(norm(vec(Q*R - A - u*v'),Inf) < norm(A)*1e1*eps) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
530 %! |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
531 %!test |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
532 %! A = [0.620405 + 0.956953i 0.480013 + 0.048806i 0.402627 + 0.338171i; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
533 %! 0.589077 + 0.658457i 0.013205 + 0.279323i 0.229284 + 0.721929i; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
534 %! 0.092758 + 0.345687i 0.928679 + 0.241052i 0.764536 + 0.832406i; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
535 %! 0.912098 + 0.721024i 0.049018 + 0.269452i 0.730029 + 0.796517i; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
536 %! 0.112849 + 0.603871i 0.486352 + 0.142337i 0.355646 + 0.151496i ]; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
537 %! |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
538 %! u = [0.20351 + 0.05401i; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
539 %! 0.13141 + 0.43708i; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
540 %! 0.29808 + 0.08789i; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
541 %! 0.69821 + 0.38844i; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
542 %! 0.74871 + 0.25821i ]; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
543 %! |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
544 %! v = [0.85839 + 0.29468i; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
545 %! 0.20820 + 0.93090i; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
546 %! 0.86184 + 0.34689i ]; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
547 %! |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
548 %! [Q,R] = qr(A); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
549 %! [Q,R] = qrupdate(Q,R,u,v); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
550 %! assert(norm(vec(Q'*Q - eye(5)),Inf) < 1e1*eps) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
551 %! assert(norm(vec(triu(R)-R),Inf) == 0) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
552 %! assert(norm(vec(Q*R - A - u*v'),Inf) < norm(A)*1e1*eps) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
553 */ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
554 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
555 DEFUN_DLD (qrinsert, args, , |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
556 "-*- texinfo -*-\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
557 @deftypefn{Loadable Function} {[@var{Q1}, @var{R1}]} = qrinsert (@var{Q}, @var{R}, @var{j}, @var{x}, @var{orient})\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
558 Given a QR@tie{}factorization of a real or complex matrix\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
559 @w{@var{A} = @var{Q}*@var{R}}, @var{Q}@tie{}unitary and\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
560 @var{R}@tie{}upper trapezoidal, return the QR@tie{}factorization of\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
561 @w{[A(:,1:j-1) x A(:,j:n)]}, where @var{u} is a column vector to be\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
562 inserted into @var{A} (if @var{orient} is @code{\"col\"}), or the\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
563 QR@tie{}factorization of @w{[A(1:j-1,:);x;A(:,j:n)]}, where @var{x}\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
564 is a row vector to be inserted into @var{A} (if @var{orient} is\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
565 @code{\"row\"}).\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
566 \n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
567 The default value of @var{orient} is @code{\"col\"}.\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
568 \n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
569 If @var{orient} is @code{\"col\"} and the matrix @var{Q} is not square,\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
570 then what gets inserted is the projection of @var{u} onto the space\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
571 spanned by columns of @var{Q}, i.e. Q*Q'*u.\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
572 \n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
573 If @var{orient} is @code{\"row\"}, @var{Q} must be square.\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
574 @seealso{qr, qrupdate, qrdelete}\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
575 @end deftypefn") |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
576 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
577 octave_idx_type nargin = args.length (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
578 octave_value_list retval; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
579 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
580 octave_value argQ,argR,argj,argx,argor; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
581 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
582 if ((nargin == 4 || nargin == 5) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
583 && (argQ = args(0), argQ.is_matrix_type ()) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
584 && (argR = args(1), argR.is_matrix_type ()) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
585 && (argj = args(2), argj.is_scalar_type ()) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
586 && (argx = args(3), argx.is_matrix_type ()) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
587 && (nargin < 5 || (argor = args (4), argor.is_string ()))) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
588 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
589 octave_idx_type m = argQ.rows (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
590 octave_idx_type n = argR.columns (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
591 octave_idx_type k = argQ.columns(); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
592 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
593 bool row = false; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
594 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
595 std::string orient = (nargin < 5) ? "col" : argor.string_value (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
596 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
597 if (nargin < 5 || (row = orient == "row") || (orient == "col")) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
598 if (argR.rows () == k |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
599 && (! row || m == k) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
600 && argx.rows () == (row ? 1 : m) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
601 && argx.columns () == (row ? n : 1)) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
602 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
603 int j = argj.int_value (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
604 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
605 if (j >= 1 && j <= (row ? n : m)+1) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
606 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
607 if (argQ.is_real_matrix () |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
608 && argR.is_real_matrix () |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
609 && argx.is_real_matrix ()) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
610 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
611 // real case |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
612 Matrix Q = argQ.matrix_value (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
613 Matrix R = argR.matrix_value (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
614 Matrix x = argx.matrix_value (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
615 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
616 QR fact (Q, R); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
617 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
618 if (row) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
619 fact.insert_row(x, j); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
620 else |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
621 fact.insert_col(x, j); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
622 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
623 retval(1) = fact.R (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
624 retval(0) = fact.Q (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
625 } |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
626 else |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
627 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
628 // complex case |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
629 ComplexMatrix Q = argQ.complex_matrix_value (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
630 ComplexMatrix R = argR.complex_matrix_value (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
631 ComplexMatrix x = argx.complex_matrix_value (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
632 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
633 ComplexQR fact (Q, R); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
634 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
635 if (row) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
636 fact.insert_row (x, j); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
637 else |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
638 fact.insert_col (x, j); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
639 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
640 retval(1) = fact.R (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
641 retval(0) = fact.Q (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
642 } |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
643 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
644 } |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
645 else |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
646 error ("qrinsert: index j out of range"); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
647 } |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
648 else |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
649 error ("qrinsert: dimension mismatch"); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
650 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
651 else |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
652 error ("qrinsert: orient must be \"col\" or \"row\""); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
653 } |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
654 else |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
655 print_usage (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
656 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
657 return retval; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
658 } |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
659 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
660 /* |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
661 %!test |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
662 %! A = [0.091364 0.613038 0.999083; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
663 %! 0.594638 0.425302 0.603537; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
664 %! 0.383594 0.291238 0.085574; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
665 %! 0.265712 0.268003 0.238409; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
666 %! 0.669966 0.743851 0.445057 ]; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
667 %! |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
668 %! x = [0.85082; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
669 %! 0.76426; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
670 %! 0.42883; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
671 %! 0.53010; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
672 %! 0.80683 ]; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
673 %! |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
674 %! [Q,R] = qr(A); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
675 %! [Q,R] = qrinsert(Q,R,3,x); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
676 %! assert(norm(vec(Q'*Q - eye(5)),Inf) < 1e1*eps) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
677 %! assert(norm(vec(triu(R)-R),Inf) == 0) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
678 %! assert(norm(vec(Q*R - [A(:,1:2) x A(:,3)]),Inf) < norm(A)*1e1*eps) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
679 %! |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
680 %!test |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
681 %! A = [0.620405 + 0.956953i 0.480013 + 0.048806i 0.402627 + 0.338171i; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
682 %! 0.589077 + 0.658457i 0.013205 + 0.279323i 0.229284 + 0.721929i; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
683 %! 0.092758 + 0.345687i 0.928679 + 0.241052i 0.764536 + 0.832406i; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
684 %! 0.912098 + 0.721024i 0.049018 + 0.269452i 0.730029 + 0.796517i; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
685 %! 0.112849 + 0.603871i 0.486352 + 0.142337i 0.355646 + 0.151496i ]; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
686 %! |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
687 %! x = [0.20351 + 0.05401i; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
688 %! 0.13141 + 0.43708i; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
689 %! 0.29808 + 0.08789i; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
690 %! 0.69821 + 0.38844i; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
691 %! 0.74871 + 0.25821i ]; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
692 %! |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
693 %! [Q,R] = qr(A); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
694 %! [Q,R] = qrinsert(Q,R,3,x); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
695 %! assert(norm(vec(Q'*Q - eye(5)),Inf) < 1e1*eps) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
696 %! assert(norm(vec(triu(R)-R),Inf) == 0) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
697 %! assert(norm(vec(Q*R - [A(:,1:2) x A(:,3)]),Inf) < norm(A)*1e1*eps) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
698 %! |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
699 %!test |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
700 %! A = [0.091364 0.613038 0.999083; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
701 %! 0.594638 0.425302 0.603537; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
702 %! 0.383594 0.291238 0.085574; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
703 %! 0.265712 0.268003 0.238409; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
704 %! 0.669966 0.743851 0.445057 ]; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
705 %! |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
706 %! x = [0.85082 0.76426 0.42883 ]; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
707 %! |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
708 %! [Q,R] = qr(A); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
709 %! [Q,R] = qrinsert(Q,R,3,x,'row'); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
710 %! assert(norm(vec(Q'*Q - eye(6)),Inf) < 1e1*eps) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
711 %! assert(norm(vec(triu(R)-R),Inf) == 0) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
712 %! assert(norm(vec(Q*R - [A(1:2,:);x;A(3:5,:)]),Inf) < norm(A)*1e1*eps) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
713 %! |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
714 %!test |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
715 %! A = [0.620405 + 0.956953i 0.480013 + 0.048806i 0.402627 + 0.338171i; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
716 %! 0.589077 + 0.658457i 0.013205 + 0.279323i 0.229284 + 0.721929i; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
717 %! 0.092758 + 0.345687i 0.928679 + 0.241052i 0.764536 + 0.832406i; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
718 %! 0.912098 + 0.721024i 0.049018 + 0.269452i 0.730029 + 0.796517i; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
719 %! 0.112849 + 0.603871i 0.486352 + 0.142337i 0.355646 + 0.151496i ]; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
720 %! |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
721 %! x = [0.20351 + 0.05401i 0.13141 + 0.43708i 0.29808 + 0.08789i ]; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
722 %! |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
723 %! [Q,R] = qr(A); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
724 %! [Q,R] = qrinsert(Q,R,3,x,'row'); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
725 %! assert(norm(vec(Q'*Q - eye(6)),Inf) < 1e1*eps) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
726 %! assert(norm(vec(triu(R)-R),Inf) == 0) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
727 %! assert(norm(vec(Q*R - [A(1:2,:);x;A(3:5,:)]),Inf) < norm(A)*1e1*eps) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
728 */ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
729 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
730 DEFUN_DLD (qrdelete, args, , |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
731 "-*- texinfo -*-\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
732 @deftypefn{Loadable Function} {[@var{Q1}, @var{R1}]} = qrdelete (@var{Q}, @var{R}, @var{j}, @var{orient})\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
733 Given a QR@tie{}factorization of a real or complex matrix\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
734 @w{@var{A} = @var{Q}*@var{R}}, @var{Q}@tie{}unitary and\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
735 @var{R}@tie{}upper trapezoidal, return the QR@tie{}factorization of\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
736 @w{[A(:,1:j-1) A(:,j+1:n)]}, i.e. @var{A} with one column deleted\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
737 (if @var{orient} is \"col\"), or the QR@tie{}factorization of\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
738 @w{[A(1:j-1,:);A(:,j+1:n)]}, i.e. @var{A} with one row deleted (if\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
739 @var{orient} is \"row\").\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
740 \n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
741 The default value of @var{orient} is \"col\".\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
742 \n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
743 If @var{orient} is \"col\", the matrix @var{Q} is not required to\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
744 be square.\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
745 \n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
746 For @sc{Matlab} compatibility, if @var{Q} is nonsquare on input, the\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
747 updated factorization is always stripped to the economical form, i.e.\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
748 @code{columns (Q) == rows (R) = columns (R)}.\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
749 \n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
750 To get the less intelligent but more natural behaviour when @var{Q}\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
751 retains it shape and @var{R} loses one column, set @var{orient} to\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
752 \"col+\" instead.\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
753 \n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
754 If @var{orient} is \"row\", @var{Q} must be square.\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
755 @seealso{qr, qrinsert, qrupdate}\n\ |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
756 @end deftypefn") |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
757 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
758 int nargin = args.length (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
759 octave_value_list retval; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
760 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
761 octave_value argQ,argR,argj,argor; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
762 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
763 if ((nargin == 3 || nargin == 4) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
764 && (argQ = args(0), argQ.is_matrix_type ()) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
765 && (argR = args(1), argR.is_matrix_type ()) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
766 && (argj = args(2), argj.is_scalar_type ()) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
767 && (nargin < 4 || (argor = args (3), argor.is_string ()))) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
768 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
769 octave_idx_type m = argQ.rows (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
770 octave_idx_type k = argQ.columns (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
771 octave_idx_type n = argR.columns (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
772 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
773 bool row = false; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
774 bool colp = false; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
775 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
776 std::string orient = (nargin < 4) ? "col" : argor.string_value (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
777 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
778 if (nargin < 4 || (row = orient == "row") |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
779 || (orient == "col") || (colp = orient == "col+")) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
780 if (argR.rows() == k) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
781 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
782 int j = argj.scalar_value (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
783 if (j >= 1 && j <= (row ? n : m)+1) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
784 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
785 if (argQ.is_real_matrix () |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
786 && argR.is_real_matrix ()) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
787 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
788 // real case |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
789 Matrix Q = argQ.matrix_value (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
790 Matrix R = argR.matrix_value (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
791 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
792 QR fact (Q, R); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
793 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
794 if (row) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
795 fact.delete_row (j); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
796 else |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
797 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
798 fact.delete_col (j); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
799 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
800 if (! colp && k < m) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
801 fact.economize (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
802 } |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
803 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
804 retval(1) = fact.R (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
805 retval(0) = fact.Q (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
806 } |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
807 else |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
808 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
809 // complex case |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
810 ComplexMatrix Q = argQ.complex_matrix_value (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
811 ComplexMatrix R = argR.complex_matrix_value (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
812 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
813 ComplexQR fact (Q, R); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
814 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
815 if (row) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
816 fact.delete_row (j); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
817 else |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
818 { |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
819 fact.delete_col (j); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
820 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
821 if (! colp && k < m) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
822 fact.economize (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
823 } |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
824 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
825 retval(1) = fact.R (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
826 retval(0) = fact.Q (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
827 } |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
828 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
829 } |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
830 else |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
831 error ("qrdelete: index j out of range"); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
832 } |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
833 else |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
834 error ("qrdelete: dimension mismatch"); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
835 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
836 else |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
837 error ("qrdelete: orient must be \"col\" or \"row\""); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
838 } |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
839 else |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
840 print_usage (); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
841 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
842 return retval; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
843 } |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
844 |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
845 /* |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
846 %!test |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
847 %! A = [0.091364 0.613038 0.027504 0.999083; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
848 %! 0.594638 0.425302 0.562834 0.603537; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
849 %! 0.383594 0.291238 0.742073 0.085574; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
850 %! 0.265712 0.268003 0.783553 0.238409; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
851 %! 0.669966 0.743851 0.457255 0.445057 ]; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
852 %! |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
853 %! [Q,R] = qr(A); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
854 %! [Q,R] = qrdelete(Q,R,3); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
855 %! assert(norm(vec(Q'*Q - eye(5)),Inf) < 1e1*eps) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
856 %! assert(norm(vec(triu(R)-R),Inf) == 0) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
857 %! assert(norm(vec(Q*R - [A(:,1:2) A(:,4)]),Inf) < norm(A)*1e1*eps) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
858 %! |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
859 %!test |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
860 %! A = [0.364554 + 0.993117i 0.669818 + 0.510234i 0.426568 + 0.041337i 0.847051 + 0.233291i; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
861 %! 0.049600 + 0.242783i 0.448946 + 0.484022i 0.141155 + 0.074420i 0.446746 + 0.392706i; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
862 %! 0.581922 + 0.657416i 0.581460 + 0.030016i 0.219909 + 0.447288i 0.201144 + 0.069132i; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
863 %! 0.694986 + 0.000571i 0.682327 + 0.841712i 0.807537 + 0.166086i 0.192767 + 0.358098i; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
864 %! 0.945002 + 0.066788i 0.350492 + 0.642638i 0.579629 + 0.048102i 0.600170 + 0.636938i ] * I; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
865 %! |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
866 %! [Q,R] = qr(A); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
867 %! [Q,R] = qrdelete(Q,R,3); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
868 %! assert(norm(vec(Q'*Q - eye(5)),Inf) < 1e1*eps) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
869 %! assert(norm(vec(triu(R)-R),Inf) == 0) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
870 %! assert(norm(vec(Q*R - [A(:,1:2) A(:,4)]),Inf) < norm(A)*1e1*eps) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
871 %! |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
872 %!test |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
873 %! A = [0.091364 0.613038 0.027504 0.999083; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
874 %! 0.594638 0.425302 0.562834 0.603537; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
875 %! 0.383594 0.291238 0.742073 0.085574; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
876 %! 0.265712 0.268003 0.783553 0.238409; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
877 %! 0.669966 0.743851 0.457255 0.445057 ]; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
878 %! |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
879 %! [Q,R] = qr(A); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
880 %! [Q,R] = qrdelete(Q,R,3,'row'); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
881 %! assert(norm(vec(Q'*Q - eye(4)),Inf) < 1e1*eps) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
882 %! assert(norm(vec(triu(R)-R),Inf) == 0) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
883 %! assert(norm(vec(Q*R - [A(1:2,:);A(4:5,:)]),Inf) < norm(A)*1e1*eps) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
884 %! |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
885 %!test |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
886 %! A = [0.364554 + 0.993117i 0.669818 + 0.510234i 0.426568 + 0.041337i 0.847051 + 0.233291i; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
887 %! 0.049600 + 0.242783i 0.448946 + 0.484022i 0.141155 + 0.074420i 0.446746 + 0.392706i; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
888 %! 0.581922 + 0.657416i 0.581460 + 0.030016i 0.219909 + 0.447288i 0.201144 + 0.069132i; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
889 %! 0.694986 + 0.000571i 0.682327 + 0.841712i 0.807537 + 0.166086i 0.192767 + 0.358098i; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
890 %! 0.945002 + 0.066788i 0.350492 + 0.642638i 0.579629 + 0.048102i 0.600170 + 0.636938i ] * I; |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
891 %! |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
892 %! [Q,R] = qr(A); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
893 %! [Q,R] = qrdelete(Q,R,3,'row'); |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
894 %! assert(norm(vec(Q'*Q - eye(4)),Inf) < 1e1*eps) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
895 %! assert(norm(vec(triu(R)-R),Inf) == 0) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
896 %! assert(norm(vec(Q*R - [A(1:2,:);A(4:5,:)]),Inf) < norm(A)*1e1*eps) |
56be6f31dd4e
implementation of QR factorization updating
Jaroslav Hajek <highegg@gmail.com>
parents:
7505
diff
changeset
|
897 */ |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
898 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
899 /* |
2928 | 900 ;;; Local Variables: *** |
901 ;;; mode: C++ *** | |
902 ;;; End: *** | |
903 */ |