GNU Octave  6.2.0
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
symbfact.cc
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////
2 //
3 // Copyright (C) 1998-2021 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 #if defined (HAVE_CONFIG_H)
27 # include "config.h"
28 #endif
29 
30 #include <cmath>
31 
32 #include <algorithm>
33 #include <string>
34 
35 #include "CSparse.h"
36 #include "boolSparse.h"
37 #include "dColVector.h"
38 #include "dSparse.h"
39 #include "oct-locbuf.h"
40 #include "oct-sparse.h"
41 #include "oct-spparms.h"
42 #include "sparse-util.h"
43 
44 #include "defun.h"
45 #include "error.h"
46 #include "errwarn.h"
47 #include "ovl.h"
48 #include "parse.h"
49 #include "utils.h"
50 
51 DEFUN (symbfact, args, nargout,
52  doc: /* -*- texinfo -*-
53 @deftypefn {} {[@var{count}, @var{h}, @var{parent}, @var{post}, @var{R}] =} symbfact (@var{S})
54 @deftypefnx {} {[@dots{}] =} symbfact (@var{S}, @var{typ})
55 @deftypefnx {} {[@dots{}] =} symbfact (@var{S}, @var{typ}, @var{mode})
56 
57 Perform a symbolic factorization analysis of the sparse matrix @var{S}.
58 
59 The input variables are
60 
61 @table @var
62 @item S
63 @var{S} is a real or complex sparse matrix.
64 
65 @item typ
66 Is the type of the factorization and can be one of
67 
68 @table @asis
69 @item @qcode{"sym"} (default)
70 Factorize @var{S}. Assumes @var{S} is symmetric and uses the upper
71 triangular portion of the matrix.
72 
73 @item @qcode{"col"}
74 Factorize @tcode{@var{S}' * @var{S}}.
75 
76 @item @qcode{"row"}
77 Factorize @tcode{@var{S} * @var{S}'}.
78 
79 @item @qcode{"lo"}
80 Factorize @tcode{@var{S}'}. Assumes @var{S} is symmetric and uses the lower
81 triangular portion of the matrix.
82 @end table
83 
84 @item mode
85 When @var{mode} is unspecified return the Cholesky@tie{}factorization for
86 @var{R}. If @var{mode} is @qcode{"lower"} or @qcode{"L"} then return
87 the conjugate transpose @tcode{@var{R}'} which is a lower triangular factor.
88 The conjugate transpose version is faster and uses less memory, but still
89 returns the same values for all other outputs: @var{count}, @var{h},
90 @var{parent}, and @var{post}.
91 @end table
92 
93 The output variables are:
94 
95 @table @var
96 @item count
97 The row counts of the Cholesky@tie{}factorization as determined by
98 @var{typ}. The computational difficulty of performing the true
99 factorization using @code{chol} is @code{sum (@var{count} .^ 2)}.
100 
101 @item h
102 The height of the elimination tree.
103 
104 @item parent
105 The elimination tree itself.
106 
107 @item post
108 A sparse boolean matrix whose structure is that of the
109 Cholesky@tie{}factorization as determined by @var{typ}.
110 @end table
111 @seealso{chol, etree, treelayout}
112 @end deftypefn */)
113 {
114 #if defined (HAVE_CHOLMOD)
115 
116  int nargin = args.length ();
117 
118  if (nargin < 1 || nargin > 3)
119  print_usage ();
120 
122 
123  double dummy;
124  cholmod_sparse Astore;
125  cholmod_sparse *A = &Astore;
126  A->packed = true;
127  A->sorted = true;
128  A->nz = nullptr;
129 #if defined (OCTAVE_ENABLE_64)
130  A->itype = CHOLMOD_LONG;
131 #else
132  A->itype = CHOLMOD_INT;
133 #endif
134  A->dtype = CHOLMOD_DOUBLE;
135  A->stype = 1;
136  A->x = &dummy;
137 
138  if (args(0).isreal ())
139  {
140  const SparseMatrix a = args(0).sparse_matrix_value ();
141  A->nrow = a.rows ();
142  A->ncol = a.cols ();
143  A->p = a.cidx ();
144  A->i = a.ridx ();
145  A->nzmax = a.nnz ();
146  A->xtype = CHOLMOD_REAL;
147 
148  if (a.rows () > 0 && a.cols () > 0)
149  A->x = a.data ();
150  }
151  else if (args(0).iscomplex ())
152  {
153  const SparseComplexMatrix a = args(0).sparse_complex_matrix_value ();
154  A->nrow = a.rows ();
155  A->ncol = a.cols ();
156  A->p = a.cidx ();
157  A->i = a.ridx ();
158  A->nzmax = a.nnz ();
159  A->xtype = CHOLMOD_COMPLEX;
160 
161  if (a.rows () > 0 && a.cols () > 0)
162  A->x = a.data ();
163  }
164  else
165  err_wrong_type_arg ("symbfact", args(0));
166 
167  bool coletree = false;
168  octave_idx_type n = A->nrow;
169 
170  if (nargin > 1)
171  {
172  std::string str = args(1).xstring_value ("TYP must be a string");
173  // FIXME: The input validation could be improved to use strncmp
174  char ch;
175  ch = tolower (str[0]);
176  if (ch == 'r') // 'row'
177  A->stype = 0;
178  else if (ch == 'c') // 'col'
179  {
180  n = A->ncol;
181  coletree = true;
182  A->stype = 0;
183  }
184  else if (ch == 's') // 'sym' (default)
185  A->stype = 1;
186  else if (ch == 'l') // 'lo'
187  A->stype = -1;
188  else
189  error (R"(symbfact: unrecognized TYP "%s")", str.c_str ());
190  }
191 
192  if (nargin == 3)
193  {
194  std::string str = args(2).xstring_value ("MODE must be a string");
195  // FIXME: The input validation could be improved to use strncmp
196  char ch;
197  ch = toupper (str[0]);
198  if (ch != 'L')
199  error (R"(symbfact: unrecognized MODE "%s")", str.c_str ());
200  }
201 
202  if (A->stype && A->nrow != A->ncol)
203  err_square_matrix_required ("symbfact", "S");
204 
210 
211  cholmod_common Common;
212  cholmod_common *cm = &Common;
213  CHOLMOD_NAME(start) (cm);
214 
215  double spu = octave_sparse_params::get_key ("spumoni");
216  if (spu == 0.0)
217  {
218  cm->print = -1;
219  SUITESPARSE_ASSIGN_FPTR (printf_func, cm->print_function, nullptr);
220  }
221  else
222  {
223  cm->print = static_cast<int> (spu) + 2;
224  SUITESPARSE_ASSIGN_FPTR (printf_func, cm->print_function, &SparseCholPrint);
225  }
226 
227  cm->error_handler = &SparseCholError;
228  SUITESPARSE_ASSIGN_FPTR2 (divcomplex_func, cm->complex_divide, divcomplex);
229  SUITESPARSE_ASSIGN_FPTR2 (hypot_func, cm->hypotenuse, hypot);
230 
231  cholmod_sparse *F = CHOLMOD_NAME(transpose) (A, 0, cm);
232  cholmod_sparse *Aup, *Alo;
233 
234  if (A->stype == 1 || coletree)
235  {
236  Aup = A;
237  Alo = F;
238  }
239  else
240  {
241  Aup = F;
242  Alo = A;
243  }
244 
245  CHOLMOD_NAME(etree) (Aup, Parent, cm);
246 
247  ColumnVector tmp (n); // Declaration must precede any goto cleanup.
248  std::string err_msg;
249 
250  if (cm->status < CHOLMOD_OK)
251  {
252  err_msg = "symbfact: matrix corrupted";
253  goto cleanup;
254  }
255 
256  if (CHOLMOD_NAME(postorder) (Parent, n, nullptr, Post, cm) != n)
257  {
258  err_msg = "symbfact: postorder failed";
259  goto cleanup;
260  }
261 
262  CHOLMOD_NAME(rowcolcounts) (Alo, nullptr, 0, Parent, Post, nullptr, ColCount,
263  First, octave::to_suitesparse_intptr (Level), cm);
264 
265  if (cm->status < CHOLMOD_OK)
266  {
267  err_msg = "symbfact: matrix corrupted";
268  goto cleanup;
269  }
270 
271  if (nargout > 4)
272  {
273  cholmod_sparse *A1, *A2;
274 
275  if (A->stype == 1)
276  {
277  A1 = A;
278  A2 = nullptr;
279  }
280  else if (A->stype == -1)
281  {
282  A1 = F;
283  A2 = nullptr;
284  }
285  else if (coletree)
286  {
287  A1 = F;
288  A2 = A;
289  }
290  else
291  {
292  A1 = A;
293  A2 = F;
294  }
295 
296  // count the total number of entries in L
297  octave_idx_type lnz = 0;
298  for (octave_idx_type j = 0 ; j < n ; j++)
299  lnz += ColCount[j];
300 
301  // allocate the output matrix L (pattern-only)
302  SparseBoolMatrix L (dim_vector (n, n), lnz);
303 
304  // initialize column pointers
305  lnz = 0;
306  for (octave_idx_type j = 0 ; j < n ; j++)
307  {
308  L.xcidx(j) = lnz;
309  lnz += ColCount[j];
310  }
311  L.xcidx(n) = lnz;
312 
313  // create a copy of the column pointers
314  octave::suitesparse_integer *W = First;
315  for (octave_idx_type j = 0 ; j < n ; j++)
316  W[j] = L.xcidx (j);
317 
318  // get workspace for computing one row of L
319  cholmod_sparse *R
320  = CHOLMOD_NAME(allocate_sparse) (n, 1, n, false, true,
321  0, CHOLMOD_PATTERN, cm);
322  octave_idx_type *Rp = static_cast<octave_idx_type *> (R->p);
323  octave_idx_type *Ri = static_cast<octave_idx_type *> (R->i);
324 
325  // compute L one row at a time
326  for (octave_idx_type k = 0 ; k < n ; k++)
327  {
328  // get the kth row of L and store in the columns of L
329  CHOLMOD_NAME(row_subtree) (A1, A2, k, Parent, R, cm);
330  for (octave_idx_type p = 0 ; p < Rp[1] ; p++)
331  L.xridx (W[Ri[p]]++) = k;
332 
333  // add the diagonal entry
334  L.xridx (W[k]++) = k;
335  }
336 
337  // free workspace
338  CHOLMOD_NAME(free_sparse) (&R, cm);
339 
340  // fill L with one's
341  std::fill_n (L.xdata (), lnz, true);
342 
343  // transpose L to get R, or leave as is
344  if (nargin < 3)
345  L = L.transpose ();
346 
347  retval(4) = L;
348  }
349 
350  if (nargout > 3)
351  {
352  for (octave_idx_type i = 0; i < n; i++)
353  tmp(i) = Post[i] + 1;
354  retval(3) = tmp;
355  }
356 
357  if (nargout > 2)
358  {
359  for (octave_idx_type i = 0; i < n; i++)
360  tmp(i) = Parent[i] + 1;
361  retval(2) = tmp;
362  }
363 
364  if (nargout > 1)
365  {
366  // compute the elimination tree height
367  octave_idx_type height = 0;
368  for (int i = 0 ; i < n ; i++)
369  height = std::max (height, Level[i]);
370  height++;
371  retval(1) = static_cast<double> (height);
372  }
373 
374  for (octave_idx_type i = 0; i < n; i++)
375  tmp(i) = ColCount[i];
376  retval(0) = tmp;
377 
378 cleanup:
379  CHOLMOD_NAME(free_sparse) (&F, cm);
380  CHOLMOD_NAME(finish) (cm);
381 
382  if (! err_msg.empty ())
383  error ("%s", err_msg.c_str ());
384 
385  return retval;
386 
387 #else
388 
389  octave_unused_parameter (args);
390  octave_unused_parameter (nargout);
391 
392  err_disabled_feature ("symbfact", "CHOLMOD");
393 
394 #endif
395 }
396 
397 /*
398 %!testif HAVE_CHOLMOD
399 %! A = sparse (magic (3));
400 %! [count, h, parent, post, r] = symbfact (A);
401 %! assert (count, [3; 2; 1]);
402 %! assert (h, 3);
403 %! assert (parent, [2; 3; 0]);
404 %! assert (r, sparse (triu (true (3))));
405 
406 %!testif HAVE_CHOLMOD
407 %! ## Test MODE "lower"
408 %! A = sparse (magic (3));
409 %! [~, ~, ~, ~, l] = symbfact (A, "sym", "lower");
410 %! assert (l, sparse (tril (true (3))));
411 
412 %!testif HAVE_CHOLMOD <*42587>
413 %! ## singular matrix
414 %! A = sparse ([1 0 8;0 1 8;8 8 1]);
415 %! [count, h, parent, post, r] = symbfact (A);
416 
417 ## Test input validation
418 %!testif HAVE_CHOLMOD
419 %! fail ("symbfact ()");
420 %! fail ("symbfact (1,2,3,4)");
421 %! fail ("symbfact ({1})", "wrong type argument 'cell'");
422 %! fail ("symbfact (sparse (1), {1})", "TYP must be a string");
423 %! fail ("symbfact (sparse (1), 'foobar')", 'unrecognized TYP "foobar"');
424 %! fail ("symbfact (sparse (1), 'sym', {'L'})", "MODE must be a string");
425 %! fail ('symbfact (sparse (1), "sym", "foobar")', 'unrecognized MODE "foobar"');
426 %! fail ("symbfact (sparse ([1, 2; 3, 4; 5, 6]))", "S must be a square matrix");
427 
428 */
charNDArray max(char d, const charNDArray &m)
Definition: chNDArray.cc:230
SparseBoolMatrix transpose(void) const
Definition: boolSparse.h:91
octave_idx_type cols(void) const
Definition: Sparse.h:251
octave_idx_type * xridx(void)
Definition: Sparse.h:485
T * data(void)
Definition: Sparse.h:470
octave_idx_type * cidx(void)
Definition: Sparse.h:492
octave_idx_type nnz(void) const
Actual number of nonzero terms.
Definition: Sparse.h:238
octave_idx_type rows(void) const
Definition: Sparse.h:250
octave_idx_type * xcidx(void)
Definition: Sparse.h:498
T * xdata(void)
Definition: Sparse.h:472
octave_idx_type * ridx(void)
Definition: Sparse.h:479
Vector representing the dimensions (size) of an Array.
Definition: dim-vector.h:95
static double get_key(const std::string &key)
Definition: oct-spparms.cc:93
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:157
OCTINTERP_API void print_usage(void)
Definition: defun.cc:53
#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:968
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
void err_wrong_type_arg(const char *name, const char *s)
Definition: errwarn.cc:166
F77_RET_T const F77_INT F77_CMPLX * A
octave_idx_type n
Definition: mx-inlines.cc:753
int suitesparse_integer
Definition: oct-sparse.h:171
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
#define CHOLMOD_NAME(name)
Definition: oct-sparse.h:128
octave_value::octave_value(const Array< char > &chm, char type) return retval
Definition: ov.cc:811
void SparseCholError(int status, char *file, int line, char *message)
Definition: sparse-util.cc:64
int SparseCholPrint(const char *fmt,...)
Definition: sparse-util.cc:76
static void transpose(octave_idx_type N, const octave_idx_type *ridx, const octave_idx_type *cidx, octave_idx_type *ridx2, octave_idx_type *cidx2)
Definition: symrcm.cc:389