view libinterp/dldfcn/symbfact.cc @ 21966:112b20240c87

move docstrings in C++ files out of C strings and into comments * __contourc__.cc, __dispatch__.cc, __dsearchn__.cc, __ichol__.cc, __ilu__.cc, __lin_interpn__.cc, __luinc__.cc, __magick_read__.cc, __pchip_deriv__.cc, __qp__.cc, balance.cc, besselj.cc, betainc.cc, bitfcns.cc, bsxfun.cc, cellfun.cc, colloc.cc, conv2.cc, daspk.cc, dasrt.cc, dassl.cc, data.cc, debug.cc, defaults.cc, det.cc, dirfns.cc, dlmread.cc, dot.cc, eig.cc, ellipj.cc, error.cc, fft.cc, fft2.cc, fftn.cc, file-io.cc, filter.cc, find.cc, gammainc.cc, gcd.cc, getgrent.cc, getpwent.cc, getrusage.cc, givens.cc, graphics.cc, hash.cc, help.cc, hess.cc, hex2num.cc, input.cc, inv.cc, kron.cc, load-path.cc, load-save.cc, lookup.cc, ls-oct-text.cc, lsode.cc, lu.cc, mappers.cc, matrix_type.cc, max.cc, mgorth.cc, nproc.cc, oct-hist.cc, octave-link.cc, ordschur.cc, pager.cc, pinv.cc, pr-output.cc, profiler.cc, psi.cc, pt-jit.cc, quad.cc, quadcc.cc, qz.cc, rand.cc, rcond.cc, regexp.cc, schur.cc, sighandlers.cc, sparse.cc, spparms.cc, sqrtm.cc, str2double.cc, strfind.cc, strfns.cc, sub2ind.cc, svd.cc, sylvester.cc, symtab.cc, syscalls.cc, sysdep.cc, time.cc, toplev.cc, tril.cc, tsearch.cc, typecast.cc, urlwrite.cc, utils.cc, variables.cc, __delaunayn__.cc, __eigs__.cc, __fltk_uigetfile__.cc, __glpk__.cc, __init_fltk__.cc, __init_gnuplot__.cc, __osmesa_print__.cc, __voronoi__.cc, amd.cc, audiodevinfo.cc, audioread.cc, ccolamd.cc, chol.cc, colamd.cc, convhulln.cc, dmperm.cc, fftw.cc, qr.cc, symbfact.cc, symrcm.cc, ov-base.cc, ov-bool-mat.cc, ov-cell.cc, ov-class.cc, ov-classdef.cc, ov-fcn-handle.cc, ov-fcn-inline.cc, ov-flt-re-mat.cc, ov-int16.cc, ov-int32.cc, ov-int64.cc, ov-int8.cc, ov-java.cc, ov-null-mat.cc, ov-oncleanup.cc, ov-range.cc, ov-re-mat.cc, ov-struct.cc, ov-typeinfo.cc, ov-uint16.cc, ov-uint32.cc, ov-uint64.cc, ov-uint8.cc, ov-usr-fcn.cc, ov.cc, octave.cc, pt-arg-list.cc, pt-binop.cc, pt-eval.cc, pt-mat.cc, lex.ll, oct-parse.in.yy: Docstrings are now comments instead of C strings. * build-aux/mk-opts.pl: Emit docstrings as comments instead of C strings. * DASPK-opts.in, LSODE-opts.in: Don't quote " in docstring fragments. * builtins.h: Include builtin-defun-decls.h unconditionally. * defun.h (DEFUN, DEFUNX, DEFCONSTFUN): Simply emit declaration. (DEFALIAS): Always expand to nothing. * defun-dld.h: No special macro expansions for MAKE_BUILTINS. (DEFUN_DLD): Use FORWARD_DECLARE_FUN. (DEFUNX_DLD): Use FORWARD_DECLARE_FUNX. * defun-int.h: No special macro expansions for MAKE_BUILTINS. (FORWARD_DECLARE_FUN, FORWARD_DECLARE_FUNX): New macros. (DEFINE_FUN_INSTALLER_FUN): If compiling an Octave source file, pass "external-doc" to DEFINE_FUNX_INSTALLER_FUN. (DEFUN_INTERNAL, DEFCONSTFUN_INTERNAL, DEFUNX_INTERNAL, DEFALIAS_INTERNAL): Delete. * common.mk (move_if_change_rule): New macro. (simple_move_if_change_rule): Define using move_if_change_rule. * find-defun-files.sh (DEFUN_PATTERN): Update. Don't transform file name extension to ".df". * libinterp/mk-pkg-add, gendoc.pl: Operate directly on source files. * mkbuiltins: New argument, SRCDIR. Operate directly on source files. * mkdefs: Delete. * libinterp/module.mk (BUILT_SOURCES): Update list to contain only files included in other source files. (GENERATED_MAKE_BUILTINS_INCS, DEF_FILES): Delete. (LIBINTERP_BUILT_DISTFILES): Include $(OPT_HANDLERS) here. (LIBINTERP_BUILT_NODISTFILES): Not here. Remove $(ALL_DEF_FILES from the list. (libinterp_EXTRA_DIST): Remove mkdefs from the list. (FOUND_DEFUN_FILES): Rename from SRC_DEF_FILES. (DLDFCN_DEFUN_FILES): Rename from DLDFCN_DEF_FILES. (SRC_DEFUN_FILES): Rename from SRC_DEF_FILES. (ALL_DEFUN_FILES): Rename from ALL_DEF_FILES. (%.df: %.cc): Delete pattern rule. (libinterp/build-env-features.cc, libinterp/builtins.cc, libinterp/dldfcn/PKG_ADD): Use mv instead of move-if-change. (libinterp/builtins.cc, libinterp/builtin-defun-decls.h): Update mkbuiltins command. ($(srcdir)/libinterp/DOCSTRINGS): Update gendoc.pl command. * liboctave/module.mk (BUILT_SOURCES): Don't include liboctave-build-info.cc in the list.
author John W. Eaton <jwe@octave.org>
date Tue, 21 Jun 2016 16:07:51 -0400
parents 55f7de37b618
children bac0d6f07a3e
line wrap: on
line source

/*

Copyright (C) 2005-2015 David Bateman
Copyright (C) 1998-2005 Andy Adler

This file is part of Octave.

Octave is free software; you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by the
Free Software Foundation; either version 3 of the License, or (at your
option) any later version.

Octave is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
for more details.

You should have received a copy of the GNU General Public License
along with Octave; see the file COPYING.  If not, see
<http://www.gnu.org/licenses/>.

*/

#if defined (HAVE_CONFIG_H)
#  include "config.h"
#endif

#include "oct-locbuf.h"
#include "oct-sparse.h"
#include "oct-spparms.h"
#include "sparse-chol.h"
#include "sparse-util.h"

#include "ov-re-sparse.h"
#include "ov-cx-sparse.h"
#include "defun-dld.h"
#include "error.h"
#include "errwarn.h"
#include "ovl.h"
#include "utils.h"

DEFUN_DLD (symbfact, args, nargout,
           doc: /* -*- texinfo -*-
@deftypefn  {} {[@var{count}, @var{h}, @var{parent}, @var{post}, @var{R}] =} symbfact (@var{S})
@deftypefnx {} {[@dots{}] =} symbfact (@var{S}, @var{typ})
@deftypefnx {} {[@dots{}] =} symbfact (@var{S}, @var{typ}, @var{mode})

Perform a symbolic factorization analysis of the sparse matrix @var{S}.

The input variables are

@table @var
@item S
@var{S} is a real or complex sparse matrix.

@item typ
Is the type of the factorization and can be one of

@table @asis
@item @qcode{"sym"} (default)
Factorize @var{S}.  Assumes @var{S} is symmetric and uses the upper
triangular portion of the matrix.

@item @qcode{"col"}
Factorize @tcode{@var{S}' * @var{S}}.

@item @qcode{"row"}
Factorize @tcode{@var{S} * @var{S}'}.

@item @qcode{"lo"}
Factorize @tcode{@var{S}'}.  Assumes @var{S} is symmetric and uses the lower
triangular portion of the matrix.
@end table

@item mode
When @var{mode} is unspecified return the Cholesky@tie{}factorization for
@var{R}.  If @var{mode} is @qcode{"lower"} or @qcode{"L"} then return
the conjugate transpose @tcode{@var{R}'} which is a lower triangular factor.
The conjugate transpose version is faster and uses less memory, but still
returns the same values for all other outputs: @var{count}, @var{h},
@var{parent}, and @var{post}.
@end table

The output variables are:

@table @var
@item count
The row counts of the Cholesky@tie{}factorization as determined by
@var{typ}.  The computational difficulty of performing the true
factorization using @code{chol} is @code{sum (@var{count} .^ 2)}.

@item h
The height of the elimination tree.

@item parent
The elimination tree itself.

@item post
A sparse boolean matrix whose structure is that of the
Cholesky@tie{}factorization as determined by @var{typ}.
@end table
@seealso{chol, etree, treelayout}
@end deftypefn */)
{
#if defined (HAVE_CHOLMOD)

  int nargin = args.length ();

  if (nargin < 1 || nargin > 3)
    print_usage ();

  octave_value_list retval;

  double dummy;
  cholmod_sparse Astore;
  cholmod_sparse *A = &Astore;
  A->packed = true;
  A->sorted = true;
  A->nz = 0;
#if defined (OCTAVE_ENABLE_64)
  A->itype = CHOLMOD_LONG;
#else
  A->itype = CHOLMOD_INT;
#endif
  A->dtype = CHOLMOD_DOUBLE;
  A->stype = 1;
  A->x = &dummy;

  if (args(0).is_real_type ())
    {
      const SparseMatrix a = args(0).sparse_matrix_value ();
      A->nrow = a.rows ();
      A->ncol = a.cols ();
      A->p = a.cidx ();
      A->i = a.ridx ();
      A->nzmax = a.nnz ();
      A->xtype = CHOLMOD_REAL;

      if (a.rows () > 0 && a.cols () > 0)
        A->x = a.data ();
    }
  else if (args(0).is_complex_type ())
    {
      const SparseComplexMatrix a = args(0).sparse_complex_matrix_value ();
      A->nrow = a.rows ();
      A->ncol = a.cols ();
      A->p = a.cidx ();
      A->i = a.ridx ();
      A->nzmax = a.nnz ();
      A->xtype = CHOLMOD_COMPLEX;

      if (a.rows () > 0 && a.cols () > 0)
        A->x = a.data ();
    }
  else
    err_wrong_type_arg ("symbfact", args(0));

  octave_idx_type coletree = false;
  octave_idx_type n = A->nrow;

  if (nargin > 1)
    {
      std::string str = args(1).xstring_value ("TYP must be a string");
      // FIXME: The input validation could be improved to use strncmp
      char ch;
      ch = tolower (str[0]);
      if (ch == 'r')          // 'row'
        A->stype = 0;
      else if (ch == 'c')     // 'col'
        {
          n = A->ncol;
          coletree = true;
          A->stype = 0;
        }
      else if (ch == 's')     // 'sym' (default)
        A->stype = 1;
      else if (ch == 'l')     // 'lo'
        A->stype = -1;
      else
        error ("symbfact: unrecognized TYP \"%s\"", str.c_str ());
    }

  if (nargin == 3)
    {
      std::string str = args(2).xstring_value ("MODE must be a string");
      // FIXME: The input validation could be improved to use strncmp
      char ch;
      ch = toupper (str[0]);
      if (ch != 'L')
        error ("symbfact: unrecognized MODE \"%s\"", str.c_str ());
    }

  if (A->stype && A->nrow != A->ncol)
    err_square_matrix_required ("symbfact", "S");

  OCTAVE_LOCAL_BUFFER (octave_idx_type, Parent, n);
  OCTAVE_LOCAL_BUFFER (octave_idx_type, Post, n);
  OCTAVE_LOCAL_BUFFER (octave_idx_type, ColCount, n);
  OCTAVE_LOCAL_BUFFER (octave_idx_type, First, n);
  OCTAVE_LOCAL_BUFFER (octave_idx_type, Level, n);

  cholmod_common Common;
  cholmod_common *cm = &Common;
  CHOLMOD_NAME(start) (cm);

  double spu = octave_sparse_params::get_key ("spumoni");
  if (spu == 0.)
    {
      cm->print = -1;
      SUITESPARSE_ASSIGN_FPTR (printf_func, cm->print_function, 0);
    }
  else
    {
      cm->print = static_cast<int> (spu) + 2;
      SUITESPARSE_ASSIGN_FPTR (printf_func, cm->print_function, &SparseCholPrint);
    }

  cm->error_handler = &SparseCholError;
  SUITESPARSE_ASSIGN_FPTR2 (divcomplex_func, cm->complex_divide, divcomplex);
  SUITESPARSE_ASSIGN_FPTR2 (hypot_func, cm->hypotenuse, hypot);

  cholmod_sparse *F = CHOLMOD_NAME(transpose) (A, 0, cm);
  cholmod_sparse *Aup, *Alo;

  if (A->stype == 1 || coletree)
    {
      Aup = A;
      Alo = F;
    }
  else
    {
      Aup = F;
      Alo = A;
    }

  CHOLMOD_NAME(etree) (Aup, Parent, cm);

  ColumnVector tmp (n);    // Declaration must precede any goto cleanup.
  std::string err_msg;

  if (cm->status < CHOLMOD_OK)
    {
      err_msg = "symbfact: matrix corrupted";
      goto cleanup;
    }

  if (CHOLMOD_NAME(postorder) (Parent, n, 0, Post, cm) != n)
    {
      err_msg = "symbfact: postorder failed";
      goto cleanup;
    }

  CHOLMOD_NAME(rowcolcounts) (Alo, 0, 0, Parent, Post, 0,
                              ColCount, First, Level, cm);

  if (cm->status < CHOLMOD_OK)
    {
      err_msg = "symbfact: matrix corrupted";
      goto cleanup;
    }

  if (nargout > 4)
    {
      cholmod_sparse *A1, *A2;

      if (A->stype == 1)
        {
          A1 = A;
          A2 = 0;
        }
      else if (A->stype == -1)
        {
          A1 = F;
          A2 = 0;
        }
      else if (coletree)
        {
          A1 = F;
          A2 = A;
        }
      else
        {
          A1 = A;
          A2 = F;
        }

      // count the total number of entries in L
      octave_idx_type lnz = 0;
      for (octave_idx_type j = 0 ; j < n ; j++)
        lnz += ColCount[j];

      // allocate the output matrix L (pattern-only)
      SparseBoolMatrix L (n, n, lnz);

      // initialize column pointers
      lnz = 0;
      for (octave_idx_type j = 0 ; j < n ; j++)
        {
          L.xcidx(j) = lnz;
          lnz += ColCount[j];
        }
      L.xcidx(n) = lnz;

      // create a copy of the column pointers
      octave_idx_type *W = First;
      for (octave_idx_type j = 0 ; j < n ; j++)
        W[j] = L.xcidx (j);

      // get workspace for computing one row of L
      cholmod_sparse *R
        = CHOLMOD_NAME(allocate_sparse) (n, 1, n, false, true,
                                         0, CHOLMOD_PATTERN, cm);
      octave_idx_type *Rp = static_cast<octave_idx_type *> (R->p);
      octave_idx_type *Ri = static_cast<octave_idx_type *> (R->i);

      // compute L one row at a time
      for (octave_idx_type k = 0 ; k < n ; k++)
        {
          // get the kth row of L and store in the columns of L
          CHOLMOD_NAME(row_subtree) (A1, A2, k, Parent, R, cm);
          for (octave_idx_type p = 0 ; p < Rp[1] ; p++)
            L.xridx (W[Ri[p]]++) = k;

          // add the diagonal entry
          L.xridx (W[k]++) = k;
        }

      // free workspace
      CHOLMOD_NAME(free_sparse) (&R, cm);

      // transpose L to get R, or leave as is
      if (nargin < 3)
        L = L.transpose ();

      // fill numerical values of L with one's
      for (octave_idx_type p = 0 ; p < lnz ; p++)
        L.xdata(p) = true;

      retval(4) = L;
    }

  if (nargout > 3)
    {
      for (octave_idx_type i = 0; i < n; i++)
        tmp(i) = Post[i] + 1;
      retval(3) = tmp;
    }

  if (nargout > 2)
    {
      for (octave_idx_type i = 0; i < n; i++)
        tmp(i) = Parent[i] + 1;
      retval(2) = tmp;
    }

  if (nargout > 1)
    {
      // compute the elimination tree height
      octave_idx_type height = 0;
      for (int i = 0 ; i < n ; i++)
        height = std::max (height, Level[i]);
      height++;
      retval(1) = static_cast<double> (height);
    }

  for (octave_idx_type i = 0; i < n; i++)
    tmp(i) = ColCount[i];
  retval(0) = tmp;

cleanup:
  CHOLMOD_NAME(free_sparse) (&F, cm);
  CHOLMOD_NAME(finish) (cm);

  if (! err_msg.empty ())
    error (err_msg.c_str ());

  return retval;

#else

  octave_unused_parameter (args);
  octave_unused_parameter (nargout);

  err_disabled_feature ("symbfact", "CHOLMOD");

#endif
}

/*
%!testif HAVE_CHOLMOD
%! A = sparse (magic (3));
%! [count, h, parent, post, r] = symbfact (A);
%! assert (count, [3; 2; 1]);
%! assert (h, 3);
%! assert (parent, [2; 3; 0]);
%! assert (r, sparse (triu (true (3))));

%!testif HAVE_CHOLMOD
%! ## Test MODE "lower"
%! A = sparse (magic (3));
%! [~, ~, ~, ~, l] = symbfact (A, "sym", "lower");
%! assert (l, sparse (tril (true (3))));

%!testif HAVE_CHOLMOD
%! ## Bug #42587, singular matrix
%! A = sparse ([1 0 8;0 1 8;8 8 1]);
%! [count, h, parent, post, r] = symbfact (A);

## Test input validation
%!testif HAVE_CHOLMOD
%! fail ("symbfact ()");
%! fail ("symbfact (1,2,3,4)");
%! fail ("symbfact ({1})", "wrong type argument 'cell'");
%! fail ("symbfact (sparse (1), {1})", "TYP must be a string");
%! fail ("symbfact (sparse (1), 'foobar')", 'unrecognized TYP "foobar"');
%! fail ("symbfact (sparse (1), 'sym', {'L'})", "MODE must be a string");
%! fail ('symbfact (sparse (1), "sym", "foobar")', 'unrecognized MODE "foobar"');
%! fail ("symbfact (sparse ([1, 2; 3, 4; 5, 6]))", "S must be a square matrix");

*/