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