Mercurial > octave-antonio
annotate src/DLD-FUNCTIONS/qr.cc @ 7505:f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
author | David Bateman <dbateman@free.fr> |
---|---|
date | Wed, 20 Feb 2008 15:52:11 -0500 |
parents | 7879ef1042a8 |
children | 56be6f31dd4e |
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 | |
23 #ifdef HAVE_CONFIG_H | |
24 #include <config.h> | |
25 #endif | |
26 | |
27 #include "CmplxQR.h" | |
28 #include "CmplxQRP.h" | |
29 #include "dbleQR.h" | |
30 #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
|
31 #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
|
32 #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
|
33 |
2928 | 34 |
35 #include "defun-dld.h" | |
36 #include "error.h" | |
37 #include "gripes.h" | |
38 #include "oct-obj.h" | |
39 #include "utils.h" | |
40 | |
3372 | 41 // [Q, R] = qr (X): form Q unitary and R upper triangular such |
42 // that Q * R = X | |
43 // | |
44 // [Q, R] = qr (X, 0): form the economy decomposition such that if X is | |
45 // m by n then only the first n columns of Q are | |
46 // computed. | |
47 // | |
48 // [Q, R, P] = qr (X): form QRP factorization of X where | |
49 // P is a permutation matrix such that | |
50 // A * P = Q * R | |
51 // | |
52 // [Q, R, P] = qr (X, 0): form the economy decomposition with | |
53 // permutation vector P such that Q * R = X (:, P) | |
54 // | |
55 // qr (X) alone returns the output of the LAPACK routine dgeqrf, such | |
56 // that R = triu (qr (X)) | |
57 | |
2928 | 58 DEFUN_DLD (qr, args, nargout, |
3548 | 59 "-*- texinfo -*-\n\ |
3372 | 60 @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
|
61 @deftypefnx {Loadable Function} {[@var{q}, @var{r}, @var{p}] =} qr (@var{a}, '0')\n\ |
3372 | 62 @cindex QR factorization\n\ |
63 Compute the QR factorization of @var{a}, using standard @sc{Lapack}\n\ | |
64 subroutines. For example, given the matrix @code{a = [1, 2; 3, 4]},\n\ | |
2928 | 65 \n\ |
3372 | 66 @example\n\ |
67 [q, r] = qr (a)\n\ | |
68 @end example\n\ | |
69 \n\ | |
70 @noindent\n\ | |
71 returns\n\ | |
72 \n\ | |
73 @example\n\ | |
74 q =\n\ | |
75 \n\ | |
76 -0.31623 -0.94868\n\ | |
77 -0.94868 0.31623\n\ | |
78 \n\ | |
79 r =\n\ | |
80 \n\ | |
81 -3.16228 -4.42719\n\ | |
82 0.00000 -0.63246\n\ | |
83 @end example\n\ | |
84 \n\ | |
85 The @code{qr} factorization has applications in the solution of least\n\ | |
86 squares problems\n\ | |
87 @iftex\n\ | |
88 @tex\n\ | |
89 $$\n\ | |
90 \\min_x \\left\\Vert A x - b \\right\\Vert_2\n\ | |
91 $$\n\ | |
92 @end tex\n\ | |
93 @end iftex\n\ | |
94 @ifinfo\n\ | |
2928 | 95 \n\ |
3372 | 96 @example\n\ |
97 @code{min norm(A x - b)}\n\ | |
98 @end example\n\ | |
99 \n\ | |
100 @end ifinfo\n\ | |
101 for overdetermined systems of equations (i.e.,\n\ | |
102 @iftex\n\ | |
103 @tex\n\ | |
104 $A$\n\ | |
105 @end tex\n\ | |
106 @end iftex\n\ | |
107 @ifinfo\n\ | |
108 @code{a}\n\ | |
109 @end ifinfo\n\ | |
110 is a tall, thin matrix). The QR factorization is\n\ | |
111 @iftex\n\ | |
112 @tex\n\ | |
113 $QR = A$ where $Q$ is an orthogonal matrix and $R$ is upper triangular.\n\ | |
114 @end tex\n\ | |
115 @end iftex\n\ | |
116 @ifinfo\n\ | |
117 @code{q * r = a} where @code{q} is an orthogonal matrix and @code{r} is\n\ | |
118 upper triangular.\n\ | |
119 @end ifinfo\n\ | |
2928 | 120 \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
|
121 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
|
122 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
|
123 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
|
124 \n\ |
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 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
|
126 @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
|
127 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
|
128 order. For example,given the matrix @code{a = [1, 2; 3, 4]},\n\ |
3372 | 129 \n\ |
130 @example\n\ | |
3600 | 131 [q, r, p] = qr(a)\n\ |
3372 | 132 @end example\n\ |
133 \n\ | |
134 @noindent\n\ | |
135 returns\n\ | |
136 \n\ | |
137 @example\n\ | |
138 q = \n\ | |
2928 | 139 \n\ |
3372 | 140 -0.44721 -0.89443\n\ |
141 -0.89443 0.44721\n\ | |
142 \n\ | |
143 r =\n\ | |
144 \n\ | |
145 -4.47214 -3.13050\n\ | |
146 0.00000 0.44721\n\ | |
147 \n\ | |
148 p =\n\ | |
149 \n\ | |
150 0 1\n\ | |
151 1 0\n\ | |
152 @end example\n\ | |
153 \n\ | |
154 The permuted @code{qr} factorization @code{[q, r, p] = qr (a)}\n\ | |
155 factorization allows the construction of an orthogonal basis of\n\ | |
156 @code{span (a)}.\n\ | |
7491 | 157 \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
|
158 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
|
159 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
|
160 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
|
161 @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
|
162 \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 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
|
164 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
|
165 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
|
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 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
|
168 @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
|
169 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
|
170 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
|
171 \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 @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
|
173 [@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
|
174 @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
|
175 @end example\n\ |
3372 | 176 @end deftypefn") |
2928 | 177 { |
178 octave_value_list retval; | |
179 | |
180 int nargin = args.length (); | |
181 | |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
182 if (nargin < 1 || nargin > (args(0).is_sparse_type() ? 3 : 2) || nargout > 3) |
2928 | 183 { |
5823 | 184 print_usage (); |
2928 | 185 return retval; |
186 } | |
187 | |
188 octave_value arg = args(0); | |
189 | |
190 int arg_is_empty = empty_arg ("qr", arg.rows (), arg.columns ()); | |
191 | |
192 if (arg_is_empty < 0) | |
193 return retval; | |
194 else if (arg_is_empty > 0) | |
195 return octave_value_list (3, Matrix ()); | |
196 | |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
197 if (arg.is_sparse_type ()) |
2928 | 198 { |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
199 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
|
200 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
|
201 int have_b = 0; |
3068 | 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 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
|
204 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
|
205 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
|
206 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
207 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
|
208 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
|
209 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
210 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
|
211 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
|
212 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
213 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
|
214 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
|
215 } |
2928 | 216 } |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
217 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
|
218 is_cmplx = true; |
2928 | 219 } |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
220 |
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 (!error_state) |
2928 | 222 { |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
223 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
|
224 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
|
225 else if (is_cmplx) |
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 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
|
228 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
|
229 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
230 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
|
231 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
232 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
|
233 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
|
234 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
|
235 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
|
236 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
237 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
|
238 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
239 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
|
240 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
|
241 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
242 else |
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(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
|
244 } |
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 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
248 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
|
249 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
|
250 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
251 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
|
252 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
253 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
|
254 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
|
255 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
|
256 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
|
257 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
258 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
|
259 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
260 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
|
261 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
|
262 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
263 else |
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(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
|
265 } |
2928 | 266 } |
267 } | |
268 } | |
269 else | |
270 { | |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
271 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
|
272 : (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
|
273 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
274 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
|
275 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
276 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
|
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 (! 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
|
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 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
|
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 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
|
283 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
|
284 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
285 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
|
286 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
|
287 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
288 break; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
289 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
290 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
|
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 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
|
293 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
|
294 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
|
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 break; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
297 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
298 default: |
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 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
|
301 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
|
302 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
|
303 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
|
304 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
305 break; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
306 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
307 } |
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 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
|
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 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
|
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 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
|
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 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
|
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 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
|
318 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
|
319 { |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
320 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
|
321 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
|
322 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
323 break; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
324 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
325 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
|
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 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
|
328 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
|
329 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
|
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 break; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
332 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
333 default: |
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 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
|
336 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
|
337 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
|
338 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
|
339 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
340 break; |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
341 } |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
342 } |
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 else |
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 gripe_wrong_type_arg ("qr", arg); |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
347 } |
2928 | 348 } |
349 | |
350 return retval; | |
351 } | |
352 | |
353 /* | |
7505
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
354 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
355 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
|
356 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
|
357 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
358 %!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
|
359 %! 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
|
360 %! 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
|
361 %! 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
|
362 %! 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
|
363 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
364 %!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
|
365 %! 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
|
366 %! 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
|
367 %! 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
|
368 %! 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
|
369 %! 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
|
370 %! 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
|
371 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
372 %!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
|
373 %! 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
|
374 %! 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
|
375 %! [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
|
376 %! 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
|
377 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
378 %!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
|
379 %! 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
|
380 %! 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
|
381 %! 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
|
382 %! [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
|
383 %! 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
|
384 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
385 %% 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
|
386 %!#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
|
387 %! 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
|
388 %! 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
|
389 %! 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
|
390 %! [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
|
391 %! 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
|
392 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
393 %!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
|
394 %! 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
|
395 %! 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
|
396 %! 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
|
397 %! 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
|
398 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
399 %!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
|
400 %! 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
|
401 %! 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
|
402 %! 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
|
403 %! 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
|
404 %! 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
|
405 %! 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
|
406 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
407 %!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
|
408 %! 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
|
409 %! 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
|
410 %! [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
|
411 %! 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
|
412 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
413 %!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
|
414 %! 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
|
415 %! 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
|
416 %! 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
|
417 %! [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
|
418 %! 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
|
419 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
420 %% 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
|
421 %!#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
|
422 %! 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
|
423 %! 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
|
424 %! 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
|
425 %! [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
|
426 %! 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
|
427 |
f5005d9510f4
Remove dispatched sparse functions and treat in the generic versions of the functions
David Bateman <dbateman@free.fr>
parents:
7491
diff
changeset
|
428 %!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
|
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 */ |
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 /* |
2928 | 434 ;;; Local Variables: *** |
435 ;;; mode: C++ *** | |
436 ;;; End: *** | |
437 */ |