5610
|
1 /* |
|
2 |
7017
|
3 Copyright (C) 2005, 2006, 2007 David Bateman |
7016
|
4 Copyright (C) 1998, 1999, 2000, 2001, 2002, 2003, 2004, 2005 Andy Adler |
|
5 |
|
6 This file is part of Octave. |
5610
|
7 |
|
8 Octave is free software; you can redistribute it and/or modify it |
|
9 under the terms of the GNU General Public License as published by the |
7016
|
10 Free Software Foundation; either version 3 of the License, or (at your |
|
11 option) any later version. |
5610
|
12 |
|
13 Octave is distributed in the hope that it will be useful, but WITHOUT |
|
14 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or |
|
15 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License |
|
16 for more details. |
|
17 |
|
18 You should have received a copy of the GNU General Public License |
7016
|
19 along with Octave; see the file COPYING. If not, see |
|
20 <http://www.gnu.org/licenses/>. |
5610
|
21 |
|
22 */ |
|
23 |
|
24 #ifdef HAVE_CONFIG_H |
|
25 #include <config.h> |
|
26 #endif |
|
27 |
|
28 #include "defun-dld.h" |
|
29 #include "error.h" |
|
30 #include "gripes.h" |
|
31 #include "oct-obj.h" |
|
32 #include "utils.h" |
|
33 |
5616
|
34 #include "oct-sparse.h" |
5610
|
35 #include "ov-re-sparse.h" |
|
36 #include "ov-cx-sparse.h" |
|
37 #include "SparseQR.h" |
|
38 #include "SparseCmplxQR.h" |
|
39 |
|
40 #ifdef IDX_TYPE_LONG |
5648
|
41 #define CXSPARSE_NAME(name) cs_dl ## name |
5610
|
42 #else |
5648
|
43 #define CXSPARSE_NAME(name) cs_di ## name |
5610
|
44 #endif |
|
45 |
|
46 // PKG_ADD: dispatch ("qr", "spqr", "sparse matrix"); |
|
47 // PKG_ADD: dispatch ("qr", "spqr", "sparse complex matrix"); |
|
48 // PKG_ADD: dispatch ("qr", "spqr", "sparse bool matrix"); |
|
49 DEFUN_DLD (spqr, args, nargout, |
|
50 "-*- texinfo -*-\n\ |
|
51 @deftypefn {Loadable Function} {@var{r} =} spqr (@var{a})\n\ |
|
52 @deftypefnx {Loadable Function} {@var{r} =} spqr (@var{a},0)\n\ |
|
53 @deftypefnx {Loadable Function} {[@var{c}, @var{r}] =} spqr (@var{a},@var{b})\n\ |
|
54 @deftypefnx {Loadable Function} {[@var{c}, @var{r}] =} spqr (@var{a},@var{b},0)\n\ |
|
55 @cindex QR factorization\n\ |
|
56 Compute the sparse QR factorization of @var{a}, using @sc{CSparse}.\n\ |
|
57 As the matrix @var{Q} is in general a full matrix, this function returns\n\ |
|
58 the @var{Q}-less factorization @var{r} of @var{a}, such that\n\ |
|
59 @code{@var{r} = chol (@var{a}' * @var{a})}.\n\ |
|
60 \n\ |
|
61 If the final argument is the scalar @code{0} and the number of rows is\n\ |
|
62 larger than the number of columns, then an economy factorization is\n\ |
|
63 returned. That is @var{r} will have only @code{size (@var{a},1)} rows.\n\ |
|
64 \n\ |
|
65 If an additional matrix @var{b} is supplied, then @code{spqr} returns\n\ |
|
66 @var{c}, where @code{@var{c} = @var{q}' * @var{b}}. This allows the\n\ |
|
67 least squares approximation of @code{@var{a} \\ @var{b}} to be calculated\n\ |
|
68 as\n\ |
|
69 \n\ |
|
70 @example\n\ |
|
71 [@var{c},@var{r}] = spqr (@var{a},@var{b})\n\ |
|
72 @var{x} = @var{r} \\ @var{c}\n\ |
|
73 @end example\n\ |
5642
|
74 @seealso{spchol, qr}\n\ |
|
75 @end deftypefn") |
5610
|
76 { |
|
77 int nargin = args.length (); |
|
78 octave_value_list retval; |
|
79 bool economy = false; |
|
80 bool is_cmplx = false; |
|
81 bool have_b = false; |
|
82 |
|
83 if (nargin < 1 || nargin > 3) |
5823
|
84 print_usage (); |
5610
|
85 else |
|
86 { |
|
87 if (args(0).is_complex_type ()) |
|
88 is_cmplx = true; |
|
89 if (nargin > 1) |
|
90 { |
|
91 have_b = true; |
|
92 if (args(nargin-1).is_scalar_type ()) |
|
93 { |
|
94 int val = args(nargin-1).int_value (); |
|
95 if (val == 0) |
|
96 { |
|
97 economy = true; |
|
98 have_b = (nargin > 2); |
|
99 } |
|
100 } |
|
101 if (have_b && args(1).is_complex_type ()) |
|
102 is_cmplx = true; |
|
103 } |
|
104 |
|
105 if (!error_state) |
|
106 { |
|
107 if (have_b && nargout < 2) |
|
108 error ("spqr: incorrect number of output arguments"); |
|
109 else if (is_cmplx) |
|
110 { |
|
111 SparseComplexQR q (args(0).sparse_complex_matrix_value ()); |
|
112 if (!error_state) |
|
113 { |
|
114 if (have_b) |
|
115 { |
|
116 retval(1) = q.R (economy); |
|
117 retval(0) = q.C (args(1).complex_matrix_value ()); |
5681
|
118 if (args(0).rows() < args(0).columns()) |
|
119 warning ("spqr: non minimum norm solution for under-determined problem"); |
5610
|
120 } |
|
121 else |
|
122 retval(0) = q.R (economy); |
|
123 } |
|
124 } |
|
125 else |
|
126 { |
|
127 SparseQR q (args(0).sparse_matrix_value ()); |
|
128 if (!error_state) |
|
129 { |
|
130 if (have_b) |
|
131 { |
|
132 retval(1) = q.R (economy); |
|
133 retval(0) = q.C (args(1).matrix_value ()); |
5681
|
134 if (args(0).rows() < args(0).columns()) |
|
135 warning ("spqr: non minimum norm solution for under-determined problem"); |
5610
|
136 } |
|
137 else |
|
138 retval(0) = q.R (economy); |
|
139 } |
|
140 } |
|
141 } |
|
142 } |
|
143 return retval; |
|
144 } |
|
145 |
|
146 /* |
|
147 |
|
148 The deactivated tests below can't be tested till rectangular back-subs is |
|
149 implemented for sparse matrices. |
|
150 |
7243
|
151 %!testif HAVE_CXSPARSE |
5610
|
152 %! n = 20; d= 0.2; |
|
153 %! a = sprandn(n,n,d)+speye(n,n); |
|
154 %! r = spqr(a); |
|
155 %! assert(r'*r,a'*a,1e-10) |
|
156 |
7243
|
157 %!testif HAVE_CXSPARSE |
5610
|
158 %! n = 20; d= 0.2; |
|
159 %! a = sprandn(n,n,d)+speye(n,n); |
|
160 %! q = symamd(a); |
|
161 %! a = a(q,q); |
|
162 %! r = spqr(a); |
|
163 %! assert(r'*r,a'*a,1e-10) |
|
164 |
7243
|
165 %!testif HAVE_CXSPARSE |
5610
|
166 %! n = 20; d= 0.2; |
|
167 %! a = sprandn(n,n,d)+speye(n,n); |
|
168 %! [c,r] = spqr(a,ones(n,1)); |
|
169 %! assert (r\c,full(a)\ones(n,1),10e-10) |
|
170 |
7243
|
171 %!testif HAVE_CXSPARSE |
5610
|
172 %! n = 20; d= 0.2; |
|
173 %! a = sprandn(n,n,d)+speye(n,n); |
|
174 %! b = randn(n,2); |
|
175 %! [c,r] = spqr(a,b); |
|
176 %! assert (r\c,full(a)\b,10e-10) |
|
177 |
|
178 %% Test under-determined systems!! |
7243
|
179 %!#testif HAVE_CXSPARSE |
5610
|
180 %! n = 20; d= 0.2; |
|
181 %! a = sprandn(n,n+1,d)+speye(n,n+1); |
|
182 %! b = randn(n,2); |
|
183 %! [c,r] = spqr(a,b); |
|
184 %! assert (r\c,full(a)\b,10e-10) |
|
185 |
7243
|
186 %!testif HAVE_CXSPARSE |
5610
|
187 %! n = 20; d= 0.2; |
|
188 %! a = 1i*sprandn(n,n,d)+speye(n,n); |
|
189 %! r = spqr(a); |
|
190 %! assert(r'*r,a'*a,1e-10) |
|
191 |
7243
|
192 %!testif HAVE_CXSPARSE |
5610
|
193 %! n = 20; d= 0.2; |
|
194 %! a = 1i*sprandn(n,n,d)+speye(n,n); |
|
195 %! q = symamd(a); |
|
196 %! a = a(q,q); |
|
197 %! r = spqr(a); |
|
198 %! assert(r'*r,a'*a,1e-10) |
|
199 |
7243
|
200 %!testif HAVE_CXSPARSE |
5610
|
201 %! n = 20; d= 0.2; |
|
202 %! a = 1i*sprandn(n,n,d)+speye(n,n); |
|
203 %! [c,r] = spqr(a,ones(n,1)); |
|
204 %! assert (r\c,full(a)\ones(n,1),10e-10) |
|
205 |
7243
|
206 %!testif HAVE_CXSPARSE |
5610
|
207 %! n = 20; d= 0.2; |
|
208 %! a = 1i*sprandn(n,n,d)+speye(n,n); |
|
209 %! b = randn(n,2); |
|
210 %! [c,r] = spqr(a,b); |
|
211 %! assert (r\c,full(a)\b,10e-10) |
|
212 |
|
213 %% Test under-determined systems!! |
7243
|
214 %!#testif HAVE_CXSPARSE |
5610
|
215 %! n = 20; d= 0.2; |
|
216 %! a = 1i*sprandn(n,n+1,d)+speye(n,n+1); |
|
217 %! b = randn(n,2); |
|
218 %! [c,r] = spqr(a,b); |
|
219 %! assert (r\c,full(a)\b,10e-10) |
|
220 |
|
221 %!error spqr(sprandn(10,10,0.2),ones(10,1)); |
|
222 |
|
223 */ |
|
224 |
|
225 static RowVector |
|
226 put_int (octave_idx_type *p, octave_idx_type n) |
|
227 { |
|
228 RowVector ret (n); |
|
229 for (octave_idx_type i = 0; i < n; i++) |
|
230 ret.xelem(i) = p[i] + 1; |
|
231 return ret; |
|
232 } |
|
233 |
6066
|
234 #if HAVE_CXSPARSE |
|
235 static octave_value_list |
|
236 dmperm_internal (bool rank, const octave_value arg, int nargout) |
5610
|
237 { |
|
238 octave_value_list retval; |
|
239 octave_idx_type nr = arg.rows (); |
|
240 octave_idx_type nc = arg.columns (); |
|
241 SparseMatrix m; |
|
242 SparseComplexMatrix cm; |
5648
|
243 CXSPARSE_NAME () csm; |
5610
|
244 csm.m = nr; |
|
245 csm.n = nc; |
|
246 csm.x = NULL; |
|
247 csm.nz = -1; |
|
248 |
|
249 if (arg.is_real_type ()) |
|
250 { |
|
251 m = arg.sparse_matrix_value (); |
|
252 csm.nzmax = m.nnz(); |
|
253 csm.p = m.xcidx (); |
|
254 csm.i = m.xridx (); |
|
255 } |
|
256 else |
|
257 { |
|
258 cm = arg.sparse_complex_matrix_value (); |
|
259 csm.nzmax = cm.nnz(); |
|
260 csm.p = cm.xcidx (); |
|
261 csm.i = cm.xridx (); |
|
262 } |
|
263 |
|
264 if (!error_state) |
|
265 { |
6066
|
266 if (nargout <= 1 || rank) |
5610
|
267 { |
5792
|
268 #if defined(CS_VER) && (CS_VER >= 2) |
|
269 octave_idx_type *jmatch = CXSPARSE_NAME (_maxtrans) (&csm, 0); |
|
270 #else |
5648
|
271 octave_idx_type *jmatch = CXSPARSE_NAME (_maxtrans) (&csm); |
5792
|
272 #endif |
6066
|
273 if (rank) |
|
274 { |
|
275 octave_idx_type r = 0; |
|
276 for (octave_idx_type i = 0; i < nc; i++) |
|
277 if (jmatch[nr+i] >= 0) |
|
278 r++; |
|
279 retval(0) = static_cast<double>(r); |
|
280 } |
|
281 else |
|
282 retval(0) = put_int (jmatch + nr, nc); |
5648
|
283 CXSPARSE_NAME (_free) (jmatch); |
5610
|
284 } |
|
285 else |
|
286 { |
5792
|
287 #if defined(CS_VER) && (CS_VER >= 2) |
|
288 CXSPARSE_NAME (d) *dm = CXSPARSE_NAME(_dmperm) (&csm, 0); |
|
289 #else |
5648
|
290 CXSPARSE_NAME (d) *dm = CXSPARSE_NAME(_dmperm) (&csm); |
5792
|
291 #endif |
6066
|
292 |
5610
|
293 //retval(5) = put_int (dm->rr, 5); |
|
294 //retval(4) = put_int (dm->cc, 5); |
5792
|
295 #if defined(CS_VER) && (CS_VER >= 2) |
|
296 retval(3) = put_int (dm->s, dm->nb+1); |
|
297 retval(2) = put_int (dm->r, dm->nb+1); |
|
298 retval(1) = put_int (dm->q, nc); |
|
299 retval(0) = put_int (dm->p, nr); |
|
300 #else |
5610
|
301 retval(3) = put_int (dm->S, dm->nb+1); |
|
302 retval(2) = put_int (dm->R, dm->nb+1); |
|
303 retval(1) = put_int (dm->Q, nc); |
|
304 retval(0) = put_int (dm->P, nr); |
5792
|
305 #endif |
5648
|
306 CXSPARSE_NAME (_dfree) (dm); |
5610
|
307 } |
|
308 } |
6066
|
309 return retval; |
|
310 } |
|
311 #endif |
|
312 |
|
313 DEFUN_DLD (dmperm, args, nargout, |
|
314 "-*- texinfo -*-\n\ |
|
315 @deftypefn {Loadable Function} {@var{p} =} dmperm (@var{s})\n\ |
|
316 @deftypefnx {Loadable Function} {[@var{p}, @var{q}, @var{r}, @var{s}] =} dmperm (@var{s})\n\ |
|
317 \n\ |
|
318 @cindex Dulmage-Mendelsohn decomposition\n\ |
|
319 Perform a Dulmage-Mendelsohn permutation on the sparse matrix @var{s}.\n\ |
|
320 With a single output argument @dfn{dmperm} performs the row permutations\n\ |
|
321 @var{p} such that @code{@var{s} (@var{p},:)} has no zero elements on the\n\ |
|
322 diagonal.\n\ |
|
323 \n\ |
|
324 Called with two or more output arguments, returns the row and column\n\ |
|
325 permutations, such that @code{@var{s} (@var{p}, @var{q})} is in block\n\ |
|
326 triangular form. The values of @var{r} and @var{s} define the boundaries\n\ |
|
327 of the blocks. If @var{s} is square then @code{@var{r} == @var{s}}.\n\ |
|
328 \n\ |
|
329 The method used is described in: A. Pothen & C.-J. Fan. Computing the block\n\ |
|
330 triangular form of a sparse matrix. ACM Trans. Math. Software,\n\ |
|
331 16(4):303-324, 1990.\n\ |
|
332 @seealso{colamd, ccolamd}\n\ |
|
333 @end deftypefn") |
|
334 { |
|
335 int nargin = args.length(); |
|
336 octave_value_list retval; |
|
337 |
|
338 if (nargin != 1) |
|
339 { |
|
340 print_usage (); |
|
341 return retval; |
|
342 } |
|
343 |
|
344 #if HAVE_CXSPARSE |
|
345 retval = dmperm_internal (false, args(0), nargout); |
5610
|
346 #else |
|
347 error ("dmperm: not available in this version of Octave"); |
|
348 #endif |
|
349 |
|
350 return retval; |
|
351 } |
|
352 |
6066
|
353 /* |
5610
|
354 |
7243
|
355 %!testif HAVE_CXSPARSE |
5610
|
356 %! n=20; |
|
357 %! a=speye(n,n);a=a(randperm(n),:); |
|
358 %! assert(a(dmperm(a),:),speye(n)) |
|
359 |
7243
|
360 %!testif HAVE_CXSPARSE |
5610
|
361 %! n=20; |
|
362 %! d=0.2; |
|
363 %! a=tril(sprandn(n,n,d),-1)+speye(n,n); |
|
364 %! a=a(randperm(n),randperm(n)); |
|
365 %! [p,q,r,s]=dmperm(a); |
|
366 %! assert(tril(a(p,q),-1),sparse(n,n)) |
|
367 |
|
368 */ |
|
369 |
6066
|
370 DEFUN_DLD (sprank, args, nargout, |
|
371 "-*- texinfo -*-\n\ |
|
372 @deftypefn {Loadable Function} {@var{p} =} sprank (@var{s})\n\ |
|
373 \n\ |
|
374 @cindex Structural Rank\n\ |
|
375 Calculates the structural rank of a sparse matrix @var{s}. Note that\n\ |
|
376 only the structure of the matrix is used in this calculation based on\n\ |
|
377 a Dulmage-Mendelsohn to block triangular form. As such the numerical\n\ |
|
378 rank of the matrix @var{s} is bounded by @code{sprank (@var{s}) >=\n\ |
|
379 rank (@var{s})}. Ignoring floating point errors @code{sprank (@var{s}) ==\n\ |
|
380 rank (@var{s})}.\n\ |
|
381 @seealso{dmperm}\n\ |
|
382 @end deftypefn") |
|
383 { |
|
384 int nargin = args.length(); |
|
385 octave_value_list retval; |
|
386 |
|
387 if (nargin != 1) |
|
388 { |
|
389 print_usage (); |
|
390 return retval; |
|
391 } |
|
392 |
|
393 #if HAVE_CXSPARSE |
|
394 retval = dmperm_internal (true, args(0), nargout); |
|
395 #else |
|
396 error ("sprank: not available in this version of Octave"); |
|
397 #endif |
|
398 |
|
399 return retval; |
|
400 } |
|
401 |
|
402 /* |
|
403 |
|
404 %!error(sprank(1,2)); |
7243
|
405 %!testif HAVE_CXSPARSE |
|
406 %! assert(sprank(speye(20)), 20) |
|
407 %!testif HAVE_CXSPARSE |
|
408 %! assert(sprank([1,0,2,0;2,0,4,0]),2) |
6066
|
409 |
|
410 */ |
5610
|
411 /* |
|
412 ;;; Local Variables: *** |
|
413 ;;; mode: C++ *** |
|
414 ;;; End: *** |
|
415 */ |