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