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 static RowVector |
|
47 put_int (octave_idx_type *p, octave_idx_type n) |
|
48 { |
|
49 RowVector ret (n); |
|
50 for (octave_idx_type i = 0; i < n; i++) |
|
51 ret.xelem(i) = p[i] + 1; |
|
52 return ret; |
|
53 } |
|
54 |
6066
|
55 #if HAVE_CXSPARSE |
|
56 static octave_value_list |
|
57 dmperm_internal (bool rank, const octave_value arg, int nargout) |
5610
|
58 { |
|
59 octave_value_list retval; |
|
60 octave_idx_type nr = arg.rows (); |
|
61 octave_idx_type nc = arg.columns (); |
|
62 SparseMatrix m; |
|
63 SparseComplexMatrix cm; |
5648
|
64 CXSPARSE_NAME () csm; |
5610
|
65 csm.m = nr; |
|
66 csm.n = nc; |
7520
|
67 csm.x = 0; |
5610
|
68 csm.nz = -1; |
|
69 |
|
70 if (arg.is_real_type ()) |
|
71 { |
|
72 m = arg.sparse_matrix_value (); |
|
73 csm.nzmax = m.nnz(); |
|
74 csm.p = m.xcidx (); |
|
75 csm.i = m.xridx (); |
|
76 } |
|
77 else |
|
78 { |
|
79 cm = arg.sparse_complex_matrix_value (); |
|
80 csm.nzmax = cm.nnz(); |
|
81 csm.p = cm.xcidx (); |
|
82 csm.i = cm.xridx (); |
|
83 } |
|
84 |
|
85 if (!error_state) |
|
86 { |
6066
|
87 if (nargout <= 1 || rank) |
5610
|
88 { |
5792
|
89 #if defined(CS_VER) && (CS_VER >= 2) |
|
90 octave_idx_type *jmatch = CXSPARSE_NAME (_maxtrans) (&csm, 0); |
|
91 #else |
5648
|
92 octave_idx_type *jmatch = CXSPARSE_NAME (_maxtrans) (&csm); |
5792
|
93 #endif |
6066
|
94 if (rank) |
|
95 { |
|
96 octave_idx_type r = 0; |
|
97 for (octave_idx_type i = 0; i < nc; i++) |
|
98 if (jmatch[nr+i] >= 0) |
|
99 r++; |
|
100 retval(0) = static_cast<double>(r); |
|
101 } |
|
102 else |
|
103 retval(0) = put_int (jmatch + nr, nc); |
5648
|
104 CXSPARSE_NAME (_free) (jmatch); |
5610
|
105 } |
|
106 else |
|
107 { |
5792
|
108 #if defined(CS_VER) && (CS_VER >= 2) |
|
109 CXSPARSE_NAME (d) *dm = CXSPARSE_NAME(_dmperm) (&csm, 0); |
|
110 #else |
5648
|
111 CXSPARSE_NAME (d) *dm = CXSPARSE_NAME(_dmperm) (&csm); |
5792
|
112 #endif |
6066
|
113 |
5610
|
114 //retval(5) = put_int (dm->rr, 5); |
|
115 //retval(4) = put_int (dm->cc, 5); |
5792
|
116 #if defined(CS_VER) && (CS_VER >= 2) |
|
117 retval(3) = put_int (dm->s, dm->nb+1); |
|
118 retval(2) = put_int (dm->r, dm->nb+1); |
|
119 retval(1) = put_int (dm->q, nc); |
|
120 retval(0) = put_int (dm->p, nr); |
|
121 #else |
5610
|
122 retval(3) = put_int (dm->S, dm->nb+1); |
|
123 retval(2) = put_int (dm->R, dm->nb+1); |
|
124 retval(1) = put_int (dm->Q, nc); |
|
125 retval(0) = put_int (dm->P, nr); |
5792
|
126 #endif |
5648
|
127 CXSPARSE_NAME (_dfree) (dm); |
5610
|
128 } |
|
129 } |
6066
|
130 return retval; |
|
131 } |
|
132 #endif |
|
133 |
|
134 DEFUN_DLD (dmperm, args, nargout, |
|
135 "-*- texinfo -*-\n\ |
|
136 @deftypefn {Loadable Function} {@var{p} =} dmperm (@var{s})\n\ |
|
137 @deftypefnx {Loadable Function} {[@var{p}, @var{q}, @var{r}, @var{s}] =} dmperm (@var{s})\n\ |
|
138 \n\ |
|
139 @cindex Dulmage-Mendelsohn decomposition\n\ |
|
140 Perform a Dulmage-Mendelsohn permutation on the sparse matrix @var{s}.\n\ |
|
141 With a single output argument @dfn{dmperm} performs the row permutations\n\ |
|
142 @var{p} such that @code{@var{s} (@var{p},:)} has no zero elements on the\n\ |
|
143 diagonal.\n\ |
|
144 \n\ |
|
145 Called with two or more output arguments, returns the row and column\n\ |
|
146 permutations, such that @code{@var{s} (@var{p}, @var{q})} is in block\n\ |
|
147 triangular form. The values of @var{r} and @var{s} define the boundaries\n\ |
|
148 of the blocks. If @var{s} is square then @code{@var{r} == @var{s}}.\n\ |
|
149 \n\ |
|
150 The method used is described in: A. Pothen & C.-J. Fan. Computing the block\n\ |
|
151 triangular form of a sparse matrix. ACM Trans. Math. Software,\n\ |
|
152 16(4):303-324, 1990.\n\ |
|
153 @seealso{colamd, ccolamd}\n\ |
|
154 @end deftypefn") |
|
155 { |
|
156 int nargin = args.length(); |
|
157 octave_value_list retval; |
|
158 |
|
159 if (nargin != 1) |
|
160 { |
|
161 print_usage (); |
|
162 return retval; |
|
163 } |
|
164 |
|
165 #if HAVE_CXSPARSE |
|
166 retval = dmperm_internal (false, args(0), nargout); |
5610
|
167 #else |
|
168 error ("dmperm: not available in this version of Octave"); |
|
169 #endif |
|
170 |
|
171 return retval; |
|
172 } |
|
173 |
6066
|
174 /* |
5610
|
175 |
7243
|
176 %!testif HAVE_CXSPARSE |
5610
|
177 %! n=20; |
|
178 %! a=speye(n,n);a=a(randperm(n),:); |
|
179 %! assert(a(dmperm(a),:),speye(n)) |
|
180 |
7243
|
181 %!testif HAVE_CXSPARSE |
5610
|
182 %! n=20; |
|
183 %! d=0.2; |
|
184 %! a=tril(sprandn(n,n,d),-1)+speye(n,n); |
|
185 %! a=a(randperm(n),randperm(n)); |
|
186 %! [p,q,r,s]=dmperm(a); |
|
187 %! assert(tril(a(p,q),-1),sparse(n,n)) |
|
188 |
|
189 */ |
|
190 |
6066
|
191 DEFUN_DLD (sprank, args, nargout, |
|
192 "-*- texinfo -*-\n\ |
|
193 @deftypefn {Loadable Function} {@var{p} =} sprank (@var{s})\n\ |
|
194 \n\ |
|
195 @cindex Structural Rank\n\ |
|
196 Calculates the structural rank of a sparse matrix @var{s}. Note that\n\ |
|
197 only the structure of the matrix is used in this calculation based on\n\ |
|
198 a Dulmage-Mendelsohn to block triangular form. As such the numerical\n\ |
|
199 rank of the matrix @var{s} is bounded by @code{sprank (@var{s}) >=\n\ |
|
200 rank (@var{s})}. Ignoring floating point errors @code{sprank (@var{s}) ==\n\ |
|
201 rank (@var{s})}.\n\ |
|
202 @seealso{dmperm}\n\ |
|
203 @end deftypefn") |
|
204 { |
|
205 int nargin = args.length(); |
|
206 octave_value_list retval; |
|
207 |
|
208 if (nargin != 1) |
|
209 { |
|
210 print_usage (); |
|
211 return retval; |
|
212 } |
|
213 |
|
214 #if HAVE_CXSPARSE |
|
215 retval = dmperm_internal (true, args(0), nargout); |
|
216 #else |
|
217 error ("sprank: not available in this version of Octave"); |
|
218 #endif |
|
219 |
|
220 return retval; |
|
221 } |
|
222 |
|
223 /* |
|
224 |
|
225 %!error(sprank(1,2)); |
7243
|
226 %!testif HAVE_CXSPARSE |
|
227 %! assert(sprank(speye(20)), 20) |
|
228 %!testif HAVE_CXSPARSE |
|
229 %! assert(sprank([1,0,2,0;2,0,4,0]),2) |
6066
|
230 |
|
231 */ |
5610
|
232 /* |
|
233 ;;; Local Variables: *** |
|
234 ;;; mode: C++ *** |
|
235 ;;; End: *** |
|
236 */ |