view libinterp/dldfcn/amd.cc @ 21200:fcac5dbbf9ed

maint: Indent #ifdef blocks in libinterp. * builtins.h, Cell.cc, __contourc__.cc, __dispatch__.cc, __dsearchn__.cc, __ichol__.cc, __ilu__.cc, __lin_interpn__.cc, __pchip_deriv__.cc, __qp__.cc, balance.cc, besselj.cc, betainc.cc, bitfcns.cc, bsxfun.cc, c-file-ptr-stream.cc, c-file-ptr-stream.h, cellfun.cc, colloc.cc, comment-list.cc, conv2.cc, daspk.cc, dasrt.cc, dassl.cc, data.cc, debug.cc, defaults.cc, defaults.in.h, defun-dld.h, defun.cc, defun.h, det.cc, dirfns.cc, display.cc, dlmread.cc, dot.cc, dynamic-ld.cc, eig.cc, ellipj.cc, error.cc, errwarn.cc, event-queue.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, gl-render.cc, gl2ps-print.cc, graphics.cc, graphics.in.h, gripes.cc, hash.cc, help.cc, hess.cc, hex2num.cc, input.cc, inv.cc, jit-ir.cc, jit-typeinfo.cc, jit-util.cc, jit-util.h, kron.cc, load-path.cc, load-save.cc, lookup.cc, ls-ascii-helper.cc, ls-hdf5.cc, ls-mat-ascii.cc, ls-mat4.cc, ls-mat5.cc, ls-oct-binary.cc, ls-oct-text.cc, ls-oct-text.h, ls-utils.cc, ls-utils.h, lsode.cc, lu.cc, luinc.cc, mappers.cc, matrix_type.cc, max.cc, mex.h, mexproto.h, mgorth.cc, nproc.cc, oct-errno.in.cc, oct-fstrm.cc, oct-hdf5-types.cc, oct-hdf5.h, oct-hist.cc, oct-iostrm.cc, oct-lvalue.cc, oct-map.cc, oct-prcstrm.cc, oct-procbuf.cc, oct-stream.cc, oct-strstrm.cc, octave-link.cc, ordschur.cc, pager.cc, pinv.cc, pr-output.cc, procstream.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-xdiv.cc, sparse-xpow.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, sysdep.h, time.cc, toplev.cc, tril.cc, tsearch.cc, txt-eng-ft.cc, txt-eng.cc, typecast.cc, urlwrite.cc, utils.cc, variables.cc, xdiv.cc, xnorm.cc, xpow.cc, zfstream.cc, __delaunayn__.cc, __eigs__.cc, __fltk_uigetfile__.cc, __glpk__.cc, __init_fltk__.cc, __init_gnuplot__.cc, __magick_read__.cc, __osmesa_print__.cc, __voronoi__.cc, amd.cc, audiodevinfo.cc, audioread.cc, ccolamd.cc, chol.cc, colamd.cc, convhulln.cc, dmperm.cc, fftw.cc, oct-qhull.h, qr.cc, symbfact.cc, symrcm.cc, oct-conf.in.cc, ov-base-diag.cc, ov-base-int.cc, ov-base-mat.cc, ov-base-scalar.cc, ov-base-sparse.cc, ov-base.cc, ov-bool-mat.cc, ov-bool-sparse.cc, ov-bool.cc, ov-builtin.cc, ov-cell.cc, ov-ch-mat.cc, ov-class.cc, ov-classdef.cc, ov-colon.cc, ov-complex.cc, ov-cs-list.cc, ov-cx-diag.cc, ov-cx-mat.cc, ov-cx-sparse.cc, ov-dld-fcn.cc, ov-fcn-handle.cc, ov-fcn-inline.cc, ov-fcn.cc, ov-float.cc, ov-flt-complex.cc, ov-flt-cx-diag.cc, ov-flt-cx-mat.cc, ov-flt-re-diag.cc, ov-flt-re-mat.cc, ov-int16.cc, ov-int32.cc, ov-int64.cc, ov-int8.cc, ov-java.cc, ov-lazy-idx.cc, ov-mex-fcn.cc, ov-null-mat.cc, ov-oncleanup.cc, ov-perm.cc, ov-range.cc, ov-re-diag.cc, ov-re-mat.cc, ov-re-sparse.cc, ov-scalar.cc, ov-str-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, ovl.cc, octave.cc, op-b-b.cc, op-b-bm.cc, op-b-sbm.cc, op-bm-b.cc, op-bm-bm.cc, op-bm-sbm.cc, op-cdm-cdm.cc, op-cell.cc, op-chm.cc, op-class.cc, op-cm-cm.cc, op-cm-cs.cc, op-cm-m.cc, op-cm-s.cc, op-cm-scm.cc, op-cm-sm.cc, op-cs-cm.cc, op-cs-cs.cc, op-cs-m.cc, op-cs-s.cc, op-cs-scm.cc, op-cs-sm.cc, op-dm-dm.cc, op-dm-scm.cc, op-dm-sm.cc, op-dm-template.cc, op-dms-template.cc, op-double-conv.cc, op-fcdm-fcdm.cc, op-fcdm-fdm.cc, op-fcm-fcm.cc, op-fcm-fcs.cc, op-fcm-fm.cc, op-fcm-fs.cc, op-fcn.cc, op-fcs-fcm.cc, op-fcs-fcs.cc, op-fcs-fm.cc, op-fcs-fs.cc, op-fdm-fdm.cc, op-float-conv.cc, op-fm-fcm.cc, op-fm-fcs.cc, op-fm-fm.cc, op-fm-fs.cc, op-fs-fcm.cc, op-fs-fcs.cc, op-fs-fm.cc, op-fs-fs.cc, op-i16-i16.cc, op-i32-i32.cc, op-i64-i64.cc, op-i8-i8.cc, op-int-concat.cc, op-int-conv.cc, op-m-cm.cc, op-m-cs.cc, op-m-m.cc, op-m-s.cc, op-m-scm.cc, op-m-sm.cc, op-pm-pm.cc, op-pm-scm.cc, op-pm-sm.cc, op-pm-template.cc, op-range.cc, op-s-cm.cc, op-s-cs.cc, op-s-m.cc, op-s-s.cc, op-s-scm.cc, op-s-sm.cc, op-sbm-b.cc, op-sbm-bm.cc, op-sbm-sbm.cc, op-scm-cm.cc, op-scm-cs.cc, op-scm-m.cc, op-scm-s.cc, op-scm-scm.cc, op-scm-sm.cc, op-sm-cm.cc, op-sm-cs.cc, op-sm-m.cc, op-sm-s.cc, op-sm-scm.cc, op-sm-sm.cc, op-str-m.cc, op-str-s.cc, op-str-str.cc, op-struct.cc, op-ui16-ui16.cc, op-ui32-ui32.cc, op-ui64-ui64.cc, op-ui8-ui8.cc, pt-arg-list.cc, pt-array-list.cc, pt-assign.cc, pt-binop.cc, pt-bp.cc, pt-cbinop.cc, pt-cell.cc, pt-check.cc, pt-classdef.cc, pt-cmd.cc, pt-colon.cc, pt-colon.h, pt-const.cc, pt-decl.cc, pt-eval.cc, pt-except.cc, pt-exp.cc, pt-fcn-handle.cc, pt-funcall.cc, pt-id.cc, pt-idx.cc, pt-jump.cc, pt-loop.cc, pt-mat.cc, pt-misc.cc, pt-pr-code.cc, pt-select.cc, pt-stmt.cc, pt-unop.cc, pt.cc, token.cc, Array-jit.cc, Array-os.cc, Array-sym.cc, Array-tc.cc, version.cc: Indent #ifdef blocks in libinterp.
author Rik <rik@octave.org>
date Fri, 05 Feb 2016 16:29:08 -0800
parents e06e600f396a
children a83e7a384ee0
line wrap: on
line source

/*

Copyright (C) 2008-2015 David Bateman

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/>.

*/

// This is the octave interface to amd, which bore the copyright given
// in the help of the functions.

#ifdef HAVE_CONFIG_H
#  include <config.h>
#endif

#include <cstdlib>

#include <string>
#include <vector>

#include "ov.h"
#include "defun-dld.h"
#include "errwarn.h"
#include "pager.h"
#include "ov-re-mat.h"

#include "ov-re-sparse.h"
#include "ov-cx-sparse.h"
#include "oct-map.h"

#include "oct-sparse.h"
#include "oct-locbuf.h"

#if defined (ENABLE_64)
#  define AMD_NAME(name) amd_l ## name
#else
#  define AMD_NAME(name) amd ## name
#endif

DEFUN_DLD (amd, args, nargout,
           "-*- texinfo -*-\n\
@deftypefn  {} {@var{p} =} amd (@var{S})\n\
@deftypefnx {} {@var{p} =} amd (@var{S}, @var{opts})\n\
\n\
Return the approximate minimum degree permutation of a matrix.\n\
\n\
This is a permutation such that the Cholesky@tie{}factorization of\n\
@code{@var{S} (@var{p}, @var{p})} tends to be sparser than the\n\
Cholesky@tie{}factorization of @var{S} itself.  @code{amd} is typically\n\
faster than @code{symamd} but serves a similar purpose.\n\
\n\
The optional parameter @var{opts} is a structure that controls the behavior\n\
of @code{amd}.  The fields of the structure are\n\
\n\
@table @asis\n\
@item @var{opts}.dense\n\
Determines what @code{amd} considers to be a dense row or column of the\n\
input matrix.  Rows or columns with more than @code{max (16, (dense *\n\
sqrt (@var{n})))} entries, where @var{n} is the order of the matrix @var{S},\n\
are ignored by @code{amd} during the calculation of the permutation.\n\
The value of dense must be a positive scalar and the default value is 10.0\n\
\n\
@item @var{opts}.aggressive\n\
If this value is a nonzero scalar, then @code{amd} performs aggressive\n\
absorption.  The default is not to perform aggressive absorption.\n\
@end table\n\
\n\
The author of the code itself is Timothy A. Davis\n\
@email{davis@@cise.ufl.edu}, University of Florida\n\
(see @url{http://www.cise.ufl.edu/research/sparse/amd}).\n\
@seealso{symamd, colamd}\n\
@end deftypefn")
{
#ifdef HAVE_AMD
  int nargin = args.length ();

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

  octave_idx_type n_row, n_col;
  const octave_idx_type *ridx, *cidx;
  SparseMatrix sm;
  SparseComplexMatrix scm;

  if (args(0).is_sparse_type ())
    {
      if (args(0).is_complex_type ())
        {
          scm = args(0).sparse_complex_matrix_value ();
          n_row = scm.rows ();
          n_col = scm.cols ();
          ridx = scm.xridx ();
          cidx = scm.xcidx ();
        }
      else
        {
          sm = args(0).sparse_matrix_value ();
          n_row = sm.rows ();
          n_col = sm.cols ();
          ridx = sm.xridx ();
          cidx = sm.xcidx ();
        }
    }
  else
    {
      if (args(0).is_complex_type ())
        sm = SparseMatrix (real (args(0).complex_matrix_value ()));
      else
        sm = SparseMatrix (args(0).matrix_value ());

      n_row = sm.rows ();
      n_col = sm.cols ();
      ridx = sm.xridx ();
      cidx = sm.xcidx ();
    }

  if (n_row != n_col)
    err_square_matrix_required ("amd", "S");

  OCTAVE_LOCAL_BUFFER (double, Control, AMD_CONTROL);
  AMD_NAME (_defaults) (Control) ;
  if (nargin > 1)
    {
      octave_scalar_map arg1 = args(1).xscalar_map_value ("amd: OPTS argument must be a scalar structure");

      octave_value tmp;

      tmp = arg1.getfield ("dense");
      if (tmp.is_defined ())
        Control[AMD_DENSE] = tmp.double_value ();

      tmp = arg1.getfield ("aggressive");
      if (tmp.is_defined ())
        Control[AMD_AGGRESSIVE] = tmp.double_value ();
    }

  OCTAVE_LOCAL_BUFFER (octave_idx_type, P, n_col);
  Matrix xinfo (AMD_INFO, 1);
  double *Info = xinfo.fortran_vec ();

  // FIXME: how can we manage the memory allocation of amd
  //        in a cleaner manner?
  SUITESPARSE_ASSIGN_FPTR (malloc_func, amd_malloc, malloc);
  SUITESPARSE_ASSIGN_FPTR (free_func, amd_free, free);
  SUITESPARSE_ASSIGN_FPTR (calloc_func, amd_calloc, calloc);
  SUITESPARSE_ASSIGN_FPTR (realloc_func, amd_realloc, realloc);
  SUITESPARSE_ASSIGN_FPTR (printf_func, amd_printf, printf);

  octave_idx_type result = AMD_NAME (_order) (n_col, cidx, ridx, P,
                                              Control, Info);

  if (result == AMD_OUT_OF_MEMORY)
    error ("amd: out of memory");
  else if (result == AMD_INVALID)
    error ("amd: matrix S is corrupted");

  Matrix Pout (1, n_col);
  for (octave_idx_type i = 0; i < n_col; i++)
    Pout.xelem (i) = P[i] + 1;

  if (nargout > 1)
    return ovl (Pout, xinfo);
  else
    return ovl (Pout);

#else
  err_disabled_feature ("amd", "AMD");
#endif
}

/*
%!shared A, A2, opts
%! A = ones (20, 30);
%! A2 = ones (30, 30);
%!
%!testif HAVE_AMD
%! assert(amd (A2), [1:30])
%! opts.dense = 25;
%! assert(amd (A2, opts), [1:30])
%! opts.aggressive = 1;
%! assert(amd (A2, opts), [1:30])

%!error <S must be a square matrix|was unavailable or disabled> amd (A)
%!error amd (A2, 2)
%!error <matrix S is corrupted|was unavailable or disabled> amd ([])
*/