GNU Octave  3.8.0
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
symbfact.cc
Go to the documentation of this file.
1 /*
2 
3 Copyright (C) 2005-2013 David Bateman
4 Copyright (C) 1998-2005 Andy Adler
5 
6 This file is part of Octave.
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
10 Free Software Foundation; either version 3 of the License, or (at your
11 option) any later version.
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
19 along with Octave; see the file COPYING. If not, see
20 <http://www.gnu.org/licenses/>.
21 
22 */
23 
24 #ifdef HAVE_CONFIG_H
25 #include <config.h>
26 #endif
27 
28 #include "SparseCmplxCHOL.h"
29 #include "SparsedbleCHOL.h"
30 #include "oct-spparms.h"
31 #include "sparse-util.h"
32 #include "oct-locbuf.h"
33 
34 #include "ov-re-sparse.h"
35 #include "ov-cx-sparse.h"
36 #include "defun-dld.h"
37 #include "error.h"
38 #include "gripes.h"
39 #include "oct-obj.h"
40 #include "utils.h"
41 
42 DEFUN_DLD (symbfact, args, nargout,
43  "-*- texinfo -*-\n\
44 @deftypefn {Loadable Function} {[@var{count}, @var{h}, @var{parent}, @var{post}, @var{r}] =} symbfact (@var{S})\n\
45 @deftypefnx {Loadable Function} {[@dots{}] =} symbfact (@var{S}, @var{typ})\n\
46 @deftypefnx {Loadable Function} {[@dots{}] =} symbfact (@var{S}, @var{typ}, @var{mode})\n\
47 \n\
48 Perform a symbolic factorization analysis on the sparse matrix @var{S}.\n\
49 Where\n\
50 \n\
51 @table @var\n\
52 @item S\n\
53 @var{S} is a complex or real sparse matrix.\n\
54 \n\
55 @item typ\n\
56 Is the type of the factorization and can be one of\n\
57 \n\
58 @table @samp\n\
59 @item sym\n\
60 Factorize @var{S}. This is the default.\n\
61 \n\
62 @item col\n\
63 Factorize @code{@var{S}' * @var{S}}.\n\
64 \n\
65 @item row\n\
66 Factorize @tcode{@var{S} * @var{S}'}.\n\
67 \n\
68 @item lo\n\
69 Factorize @tcode{@var{S}'}\n\
70 @end table\n\
71 \n\
72 @item mode\n\
73 The default is to return the Cholesky@tie{}factorization for @var{r}, and if\n\
74 @var{mode} is @qcode{'L'}, the conjugate transpose of the\n\
75 Cholesky@tie{}factorization is returned. The conjugate transpose version is\n\
76 faster and uses less memory, but returns the same values for @var{count},\n\
77 @var{h}, @var{parent} and @var{post} outputs.\n\
78 @end table\n\
79 \n\
80 The output variables are\n\
81 \n\
82 @table @var\n\
83 @item count\n\
84 The row counts of the Cholesky@tie{}factorization as determined by @var{typ}.\n\
85 \n\
86 @item h\n\
87 The height of the elimination tree.\n\
88 \n\
89 @item parent\n\
90 The elimination tree itself.\n\
91 \n\
92 @item post\n\
93 A sparse boolean matrix whose structure is that of the Cholesky\n\
94 factorization as determined by @var{typ}.\n\
95 @end table\n\
96 @end deftypefn")
97 {
98  octave_value_list retval;
99  int nargin = args.length ();
100 
101  if (nargin < 1 || nargin > 3 || nargout > 5)
102  {
103  print_usage ();
104  return retval;
105  }
106 
107 #ifdef HAVE_CHOLMOD
108 
109  cholmod_common Common;
110  cholmod_common *cm = &Common;
111  CHOLMOD_NAME(start) (cm);
112 
113  double spu = octave_sparse_params::get_key ("spumoni");
114  if (spu == 0.)
115  {
116  cm->print = -1;
117  cm->print_function = 0;
118  }
119  else
120  {
121  cm->print = static_cast<int> (spu) + 2;
122  cm->print_function =&SparseCholPrint;
123  }
124 
125  cm->error_handler = &SparseCholError;
126  cm->complex_divide = CHOLMOD_NAME(divcomplex);
127  cm->hypotenuse = CHOLMOD_NAME(hypot);
128 
129  double dummy;
130  cholmod_sparse Astore;
131  cholmod_sparse *A = &Astore;
132  A->packed = true;
133  A->sorted = true;
134  A->nz = 0;
135 #ifdef USE_64_BIT_IDX_T
136  A->itype = CHOLMOD_LONG;
137 #else
138  A->itype = CHOLMOD_INT;
139 #endif
140  A->dtype = CHOLMOD_DOUBLE;
141  A->stype = 1;
142  A->x = &dummy;
143 
144  if (args(0).is_real_type ())
145  {
146  const SparseMatrix a = args(0).sparse_matrix_value ();
147  A->nrow = a.rows ();
148  A->ncol = a.cols ();
149  A->p = a.cidx ();
150  A->i = a.ridx ();
151  A->nzmax = a.nnz ();
152  A->xtype = CHOLMOD_REAL;
153 
154  if (a.rows () > 0 && a.cols () > 0)
155  A->x = a.data ();
156  }
157  else if (args(0).is_complex_type ())
158  {
159  const SparseComplexMatrix a = args(0).sparse_complex_matrix_value ();
160  A->nrow = a.rows ();
161  A->ncol = a.cols ();
162  A->p = a.cidx ();
163  A->i = a.ridx ();
164  A->nzmax = a.nnz ();
165  A->xtype = CHOLMOD_COMPLEX;
166 
167  if (a.rows () > 0 && a.cols () > 0)
168  A->x = a.data ();
169  }
170  else
171  gripe_wrong_type_arg ("symbfact", args(0));
172 
173  octave_idx_type coletree = false;
174  octave_idx_type n = A->nrow;
175 
176  if (nargin > 1)
177  {
178  char ch;
179  std::string str = args(1).string_value ();
180  ch = tolower (str.c_str ()[0]);
181  if (ch == 'r')
182  A->stype = 0;
183  else if (ch == 'c')
184  {
185  n = A->ncol;
186  coletree = true;
187  A->stype = 0;
188  }
189  else if (ch == 's')
190  A->stype = 1;
191  else if (ch == 's')
192  A->stype = -1;
193  else
194  error ("symbfact: unrecognized TYP in symbolic factorization");
195  }
196 
197  if (A->stype && A->nrow != A->ncol)
198  error ("symbfact: S must be a square matrix");
199 
200  if (!error_state)
201  {
204  OCTAVE_LOCAL_BUFFER (octave_idx_type, ColCount, n);
207 
208  cholmod_sparse *F = CHOLMOD_NAME(transpose) (A, 0, cm);
209  cholmod_sparse *Aup, *Alo;
210 
211  if (A->stype == 1 || coletree)
212  {
213  Aup = A ;
214  Alo = F ;
215  }
216  else
217  {
218  Aup = F ;
219  Alo = A ;
220  }
221 
222  CHOLMOD_NAME(etree) (Aup, Parent, cm);
223 
224  if (cm->status < CHOLMOD_OK)
225  {
226  error ("matrix corrupted");
227  goto symbfact_error;
228  }
229 
230  if (CHOLMOD_NAME(postorder) (Parent, n, 0, Post, cm) != n)
231  {
232  error ("postorder failed");
233  goto symbfact_error;
234  }
235 
236  CHOLMOD_NAME(rowcolcounts) (Alo, 0, 0, Parent, Post, 0,
237  ColCount, First, Level, cm);
238 
239  if (cm->status < CHOLMOD_OK)
240  {
241  error ("matrix corrupted");
242  goto symbfact_error;
243  }
244 
245  if (nargout > 4)
246  {
247  cholmod_sparse *A1, *A2;
248 
249  if (A->stype == 1)
250  {
251  A1 = A;
252  A2 = 0;
253  }
254  else if (A->stype == -1)
255  {
256  A1 = F;
257  A2 = 0;
258  }
259  else if (coletree)
260  {
261  A1 = F;
262  A2 = A;
263  }
264  else
265  {
266  A1 = A;
267  A2 = F;
268  }
269 
270  // count the total number of entries in L
271  octave_idx_type lnz = 0 ;
272  for (octave_idx_type j = 0 ; j < n ; j++)
273  lnz += ColCount[j];
274 
275 
276  // allocate the output matrix L (pattern-only)
277  SparseBoolMatrix L (n, n, lnz);
278 
279  // initialize column pointers
280  lnz = 0;
281  for (octave_idx_type j = 0 ; j < n ; j++)
282  {
283  L.xcidx(j) = lnz;
284  lnz += ColCount[j];
285  }
286  L.xcidx(n) = lnz;
287 
288 
289  /* create a copy of the column pointers */
290  octave_idx_type *W = First;
291  for (octave_idx_type j = 0 ; j < n ; j++)
292  W[j] = L.xcidx (j);
293 
294  // get workspace for computing one row of L
295  cholmod_sparse *R
296  = CHOLMOD_NAME (allocate_sparse) (n, 1, n, false, true,
297  0, CHOLMOD_PATTERN, cm);
298  octave_idx_type *Rp = static_cast<octave_idx_type *>(R->p);
299  octave_idx_type *Ri = static_cast<octave_idx_type *>(R->i);
300 
301  // compute L one row at a time
302  for (octave_idx_type k = 0 ; k < n ; k++)
303  {
304  // get the kth row of L and store in the columns of L
305  CHOLMOD_NAME (row_subtree) (A1, A2, k, Parent, R, cm) ;
306  for (octave_idx_type p = 0 ; p < Rp[1] ; p++)
307  L.xridx (W[Ri[p]]++) = k ;
308 
309  // add the diagonal entry
310  L.xridx (W[k]++) = k ;
311  }
312 
313  // free workspace
314  CHOLMOD_NAME (free_sparse) (&R, cm) ;
315 
316 
317  // transpose L to get R, or leave as is
318  if (nargin < 3)
319  L = L.transpose ();
320 
321  // fill numerical values of L with one's
322  for (octave_idx_type p = 0 ; p < lnz ; p++)
323  L.xdata(p) = true;
324 
325  retval(4) = L;
326  }
327 
328  ColumnVector tmp (n);
329  if (nargout > 3)
330  {
331  for (octave_idx_type i = 0; i < n; i++)
332  tmp(i) = Post[i] + 1;
333  retval(3) = tmp;
334  }
335 
336  if (nargout > 2)
337  {
338  for (octave_idx_type i = 0; i < n; i++)
339  tmp(i) = Parent[i] + 1;
340  retval(2) = tmp;
341  }
342 
343  if (nargout > 1)
344  {
345  /* compute the elimination tree height */
346  octave_idx_type height = 0 ;
347  for (int i = 0 ; i < n ; i++)
348  height = (height > Level[i] ? height : Level[i]);
349  height++ ;
350  retval(1) = static_cast<double> (height);
351  }
352 
353  for (octave_idx_type i = 0; i < n; i++)
354  tmp(i) = ColCount[i];
355  retval(0) = tmp;
356  }
357 
358 symbfact_error:
359 #else
360  error ("symbfact: not available in this version of Octave");
361 #endif
362 
363  return retval;
364 }