view libinterp/corefcn/chol.cc @ 31605:e88a07dec498 stable

maint: Use macros to begin/end C++ namespaces. * oct-conf-post-public.in.h: Define two macros (OCTAVE_BEGIN_NAMESPACE, OCTAVE_END_NAMESPACE) that can be used to start/end a namespace. * mk-opts.pl, build-env.h, build-env.in.cc, __betainc__.cc, __contourc__.cc, __dsearchn__.cc, __eigs__.cc, __expint__.cc, __ftp__.cc, __gammainc__.cc, __ichol__.cc, __ilu__.cc, __isprimelarge__.cc, __lin_interpn__.cc, __magick_read__.cc, __pchip_deriv__.cc, __qp__.cc, amd.cc, auto-shlib.cc, auto-shlib.h, balance.cc, base-text-renderer.cc, base-text-renderer.h, besselj.cc, bitfcns.cc, bsxfun.cc, c-file-ptr-stream.cc, c-file-ptr-stream.h, call-stack.cc, call-stack.h, ccolamd.cc, cellfun.cc, chol.cc, colamd.cc, colloc.cc, conv2.cc, daspk.cc, dasrt.cc, dassl.cc, data.cc, data.h, debug.cc, defaults.cc, defaults.h, defun-int.h, defun.cc, det.cc, dirfns.cc, display.cc, display.h, dlmread.cc, dmperm.cc, dot.cc, dynamic-ld.cc, dynamic-ld.h, eig.cc, ellipj.cc, environment.cc, environment.h, error.cc, error.h, errwarn.h, event-manager.cc, event-manager.h, event-queue.cc, event-queue.h, fcn-info.cc, fcn-info.h, fft.cc, fft2.cc, fftn.cc, file-io.cc, filter.cc, find.cc, ft-text-renderer.cc, ft-text-renderer.h, gcd.cc, getgrent.cc, getpwent.cc, getrusage.cc, givens.cc, gl-render.cc, gl-render.h, gl2ps-print.cc, gl2ps-print.h, graphics-toolkit.cc, graphics-toolkit.h, graphics.cc, graphics.in.h, gsvd.cc, gtk-manager.cc, gtk-manager.h, hash.cc, help.cc, help.h, hess.cc, hex2num.cc, hook-fcn.cc, hook-fcn.h, input.cc, input.h, interpreter-private.cc, interpreter-private.h, interpreter.cc, interpreter.h, inv.cc, jsondecode.cc, jsonencode.cc, kron.cc, latex-text-renderer.cc, latex-text-renderer.h, load-path.cc, load-path.h, load-save.cc, load-save.h, lookup.cc, ls-ascii-helper.cc, ls-ascii-helper.h, ls-oct-text.cc, ls-utils.cc, ls-utils.h, lsode.cc, lu.cc, mappers.cc, matrix_type.cc, max.cc, mex-private.h, mex.cc, mgorth.cc, nproc.cc, oct-fstrm.cc, oct-fstrm.h, oct-hdf5-types.cc, oct-hdf5-types.h, oct-hist.cc, oct-hist.h, oct-iostrm.cc, oct-iostrm.h, oct-opengl.h, oct-prcstrm.cc, oct-prcstrm.h, oct-procbuf.cc, oct-procbuf.h, oct-process.cc, oct-process.h, oct-stdstrm.h, oct-stream.cc, oct-stream.h, oct-strstrm.cc, oct-strstrm.h, oct-tex-lexer.in.ll, oct-tex-parser.yy, ordqz.cc, ordschur.cc, pager.cc, pager.h, pinv.cc, pow2.cc, pr-flt-fmt.cc, pr-output.cc, procstream.cc, procstream.h, psi.cc, qr.cc, quad.cc, quadcc.cc, qz.cc, rand.cc, rcond.cc, regexp.cc, schur.cc, settings.cc, settings.h, sighandlers.cc, sighandlers.h, sparse-xdiv.cc, sparse-xdiv.h, sparse-xpow.cc, sparse-xpow.h, sparse.cc, spparms.cc, sqrtm.cc, stack-frame.cc, stack-frame.h, stream-euler.cc, strfind.cc, strfns.cc, sub2ind.cc, svd.cc, sylvester.cc, symbfact.cc, syminfo.cc, syminfo.h, symrcm.cc, symrec.cc, symrec.h, symscope.cc, symscope.h, symtab.cc, symtab.h, syscalls.cc, sysdep.cc, sysdep.h, text-engine.cc, text-engine.h, text-renderer.cc, text-renderer.h, time.cc, toplev.cc, tril.cc, tsearch.cc, typecast.cc, url-handle-manager.cc, url-handle-manager.h, urlwrite.cc, utils.cc, utils.h, variables.cc, variables.h, xdiv.cc, xdiv.h, xnorm.cc, xnorm.h, xpow.cc, xpow.h, __delaunayn__.cc, __fltk_uigetfile__.cc, __glpk__.cc, __init_fltk__.cc, __init_gnuplot__.cc, __ode15__.cc, __voronoi__.cc, audiodevinfo.cc, audioread.cc, convhulln.cc, fftw.cc, gzip.cc, mk-build-env-features.sh, mk-builtins.pl, cdef-class.cc, cdef-class.h, cdef-fwd.h, cdef-manager.cc, cdef-manager.h, cdef-method.cc, cdef-method.h, cdef-object.cc, cdef-object.h, cdef-package.cc, cdef-package.h, cdef-property.cc, cdef-property.h, cdef-utils.cc, cdef-utils.h, ov-base.cc, ov-base.h, ov-bool-mat.cc, ov-builtin.h, ov-cell.cc, ov-class.cc, ov-class.h, ov-classdef.cc, ov-classdef.h, ov-complex.cc, ov-fcn-handle.cc, ov-fcn-handle.h, ov-fcn.h, ov-java.cc, ov-java.h, ov-mex-fcn.h, ov-null-mat.cc, ov-oncleanup.cc, ov-struct.cc, ov-typeinfo.cc, ov-typeinfo.h, ov-usr-fcn.cc, ov-usr-fcn.h, ov.cc, ov.h, octave.cc, octave.h, mk-ops.sh, 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-fcdm-fcdm.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-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-m-cm.cc, op-m-cs.cc, op-m-m.cc, op-m-s.cc, op-m-scm.cc, op-m-sm.cc, op-mi.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, ops.h, anon-fcn-validator.cc, anon-fcn-validator.h, bp-table.cc, bp-table.h, comment-list.cc, comment-list.h, filepos.h, lex.h, lex.ll, oct-lvalue.cc, oct-lvalue.h, oct-parse.yy, parse.h, profiler.cc, profiler.h, pt-anon-scopes.cc, pt-anon-scopes.h, pt-arg-list.cc, pt-arg-list.h, pt-args-block.cc, pt-args-block.h, pt-array-list.cc, pt-array-list.h, pt-assign.cc, pt-assign.h, pt-binop.cc, pt-binop.h, pt-bp.cc, pt-bp.h, pt-cbinop.cc, pt-cbinop.h, pt-cell.cc, pt-cell.h, pt-check.cc, pt-check.h, pt-classdef.cc, pt-classdef.h, pt-cmd.h, pt-colon.cc, pt-colon.h, pt-const.cc, pt-const.h, pt-decl.cc, pt-decl.h, pt-eval.cc, pt-eval.h, pt-except.cc, pt-except.h, pt-exp.cc, pt-exp.h, pt-fcn-handle.cc, pt-fcn-handle.h, pt-id.cc, pt-id.h, pt-idx.cc, pt-idx.h, pt-jump.h, pt-loop.cc, pt-loop.h, pt-mat.cc, pt-mat.h, pt-misc.cc, pt-misc.h, pt-pr-code.cc, pt-pr-code.h, pt-select.cc, pt-select.h, pt-spmd.cc, pt-spmd.h, pt-stmt.cc, pt-stmt.h, pt-tm-const.cc, pt-tm-const.h, pt-unop.cc, pt-unop.h, pt-vm-eval.cc, pt-walk.cc, pt-walk.h, pt.cc, pt.h, token.cc, token.h, Range.cc, Range.h, idx-vector.cc, idx-vector.h, range-fwd.h, CollocWt.cc, CollocWt.h, aepbalance.cc, aepbalance.h, chol.cc, chol.h, gepbalance.cc, gepbalance.h, gsvd.cc, gsvd.h, hess.cc, hess.h, lo-mappers.cc, lo-mappers.h, lo-specfun.cc, lo-specfun.h, lu.cc, lu.h, oct-convn.cc, oct-convn.h, oct-fftw.cc, oct-fftw.h, oct-norm.cc, oct-norm.h, oct-rand.cc, oct-rand.h, oct-spparms.cc, oct-spparms.h, qr.cc, qr.h, qrp.cc, qrp.h, randgamma.cc, randgamma.h, randmtzig.cc, randmtzig.h, randpoisson.cc, randpoisson.h, schur.cc, schur.h, sparse-chol.cc, sparse-chol.h, sparse-lu.cc, sparse-lu.h, sparse-qr.cc, sparse-qr.h, svd.cc, svd.h, child-list.cc, child-list.h, dir-ops.cc, dir-ops.h, file-ops.cc, file-ops.h, file-stat.cc, file-stat.h, lo-sysdep.cc, lo-sysdep.h, lo-sysinfo.cc, lo-sysinfo.h, mach-info.cc, mach-info.h, oct-env.cc, oct-env.h, oct-group.cc, oct-group.h, oct-password.cc, oct-password.h, oct-syscalls.cc, oct-syscalls.h, oct-time.cc, oct-time.h, oct-uname.cc, oct-uname.h, action-container.cc, action-container.h, base-list.h, cmd-edit.cc, cmd-edit.h, cmd-hist.cc, cmd-hist.h, f77-fcn.h, file-info.cc, file-info.h, lo-array-errwarn.cc, lo-array-errwarn.h, lo-hash.cc, lo-hash.h, lo-ieee.h, lo-regexp.cc, lo-regexp.h, lo-utils.cc, lo-utils.h, oct-base64.cc, oct-base64.h, oct-glob.cc, oct-glob.h, oct-inttypes.h, oct-mutex.cc, oct-mutex.h, oct-refcount.h, oct-shlib.cc, oct-shlib.h, oct-sparse.cc, oct-sparse.h, oct-string.h, octave-preserve-stream-state.h, pathsearch.cc, pathsearch.h, quit.cc, quit.h, unwind-prot.cc, unwind-prot.h, url-transfer.cc, url-transfer.h : Use new macros to begin/end C++ namespaces.
author Rik <rik@octave.org>
date Thu, 01 Dec 2022 14:23:45 -0800
parents 32d2b6604a9f
children aac27ad79be6
line wrap: on
line source

////////////////////////////////////////////////////////////////////////
//
// Copyright (C) 1996-2022 The Octave Project Developers
//
// See the file COPYRIGHT.md in the top-level directory of this
// distribution or <https://octave.org/copyright/>.
//
// 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
// <https://www.gnu.org/licenses/>.
//
////////////////////////////////////////////////////////////////////////

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

#include <string>

#include "Matrix.h"
#include "chol.h"
#include "oct-string.h"
#include "sparse-chol.h"
#include "sparse-util.h"

#include "defun.h"
#include "error.h"
#include "errwarn.h"
#include "ov.h"
#include "ovl.h"

OCTAVE_BEGIN_NAMESPACE(octave)

template <typename CHOLT>
static octave_value
get_chol (const CHOLT& fact)
{
  return octave_value (fact.chol_matrix());
}

template <typename CHOLT>
static octave_value
get_chol_r (const CHOLT& fact)
{
  return octave_value (fact.chol_matrix (),
                       MatrixType (MatrixType::Upper));
}

template <typename CHOLT>
static octave_value
get_chol_l (const CHOLT& fact)
{
  return octave_value (fact.chol_matrix ().transpose (),
                       MatrixType (MatrixType::Lower));
}

DEFUN (chol, args, nargout,
       doc: /* -*- texinfo -*-
@deftypefn  {} {@var{R} =} chol (@var{A})
@deftypefnx {} {[@var{R}, @var{p}] =} chol (@var{A})
@deftypefnx {} {[@var{R}, @var{p}, @var{Q}] =} chol (@var{A})
@deftypefnx {} {[@var{R}, @var{p}, @var{Q}] =} chol (@var{A}, "vector")
@deftypefnx {} {[@var{L}, @dots{}] =} chol (@dots{}, "lower")
@deftypefnx {} {[@var{R}, @dots{}] =} chol (@dots{}, "upper")
@cindex Cholesky factorization
Compute the upper Cholesky@tie{}factor, @var{R}, of the real symmetric
or complex Hermitian positive definite matrix @var{A}.

The upper Cholesky@tie{}factor @var{R} is computed by using the upper
triangular part of matrix @var{A} and is defined by
@tex
$ R^T R = A $.
@end tex
@ifnottex

@example
@var{R}' * @var{R} = @var{A}.
@end example

@end ifnottex

Calling @code{chol} using the optional @qcode{"upper"} flag has the
same behavior.  In contrast, using the optional @qcode{"lower"} flag,
@code{chol} returns the lower triangular factorization, computed by using
the lower triangular part of matrix @var{A}, such that
@tex
$ L L^T = A $.
@end tex
@ifnottex

@example
@var{L} * @var{L}' = @var{A}.
@end example

@end ifnottex

Called with one output argument @code{chol} fails if matrix @var{A} is
not positive definite.  Note that if matrix @var{A} is not real symmetric
or complex Hermitian then the lower triangular part is considered to be
the (complex conjugate) transpose of the upper triangular part, or vice
versa, given the @qcode{"lower"} flag.

Called with two or more output arguments @var{p} flags whether the matrix
@var{A} was positive definite and @code{chol} does not fail.  A zero value
of @var{p} indicates that matrix @var{A} is positive definite and @var{R}
gives the factorization.  Otherwise, @var{p} will have a positive value.

If called with three output arguments matrix @var{A} must be sparse and
a sparsity preserving row/column permutation is applied to matrix @var{A}
prior to the factorization.  That is @var{R} is the factorization of
@code{@var{A}(@var{Q},@var{Q})} such that
@tex
$ R^T R = Q^T A Q$.
@end tex
@ifnottex

@example
@var{R}' * @var{R} = @var{Q}' * @var{A} * @var{Q}.
@end example

@end ifnottex

The sparsity preserving permutation is generally returned as a matrix.
However, given the optional flag @qcode{"vector"}, @var{Q} will be
returned as a vector such that
@tex
$ R^T R = A (Q, Q)$.
@end tex
@ifnottex

@example
@var{R}' * @var{R} = @var{A}(@var{Q}, @var{Q}).
@end example

@end ifnottex

In general the lower triangular factorization is significantly faster for
sparse matrices.
@seealso{hess, lu, qr, qz, schur, svd, ichol, cholinv, chol2inv, cholupdate,
cholinsert, choldelete, cholshift}
@end deftypefn */)
{
  int nargin = args.length ();

  if (nargin < 1 || nargin > 3 || nargout > 3)
    print_usage ();
  if (nargout > 2 && ! args(0).issparse ())
    error ("chol: using three output arguments, matrix A must be sparse");

  bool LLt = false;
  bool vecout = false;

  int n = 1;
  while (n < nargin)
    {
      std::string tmp = args(n++).xstring_value ("chol: optional arguments must be strings");

      if (string::strcmpi (tmp, "vector"))
        vecout = true;
      else if (string::strcmpi (tmp, "lower"))
        LLt = true;
      else if (string::strcmpi (tmp, "upper"))
        LLt = false;
      else
        error (R"(chol: optional argument must be one of "vector", "lower", or "upper")");
    }

  octave_value_list retval;
  octave_value arg = args(0);

  if (arg.isempty ())
    return ovl (Matrix ());

  if (arg.issparse ())
    {
      octave_idx_type info;
      bool natural = (nargout != 3);
      bool force = nargout > 1;

      if (arg.isreal ())
        {
          SparseMatrix m = arg.sparse_matrix_value ();

          math::sparse_chol<SparseMatrix> fact (m, info, natural, force);

          if (nargout == 3)
            {
              if (vecout)
                retval(2) = fact.perm ();
              else
                retval(2) = fact.Q ();
            }

          if (nargout >= 2 || info == 0)
            {
              retval(1) = info;
              if (LLt)
                retval(0) = fact.L ();
              else
                retval(0) = fact.R ();
            }
          else
            error ("chol: input matrix must be positive definite");
        }
      else if (arg.iscomplex ())
        {
          SparseComplexMatrix m = arg.sparse_complex_matrix_value ();

          math::sparse_chol<SparseComplexMatrix> fact (m, info, natural, force);

          if (nargout == 3)
            {
              if (vecout)
                retval(2) = fact.perm ();
              else
                retval(2) = fact.Q ();
            }

          if (nargout >= 2 || info == 0)
            {
              retval(1) = info;
              if (LLt)
                retval(0) = fact.L ();
              else
                retval(0) = fact.R ();
            }
          else
            error ("chol: input matrix must be positive definite");
        }
      else
        err_wrong_type_arg ("chol", arg);
    }
  else if (arg.is_single_type ())
    {
      if (vecout)
        error (R"(chol: A must be sparse for the "vector" option)");
      if (arg.isreal ())
        {
          FloatMatrix m = arg.float_matrix_value ();

          octave_idx_type info;

          math::chol<FloatMatrix> fact (m, info, LLt != true);

          if (nargout == 2 || info == 0)
            retval = ovl (get_chol (fact), info);
          else
            error ("chol: input matrix must be positive definite");
        }
      else if (arg.iscomplex ())
        {
          FloatComplexMatrix m = arg.float_complex_matrix_value ();

          octave_idx_type info;

          math::chol<FloatComplexMatrix> fact (m, info, LLt != true);

          if (nargout == 2 || info == 0)
            retval = ovl (get_chol (fact), info);
          else
            error ("chol: input matrix must be positive definite");
        }
      else
        err_wrong_type_arg ("chol", arg);
    }
  else
    {
      if (vecout)
        error (R"(chol: A must be sparse for the "vector" option)");
      if (arg.isreal ())
        {
          Matrix m = arg.matrix_value ();

          octave_idx_type info;

          math::chol<Matrix> fact (m, info, LLt != true);

          if (nargout == 2 || info == 0)
            retval = ovl (get_chol (fact), info);
          else
            error ("chol: input matrix must be positive definite");
        }
      else if (arg.iscomplex ())
        {
          ComplexMatrix m = arg.complex_matrix_value ();

          octave_idx_type info;

          math::chol<ComplexMatrix> fact (m, info, LLt != true);

          if (nargout == 2 || info == 0)
            retval = ovl (get_chol (fact), info);
          else
            error ("chol: input matrix must be positive definite");
        }
      else
        err_wrong_type_arg ("chol", arg);
    }

  return retval;
}

/*
%!assert (chol ([2, 1; 1, 1]), [sqrt(2), 1/sqrt(2); 0, 1/sqrt(2)], sqrt (eps))
%!assert (chol (single ([2, 1; 1, 1])),
%!        single ([sqrt(2), 1/sqrt(2); 0, 1/sqrt(2)]), sqrt (eps ("single")))

%!assert (chol ([2, 1; 1, 1], "upper"), [sqrt(2), 1/sqrt(2); 0, 1/sqrt(2)],
%!        sqrt (eps))
%!assert (chol ([2, 1; 1, 1], "lower"), [sqrt(2), 0; 1/sqrt(2), 1/sqrt(2)],
%!        sqrt (eps))

%!assert (chol ([2, 1; 1, 1], "lower"), chol ([2, 1; 1, 1], "LoweR"))
%!assert (chol ([2, 1; 1, 1], "upper"), chol ([2, 1; 1, 1], "Upper"))

## Check the "vector" option which only affects the 3rd argument and
## is only valid for sparse input.
%!testif HAVE_CHOLMOD
%! a = sparse ([2 1; 1 1]);
%! r = sparse ([sqrt(2), 1/sqrt(2); 0, 1/sqrt(2)]);
%! [rd, pd, qd] = chol (a);
%! [rv, pv, qv] = chol (a, "vector");
%! assert (r, rd, eps)
%! assert (r, rv, eps)
%! assert (pd, 0)
%! assert (pd, pv)
%! assert (qd, sparse (eye (2)))
%! assert (qv, [1 2])
%!
%! [rv, pv, qv] = chol (a, "Vector"); # check case sensitivity
%! assert (r, rv, eps)
%! assert (pd, pv)
%! assert (qv, [1 2])

%!testif HAVE_CHOLMOD <*42587>
%! A = sparse ([1 0 8;0 1 8;8 8 1]);
%! [Q, p] = chol (A);
%! assert (p != 0);

%!error chol ()
%!error <matrix must be positive definite> chol ([1, 2; 3, 4])
%!error <requires square matrix> chol ([1, 2; 3, 4; 5, 6])
%!error <optional arguments must be strings> chol (1, 2)
%!error <optional argument must be one of "vector", "lower"> chol (1, "foobar")
%!error <matrix A must be sparse> [L,p,Q] = chol ([1, 2; 3, 4])
%!error <A must be sparse> [L, p] = chol ([1, 2; 3, 4], "vector")
*/

DEFUN (cholinv, args, ,
       doc: /* -*- texinfo -*-
@deftypefn {} {@var{Ainv} =} cholinv (@var{A})
Compute the inverse of the symmetric positive definite matrix @var{A} using
the Cholesky@tie{}factorization.
@seealso{chol, chol2inv, inv}
@end deftypefn */)
{
  if (args.length () != 1)
    print_usage ();

  octave_value retval;
  octave_value arg = args(0);

  octave_idx_type nr = arg.rows ();
  octave_idx_type nc = arg.columns ();

  if (nr == 0 || nc == 0)
    retval = Matrix ();
  else
    {
      if (arg.issparse ())
        {
          octave_idx_type info;

          if (arg.isreal ())
            {
              SparseMatrix m = arg.sparse_matrix_value ();

              math::sparse_chol<SparseMatrix> chol (m, info);

              if (info == 0)
                retval = chol.inverse ();
              else
                error ("cholinv: A must be positive definite");
            }
          else if (arg.iscomplex ())
            {
              SparseComplexMatrix m = arg.sparse_complex_matrix_value ();

              math::sparse_chol<SparseComplexMatrix> chol (m, info);

              if (info == 0)
                retval = chol.inverse ();
              else
                error ("cholinv: A must be positive definite");
            }
          else
            err_wrong_type_arg ("cholinv", arg);
        }
      else if (arg.is_single_type ())
        {
          if (arg.isreal ())
            {
              FloatMatrix m = arg.float_matrix_value ();

              octave_idx_type info;
              math::chol<FloatMatrix> chol (m, info);
              if (info == 0)
                retval = chol.inverse ();
              else
                error ("cholinv: A must be positive definite");
            }
          else if (arg.iscomplex ())
            {
              FloatComplexMatrix m = arg.float_complex_matrix_value ();

              octave_idx_type info;
              math::chol<FloatComplexMatrix> chol (m, info);
              if (info == 0)
                retval = chol.inverse ();
              else
                error ("cholinv: A must be positive definite");
            }
          else
            err_wrong_type_arg ("chol", arg);
        }
      else
        {
          if (arg.isreal ())
            {
              Matrix m = arg.matrix_value ();

              octave_idx_type info;
              math::chol<Matrix> chol (m, info);
              if (info == 0)
                retval = chol.inverse ();
              else
                error ("cholinv: A must be positive definite");
            }
          else if (arg.iscomplex ())
            {
              ComplexMatrix m = arg.complex_matrix_value ();

              octave_idx_type info;
              math::chol<ComplexMatrix> chol (m, info);
              if (info == 0)
                retval = chol.inverse ();
              else
                error ("cholinv: A must be positive definite");
            }
          else
            err_wrong_type_arg ("chol", arg);
        }
    }

  return retval;
}

/*
%!shared A, Ainv
%! A = [2,0.2;0.2,1];
%! Ainv = inv (A);
%!test
%! Ainv1 = cholinv (A);
%! assert (norm (Ainv-Ainv1), 0, 1e-10);
%!testif HAVE_CHOLMOD
%! Ainv2 = inv (sparse (A));
%! assert (norm (Ainv-Ainv2), 0, 1e-10);
%!testif HAVE_CHOLMOD
%! Ainv3 = cholinv (sparse (A));
%! assert (norm (Ainv-Ainv3), 0, 1e-10);
*/

DEFUN (chol2inv, args, ,
       doc: /* -*- texinfo -*-
@deftypefn {} {@var{Ainv} =} chol2inv (@var{R})
Invert a symmetric, positive definite square matrix from its Cholesky
decomposition, @var{R}.

Note that @var{R} should be an upper-triangular matrix with positive diagonal
elements.  @code{chol2inv (@var{U})} provides @code{inv (@var{R}'*@var{R})} but
is much faster than using @code{inv}.
@seealso{chol, cholinv, inv}
@end deftypefn */)
{
  if (args.length () != 1)
    print_usage ();

  octave_value retval;

  octave_value arg = args(0);

  octave_idx_type nr = arg.rows ();
  octave_idx_type nc = arg.columns ();

  if (nr == 0 || nc == 0)
    retval = Matrix ();
  else
    {
      if (arg.issparse ())
        {
          if (arg.isreal ())
            {
              SparseMatrix r = arg.sparse_matrix_value ();

              retval = math::chol2inv (r);
            }
          else if (arg.iscomplex ())
            {
              SparseComplexMatrix r = arg.sparse_complex_matrix_value ();

              retval = math::chol2inv (r);
            }
          else
            err_wrong_type_arg ("chol2inv", arg);
        }
      else if (arg.is_single_type ())
        {
          if (arg.isreal ())
            {
              FloatMatrix r = arg.float_matrix_value ();

              retval = math::chol2inv (r);
            }
          else if (arg.iscomplex ())
            {
              FloatComplexMatrix r = arg.float_complex_matrix_value ();

              retval = math::chol2inv (r);
            }
          else
            err_wrong_type_arg ("chol2inv", arg);

        }
      else
        {
          if (arg.isreal ())
            {
              Matrix r = arg.matrix_value ();

              retval = math::chol2inv (r);
            }
          else if (arg.iscomplex ())
            {
              ComplexMatrix r = arg.complex_matrix_value ();

              retval = math::chol2inv (r);
            }
          else
            err_wrong_type_arg ("chol2inv", arg);
        }
    }

  return retval;
}

/*

## Test for bug #36437
%!function sparse_chol2inv (T, tol)
%!  iT = inv (T);
%!  ciT = chol2inv (chol (T));
%!  assert (ciT, iT, tol);
%!  assert (chol2inv (chol ( full (T))), ciT, tol*2);
%!endfunction

%!testif HAVE_CHOLMOD
%! A = gallery ("poisson", 3);
%! sparse_chol2inv (A, eps);

%!testif HAVE_CHOLMOD
%! n = 10;
%! B = spdiags (ones (n, 1) * [1 2 1], [-1 0 1], n, n);
%! sparse_chol2inv (B, eps*100);

%!testif HAVE_CHOLMOD
%! C = gallery ("tridiag", 5);
%! sparse_chol2inv (C, eps*10);

%!testif HAVE_CHOLMOD
%! D = gallery ("wathen", 1, 1);
%! sparse_chol2inv (D, eps*10^4);

*/

DEFUN (cholupdate, args, nargout,
       doc: /* -*- texinfo -*-
@deftypefn {} {[@var{R1}, @var{info}] =} cholupdate (@var{R}, @var{u}, @var{op})
Update or downdate a Cholesky@tie{}factorization.

Given an upper triangular matrix @var{R} and a column vector @var{u},
attempt to determine another upper triangular matrix @var{R1} such that

@itemize @bullet
@item
@var{R1}'*@var{R1} = @var{R}'*@var{R} + @var{u}*@var{u}'
if @var{op} is @qcode{"+"}

@item
@var{R1}'*@var{R1} = @var{R}'*@var{R} - @var{u}*@var{u}'
if @var{op} is @qcode{"-"}
@end itemize

If @var{op} is @qcode{"-"}, @var{info} is set to

@itemize
@item 0 if the downdate was successful,

@item 1 if @var{R}'*@var{R} - @var{u}*@var{u}' is not positive definite,

@item 2 if @var{R} is singular.
@end itemize

If @var{info} is not present, an error message is printed in cases 1 and 2.
@seealso{chol, cholinsert, choldelete, cholshift}
@end deftypefn */)
{
  int nargin = args.length ();

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

  octave_value argr = args(0);
  octave_value argu = args(1);

  if (! argr.isnumeric () || ! argu.isnumeric ()
      || (nargin > 2 && ! args(2).is_string ()))
    print_usage ();

  octave_value_list retval (nargout == 2 ? 2 : 1);

  octave_idx_type n = argr.rows ();

  std::string op = (nargin < 3) ? "+" : args(2).string_value ();

  bool down = (op == "-");

  if (! down && op != "+")
    error (R"(cholupdate: OP must be "+" or "-")");

  if (argr.columns () != n || argu.rows () != n || argu.columns () != 1)
    error ("cholupdate: dimension mismatch between R and U");

  int err = 0;
  if (argr.is_single_type () || argu.is_single_type ())
    {
      if (argr.isreal () && argu.isreal ())
        {
          // real case
          FloatMatrix R = argr.float_matrix_value ();
          FloatColumnVector u = argu.float_column_vector_value ();

          math::chol<FloatMatrix> fact;
          fact.set (R);

          if (down)
            err = fact.downdate (u);
          else
            fact.update (u);

          retval = ovl (get_chol_r (fact));
        }
      else
        {
          // complex case
          FloatComplexMatrix R = argr.float_complex_matrix_value ();
          FloatComplexColumnVector u
            = argu.float_complex_column_vector_value ();

          math::chol<FloatComplexMatrix> fact;
          fact.set (R);

          if (down)
            err = fact.downdate (u);
          else
            fact.update (u);

          retval = ovl (get_chol_r (fact));
        }
    }
  else
    {
      if (argr.isreal () && argu.isreal ())
        {
          // real case
          Matrix R = argr.matrix_value ();
          ColumnVector u = argu.column_vector_value ();

          math::chol<Matrix> fact;
          fact.set (R);

          if (down)
            err = fact.downdate (u);
          else
            fact.update (u);

          retval = ovl (get_chol_r (fact));
        }
      else
        {
          // complex case
          ComplexMatrix R = argr.complex_matrix_value ();
          ComplexColumnVector u = argu.complex_column_vector_value ();

          math::chol<ComplexMatrix> fact;
          fact.set (R);

          if (down)
            err = fact.downdate (u);
          else
            fact.update (u);

          retval = ovl (get_chol_r (fact));
        }
    }

  if (nargout > 1)
    retval(1) = err;
  else if (err == 1)
    error ("cholupdate: downdate violates positiveness");
  else if (err == 2)
    error ("cholupdate: singular matrix");

  return retval;
}

/*
%!shared A, u, Ac, uc
%! A = [  0.436997  -0.131721   0.124120  -0.061673 ;
%!       -0.131721   0.738529   0.019851  -0.140295 ;
%!        0.124120   0.019851   0.354879  -0.059472 ;
%!       -0.061673  -0.140295  -0.059472   0.600939 ];
%!
%! u = [  0.98950 ;
%!        0.39844 ;
%!        0.63484 ;
%!        0.13351 ];
%! Ac = [  0.5585528 + 0.0000000i  -0.1662088 - 0.0315341i   0.0107873 + 0.0236411i  -0.0276775 - 0.0186073i ;
%!        -0.1662088 + 0.0315341i   0.6760061 + 0.0000000i   0.0011452 - 0.0475528i   0.0145967 + 0.0247641i ;
%!         0.0107873 - 0.0236411i   0.0011452 + 0.0475528i   0.6263149 - 0.0000000i  -0.1585837 - 0.0719763i ;
%!        -0.0276775 + 0.0186073i   0.0145967 - 0.0247641i  -0.1585837 + 0.0719763i   0.6034234 - 0.0000000i ];
%!
%! uc = [ 0.54267 + 0.91519i ;
%!        0.99647 + 0.43141i ;
%!        0.83760 + 0.68977i ;
%!        0.39160 + 0.90378i ];

%!test
%! R = chol (A);
%! R1 = cholupdate (R, u);
%! assert (norm (triu (R1)-R1, Inf), 0);
%! assert (norm (R1'*R1 - R'*R - u*u', Inf) < 1e1*eps);
%!
%! R1 = cholupdate (R1, u, "-");
%! assert (norm (triu (R1)-R1, Inf), 0);
%! assert (norm (R1 - R, Inf) < 1e1*eps);

%!test
%! R = chol (Ac);
%! R1 = cholupdate (R, uc);
%! assert (norm (triu (R1)-R1, Inf), 0);
%! assert (norm (R1'*R1 - R'*R - uc*uc', Inf) < 1e1*eps);
%!
%! R1 = cholupdate (R1, uc, "-");
%! assert (norm (triu (R1)-R1, Inf), 0);
%! assert (norm (R1 - R, Inf) < 1e1*eps);

%!test
%! R = chol (single (A));
%! R1 = cholupdate (R, single (u));
%! assert (norm (triu (R1)-R1, Inf), single (0));
%! assert (norm (R1'*R1 - R'*R - single (u*u'), Inf) < 1e1*eps ("single"));
%!
%! R1 = cholupdate (R1, single (u), "-");
%! assert (norm (triu (R1)-R1, Inf), single (0));
%! assert (norm (R1 - R, Inf) < 2e1*eps ("single"));

%!test
%! R = chol (single (Ac));
%! R1 = cholupdate (R, single (uc));
%! assert (norm (triu (R1)-R1, Inf), single (0));
%! assert (norm (R1'*R1 - R'*R - single (uc*uc'), Inf) < 1e1*eps ("single"));
%!
%! R1 = cholupdate (R1, single (uc), "-");
%! assert (norm (triu (R1)-R1, Inf), single (0));
%! assert (norm (R1 - R, Inf) < 2e1*eps ("single"));
*/

DEFUN (cholinsert, args, nargout,
       doc: /* -*- texinfo -*-
@deftypefn  {} {@var{R1} =} cholinsert (@var{R}, @var{j}, @var{u})
@deftypefnx {} {[@var{R1}, @var{info}] =} cholinsert (@var{R}, @var{j}, @var{u})
Update a Cholesky factorization given a row or column to insert in the
original factored matrix.

Given a Cholesky@tie{}factorization of a real symmetric or complex Hermitian
positive definite matrix @w{@var{A} = @var{R}'*@var{R}}, @var{R}@tie{}upper
triangular, return the Cholesky@tie{}factorization of
@var{A1}, where @w{A1(p,p) = A}, @w{A1(:,j) = A1(j,:)' = u} and
@w{p = [1:j-1,j+1:n+1]}.  @w{u(j)} should be positive.

On return, @var{info} is set to

@itemize
@item 0 if the insertion was successful,

@item 1 if @var{A1} is not positive definite,

@item 2 if @var{R} is singular.
@end itemize

If @var{info} is not present, an error message is printed in cases 1 and 2.
@seealso{chol, cholupdate, choldelete, cholshift}
@end deftypefn */)
{
  if (args.length () != 3)
    print_usage ();

  octave_value argr = args(0);
  octave_value argj = args(1);
  octave_value argu = args(2);

  if (! argr.isnumeric () || ! argu.isnumeric ()
      || ! argj.is_real_scalar ())
    print_usage ();

  octave_idx_type n = argr.rows ();
  octave_idx_type j = argj.scalar_value ();

  if (argr.columns () != n || argu.rows () != n+1 || argu.columns () != 1)
    error ("cholinsert: dimension mismatch between R and U");

  if (j < 1 || j > n+1)
    error ("cholinsert: index J out of range");

  octave_value_list retval (nargout == 2 ? 2 : 1);

  int err = 0;
  if (argr.is_single_type () || argu.is_single_type ())
    {
      if (argr.isreal () && argu.isreal ())
        {
          // real case
          FloatMatrix R = argr.float_matrix_value ();
          FloatColumnVector u = argu.float_column_vector_value ();

          math::chol<FloatMatrix> fact;
          fact.set (R);
          err = fact.insert_sym (u, j-1);

          retval = ovl (get_chol_r (fact));
        }
      else
        {
          // complex case
          FloatComplexMatrix R = argr.float_complex_matrix_value ();
          FloatComplexColumnVector u
            = argu.float_complex_column_vector_value ();

          math::chol<FloatComplexMatrix> fact;
          fact.set (R);
          err = fact.insert_sym (u, j-1);

          retval = ovl (get_chol_r (fact));
        }
    }
  else
    {
      if (argr.isreal () && argu.isreal ())
        {
          // real case
          Matrix R = argr.matrix_value ();
          ColumnVector u = argu.column_vector_value ();

          math::chol<Matrix> fact;
          fact.set (R);
          err = fact.insert_sym (u, j-1);

          retval = ovl (get_chol_r (fact));
        }
      else
        {
          // complex case
          ComplexMatrix R = argr.complex_matrix_value ();
          ComplexColumnVector u = argu.complex_column_vector_value ();

          math::chol<ComplexMatrix> fact;
          fact.set (R);
          err = fact.insert_sym (u, j-1);

          retval = ovl (get_chol_r (fact));
        }
    }

  if (nargout > 1)
    retval(1) = err;
  else if (err == 1)
    error ("cholinsert: insertion violates positiveness");
  else if (err == 2)
    error ("cholinsert: singular matrix");
  else if (err == 3)
    error ("cholinsert: diagonal element must be real");

  return retval;
}

/*
%!test
%! u2 = [  0.35080 ;
%!         0.63930 ;
%!         3.31057 ;
%!        -0.13825 ;
%!         0.45266 ];
%!
%! R = chol (A);
%!
%! j = 3;  p = [1:j-1, j+1:5];
%! R1 = cholinsert (R, j, u2);
%! A1 = R1'*R1;
%!
%! assert (norm (triu (R1)-R1, Inf), 0);
%! assert (norm (A1(p,p) - A, Inf) < 1e1*eps);

%!test
%! u2 = [  0.35080  + 0.04298i;
%!         0.63930  + 0.23778i;
%!         3.31057  + 0.00000i;
%!        -0.13825  + 0.19879i;
%!         0.45266  + 0.50020i];
%!
%! R = chol (Ac);
%!
%! j = 3;  p = [1:j-1, j+1:5];
%! R1 = cholinsert (R, j, u2);
%! A1 = R1'*R1;
%!
%! assert (norm (triu (R1)-R1, Inf), 0);
%! assert (norm (A1(p,p) - Ac, Inf) < 1e1*eps);

%!test
%! u2 = single ([  0.35080 ;
%!                 0.63930 ;
%!                 3.31057 ;
%!                -0.13825 ;
%!                 0.45266 ]);
%!
%! R = chol (single (A));
%!
%! j = 3;  p = [1:j-1, j+1:5];
%! R1 = cholinsert (R, j, u2);
%! A1 = R1'*R1;
%!
%! assert (norm (triu (R1)-R1, Inf), single (0));
%! assert (norm (A1(p,p) - A, Inf) < 1e1*eps ("single"));

%!test
%! u2 = single ([  0.35080  + 0.04298i;
%!                 0.63930  + 0.23778i;
%!                 3.31057  + 0.00000i;
%!                -0.13825  + 0.19879i;
%!                 0.45266  + 0.50020i]);
%!
%! R = chol (single (Ac));
%!
%! j = 3;  p = [1:j-1, j+1:5];
%! R1 = cholinsert (R, j, u2);
%! A1 = R1'*R1;
%!
%! assert (norm (triu (R1)-R1, Inf), single (0));
%! assert (norm (A1(p,p) - single (Ac), Inf) < 2e1*eps ("single"));

%!test
%! cu = chol (triu (A), "upper");
%! cl = chol (tril (A), "lower");
%! assert (cu, cl', eps);

%!test
%! cca   = chol (Ac);
%!
%! ccal  = chol (Ac, "lower");
%! ccal2 = chol (tril (Ac), "lower");
%!
%! ccau  = chol (Ac, "upper");
%! ccau2 = chol (triu (Ac), "upper");
%!
%! assert (cca'*cca,     Ac, eps);
%! assert (ccau'*ccau,   Ac, eps);
%! assert (ccau2'*ccau2, Ac, eps);
%!
%! assert (cca, ccal',  eps);
%! assert (cca, ccau,   eps);
%! assert (cca, ccal2', eps);
%! assert (cca, ccau2,  eps);

%!test
%! cca   = chol (single (Ac));
%!
%! ccal  = chol (single (Ac), "lower");
%! ccal2 = chol (tril (single (Ac)), "lower");
%!
%! ccau  = chol (single (Ac), "upper");
%! ccau2 = chol (triu (single (Ac)), "upper");
%!
%! assert (cca'*cca,     single (Ac), eps ("single"));
%! assert (ccau'*ccau,   single (Ac), eps ("single"));
%! assert (ccau2'*ccau2, single (Ac), eps ("single"));
%!
%! assert (cca, ccal',  eps ("single"));
%! assert (cca, ccau,   eps ("single"));
%! assert (cca, ccal2', eps ("single"));
%! assert (cca, ccau2,  eps ("single"));

%!test
%! a = [12,  2,  3,  4;
%!       2, 14,  5,  3;
%!       3,  5, 16,  6;
%!       4,  3,  6, 16];
%!
%! b = [0,  1,  2,  3;
%!     -1,  0,  1,  2;
%!     -2, -1,  0,  1;
%!     -3, -2, -1,  0];
%!
%! ca = a + i*b;
%!
%! cca   = chol (ca);
%!
%! ccal  = chol (ca, "lower");
%! ccal2 = chol (tril (ca), "lower");
%!
%! ccau  = chol (ca, "upper");
%! ccau2 = chol (triu (ca), "upper");
%!
%! assert (cca'*cca,     ca, 16*eps);
%! assert (ccau'*ccau,   ca, 16*eps);
%! assert (ccau2'*ccau2, ca, 16*eps);
%!
%! assert (cca, ccal',  16*eps);
%! assert (cca, ccau,   16*eps);
%! assert (cca, ccal2', 16*eps);
%! assert (cca, ccau2,  16*eps);
*/

DEFUN (choldelete, args, ,
       doc: /* -*- texinfo -*-
@deftypefn {} {@var{R1} =} choldelete (@var{R}, @var{j})
Update a Cholesky factorization given a row or column to delete from the
original factored matrix.

Given a Cholesky@tie{}factorization of a real symmetric or complex Hermitian
positive definite matrix @w{@var{A} = @var{R}'*@var{R}}, @var{R}@tie{}upper
triangular, return the Cholesky@tie{}factorization of @w{A(p,p)}, where
@w{p = [1:j-1,j+1:n+1]}.
@seealso{chol, cholupdate, cholinsert, cholshift}
@end deftypefn */)
{
  if (args.length () != 2)
    print_usage ();

  octave_value argr = args(0);
  octave_value argj = args(1);

  if (! argr.isnumeric () || ! argj.is_real_scalar ())
    print_usage ();

  octave_idx_type n = argr.rows ();
  octave_idx_type j = argj.scalar_value ();

  if (argr.columns () != n)
    err_square_matrix_required ("choldelete", "R");

  if (j < 0 && j > n)
    error ("choldelete: index J out of range");

  octave_value_list retval;

  if (argr.is_single_type ())
    {
      if (argr.isreal ())
        {
          // real case
          FloatMatrix R = argr.float_matrix_value ();

          math::chol<FloatMatrix> fact;
          fact.set (R);
          fact.delete_sym (j-1);

          retval = ovl (get_chol_r (fact));
        }
      else
        {
          // complex case
          FloatComplexMatrix R = argr.float_complex_matrix_value ();

          math::chol<FloatComplexMatrix> fact;
          fact.set (R);
          fact.delete_sym (j-1);

          retval = ovl (get_chol_r (fact));
        }
    }
  else
    {
      if (argr.isreal ())
        {
          // real case
          Matrix R = argr.matrix_value ();

          math::chol<Matrix> fact;
          fact.set (R);
          fact.delete_sym (j-1);

          retval = ovl (get_chol_r (fact));
        }
      else
        {
          // complex case
          ComplexMatrix R = argr.complex_matrix_value ();

          math::chol<ComplexMatrix> fact;
          fact.set (R);
          fact.delete_sym (j-1);

          retval = ovl (get_chol_r (fact));
        }
    }

  return retval;
}

/*
%!test
%! R = chol (A);
%!
%! j = 3;  p = [1:j-1,j+1:4];
%! R1 = choldelete (R, j);
%!
%! assert (norm (triu (R1)-R1, Inf), 0);
%! assert (norm (R1'*R1 - A(p,p), Inf) < 1e1*eps);

%!test
%! R = chol (Ac);
%!
%! j = 3;  p = [1:j-1,j+1:4];
%! R1 = choldelete (R, j);
%!
%! assert (norm (triu (R1)-R1, Inf), 0);
%! assert (norm (R1'*R1 - Ac(p,p), Inf) < 1e1*eps);

%!test
%! R = chol (single (A));
%!
%! j = 3;  p = [1:j-1,j+1:4];
%! R1 = choldelete (R, j);
%!
%! assert (norm (triu (R1)-R1, Inf), single (0));
%! assert (norm (R1'*R1 - single (A(p,p)), Inf) < 1e1*eps ("single"));

%!test
%! R = chol (single (Ac));
%!
%! j = 3;  p = [1:j-1,j+1:4];
%! R1 = choldelete (R,j);
%!
%! assert (norm (triu (R1)-R1, Inf), single (0));
%! assert (norm (R1'*R1 - single (Ac(p,p)), Inf) < 1e1*eps ("single"));
*/

DEFUN (cholshift, args, ,
       doc: /* -*- texinfo -*-
@deftypefn {} {@var{R1} =} cholshift (@var{R}, @var{i}, @var{j})
Update a Cholesky factorization given a range of columns to shift in the
original factored matrix.

Given a Cholesky@tie{}factorization of a real symmetric or complex Hermitian
positive definite matrix @w{@var{A} = @var{R}'*@var{R}}, @var{R}@tie{}upper
triangular, return the Cholesky@tie{}factorization of
@w{@var{A}(p,p)}, where @w{p} is the permutation @*
@code{p = [1:i-1, shift(i:j, 1), j+1:n]} if @w{@var{i} < @var{j}} @*
 or @*
@code{p = [1:j-1, shift(j:i,-1), i+1:n]} if @w{@var{j} < @var{i}}.  @*

@seealso{chol, cholupdate, cholinsert, choldelete}
@end deftypefn */)
{
  if (args.length () != 3)
    print_usage ();

  octave_value argr = args(0);
  octave_value argi = args(1);
  octave_value argj = args(2);

  if (! argr.isnumeric () || ! argi.is_real_scalar ()
      || ! argj.is_real_scalar ())
    print_usage ();

  octave_idx_type n = argr.rows ();
  octave_idx_type i = argi.scalar_value ();
  octave_idx_type j = argj.scalar_value ();

  if (argr.columns () != n)
    err_square_matrix_required ("cholshift", "R");

  if (j < 0 || j > n+1 || i < 0 || i > n+1)
    error ("cholshift: index I or J is out of range");

  octave_value_list retval;

  if (argr.is_single_type () && argi.is_single_type ()
      && argj.is_single_type ())
    {
      if (argr.isreal ())
        {
          // real case
          FloatMatrix R = argr.float_matrix_value ();

          math::chol<FloatMatrix> fact;
          fact.set (R);
          fact.shift_sym (i-1, j-1);

          retval = ovl (get_chol_r (fact));
        }
      else
        {
          // complex case
          FloatComplexMatrix R = argr.float_complex_matrix_value ();

          math::chol<FloatComplexMatrix> fact;
          fact.set (R);
          fact.shift_sym (i-1, j-1);

          retval = ovl (get_chol_r (fact));
        }
    }
  else
    {
      if (argr.isreal ())
        {
          // real case
          Matrix R = argr.matrix_value ();

          math::chol<Matrix> fact;
          fact.set (R);
          fact.shift_sym (i-1, j-1);

          retval = ovl (get_chol_r (fact));
        }
      else
        {
          // complex case
          ComplexMatrix R = argr.complex_matrix_value ();

          math::chol<ComplexMatrix> fact;
          fact.set (R);
          fact.shift_sym (i-1, j-1);

          retval = ovl (get_chol_r (fact));
        }
    }

  return retval;
}

/*
%!test
%! R = chol (A);
%!
%! i = 1;  j = 3;  p = [1:i-1, shift(i:j,-1), j+1:4];
%! R1 = cholshift (R, i, j);
%!
%! assert (norm (triu (R1)-R1, Inf), 0);
%! assert (norm (R1'*R1 - A(p,p), Inf) < 1e1*eps);
%!
%! j = 1;  i = 3;  p = [1:j-1, shift(j:i,+1), i+1:4];
%! R1 = cholshift (R, i, j);
%!
%! assert (norm (triu (R1) - R1, Inf), 0);
%! assert (norm (R1'*R1 - A(p,p), Inf) < 1e1*eps);

%!test
%! R = chol (Ac);
%!
%! i = 1;  j = 3;  p = [1:i-1, shift(i:j,-1), j+1:4];
%! R1 = cholshift (R, i, j);
%!
%! assert (norm (triu (R1)-R1, Inf), 0);
%! assert (norm (R1'*R1 - Ac(p,p), Inf) < 1e1*eps);
%!
%! j = 1;  i = 3;  p = [1:j-1, shift(j:i,+1), i+1:4];
%! R1 = cholshift (R, i, j);
%!
%! assert (norm (triu (R1)-R1, Inf), 0);
%! assert (norm (R1'*R1 - Ac(p,p), Inf) < 1e1*eps);

%!test
%! R = chol (single (A));
%!
%! i = 1;  j = 3;  p = [1:i-1, shift(i:j,-1), j+1:4];
%! R1 = cholshift (R, i, j);
%!
%! assert (norm (triu (R1)-R1, Inf), 0);
%! assert (norm (R1'*R1 - single (A(p,p)), Inf) < 1e1*eps ("single"));
%!
%! j = 1;  i = 3;  p = [1:j-1, shift(j:i,+1), i+1:4];
%! R1 = cholshift (R, i, j);
%!
%! assert (norm (triu (R1)-R1, Inf), 0);
%! assert (norm (R1'*R1 - single (A(p,p)), Inf) < 1e1*eps ("single"));

%!test
%! R = chol (single (Ac));
%!
%! i = 1;  j = 3;  p = [1:i-1, shift(i:j,-1), j+1:4];
%! R1 = cholshift (R, i, j);
%!
%! assert (norm (triu (R1)-R1, Inf), 0);
%! assert (norm (R1'*R1 - single (Ac(p,p)), Inf) < 1e1*eps ("single"));
%!
%! j = 1; i = 3; p = [1:j-1, shift(j:i,+1), i+1:4];
%! R1 = cholshift (R, i, j);
%!
%! assert (norm (triu (R1)-R1, Inf), 0);
%! assert (norm (R1'*R1 - single (Ac(p,p)), Inf) < 1e1*eps ("single"));
*/

OCTAVE_END_NAMESPACE(octave)