5164
|
1 /* |
|
2 |
|
3 Copyright (C) 2004 David Bateman |
|
4 Copyright (C) 1998-2004 Andy Adler |
|
5 |
|
6 Octave is free software; you can redistribute it and/or modify it |
|
7 under the terms of the GNU General Public License as published by the |
|
8 Free Software Foundation; either version 2, or (at your option) any |
|
9 later version. |
|
10 |
|
11 Octave is distributed in the hope that it will be useful, but WITHOUT |
|
12 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or |
|
13 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License |
|
14 for more details. |
|
15 |
|
16 You should have received a copy of the GNU General Public License |
5307
|
17 along with this program; see the file COPYING. If not, write to the |
|
18 Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, |
|
19 Boston, MA 02110-1301, USA. |
5164
|
20 |
|
21 */ |
|
22 |
|
23 // This is the octave interface to colamd, which bore the copyright given |
|
24 // in the help of the functions. |
|
25 |
|
26 #ifdef HAVE_CONFIG_H |
|
27 #include <config.h> |
|
28 #endif |
|
29 |
|
30 #include <cstdlib> |
|
31 |
|
32 #include <string> |
|
33 #include <vector> |
|
34 |
|
35 #include "ov.h" |
|
36 #include "defun-dld.h" |
|
37 #include "pager.h" |
|
38 #include "ov-re-mat.h" |
|
39 |
|
40 #include "ov-re-sparse.h" |
|
41 #include "ov-cx-sparse.h" |
|
42 |
|
43 // External COLAMD functions in C |
|
44 extern "C" { |
|
45 #include "COLAMD/colamd.h" |
|
46 } |
|
47 |
5440
|
48 #ifdef IDX_TYPE_LONG |
|
49 #define COLAMD_NAME(name) colamd_l ## name |
|
50 #define SYMAMD_NAME(name) symamd_l ## name |
|
51 #else |
|
52 #define COLAMD_NAME(name) colamd ## name |
|
53 #define SYMAMD_NAME(name) symamd ## name |
|
54 #endif |
|
55 |
5164
|
56 // The symmetric column elimination tree code take from the Davis LDL code. |
|
57 // Copyright given elsewhere in this file. |
|
58 static |
5440
|
59 void symetree (const octave_idx_type *ridx, const octave_idx_type *cidx, |
|
60 octave_idx_type *Parent, octave_idx_type *P, octave_idx_type n) |
5164
|
61 { |
5440
|
62 OCTAVE_LOCAL_BUFFER (octave_idx_type, Flag, n); |
|
63 OCTAVE_LOCAL_BUFFER (octave_idx_type, Pinv, (P ? n : 0)); |
5164
|
64 if (P) |
|
65 // If P is present then compute Pinv, the inverse of P |
5440
|
66 for (octave_idx_type k = 0 ; k < n ; k++) |
5164
|
67 Pinv [P [k]] = k ; |
|
68 |
5440
|
69 for (octave_idx_type k = 0 ; k < n ; k++) |
5164
|
70 { |
|
71 // L(k,:) pattern: all nodes reachable in etree from nz in A(0:k-1,k) |
|
72 Parent [k] = n ; // parent of k is not yet known |
|
73 Flag [k] = k ; // mark node k as visited |
5440
|
74 octave_idx_type kk = (P) ? (P [k]) : (k) ; // kth original, or permuted, column |
|
75 octave_idx_type p2 = cidx [kk+1] ; |
|
76 for (octave_idx_type p = cidx [kk] ; p < p2 ; p++) |
5164
|
77 { |
|
78 // A (i,k) is nonzero (original or permuted A) |
5440
|
79 octave_idx_type i = (Pinv) ? (Pinv [ridx [p]]) : (ridx [p]) ; |
5164
|
80 if (i < k) |
|
81 { |
|
82 // follow path from i to root of etree, stop at flagged node |
|
83 for ( ; Flag [i] != k ; i = Parent [i]) |
|
84 { |
|
85 // find parent of i if not yet determined |
|
86 if (Parent [i] == n) |
|
87 Parent [i] = k ; |
|
88 Flag [i] = k ; // mark i as visited |
|
89 } |
|
90 } |
|
91 } |
|
92 } |
|
93 } |
|
94 |
|
95 // The elimination tree post-ordering code below is taken from SuperLU |
|
96 static inline |
5440
|
97 octave_idx_type make_set (octave_idx_type i, octave_idx_type *pp) |
5164
|
98 { |
|
99 pp[i] = i; |
|
100 return i; |
|
101 } |
|
102 |
|
103 static inline |
5440
|
104 octave_idx_type link (octave_idx_type s, octave_idx_type t, octave_idx_type *pp) |
5164
|
105 { |
|
106 pp[s] = t; |
|
107 return t; |
|
108 } |
|
109 |
|
110 static inline |
5440
|
111 octave_idx_type find (octave_idx_type i, octave_idx_type *pp) |
5164
|
112 { |
5440
|
113 register octave_idx_type p, gp; |
5164
|
114 |
|
115 p = pp[i]; |
|
116 gp = pp[p]; |
|
117 while (gp != p) { |
|
118 pp[i] = gp; |
|
119 i = gp; |
|
120 p = pp[i]; |
|
121 gp = pp[p]; |
|
122 } |
|
123 return (p); |
|
124 } |
|
125 |
|
126 static |
5440
|
127 octave_idx_type etdfs (octave_idx_type v, octave_idx_type *first_kid, |
|
128 octave_idx_type *next_kid, octave_idx_type *post, |
|
129 octave_idx_type postnum) |
5164
|
130 { |
5440
|
131 for (octave_idx_type w = first_kid[v]; w != -1; w = next_kid[w]) { |
5164
|
132 postnum = etdfs (w, first_kid, next_kid, post, postnum); |
|
133 } |
|
134 post[postnum++] = v; |
|
135 |
|
136 return postnum; |
|
137 } |
|
138 |
|
139 static |
5440
|
140 void TreePostorder(octave_idx_type n, octave_idx_type *parent, octave_idx_type *post) |
5164
|
141 { |
|
142 // Allocate storage for working arrays and results |
5440
|
143 OCTAVE_LOCAL_BUFFER (octave_idx_type, first_kid, n+1); |
|
144 OCTAVE_LOCAL_BUFFER (octave_idx_type, next_kid, n+1); |
5164
|
145 |
|
146 // Set up structure describing children |
5440
|
147 for (octave_idx_type v = 0; v <= n; first_kid[v++] = -1); |
|
148 for (octave_idx_type v = n-1; v >= 0; v--) |
5164
|
149 { |
5440
|
150 octave_idx_type dad = parent[v]; |
5164
|
151 next_kid[v] = first_kid[dad]; |
|
152 first_kid[dad] = v; |
|
153 } |
|
154 |
|
155 // Depth-first search from dummy root vertex #n |
|
156 etdfs (n, first_kid, next_kid, post, 0); |
|
157 } |
|
158 |
|
159 static |
5440
|
160 void coletree (const octave_idx_type *ridx, const octave_idx_type *colbeg, |
|
161 octave_idx_type *colend, octave_idx_type *parent, |
|
162 octave_idx_type nr, octave_idx_type nc) |
5164
|
163 { |
5440
|
164 OCTAVE_LOCAL_BUFFER (octave_idx_type, root, nc); |
|
165 OCTAVE_LOCAL_BUFFER (octave_idx_type, pp, nc); |
|
166 OCTAVE_LOCAL_BUFFER (octave_idx_type, firstcol, nr); |
5164
|
167 |
|
168 // Compute firstcol[row] = first nonzero column in row |
5440
|
169 for (octave_idx_type row = 0; row < nr; firstcol[row++] = nc); |
|
170 for (octave_idx_type col = 0; col < nc; col++) |
|
171 for (octave_idx_type p = colbeg[col]; p < colend[col]; p++) |
5164
|
172 { |
5440
|
173 octave_idx_type row = ridx[p]; |
5164
|
174 if (firstcol[row] > col) |
|
175 firstcol[row] = col; |
|
176 } |
|
177 |
|
178 // Compute etree by Liu's algorithm for symmetric matrices, |
|
179 // except use (firstcol[r],c) in place of an edge (r,c) of A. |
|
180 // Thus each row clique in A'*A is replaced by a star |
|
181 // centered at its first vertex, which has the same fill. |
5440
|
182 for (octave_idx_type col = 0; col < nc; col++) |
5164
|
183 { |
5440
|
184 octave_idx_type cset = make_set (col, pp); |
5164
|
185 root[cset] = col; |
|
186 parent[col] = nc; |
5440
|
187 for (octave_idx_type p = colbeg[col]; p < colend[col]; p++) |
5164
|
188 { |
5440
|
189 octave_idx_type row = firstcol[ridx[p]]; |
5164
|
190 if (row >= col) |
|
191 continue; |
5440
|
192 octave_idx_type rset = find (row, pp); |
|
193 octave_idx_type rroot = root[rset]; |
5164
|
194 if (rroot != col) |
|
195 { |
|
196 parent[rroot] = col; |
|
197 cset = link (cset, rset, pp); |
|
198 root[cset] = col; |
|
199 } |
|
200 } |
|
201 } |
|
202 } |
|
203 |
|
204 DEFUN_DLD (colamd, args, nargout, |
|
205 "-*- texinfo -*-\n\ |
|
206 @deftypefn {Loadable Function} {@var{p} =} colamd (@var{s})\n\ |
|
207 @deftypefnx {Loadable Function} {@var{p} =} colamd (@var{s}, @var{knobs})\n\ |
|
208 @deftypefnx {Loadable Function} {[@var{p}, @var{stats}] =} colamd (@var{s})\n\ |
|
209 @deftypefnx {Loadable Function} {[@var{p}, @var{stats}] =} colamd (@var{s}, @var{knobs})\n\ |
|
210 \n\ |
|
211 Column approximate minimum degree permutation. @code{@var{p} = colamd\n\ |
|
212 (@var{s})} returns the column approximate minimum degree permutation\n\ |
|
213 vector for the sparse matrix @var{s}. For a non-symmetric matrix @var{s},\n\ |
|
214 @code{@var{s} (:,@var{p})} tends to have sparser LU factors than @var{s}.\n\ |
|
215 The Cholesky factorization of @code{@var{s} (:,@var{p})' * @var{s}\n\ |
|
216 (:,@var{p})} also tends to be sparser than that of @code{@var{s}' *\n\ |
|
217 @var{s}}.\n\ |
|
218 \n\ |
5440
|
219 @var{knobs} is an optional one- to three-element input vector. If @var{s} is\n\ |
|
220 m-by-n, then rows with more than @code{max(16,@var{knobs}(1)*sqrt(n))} entries\n\ |
|
221 are ignored. Columns with more than @code{max(16,knobs(2)*sqrt(min(m,n)))}\n\ |
|
222 entries are removed prior to ordering, and ordered last in the output\n\ |
|
223 permutation @var{p}. Only completely dense rows or columns are removed\n\ |
|
224 if @code{@var{knobs} (1)} and @code{@var{knobs} (2)} are < 0, respectively.\n\ |
|
225 If @code{@var{knobs} (3)} is nonzero, @var{stats} and @var{knobs} are\n\ |
|
226 printed. The default is @code{@var{knobs} = [10 10 0]}. Note that\n\ |
|
227 @var{knobs} differs from earlier versions of colamd\n\ |
5164
|
228 \n\ |
|
229 @var{stats} is an optional 20-element output vector that provides data\n\ |
|
230 about the ordering and the validity of the input matrix @var{s}. Ordering\n\ |
|
231 statistics are in @code{@var{stats} (1:3)}. @code{@var{stats} (1)} and\n\ |
|
232 @code{@var{stats} (2)} are the number of dense or empty rows and columns\n\ |
|
233 ignored by COLAMD and @code{@var{stats} (3)} is the number of garbage\n\ |
|
234 collections performed on the internal data structure used by COLAMD\n\ |
|
235 (roughly of size @code{2.2 * nnz(@var{s}) + 4 * @var{m} + 7 * @var{n}}\n\ |
|
236 integers).\n\ |
|
237 \n\ |
|
238 Octave built-in functions are intended to generate valid sparse matrices,\n\ |
|
239 with no duplicate entries, with ascending row indices of the nonzeros\n\ |
|
240 in each column, with a non-negative number of entries in each column (!)\n\ |
|
241 and so on. If a matrix is invalid, then COLAMD may or may not be able\n\ |
|
242 to continue. If there are duplicate entries (a row index appears two or\n\ |
|
243 more times in the same column) or if the row indices in a column are out\n\ |
|
244 of order, then COLAMD can correct these errors by ignoring the duplicate\n\ |
|
245 entries and sorting each column of its internal copy of the matrix\n\ |
|
246 @var{s} (the input matrix @var{s} is not repaired, however). If a matrix\n\ |
|
247 is invalid in other ways then COLAMD cannot continue, an error message is\n\ |
|
248 printed, and no output arguments (@var{p} or @var{stats}) are returned.\n\ |
|
249 COLAMD is thus a simple way to check a sparse matrix to see if it's\n\ |
|
250 valid.\n\ |
|
251 \n\ |
|
252 @code{@var{stats} (4:7)} provide information if COLAMD was able to\n\ |
|
253 continue. The matrix is OK if @code{@var{stats} (4)} is zero, or 1 if\n\ |
|
254 invalid. @code{@var{stats} (5)} is the rightmost column index that is\n\ |
|
255 unsorted or contains duplicate entries, or zero if no such column exists.\n\ |
|
256 @code{@var{stats} (6)} is the last seen duplicate or out-of-order row\n\ |
|
257 index in the column index given by @code{@var{stats} (5)}, or zero if no\n\ |
|
258 such row index exists. @code{@var{stats} (7)} is the number of duplicate\n\ |
|
259 or out-of-order row indices. @code{@var{stats} (8:20)} is always zero in\n\ |
|
260 the current version of COLAMD (reserved for future use).\n\ |
|
261 \n\ |
|
262 The ordering is followed by a column elimination tree post-ordering.\n\ |
|
263 \n\ |
|
264 The authors of the code itself are Stefan I. Larimore and Timothy A.\n\ |
|
265 Davis (davis@@cise.ufl.edu), University of Florida. The algorithm was\n\ |
|
266 developed in collaboration with John Gilbert, Xerox PARC, and Esmond\n\ |
|
267 Ng, Oak Ridge National Laboratory. (see\n\ |
|
268 @url{http://www.cise.ufl.edu/research/sparse/colamd})\n\ |
|
269 @end deftypefn\n\ |
|
270 @seealso{colperm, symamd}") |
|
271 { |
|
272 octave_value_list retval; |
|
273 int nargin = args.length (); |
|
274 int spumoni = 0; |
|
275 |
|
276 if (nargout < 0 || nargout > 2 || nargin < 0 || nargin > 2) |
|
277 usage ("colamd: incorrect number of input and/or output arguments"); |
|
278 else |
|
279 { |
|
280 // Get knobs |
5440
|
281 OCTAVE_LOCAL_BUFFER (double, knobs, COLAMD_KNOBS); |
|
282 COLAMD_NAME (_set_defaults) (knobs); |
5164
|
283 |
|
284 // Check for user-passed knobs |
|
285 if (nargin == 2) |
|
286 { |
|
287 NDArray User_knobs = args(1).array_value (); |
|
288 int nel_User_knobs = User_knobs.length (); |
|
289 |
|
290 if (nel_User_knobs > 0) |
5440
|
291 knobs [COLAMD_DENSE_ROW] = User_knobs (0); |
5164
|
292 if (nel_User_knobs > 1) |
5440
|
293 knobs [COLAMD_DENSE_COL] = User_knobs (1) ; |
5164
|
294 if (nel_User_knobs > 2) |
|
295 spumoni = (int) User_knobs (2); |
5440
|
296 |
|
297 // print knob settings if spumoni is set |
|
298 if (spumoni) |
|
299 { |
|
300 |
|
301 octave_stdout << "\ncolamd version " << COLAMD_MAIN_VERSION << "." |
|
302 << COLAMD_SUB_VERSION << ", " << COLAMD_DATE << ":\n"; |
5164
|
303 |
5440
|
304 if (knobs [COLAMD_DENSE_ROW] >= 0) |
|
305 octave_stdout << "knobs(1): " << User_knobs (0) |
|
306 << ", rows with > max(16," |
|
307 << knobs [COLAMD_DENSE_ROW] << "*sqrt(size(A,2)))" |
|
308 << " entries removed\n"; |
|
309 else |
|
310 octave_stdout << "knobs(1): " << User_knobs (0) |
|
311 << ", only completely dense rows removed\n"; |
|
312 |
|
313 if (knobs [COLAMD_DENSE_COL] >= 0) |
|
314 octave_stdout << "knobs(2): " << User_knobs (1) |
|
315 << ", cols with > max(16," |
|
316 << knobs [COLAMD_DENSE_COL] << "*sqrt(size(A)))" |
|
317 << " entries removed\n"; |
|
318 else |
|
319 octave_stdout << "knobs(2): " << User_knobs (1) |
|
320 << ", only completely dense columns removed\n"; |
|
321 |
|
322 octave_stdout << "knobs(3): " << User_knobs (2) |
|
323 << ", statistics and knobs printed\n"; |
|
324 |
|
325 } |
5164
|
326 } |
|
327 |
5440
|
328 octave_idx_type n_row, n_col, nnz; |
|
329 octave_idx_type *ridx, *cidx; |
5164
|
330 SparseComplexMatrix scm; |
|
331 SparseMatrix sm; |
|
332 |
|
333 if (args(0).class_name () == "sparse") |
|
334 { |
|
335 if (args(0).is_complex_type ()) |
|
336 { |
|
337 scm = args(0). sparse_complex_matrix_value (); |
|
338 n_row = scm.rows (); |
|
339 n_col = scm.cols (); |
|
340 nnz = scm.nnz (); |
|
341 ridx = scm.xridx (); |
|
342 cidx = scm.xcidx (); |
|
343 } |
|
344 else |
|
345 { |
|
346 sm = args(0).sparse_matrix_value (); |
|
347 |
|
348 n_row = sm.rows (); |
|
349 n_col = sm.cols (); |
|
350 nnz = sm.nnz (); |
|
351 ridx = sm.xridx (); |
|
352 cidx = sm.xcidx (); |
|
353 } |
|
354 } |
|
355 else |
|
356 { |
|
357 if (args(0).is_complex_type ()) |
|
358 sm = SparseMatrix (real (args(0).complex_matrix_value ())); |
|
359 else |
|
360 sm = SparseMatrix (args(0).matrix_value ()); |
|
361 |
|
362 n_row = sm.rows (); |
|
363 n_col = sm.cols (); |
|
364 nnz = sm.nnz (); |
|
365 ridx = sm.xridx (); |
|
366 cidx = sm.xcidx (); |
|
367 } |
|
368 |
|
369 // Allocate workspace for colamd |
5440
|
370 OCTAVE_LOCAL_BUFFER (octave_idx_type, p, n_col+1); |
|
371 for (octave_idx_type i = 0; i < n_col+1; i++) |
5164
|
372 p[i] = cidx [i]; |
|
373 |
5440
|
374 octave_idx_type Alen = COLAMD_NAME (_recommended) (nnz, n_row, n_col); |
|
375 OCTAVE_LOCAL_BUFFER (octave_idx_type, A, Alen); |
|
376 for (octave_idx_type i = 0; i < nnz; i++) |
5164
|
377 A[i] = ridx [i]; |
|
378 |
|
379 // Order the columns (destroys A) |
5440
|
380 OCTAVE_LOCAL_BUFFER (octave_idx_type, stats, COLAMD_STATS); |
|
381 if (! COLAMD_NAME () (n_row, n_col, Alen, A, p, knobs, stats)) |
5164
|
382 { |
5440
|
383 COLAMD_NAME (_report) (stats) ; |
5164
|
384 error ("colamd: internal error!"); |
|
385 return retval; |
|
386 } |
|
387 |
|
388 // column elimination tree post-ordering (reuse variables) |
5440
|
389 OCTAVE_LOCAL_BUFFER (octave_idx_type, colbeg, n_col + 1); |
|
390 OCTAVE_LOCAL_BUFFER (octave_idx_type, colend, n_col + 1); |
|
391 OCTAVE_LOCAL_BUFFER (octave_idx_type, etree, n_col + 1); |
5164
|
392 |
5440
|
393 for (octave_idx_type i = 0; i < n_col; i++) |
5164
|
394 { |
|
395 colbeg[i] = cidx[p[i]]; |
|
396 colend[i] = cidx[p[i]+1]; |
|
397 } |
|
398 |
|
399 coletree (ridx, colbeg, colend, etree, n_row, n_col); |
|
400 |
|
401 // Calculate the tree post-ordering |
|
402 TreePostorder (n_col, etree, colbeg); |
|
403 |
|
404 // return the permutation vector |
|
405 NDArray out_perm (dim_vector (1, n_col)); |
5440
|
406 for (octave_idx_type i = 0; i < n_col; i++) |
5164
|
407 out_perm(i) = p [colbeg [i]] + 1; |
|
408 |
|
409 retval (0) = out_perm; |
|
410 |
|
411 // print stats if spumoni > 0 |
|
412 if (spumoni > 0) |
5440
|
413 COLAMD_NAME (_report) (stats) ; |
5164
|
414 |
|
415 // Return the stats vector |
|
416 if (nargout == 2) |
|
417 { |
|
418 NDArray out_stats (dim_vector (1, COLAMD_STATS)); |
5440
|
419 for (octave_idx_type i = 0 ; i < COLAMD_STATS ; i++) |
5164
|
420 out_stats (i) = stats [i] ; |
|
421 retval(1) = out_stats; |
|
422 |
|
423 // fix stats (5) and (6), for 1-based information on |
|
424 // jumbled matrix. note that this correction doesn't |
|
425 // occur if symamd returns FALSE |
|
426 out_stats (COLAMD_INFO1) ++ ; |
|
427 out_stats (COLAMD_INFO2) ++ ; |
|
428 } |
|
429 } |
|
430 |
|
431 return retval; |
|
432 } |
|
433 |
|
434 DEFUN_DLD (symamd, args, nargout, |
|
435 "-*- texinfo -*-\n\ |
|
436 @deftypefn {Loadable Function} {@var{p} =} symamd (@var{s})\n\ |
|
437 @deftypefnx {Loadable Function} {@var{p} =} symamd (@var{s}, @var{knobs})\n\ |
|
438 @deftypefnx {Loadable Function} {[@var{p}, @var{stats}] =} symamd (@var{s})\n\ |
|
439 @deftypefnx {Loadable Function} {[@var{p}, @var{stats}] =} symamd (@var{s}, @var{knobs})\n\ |
|
440 \n\ |
|
441 For a symmetric positive definite matrix @var{s}, returns the permutation\n\ |
|
442 vector p such that @code{@var{s} (@var{p}, @var{p})} tends to have a\n\ |
|
443 sparser Cholesky factor than @var{s}. Sometimes SYMAMD works well for\n\ |
|
444 symmetric indefinite matrices too. The matrix @var{s} is assumed to be\n\ |
|
445 symmetric; only the strictly lower triangular part is referenced. @var{s}\n\ |
|
446 must be square.\n\ |
|
447 \n\ |
5440
|
448 @var{knobs} is an optional one- to two-element input vector. If @var{s} is\n\ |
|
449 n-by-n, then rows and columns with more than\n\ |
|
450 @code{max(16,@var{knobs}(1)*sqrt(n))} entries are removed prior to ordering,\n\ |
|
451 and ordered last in the output permutation @var{p}. No rows/columns are\n\ |
|
452 removed if @code{@var{knobs}(1) < 0}. If @code{@var{knobs} (2)} is nonzero,\n\ |
|
453 @code{stats} and @var{knobs} are printed. The default is @code{@var{knobs} \n\ |
|
454 = [10 0]}. Note that @var{knobs} differs from earlier versions of symamd.\n\ |
5164
|
455 \n\ |
|
456 @var{stats} is an optional 20-element output vector that provides data\n\ |
|
457 about the ordering and the validity of the input matrix @var{s}. Ordering\n\ |
|
458 statistics are in @code{@var{stats} (1:3)}. @code{@var{stats} (1) =\n\ |
|
459 @var{stats} (2)} is the number of dense or empty rows and columns\n\ |
|
460 ignored by SYMAMD and @code{@var{stats} (3)} is the number of garbage\n\ |
|
461 collections performed on the internal data structure used by SYMAMD\n\ |
|
462 (roughly of size @code{8.4 * nnz (tril (@var{s}, -1)) + 9 * @var{n}}\n\ |
|
463 integers).\n\ |
|
464 \n\ |
|
465 Octave built-in functions are intended to generate valid sparse matrices,\n\ |
|
466 with no duplicate entries, with ascending row indices of the nonzeros\n\ |
|
467 in each column, with a non-negative number of entries in each column (!)\n\ |
|
468 and so on. If a matrix is invalid, then SYMAMD may or may not be able\n\ |
|
469 to continue. If there are duplicate entries (a row index appears two or\n\ |
|
470 more times in the same column) or if the row indices in a column are out\n\ |
|
471 of order, then SYMAMD can correct these errors by ignoring the duplicate\n\ |
|
472 entries and sorting each column of its internal copy of the matrix S (the\n\ |
|
473 input matrix S is not repaired, however). If a matrix is invalid in\n\ |
|
474 other ways then SYMAMD cannot continue, an error message is printed, and\n\ |
|
475 no output arguments (@var{p} or @var{stats}) are returned. SYMAMD is\n\ |
|
476 thus a simple way to check a sparse matrix to see if it's valid.\n\ |
|
477 \n\ |
|
478 @code{@var{stats} (4:7)} provide information if SYMAMD was able to\n\ |
|
479 continue. The matrix is OK if @code{@var{stats} (4)} is zero, or 1\n\ |
|
480 if invalid. @code{@var{stats} (5)} is the rightmost column index that\n\ |
|
481 is unsorted or contains duplicate entries, or zero if no such column\n\ |
|
482 exists. @code{@var{stats} (6)} is the last seen duplicate or out-of-order\n\ |
|
483 row index in the column index given by @code{@var{stats} (5)}, or zero\n\ |
|
484 if no such row index exists. @code{@var{stats} (7)} is the number of\n\ |
|
485 duplicate or out-of-order row indices. @code{@var{stats} (8:20)} is\n\ |
|
486 always zero in the current version of SYMAMD (reserved for future use).\n\ |
|
487 \n\ |
|
488 The ordering is followed by a column elimination tree post-ordering.\n\ |
|
489 \n\ |
|
490 \n\ |
|
491 The authors of the code itself are Stefan I. Larimore and Timothy A.\n\ |
|
492 Davis (davis@@cise.ufl.edu), University of Florida. The algorithm was\n\ |
|
493 developed in collaboration with John Gilbert, Xerox PARC, and Esmond\n\ |
|
494 Ng, Oak Ridge National Laboratory. (see\n\ |
|
495 @url{http://www.cise.ufl.edu/research/sparse/colamd})\n\ |
|
496 @end deftypefn\n\ |
|
497 @seealso{colperm, colamd}") |
|
498 { |
|
499 octave_value_list retval; |
|
500 int nargin = args.length (); |
|
501 int spumoni = 0; |
|
502 |
|
503 if (nargout < 0 || nargout > 2 || nargin < 0 || nargin > 2) |
|
504 usage ("symamd: incorrect number of input and/or output arguments"); |
|
505 else |
|
506 { |
|
507 // Get knobs |
|
508 OCTAVE_LOCAL_BUFFER (double, knobs, COLAMD_KNOBS); |
5440
|
509 COLAMD_NAME (_set_defaults) (knobs); |
5164
|
510 |
|
511 // Check for user-passed knobs |
|
512 if (nargin == 2) |
|
513 { |
|
514 NDArray User_knobs = args(1).array_value (); |
|
515 int nel_User_knobs = User_knobs.length (); |
|
516 |
|
517 if (nel_User_knobs > 0) |
|
518 knobs [COLAMD_DENSE_ROW] = User_knobs (COLAMD_DENSE_ROW); |
|
519 if (nel_User_knobs > 1) |
|
520 spumoni = (int) User_knobs (1); |
|
521 } |
|
522 |
|
523 // print knob settings if spumoni is set |
|
524 if (spumoni > 0) |
|
525 octave_stdout << "symamd: dense row/col fraction: " |
|
526 << knobs [COLAMD_DENSE_ROW] << std::endl; |
|
527 |
5440
|
528 octave_idx_type n_row, n_col, nnz; |
|
529 octave_idx_type *ridx, *cidx; |
5164
|
530 SparseMatrix sm; |
|
531 SparseComplexMatrix scm; |
|
532 |
|
533 if (args(0).class_name () == "sparse") |
|
534 { |
|
535 if (args(0).is_complex_type ()) |
|
536 { |
|
537 scm = args(0).sparse_complex_matrix_value (); |
|
538 n_row = scm.rows (); |
|
539 n_col = scm.cols (); |
|
540 nnz = scm.nnz (); |
|
541 ridx = scm.xridx (); |
|
542 cidx = scm.xcidx (); |
|
543 } |
|
544 else |
|
545 { |
|
546 sm = args(0).sparse_matrix_value (); |
|
547 n_row = sm.rows (); |
|
548 n_col = sm.cols (); |
|
549 nnz = sm.nnz (); |
|
550 ridx = sm.xridx (); |
|
551 cidx = sm.xcidx (); |
|
552 } |
|
553 } |
|
554 else |
|
555 { |
|
556 if (args(0).is_complex_type ()) |
|
557 sm = SparseMatrix (real (args(0).complex_matrix_value ())); |
|
558 else |
|
559 sm = SparseMatrix (args(0).matrix_value ()); |
|
560 |
|
561 n_row = sm.rows (); |
|
562 n_col = sm.cols (); |
|
563 nnz = sm.nnz (); |
|
564 ridx = sm.xridx (); |
|
565 cidx = sm.xcidx (); |
|
566 } |
|
567 |
|
568 if (n_row != n_col) |
|
569 { |
|
570 error ("symamd: matrix must be square"); |
|
571 return retval; |
|
572 } |
|
573 |
|
574 // Allocate workspace for symamd |
5440
|
575 OCTAVE_LOCAL_BUFFER (octave_idx_type, perm, n_col+1); |
|
576 OCTAVE_LOCAL_BUFFER (octave_idx_type, stats, COLAMD_STATS); |
|
577 if (!SYMAMD_NAME () (n_col, ridx, cidx, perm, knobs, stats, &calloc, &free)) |
5164
|
578 { |
5440
|
579 SYMAMD_NAME (_report) (stats) ; |
5164
|
580 error ("symamd: internal error!") ; |
|
581 return retval; |
|
582 } |
|
583 |
|
584 // column elimination tree post-ordering |
5440
|
585 OCTAVE_LOCAL_BUFFER (octave_idx_type, etree, n_col + 1); |
5164
|
586 symetree (ridx, cidx, etree, perm, n_col); |
|
587 |
|
588 // Calculate the tree post-ordering |
5440
|
589 OCTAVE_LOCAL_BUFFER (octave_idx_type, post, n_col + 1); |
5164
|
590 TreePostorder (n_col, etree, post); |
|
591 |
|
592 // return the permutation vector |
|
593 NDArray out_perm (dim_vector (1, n_col)); |
5440
|
594 for (octave_idx_type i = 0; i < n_col; i++) |
5164
|
595 out_perm(i) = perm [post [i]] + 1; |
|
596 |
|
597 retval (0) = out_perm; |
|
598 |
|
599 // print stats if spumoni > 0 |
|
600 if (spumoni > 0) |
5440
|
601 SYMAMD_NAME (_report) (stats) ; |
5164
|
602 |
|
603 // Return the stats vector |
|
604 if (nargout == 2) |
|
605 { |
|
606 NDArray out_stats (dim_vector (1, COLAMD_STATS)); |
5440
|
607 for (octave_idx_type i = 0 ; i < COLAMD_STATS ; i++) |
5164
|
608 out_stats (i) = stats [i] ; |
|
609 retval(1) = out_stats; |
|
610 |
|
611 // fix stats (5) and (6), for 1-based information on |
|
612 // jumbled matrix. note that this correction doesn't |
|
613 // occur if symamd returns FALSE |
|
614 out_stats (COLAMD_INFO1) ++ ; |
|
615 out_stats (COLAMD_INFO2) ++ ; |
|
616 } |
|
617 } |
|
618 |
|
619 return retval; |
|
620 } |
|
621 |
|
622 DEFUN_DLD (etree, args, nargout, |
|
623 "-*- texinfo -*-\n\ |
|
624 @deftypefn {Loadable Function} {@var{p} =} etree (@var{s})\n\ |
|
625 @deftypefnx {Loadable Function} {@var{p} =} etree (@var{s}, @var{typ})\n\ |
|
626 @deftypefnx {Loadable Function} {[@var{p}, @var{q}] =} etree (@var{s}, @var{typ})\n\ |
|
627 \n\ |
|
628 Returns the elimination tree for the matrix @var{s}. By default @var{s}\n\ |
|
629 is assumed to be symmetric and the symmetric elimination tree is\n\ |
|
630 returned. The argument @var{typ} controls whether a symmetric or\n\ |
|
631 column elimination tree is returned. Valid values of @var{typ} are\n\ |
|
632 'sym' or 'col', for symmetric or column elimination tree respectively\n\ |
|
633 \n\ |
|
634 Called with a second argument, @dfn{etree} also returns the postorder\n\ |
|
635 permutations on the tree.\n\ |
|
636 @end deftypefn") |
|
637 { |
|
638 octave_value_list retval; |
5297
|
639 |
5164
|
640 int nargin = args.length (); |
|
641 |
|
642 if (nargout < 0 || nargout > 2 || nargin < 0 || nargin > 2) |
|
643 usage ("etree: incorrect number of input and/or output arguments"); |
|
644 else |
|
645 { |
5440
|
646 octave_idx_type n_row, n_col, nnz; |
|
647 octave_idx_type *ridx, *cidx; |
5164
|
648 bool is_sym = true; |
|
649 SparseMatrix sm; |
|
650 SparseComplexMatrix scm; |
|
651 |
|
652 if (args(0).class_name () == "sparse") |
|
653 { |
|
654 if (args(0).is_complex_type ()) |
|
655 { |
|
656 scm = args(0).sparse_complex_matrix_value (); |
|
657 n_row = scm.rows (); |
|
658 n_col = scm.cols (); |
|
659 nnz = scm.nnz (); |
|
660 ridx = scm.xridx (); |
|
661 cidx = scm.xcidx (); |
|
662 } |
|
663 else |
|
664 { |
|
665 sm = args(0).sparse_matrix_value (); |
|
666 n_row = sm.rows (); |
|
667 n_col = sm.cols (); |
|
668 nnz = sm.nnz (); |
|
669 ridx = sm.xridx (); |
|
670 cidx = sm.xcidx (); |
|
671 } |
|
672 |
|
673 } |
|
674 else |
|
675 { |
|
676 error ("etree: must be called with a sparse matrix"); |
|
677 return retval; |
|
678 } |
|
679 |
|
680 if (nargin == 2) |
|
681 if (args(1).is_string ()) |
|
682 { |
|
683 std::string str = args(1).string_value (); |
|
684 if (str.find("C") == 0 || str.find("c") == 0) |
|
685 is_sym = false; |
|
686 } |
|
687 else |
|
688 { |
|
689 error ("etree: second argument must be a string"); |
|
690 return retval; |
|
691 } |
|
692 |
|
693 // column elimination tree post-ordering (reuse variables) |
5440
|
694 OCTAVE_LOCAL_BUFFER (octave_idx_type, etree, n_col + 1); |
5164
|
695 |
|
696 |
|
697 if (is_sym) |
|
698 { |
|
699 if (n_row != n_col) |
|
700 { |
|
701 error ("etree: matrix is marked as symmetric, but not square"); |
|
702 return retval; |
|
703 } |
|
704 symetree (ridx, cidx, etree, NULL, n_col); |
|
705 } |
|
706 else |
|
707 { |
5440
|
708 OCTAVE_LOCAL_BUFFER (octave_idx_type, colbeg, n_col); |
|
709 OCTAVE_LOCAL_BUFFER (octave_idx_type, colend, n_col); |
5164
|
710 |
5440
|
711 for (octave_idx_type i = 0; i < n_col; i++) |
5164
|
712 { |
|
713 colbeg[i] = cidx[i]; |
|
714 colend[i] = cidx[i+1]; |
|
715 } |
|
716 |
|
717 coletree (ridx, colbeg, colend, etree, n_row, n_col); |
|
718 } |
|
719 |
|
720 NDArray tree (dim_vector (1, n_col)); |
5440
|
721 for (octave_idx_type i = 0; i < n_col; i++) |
5164
|
722 // We flag a root with n_col while Matlab does it with zero |
|
723 // Convert for matlab compatiable output |
|
724 if (etree[i] == n_col) |
|
725 tree (i) = 0; |
|
726 else |
|
727 tree (i) = etree[i] + 1; |
|
728 |
|
729 retval (0) = tree; |
|
730 |
|
731 if (nargout == 2) |
|
732 { |
|
733 // Calculate the tree post-ordering |
5440
|
734 OCTAVE_LOCAL_BUFFER (octave_idx_type, post, n_col + 1); |
5164
|
735 TreePostorder (n_col, etree, post); |
|
736 |
|
737 NDArray postorder (dim_vector (1, n_col)); |
5440
|
738 for (octave_idx_type i = 0; i < n_col; i++) |
5164
|
739 postorder (i) = post[i] + 1; |
|
740 |
|
741 retval (1) = postorder; |
|
742 } |
|
743 } |
|
744 |
|
745 return retval; |
|
746 } |
|
747 |
|
748 DEFUN_DLD (symbfact, args, nargout, |
|
749 "-*- texinfo -*-\n\ |
|
750 @deftypefn {Loadable Function} {[@var{count}, @var{h}, @var{parent}, @var{post}, @var{r}]} = symbfact (@var{s}, @var{typ})\n\ |
|
751 \n\ |
|
752 Performs a symbolic factorization analysis on the sparse matrix @var{s}.\n\ |
|
753 @end deftypefn") |
|
754 { |
|
755 error ("symbfact: not implemented yet"); |
|
756 return octave_value (); |
|
757 } |
|
758 |
|
759 /* |
|
760 ;;; Local Variables: *** |
|
761 ;;; mode: C++ *** |
|
762 ;;; End: *** |
|
763 */ |