colamd.cc

Go to the documentation of this file.
00001 /*
00002 
00003 Copyright (C) 2004-2012 David Bateman
00004 Copyright (C) 1998-2004 Andy Adler
00005 
00006 This file is part of Octave.
00007 
00008 Octave is free software; you can redistribute it and/or modify it
00009 under the terms of the GNU General Public License as published by the
00010 Free Software Foundation; either version 3 of the License, or (at your
00011 option) any later version.
00012 
00013 Octave is distributed in the hope that it will be useful, but WITHOUT
00014 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00015 FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00016 for more details.
00017 
00018 You should have received a copy of the GNU General Public License
00019 along with Octave; see the file COPYING.  If not, see
00020 <http://www.gnu.org/licenses/>.
00021 
00022 */
00023 
00024 // This is the octave interface to colamd, which bore the copyright given
00025 // in the help of the functions.
00026 
00027 #ifdef HAVE_CONFIG_H
00028 #include <config.h>
00029 #endif
00030 
00031 #include <cstdlib>
00032 
00033 #include <string>
00034 #include <vector>
00035 
00036 #include "ov.h"
00037 #include "defun-dld.h"
00038 #include "pager.h"
00039 #include "ov-re-mat.h"
00040 
00041 #include "ov-re-sparse.h"
00042 #include "ov-cx-sparse.h"
00043 
00044 #include "oct-sparse.h"
00045 #include "oct-locbuf.h"
00046 
00047 #ifdef IDX_TYPE_LONG
00048 #define COLAMD_NAME(name) colamd_l ## name
00049 #define SYMAMD_NAME(name) symamd_l ## name
00050 #else
00051 #define COLAMD_NAME(name) colamd ## name
00052 #define SYMAMD_NAME(name) symamd ## name
00053 #endif
00054 
00055 // The symmetric column elimination tree code take from the Davis LDL code.
00056 // Copyright given elsewhere in this file.
00057 static void
00058 symetree (const octave_idx_type *ridx, const octave_idx_type *cidx,
00059           octave_idx_type *Parent, octave_idx_type *P, octave_idx_type n)
00060 {
00061   OCTAVE_LOCAL_BUFFER (octave_idx_type, Flag, n);
00062   OCTAVE_LOCAL_BUFFER (octave_idx_type, Pinv, (P ? n : 0));
00063   if (P)
00064     // If P is present then compute Pinv, the inverse of P
00065     for (octave_idx_type k = 0 ; k < n ; k++)
00066       Pinv [P [k]] = k ;
00067 
00068   for (octave_idx_type k = 0 ; k < n ; k++)
00069     {
00070       // L(k,:) pattern: all nodes reachable in etree from nz in A(0:k-1,k)
00071       Parent [k] = n ;                // parent of k is not yet known
00072       Flag [k] = k ;                  // mark node k as visited
00073       octave_idx_type kk = (P) ? (P [k]) : (k) ;  // kth original, or permuted, column
00074       octave_idx_type p2 = cidx [kk+1] ;
00075       for (octave_idx_type p = cidx [kk] ; p < p2 ; p++)
00076         {
00077           // A (i,k) is nonzero (original or permuted A)
00078           octave_idx_type i = (Pinv) ? (Pinv [ridx [p]]) : (ridx [p]) ;
00079           if (i < k)
00080             {
00081               // follow path from i to root of etree, stop at flagged node
00082               for ( ; Flag [i] != k ; i = Parent [i])
00083                 {
00084                   // find parent of i if not yet determined
00085                   if (Parent [i] == n)
00086                     Parent [i] = k ;
00087                   Flag [i] = k ;        // mark i as visited
00088                 }
00089             }
00090         }
00091     }
00092 }
00093 
00094 // The elimination tree post-ordering code below is taken from SuperLU
00095 static inline octave_idx_type
00096 make_set (octave_idx_type i, octave_idx_type *pp)
00097 {
00098   pp[i] = i;
00099   return i;
00100 }
00101 
00102 static inline octave_idx_type
00103 link (octave_idx_type s, octave_idx_type t, octave_idx_type *pp)
00104 {
00105   pp[s] = t;
00106   return t;
00107 }
00108 
00109 static inline octave_idx_type
00110 find (octave_idx_type i, octave_idx_type *pp)
00111 {
00112   register octave_idx_type p, gp;
00113 
00114   p = pp[i];
00115   gp = pp[p];
00116 
00117   while (gp != p)
00118     {
00119       pp[i] = gp;
00120       i = gp;
00121       p = pp[i];
00122       gp = pp[p];
00123     }
00124 
00125   return p;
00126 }
00127 
00128 static octave_idx_type
00129 etdfs (octave_idx_type v, octave_idx_type *first_kid,
00130        octave_idx_type *next_kid, octave_idx_type *post,
00131        octave_idx_type postnum)
00132 {
00133   for (octave_idx_type w = first_kid[v]; w != -1; w = next_kid[w])
00134     postnum = etdfs (w, first_kid, next_kid, post, postnum);
00135 
00136   post[postnum++] = v;
00137 
00138   return postnum;
00139 }
00140 
00141 static void
00142 tree_postorder (octave_idx_type n, octave_idx_type *parent,
00143                 octave_idx_type *post)
00144 {
00145   // Allocate storage for working arrays and results
00146   OCTAVE_LOCAL_BUFFER (octave_idx_type, first_kid, n+1);
00147   OCTAVE_LOCAL_BUFFER (octave_idx_type, next_kid, n+1);
00148 
00149   // Set up structure describing children
00150   for (octave_idx_type v = 0; v <= n; first_kid[v++] = -1)
00151     /* do nothing */;
00152 
00153   for (octave_idx_type v = n-1; v >= 0; v--)
00154     {
00155       octave_idx_type dad = parent[v];
00156       next_kid[v] = first_kid[dad];
00157       first_kid[dad] = v;
00158     }
00159 
00160   // Depth-first search from dummy root vertex #n
00161   etdfs (n, first_kid, next_kid, post, 0);
00162 }
00163 
00164 static void
00165 coletree (const octave_idx_type *ridx, const octave_idx_type *colbeg,
00166           octave_idx_type *colend, octave_idx_type *parent,
00167           octave_idx_type nr, octave_idx_type nc)
00168 {
00169   OCTAVE_LOCAL_BUFFER (octave_idx_type, root, nc);
00170   OCTAVE_LOCAL_BUFFER (octave_idx_type, pp, nc);
00171   OCTAVE_LOCAL_BUFFER (octave_idx_type, firstcol, nr);
00172 
00173   // Compute firstcol[row] = first nonzero column in row
00174   for (octave_idx_type row = 0; row < nr; firstcol[row++] = nc)
00175     /* do nothing */;
00176 
00177   for (octave_idx_type col = 0; col < nc; col++)
00178     for (octave_idx_type p = colbeg[col]; p < colend[col]; p++)
00179       {
00180         octave_idx_type row = ridx[p];
00181         if (firstcol[row] > col)
00182           firstcol[row] = col;
00183       }
00184 
00185   // Compute etree by Liu's algorithm for symmetric matrices,
00186   // except use (firstcol[r],c) in place of an edge (r,c) of A.
00187   // Thus each row clique in A'*A is replaced by a star
00188   // centered at its first vertex, which has the same fill.
00189   for (octave_idx_type col = 0; col < nc; col++)
00190     {
00191       octave_idx_type cset = make_set (col, pp);
00192       root[cset] = col;
00193       parent[col] = nc;
00194       for (octave_idx_type p = colbeg[col]; p < colend[col]; p++)
00195         {
00196           octave_idx_type row = firstcol[ridx[p]];
00197           if (row >= col)
00198             continue;
00199           octave_idx_type rset = find (row, pp);
00200           octave_idx_type rroot = root[rset];
00201           if (rroot != col)
00202             {
00203               parent[rroot] = col;
00204               cset = link (cset, rset, pp);
00205               root[cset] = col;
00206             }
00207         }
00208     }
00209 }
00210 
00211 DEFUN_DLD (colamd, args, nargout,
00212     "-*- texinfo -*-\n\
00213 @deftypefn  {Loadable Function} {@var{p} =} colamd (@var{S})\n\
00214 @deftypefnx {Loadable Function} {@var{p} =} colamd (@var{S}, @var{knobs})\n\
00215 @deftypefnx {Loadable Function} {[@var{p}, @var{stats}] =} colamd (@var{S})\n\
00216 @deftypefnx {Loadable Function} {[@var{p}, @var{stats}] =} colamd (@var{S}, @var{knobs})\n\
00217 \n\
00218 Column approximate minimum degree permutation.\n\
00219 @code{@var{p} = colamd (@var{S})} returns the column approximate minimum\n\
00220 degree permutation vector for the sparse matrix @var{S}.  For a\n\
00221 non-symmetric matrix @var{S}, @code{@var{S}(:,@var{p})} tends to have\n\
00222 sparser LU@tie{}factors than @var{S}.  The Cholesky@tie{}factorization of\n\
00223 @code{@var{S}(:,@var{p})' * @var{S}(:,@var{p})} also tends to be sparser\n\
00224 than that of @code{@var{S}' * @var{S}}.\n\
00225 \n\
00226 @var{knobs} is an optional one- to three-element input vector.  If @var{S} is\n\
00227 m-by-n, then rows with more than @code{max(16,@var{knobs}(1)*sqrt(n))}\n\
00228 entries are ignored.  Columns with more than\n\
00229 @code{max(16,@var{knobs}(2)*sqrt(min(m,n)))} entries are removed prior to\n\
00230 ordering, and ordered last in the output permutation @var{p}.  Only\n\
00231 completely dense rows or columns are removed if @code{@var{knobs}(1)} and\n\
00232 @code{@var{knobs}(2)} are < 0, respectively.  If @code{@var{knobs}(3)} is\n\
00233 nonzero, @var{stats} and @var{knobs} are printed.  The default is\n\
00234 @code{@var{knobs} = [10 10 0]}.  Note that @var{knobs} differs from earlier\n\
00235 versions of colamd.\n\
00236 \n\
00237 @var{stats} is an optional 20-element output vector that provides data\n\
00238 about the ordering and the validity of the input matrix @var{S}.  Ordering\n\
00239 statistics are in @code{@var{stats}(1:3)}.  @code{@var{stats}(1)} and\n\
00240 @code{@var{stats}(2)} are the number of dense or empty rows and columns\n\
00241 ignored by @sc{colamd} and @code{@var{stats}(3)} is the number of garbage\n\
00242 collections performed on the internal data structure used by @sc{colamd}\n\
00243 (roughly of size @code{2.2 * nnz(@var{S}) + 4 * @var{m} + 7 * @var{n}}\n\
00244 integers).\n\
00245 \n\
00246 Octave built-in functions are intended to generate valid sparse matrices,\n\
00247 with no duplicate entries, with ascending row indices of the nonzeros\n\
00248 in each column, with a non-negative number of entries in each column (!)\n\
00249 and so on.  If a matrix is invalid, then @sc{colamd} may or may not be able\n\
00250 to continue.  If there are duplicate entries (a row index appears two or\n\
00251 more times in the same column) or if the row indices in a column are out\n\
00252 of order, then @sc{colamd} can correct these errors by ignoring the duplicate\n\
00253 entries and sorting each column of its internal copy of the matrix\n\
00254 @var{S} (the input matrix @var{S} is not repaired, however).  If a matrix\n\
00255 is invalid in other ways then @sc{colamd} cannot continue, an error message\n\
00256 is printed, and no output arguments (@var{p} or @var{stats}) are returned.\n\
00257 @sc{colamd} is thus a simple way to check a sparse matrix to see if it's\n\
00258 valid.\n\
00259 \n\
00260 @code{@var{stats}(4:7)} provide information if COLAMD was able to\n\
00261 continue.  The matrix is OK if @code{@var{stats}(4)} is zero, or 1 if\n\
00262 invalid.  @code{@var{stats}(5)} is the rightmost column index that is\n\
00263 unsorted or contains duplicate entries, or zero if no such column exists.\n\
00264 @code{@var{stats}(6)} is the last seen duplicate or out-of-order row\n\
00265 index in the column index given by @code{@var{stats}(5)}, or zero if no\n\
00266 such row index exists.  @code{@var{stats}(7)} is the number of duplicate\n\
00267 or out-of-order row indices.  @code{@var{stats}(8:20)} is always zero in\n\
00268 the current version of @sc{colamd} (reserved for future use).\n\
00269 \n\
00270 The ordering is followed by a column elimination tree post-ordering.\n\
00271 \n\
00272 The authors of the code itself are Stefan I. Larimore and Timothy A.\n\
00273 Davis @email{davis@@cise.ufl.edu}, University of Florida.  The algorithm was\n\
00274 developed in collaboration with John Gilbert, Xerox PARC, and Esmond\n\
00275 Ng, Oak Ridge National Laboratory.  (see\n\
00276 @url{http://www.cise.ufl.edu/research/sparse/colamd})\n\
00277 @seealso{colperm, symamd, ccolamd}\n\
00278 @end deftypefn")
00279 {
00280   octave_value_list retval;
00281 
00282 #ifdef HAVE_COLAMD
00283 
00284   int nargin = args.length ();
00285   int spumoni = 0;
00286 
00287   if (nargout > 2 || nargin < 1 || nargin > 2)
00288     print_usage ();
00289   else
00290     {
00291       // Get knobs
00292       OCTAVE_LOCAL_BUFFER (double, knobs, COLAMD_KNOBS);
00293       COLAMD_NAME (_set_defaults) (knobs);
00294 
00295       // Check for user-passed knobs
00296       if (nargin == 2)
00297         {
00298           NDArray User_knobs = args(1).array_value ();
00299           int nel_User_knobs = User_knobs.length ();
00300 
00301           if (nel_User_knobs > 0)
00302             knobs [COLAMD_DENSE_ROW] = User_knobs (0);
00303           if (nel_User_knobs > 1)
00304             knobs [COLAMD_DENSE_COL] = User_knobs (1) ;
00305           if (nel_User_knobs > 2)
00306             spumoni = static_cast<int> (User_knobs (2));
00307 
00308           // print knob settings if spumoni is set
00309           if (spumoni)
00310             {
00311 
00312               octave_stdout << "\ncolamd version " << COLAMD_MAIN_VERSION << "."
00313                             <<  COLAMD_SUB_VERSION << ", " << COLAMD_DATE << ":\n";
00314 
00315               if (knobs [COLAMD_DENSE_ROW] >= 0)
00316                 octave_stdout << "knobs(1): " << User_knobs (0)
00317                               << ", rows with > max(16,"
00318                               << knobs [COLAMD_DENSE_ROW] << "*sqrt(size(A,2)))"
00319                               << " entries removed\n";
00320               else
00321                 octave_stdout << "knobs(1): " << User_knobs (0)
00322                               << ", only completely dense rows removed\n";
00323 
00324               if (knobs [COLAMD_DENSE_COL] >= 0)
00325                 octave_stdout << "knobs(2): " << User_knobs (1)
00326                               << ", cols with > max(16,"
00327                               << knobs [COLAMD_DENSE_COL] << "*sqrt(size(A)))"
00328                               << " entries removed\n";
00329               else
00330                 octave_stdout << "knobs(2): " << User_knobs (1)
00331                               << ", only completely dense columns removed\n";
00332 
00333               octave_stdout << "knobs(3): " << User_knobs (2)
00334                             << ", statistics and knobs printed\n";
00335 
00336             }
00337         }
00338 
00339       octave_idx_type n_row, n_col, nnz;
00340       octave_idx_type *ridx, *cidx;
00341       SparseComplexMatrix scm;
00342       SparseMatrix sm;
00343 
00344       if (args(0).is_sparse_type ())
00345         {
00346           if (args(0).is_complex_type ())
00347             {
00348               scm = args(0). sparse_complex_matrix_value ();
00349               n_row = scm.rows ();
00350               n_col = scm.cols ();
00351               nnz = scm.nnz ();
00352               ridx = scm.xridx ();
00353               cidx = scm.xcidx ();
00354             }
00355           else
00356             {
00357               sm = args(0).sparse_matrix_value ();
00358 
00359               n_row = sm.rows ();
00360               n_col = sm.cols ();
00361               nnz = sm.nnz ();
00362               ridx = sm.xridx ();
00363               cidx = sm.xcidx ();
00364             }
00365         }
00366       else
00367         {
00368           if (args(0).is_complex_type ())
00369             sm = SparseMatrix (real (args(0).complex_matrix_value ()));
00370           else
00371             sm = SparseMatrix (args(0).matrix_value ());
00372 
00373           n_row = sm.rows ();
00374           n_col = sm.cols ();
00375           nnz = sm.nnz ();
00376           ridx = sm.xridx ();
00377           cidx = sm.xcidx ();
00378         }
00379 
00380       // Allocate workspace for colamd
00381       OCTAVE_LOCAL_BUFFER (octave_idx_type, p, n_col+1);
00382       for (octave_idx_type i = 0; i < n_col+1; i++)
00383         p[i] = cidx [i];
00384 
00385       octave_idx_type Alen = COLAMD_NAME (_recommended) (nnz, n_row, n_col);
00386       OCTAVE_LOCAL_BUFFER (octave_idx_type, A, Alen);
00387       for (octave_idx_type i = 0; i < nnz; i++)
00388         A[i] = ridx [i];
00389 
00390       // Order the columns (destroys A)
00391       OCTAVE_LOCAL_BUFFER (octave_idx_type, stats, COLAMD_STATS);
00392       if (! COLAMD_NAME () (n_row, n_col, Alen, A, p, knobs, stats))
00393         {
00394           COLAMD_NAME (_report) (stats) ;
00395           error ("colamd: internal error!");
00396           return retval;
00397         }
00398 
00399       // column elimination tree post-ordering (reuse variables)
00400       OCTAVE_LOCAL_BUFFER (octave_idx_type, colbeg, n_col + 1);
00401       OCTAVE_LOCAL_BUFFER (octave_idx_type, colend, n_col + 1);
00402       OCTAVE_LOCAL_BUFFER (octave_idx_type, etree, n_col + 1);
00403 
00404       for (octave_idx_type i = 0; i < n_col; i++)
00405         {
00406           colbeg[i] = cidx[p[i]];
00407           colend[i] = cidx[p[i]+1];
00408         }
00409 
00410       coletree (ridx, colbeg, colend, etree, n_row, n_col);
00411 
00412       // Calculate the tree post-ordering
00413       tree_postorder (n_col, etree, colbeg);
00414 
00415       // return the permutation vector
00416       NDArray out_perm (dim_vector (1, n_col));
00417       for (octave_idx_type i = 0; i < n_col; i++)
00418         out_perm(i) = p [colbeg [i]] + 1;
00419 
00420       retval(0) = out_perm;
00421 
00422       // print stats if spumoni > 0
00423       if (spumoni > 0)
00424         COLAMD_NAME (_report) (stats) ;
00425 
00426       // Return the stats vector
00427       if (nargout == 2)
00428         {
00429           NDArray out_stats (dim_vector (1, COLAMD_STATS));
00430           for (octave_idx_type i = 0 ; i < COLAMD_STATS ; i++)
00431             out_stats (i) = stats [i] ;
00432           retval(1) = out_stats;
00433 
00434           // fix stats (5) and (6), for 1-based information on
00435           // jumbled matrix.  note that this correction doesn't
00436           // occur if symamd returns FALSE
00437           out_stats (COLAMD_INFO1) ++ ;
00438           out_stats (COLAMD_INFO2) ++ ;
00439         }
00440     }
00441 
00442 #else
00443 
00444   error ("colamd: not available in this version of Octave");
00445 
00446 #endif
00447 
00448   return retval;
00449 }
00450 
00451 DEFUN_DLD (symamd, args, nargout,
00452     "-*- texinfo -*-\n\
00453 @deftypefn  {Loadable Function} {@var{p} =} symamd (@var{S})\n\
00454 @deftypefnx {Loadable Function} {@var{p} =} symamd (@var{S}, @var{knobs})\n\
00455 @deftypefnx {Loadable Function} {[@var{p}, @var{stats}] =} symamd (@var{S})\n\
00456 @deftypefnx {Loadable Function} {[@var{p}, @var{stats}] =} symamd (@var{S}, @var{knobs})\n\
00457 \n\
00458 For a symmetric positive definite matrix @var{S}, returns the permutation\n\
00459 vector p such that @code{@var{S}(@var{p}, @var{p})} tends to have a\n\
00460 sparser Cholesky@tie{}factor than @var{S}.  Sometimes @code{symamd} works\n\
00461 well for symmetric indefinite matrices too.  The matrix @var{S} is assumed\n\
00462 to be symmetric; only the strictly lower triangular part is referenced.\n\
00463 @var{S} must be square.\n\
00464 \n\
00465 @var{knobs} is an optional one- to two-element input vector.  If @var{S} is\n\
00466 n-by-n, then rows and columns with more than\n\
00467 @code{max(16,@var{knobs}(1)*sqrt(n))} entries are removed prior to ordering,\n\
00468 and ordered last in the output permutation @var{p}.  No rows/columns are\n\
00469 removed if @code{@var{knobs}(1) < 0}.  If @code{@var{knobs} (2)} is nonzero,\n\
00470 @code{stats} and @var{knobs} are printed.  The default is @code{@var{knobs}\n\
00471 = [10 0]}.  Note that @var{knobs} differs from earlier versions of symamd.\n\
00472 \n\
00473 @var{stats} is an optional 20-element output vector that provides data\n\
00474 about the ordering and the validity of the input matrix @var{S}.  Ordering\n\
00475 statistics are in @code{@var{stats}(1:3)}.  @code{@var{stats}(1) =\n\
00476 @var{stats}(2)} is the number of dense or empty rows and columns\n\
00477 ignored by SYMAMD and @code{@var{stats}(3)} is the number of garbage\n\
00478 collections performed on the internal data structure used by SYMAMD\n\
00479 (roughly of size @code{8.4 * nnz (tril (@var{S}, -1)) + 9 * @var{n}}\n\
00480 integers).\n\
00481 \n\
00482 Octave built-in functions are intended to generate valid sparse matrices,\n\
00483 with no duplicate entries, with ascending row indices of the nonzeros\n\
00484 in each column, with a non-negative number of entries in each column (!)\n\
00485 and so on.  If a matrix is invalid, then SYMAMD may or may not be able\n\
00486 to continue.  If there are duplicate entries (a row index appears two or\n\
00487 more times in the same column) or if the row indices in a column are out\n\
00488 of order, then SYMAMD can correct these errors by ignoring the duplicate\n\
00489 entries and sorting each column of its internal copy of the matrix S (the\n\
00490 input matrix S is not repaired, however).  If a matrix is invalid in\n\
00491 other ways then SYMAMD cannot continue, an error message is printed, and\n\
00492 no output arguments (@var{p} or @var{stats}) are returned.  SYMAMD is\n\
00493 thus a simple way to check a sparse matrix to see if it's valid.\n\
00494 \n\
00495 @code{@var{stats}(4:7)} provide information if SYMAMD was able to\n\
00496 continue.  The matrix is OK if @code{@var{stats} (4)} is zero, or 1\n\
00497 if invalid.  @code{@var{stats}(5)} is the rightmost column index that\n\
00498 is unsorted or contains duplicate entries, or zero if no such column\n\
00499 exists.  @code{@var{stats}(6)} is the last seen duplicate or out-of-order\n\
00500 row index in the column index given by @code{@var{stats}(5)}, or zero\n\
00501 if no such row index exists.  @code{@var{stats}(7)} is the number of\n\
00502 duplicate or out-of-order row indices.  @code{@var{stats}(8:20)} is\n\
00503 always zero in the current version of SYMAMD (reserved for future use).\n\
00504 \n\
00505 The ordering is followed by a column elimination tree post-ordering.\n\
00506 \n\
00507 The authors of the code itself are Stefan I. Larimore and Timothy A.\n\
00508 Davis @email{davis@@cise.ufl.edu}, University of Florida.  The algorithm was\n\
00509 developed in collaboration with John Gilbert, Xerox PARC, and Esmond\n\
00510 Ng, Oak Ridge National Laboratory.  (see\n\
00511 @url{http://www.cise.ufl.edu/research/sparse/colamd})\n\
00512 @seealso{colperm, colamd}\n\
00513 @end deftypefn")
00514 {
00515   octave_value_list retval;
00516 
00517 #ifdef HAVE_COLAMD
00518 
00519   int nargin = args.length ();
00520   int spumoni = 0;
00521 
00522   if (nargout > 2 || nargin < 1 || nargin > 2)
00523     print_usage ();
00524   else
00525     {
00526       // Get knobs
00527       OCTAVE_LOCAL_BUFFER (double, knobs, COLAMD_KNOBS);
00528       COLAMD_NAME (_set_defaults) (knobs);
00529 
00530       // Check for user-passed knobs
00531       if (nargin == 2)
00532         {
00533           NDArray User_knobs = args(1).array_value ();
00534           int nel_User_knobs = User_knobs.length ();
00535 
00536           if (nel_User_knobs > 0)
00537             knobs [COLAMD_DENSE_ROW] = User_knobs (COLAMD_DENSE_ROW);
00538           if (nel_User_knobs > 1)
00539             spumoni = static_cast<int> (User_knobs (1));
00540         }
00541 
00542       // print knob settings if spumoni is set
00543       if (spumoni > 0)
00544         octave_stdout << "symamd: dense row/col fraction: "
00545                       << knobs [COLAMD_DENSE_ROW] << std::endl;
00546 
00547       octave_idx_type n_row, n_col;
00548       octave_idx_type *ridx, *cidx;
00549       SparseMatrix sm;
00550       SparseComplexMatrix scm;
00551 
00552       if (args(0).is_sparse_type ())
00553         {
00554           if (args(0).is_complex_type ())
00555             {
00556               scm = args(0).sparse_complex_matrix_value ();
00557               n_row = scm.rows ();
00558               n_col = scm.cols ();
00559               ridx = scm.xridx ();
00560               cidx = scm.xcidx ();
00561             }
00562           else
00563             {
00564               sm = args(0).sparse_matrix_value ();
00565               n_row = sm.rows ();
00566               n_col = sm.cols ();
00567               ridx = sm.xridx ();
00568               cidx = sm.xcidx ();
00569             }
00570         }
00571       else
00572         {
00573           if (args(0).is_complex_type ())
00574             sm = SparseMatrix (real (args(0).complex_matrix_value ()));
00575           else
00576             sm = SparseMatrix (args(0).matrix_value ());
00577 
00578           n_row = sm.rows ();
00579           n_col = sm.cols ();
00580           ridx = sm.xridx ();
00581           cidx = sm.xcidx ();
00582         }
00583 
00584       if (n_row != n_col)
00585         {
00586           error ("symamd: matrix S must be square");
00587           return retval;
00588         }
00589 
00590       // Allocate workspace for symamd
00591       OCTAVE_LOCAL_BUFFER (octave_idx_type, perm, n_col+1);
00592       OCTAVE_LOCAL_BUFFER (octave_idx_type, stats, COLAMD_STATS);
00593       if (!SYMAMD_NAME () (n_col, ridx, cidx, perm, knobs, stats, &calloc, &free))
00594         {
00595           SYMAMD_NAME (_report) (stats) ;
00596           error ("symamd: internal error!") ;
00597           return retval;
00598         }
00599 
00600       // column elimination tree post-ordering
00601       OCTAVE_LOCAL_BUFFER (octave_idx_type, etree, n_col + 1);
00602       symetree (ridx, cidx, etree, perm, n_col);
00603 
00604       // Calculate the tree post-ordering
00605       OCTAVE_LOCAL_BUFFER (octave_idx_type, post, n_col + 1);
00606       tree_postorder (n_col, etree, post);
00607 
00608       // return the permutation vector
00609       NDArray out_perm (dim_vector (1, n_col));
00610       for (octave_idx_type i = 0; i < n_col; i++)
00611         out_perm(i) = perm [post [i]] + 1;
00612 
00613       retval(0) = out_perm;
00614 
00615       // print stats if spumoni > 0
00616       if (spumoni > 0)
00617         SYMAMD_NAME (_report) (stats) ;
00618 
00619       // Return the stats vector
00620       if (nargout == 2)
00621         {
00622           NDArray out_stats (dim_vector (1, COLAMD_STATS));
00623           for (octave_idx_type i = 0 ; i < COLAMD_STATS ; i++)
00624             out_stats (i) = stats [i] ;
00625           retval(1) = out_stats;
00626 
00627           // fix stats (5) and (6), for 1-based information on
00628           // jumbled matrix.  note that this correction doesn't
00629           // occur if symamd returns FALSE
00630           out_stats (COLAMD_INFO1) ++ ;
00631           out_stats (COLAMD_INFO2) ++ ;
00632         }
00633     }
00634 
00635 #else
00636 
00637   error ("symamd: not available in this version of Octave");
00638 
00639 #endif
00640 
00641   return retval;
00642 }
00643 
00644 DEFUN_DLD (etree, args, nargout,
00645     "-*- texinfo -*-\n\
00646 @deftypefn  {Loadable Function} {@var{p} =} etree (@var{S})\n\
00647 @deftypefnx {Loadable Function} {@var{p} =} etree (@var{S}, @var{typ})\n\
00648 @deftypefnx {Loadable Function} {[@var{p}, @var{q}] =} etree (@var{S}, @var{typ})\n\
00649 \n\
00650 Return the elimination tree for the matrix @var{S}.  By default @var{S}\n\
00651 is assumed to be symmetric and the symmetric elimination tree is\n\
00652 returned.  The argument @var{typ} controls whether a symmetric or\n\
00653 column elimination tree is returned.  Valid values of @var{typ} are\n\
00654 'sym' or 'col', for symmetric or column elimination tree respectively\n\
00655 \n\
00656 Called with a second argument, @code{etree} also returns the postorder\n\
00657 permutations on the tree.\n\
00658 @end deftypefn")
00659 {
00660   octave_value_list retval;
00661 
00662   int nargin = args.length ();
00663 
00664   if (nargout > 2 || nargin < 1 || nargin > 2)
00665     print_usage ();
00666   else
00667     {
00668       octave_idx_type n_row, n_col;
00669       octave_idx_type *ridx, *cidx;
00670       bool is_sym = true;
00671       SparseMatrix sm;
00672       SparseComplexMatrix scm;
00673 
00674       if (args(0).is_sparse_type ())
00675         {
00676           if (args(0).is_complex_type ())
00677             {
00678               scm = args(0).sparse_complex_matrix_value ();
00679               n_row = scm.rows ();
00680               n_col = scm.cols ();
00681               ridx = scm.xridx ();
00682               cidx = scm.xcidx ();
00683             }
00684           else
00685             {
00686               sm = args(0).sparse_matrix_value ();
00687               n_row = sm.rows ();
00688               n_col = sm.cols ();
00689               ridx = sm.xridx ();
00690               cidx = sm.xcidx ();
00691             }
00692 
00693         }
00694       else
00695         {
00696           error ("etree: S must be a sparse matrix");
00697           return retval;
00698         }
00699 
00700       if (nargin == 2)
00701         {
00702           if (args(1).is_string ())
00703             {
00704               std::string str = args(1).string_value ();
00705               if (str.find ("C") == 0 || str.find ("c") == 0)
00706                 is_sym = false;
00707             }
00708           else
00709             {
00710               error ("etree: TYP must be a string");
00711               return retval;
00712             }
00713         }
00714 
00715       // column elimination tree post-ordering (reuse variables)
00716       OCTAVE_LOCAL_BUFFER (octave_idx_type, etree, n_col + 1);
00717 
00718       if (is_sym)
00719         {
00720           if (n_row != n_col)
00721             {
00722               error ("etree: S is marked as symmetric, but is not square");
00723               return retval;
00724             }
00725 
00726           symetree (ridx, cidx, etree, 0, n_col);
00727         }
00728       else
00729         {
00730           OCTAVE_LOCAL_BUFFER (octave_idx_type, colbeg, n_col);
00731           OCTAVE_LOCAL_BUFFER (octave_idx_type, colend, n_col);
00732 
00733           for (octave_idx_type i = 0; i < n_col; i++)
00734             {
00735               colbeg[i] = cidx[i];
00736               colend[i] = cidx[i+1];
00737             }
00738 
00739           coletree (ridx, colbeg, colend, etree, n_row, n_col);
00740         }
00741 
00742       NDArray tree (dim_vector (1, n_col));
00743       for (octave_idx_type i = 0; i < n_col; i++)
00744         // We flag a root with n_col while Matlab does it with zero
00745         // Convert for matlab compatiable output
00746         if (etree[i] == n_col)
00747           tree(i) = 0;
00748         else
00749           tree(i) = etree[i] + 1;
00750 
00751       retval(0) = tree;
00752 
00753       if (nargout == 2)
00754         {
00755           // Calculate the tree post-ordering
00756           OCTAVE_LOCAL_BUFFER (octave_idx_type, post, n_col + 1);
00757           tree_postorder (n_col, etree, post);
00758 
00759           NDArray postorder (dim_vector (1, n_col));
00760           for (octave_idx_type i = 0; i < n_col; i++)
00761             postorder(i) = post[i] + 1;
00762 
00763           retval(1) = postorder;
00764         }
00765     }
00766 
00767   return retval;
00768 }
 All Classes Files Functions Variables Typedefs Enumerations Enumerator Friends Defines