view liboctave/numeric/chol.cc @ 21301:40de9f8f23a6

Use '#include "config.h"' rather than <config.h>. * mk-octave-config-h.sh, mk-opts.pl, Backend.cc, BaseControl.cc, ButtonControl.cc, Canvas.cc, CheckBoxControl.cc, Container.cc, ContextMenu.cc, EditControl.cc, Figure.cc, FigureWindow.cc, GLCanvas.cc, KeyMap.cc, ListBoxControl.cc, Logger.cc, Menu.cc, MouseModeActionGroup.cc, Object.cc, ObjectFactory.cc, ObjectProxy.cc, Panel.cc, PopupMenuControl.cc, PushButtonControl.cc, PushTool.cc, QtHandlesUtils.cc, RadioButtonControl.cc, SliderControl.cc, TextControl.cc, TextEdit.cc, ToggleButtonControl.cc, ToggleTool.cc, ToolBar.cc, ToolBarButton.cc, __init_qt__.cc, annotation-dialog.cc, gl-select.cc, module.mk, kpty.cpp, color-picker.cc, dialog.cc, documentation-dock-widget.cc, files-dock-widget.cc, find-files-dialog.cc, find-files-model.cc, history-dock-widget.cc, file-editor-tab.cc, file-editor-tab.h, file-editor.cc, find-dialog.cc, marker.cc, octave-qscintilla.cc, octave-txt-lexer.cc, main-window.cc, octave-cmd.cc, octave-dock-widget.cc, octave-gui.cc, octave-interpreter.cc, octave-qt-link.cc, parser.cc, webinfo.cc, resource-manager.cc, settings-dialog.cc, shortcut-manager.cc, terminal-dock-widget.cc, thread-manager.cc, welcome-wizard.cc, workspace-model.cc, workspace-view.cc, build-env-features.sh, build-env.in.cc, 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, cdisplay.c, cellfun.cc, coct-hdf5-types.c, colloc.cc, comment-list.cc, conv2.cc, daspk.cc, dasrt.cc, dassl.cc, data.cc, debug.cc, defaults.cc, defun.cc, 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, ft-text-renderer.cc, gammainc.cc, gcd.cc, getgrent.cc, getpwent.cc, getrusage.cc, givens.cc, gl-render.cc, gl2ps-print.cc, graphics.cc, gripes.cc, hash.cc, help.cc, hess.cc, hex2num.cc, hook-fcn.cc, input.cc, inv.cc, jit-ir.cc, jit-typeinfo.cc, jit-util.cc, 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-utils.cc, lsode.cc, lu.cc, luinc.cc, mappers.cc, matrix_type.cc, max.cc, mex.cc, mgorth.cc, nproc.cc, oct-errno.in.cc, oct-fstrm.cc, oct-hdf5-types.cc, oct-hist.cc, oct-iostrm.cc, oct-lvalue.cc, oct-map.cc, oct-prcstrm.cc, oct-procbuf.cc, oct-stream.cc, oct-strstrm.cc, oct-tex-lexer.in.ll, oct-tex-parser.in.yy, 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, siglist.c, 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, text-renderer.cc, time.cc, toplev.cc, tril.cc, tsearch.cc, txt-eng.cc, typecast.cc, urlwrite.cc, utils.cc, variables.cc, xdiv.cc, xgl2ps.c, 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, qr.cc, symbfact.cc, symrcm.cc, mkbuiltins, mkops, 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, lex.ll, oct-parse.in.yy, 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-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, Array-C.cc, Array-b.cc, Array-ch.cc, Array-d.cc, Array-f.cc, Array-fC.cc, Array-i.cc, Array-idx-vec.cc, Array-s.cc, Array-str.cc, Array-util.cc, Array-voidp.cc, Array.cc, CColVector.cc, CDiagMatrix.cc, CMatrix.cc, CNDArray.cc, CRowVector.cc, CSparse.cc, DiagArray2.cc, MArray-C.cc, MArray-d.cc, MArray-f.cc, MArray-fC.cc, MArray-i.cc, MArray-s.cc, MArray.cc, MDiagArray2.cc, MSparse-C.cc, MSparse-d.cc, MatrixType.cc, PermMatrix.cc, Range.cc, Sparse-C.cc, Sparse-b.cc, Sparse-d.cc, Sparse.cc, boolMatrix.cc, boolNDArray.cc, boolSparse.cc, chMatrix.cc, chNDArray.cc, dColVector.cc, dDiagMatrix.cc, dMatrix.cc, dNDArray.cc, dRowVector.cc, dSparse.cc, dim-vector.cc, fCColVector.cc, fCDiagMatrix.cc, fCMatrix.cc, fCNDArray.cc, fCRowVector.cc, fColVector.cc, fDiagMatrix.cc, fMatrix.cc, fNDArray.cc, fRowVector.cc, idx-vector.cc, int16NDArray.cc, int32NDArray.cc, int64NDArray.cc, int8NDArray.cc, intNDArray.cc, uint16NDArray.cc, uint32NDArray.cc, uint64NDArray.cc, uint8NDArray.cc, blaswrap.c, cquit.c, f77-extern.cc, f77-fcn.c, lo-error.c, quit.cc, CollocWt.cc, DASPK.cc, DASRT.cc, DASSL.cc, EIG.cc, LSODE.cc, ODES.cc, Quad.cc, aepbalance.cc, chol.cc, eigs-base.cc, fEIG.cc, gepbalance.cc, hess.cc, lo-mappers.cc, lo-specfun.cc, lu.cc, oct-convn.cc, oct-fftw.cc, oct-norm.cc, oct-rand.cc, oct-spparms.cc, qr.cc, qrp.cc, randgamma.c, randmtzig.c, randpoisson.c, schur.cc, sparse-chol.cc, sparse-dmsolve.cc, sparse-lu.cc, sparse-qr.cc, svd.cc, mk-ops.awk, dir-ops.cc, file-ops.cc, file-stat.cc, lo-sysdep.cc, mach-info.cc, oct-env.cc, oct-group.cc, oct-passwd.cc, oct-syscalls.cc, oct-time.cc, oct-uname.cc, cmd-edit.cc, cmd-hist.cc, data-conv.cc, f2c-main.c, glob-match.cc, kpse.cc, lo-array-errwarn.cc, lo-array-gripes.cc, lo-cutils.c, lo-ieee.cc, lo-regexp.cc, lo-utils.cc, oct-base64.cc, oct-glob.cc, oct-inttypes.cc, oct-locbuf.cc, oct-mutex.cc, oct-rl-edit.c, oct-rl-hist.c, oct-shlib.cc, oct-sort.cc, pathsearch.cc, singleton-cleanup.cc, sparse-sort.cc, sparse-util.cc, str-vec.cc, unwind-prot.cc, url-transfer.cc, display-available.c, main-cli.cc, main-gui.cc, main.in.cc, mkoctfile.in.cc, octave-config.in.cc: Use '#include "config.h"' rather than <config.h>.
author Rik <rik@octave.org>
date Thu, 18 Feb 2016 13:34:50 -0800
parents d3b265a83adc
children aba2e6293dd8
line wrap: on
line source

/*

Copyright (C) 1994-2015 John W. Eaton
Copyright (C) 2008-2009 Jaroslav Hajek

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

*/

#ifdef HAVE_CONFIG_H
#  include "config.h"
#endif

#include <vector>

#include "CColVector.h"
#include "CMatrix.h"
#include "CRowVector.h"
#include "chol.h"
#include "dColVector.h"
#include "dMatrix.h"
#include "dRowVector.h"
#include "f77-fcn.h"
#include "fCColVector.h"
#include "fCMatrix.h"
#include "fCRowVector.h"
#include "fColVector.h"
#include "fMatrix.h"
#include "fRowVector.h"
#include "lo-error.h"
#include "oct-locbuf.h"
#include "oct-norm.h"

#if ! defined (HAVE_QRUPDATE)
#  include "qr.h"
#endif

extern "C"
{
  F77_RET_T
  F77_FUNC (dpotrf, DPOTRF) (F77_CONST_CHAR_ARG_DECL,
                             const octave_idx_type&, double*,
                             const octave_idx_type&, octave_idx_type&
                             F77_CHAR_ARG_LEN_DECL);

  F77_RET_T
  F77_FUNC (dpotri, DPOTRI) (F77_CONST_CHAR_ARG_DECL,
                             const octave_idx_type&, double*,
                             const octave_idx_type&, octave_idx_type&
                             F77_CHAR_ARG_LEN_DECL);

  F77_RET_T
  F77_FUNC (dpocon, DPOCON) (F77_CONST_CHAR_ARG_DECL,
                             const octave_idx_type&, double*,
                             const octave_idx_type&, const double&,
                             double&, double*, octave_idx_type*,
                             octave_idx_type&
                             F77_CHAR_ARG_LEN_DECL);
#ifdef HAVE_QRUPDATE

  F77_RET_T
  F77_FUNC (dch1up, DCH1UP) (const octave_idx_type&, double*,
                             const octave_idx_type&, double*, double*);

  F77_RET_T
  F77_FUNC (dch1dn, DCH1DN) (const octave_idx_type&, double*,
                             const octave_idx_type&, double*, double*,
                             octave_idx_type&);

  F77_RET_T
  F77_FUNC (dchinx, DCHINX) (const octave_idx_type&, double*,
                             const octave_idx_type&, const octave_idx_type&,
                             double*, double*, octave_idx_type&);

  F77_RET_T
  F77_FUNC (dchdex, DCHDEX) (const octave_idx_type&, double*,
                             const octave_idx_type&, const octave_idx_type&,
                             double*);

  F77_RET_T
  F77_FUNC (dchshx, DCHSHX) (const octave_idx_type&, double*,
                             const octave_idx_type&, const octave_idx_type&,
                             const octave_idx_type&, double*);
#endif

  F77_RET_T
  F77_FUNC (spotrf, SPOTRF) (F77_CONST_CHAR_ARG_DECL,
                             const octave_idx_type&, float*,
                             const octave_idx_type&, octave_idx_type&
                             F77_CHAR_ARG_LEN_DECL);

  F77_RET_T
  F77_FUNC (spotri, SPOTRI) (F77_CONST_CHAR_ARG_DECL,
                             const octave_idx_type&, float*,
                             const octave_idx_type&, octave_idx_type&
                             F77_CHAR_ARG_LEN_DECL);

  F77_RET_T
  F77_FUNC (spocon, SPOCON) (F77_CONST_CHAR_ARG_DECL,
                             const octave_idx_type&, float*,
                             const octave_idx_type&, const float&,
                             float&, float*, octave_idx_type*,
                             octave_idx_type&
                             F77_CHAR_ARG_LEN_DECL);
#ifdef HAVE_QRUPDATE

  F77_RET_T
  F77_FUNC (sch1up, SCH1UP) (const octave_idx_type&, float*,
                             const octave_idx_type&, float*, float*);

  F77_RET_T
  F77_FUNC (sch1dn, SCH1DN) (const octave_idx_type&, float*,
                             const octave_idx_type&, float*, float*,
                             octave_idx_type&);

  F77_RET_T
  F77_FUNC (schinx, SCHINX) (const octave_idx_type&, float*,
                             const octave_idx_type&, const octave_idx_type&,
                             float*, float*, octave_idx_type&);

  F77_RET_T
  F77_FUNC (schdex, SCHDEX) (const octave_idx_type&, float*,
                             const octave_idx_type&, const octave_idx_type&,
                             float*);

  F77_RET_T
  F77_FUNC (schshx, SCHSHX) (const octave_idx_type&, float*,
                             const octave_idx_type&, const octave_idx_type&,
                             const octave_idx_type&, float*);
#endif

  F77_RET_T
  F77_FUNC (zpotrf, ZPOTRF) (F77_CONST_CHAR_ARG_DECL,
                             const octave_idx_type&, Complex*,
                             const octave_idx_type&, octave_idx_type&
                             F77_CHAR_ARG_LEN_DECL);
  F77_RET_T
  F77_FUNC (zpotri, ZPOTRI) (F77_CONST_CHAR_ARG_DECL,
                             const octave_idx_type&, Complex*,
                             const octave_idx_type&, octave_idx_type&
                             F77_CHAR_ARG_LEN_DECL);

  F77_RET_T
  F77_FUNC (zpocon, ZPOCON) (F77_CONST_CHAR_ARG_DECL,
                             const octave_idx_type&, Complex*,
                             const octave_idx_type&, const double&,
                             double&, Complex*, double*, octave_idx_type&
                             F77_CHAR_ARG_LEN_DECL);
#ifdef HAVE_QRUPDATE

  F77_RET_T
  F77_FUNC (zch1up, ZCH1UP) (const octave_idx_type&, Complex*,
                             const octave_idx_type&, Complex*, double*);

  F77_RET_T
  F77_FUNC (zch1dn, ZCH1DN) (const octave_idx_type&, Complex*,
                             const octave_idx_type&, Complex*, double*,
                             octave_idx_type&);

  F77_RET_T
  F77_FUNC (zchinx, ZCHINX) (const octave_idx_type&, Complex*,
                             const octave_idx_type&, const octave_idx_type&,
                             Complex*, double*, octave_idx_type&);

  F77_RET_T
  F77_FUNC (zchdex, ZCHDEX) (const octave_idx_type&, Complex*,
                             const octave_idx_type&, const octave_idx_type&,
                             double*);

  F77_RET_T
  F77_FUNC (zchshx, ZCHSHX) (const octave_idx_type&, Complex*,
                             const octave_idx_type&, const octave_idx_type&,
                             const octave_idx_type&, Complex*, double*);
#endif

  F77_RET_T
  F77_FUNC (cpotrf, CPOTRF) (F77_CONST_CHAR_ARG_DECL,
                             const octave_idx_type&, FloatComplex*,
                             const octave_idx_type&, octave_idx_type&
                             F77_CHAR_ARG_LEN_DECL);
  F77_RET_T
  F77_FUNC (cpotri, CPOTRI) (F77_CONST_CHAR_ARG_DECL,
                             const octave_idx_type&, FloatComplex*,
                             const octave_idx_type&, octave_idx_type&
                             F77_CHAR_ARG_LEN_DECL);

  F77_RET_T
  F77_FUNC (cpocon, CPOCON) (F77_CONST_CHAR_ARG_DECL,
                             const octave_idx_type&, FloatComplex*,
                             const octave_idx_type&, const float&,
                             float&, FloatComplex*, float*, octave_idx_type&
                             F77_CHAR_ARG_LEN_DECL);
#ifdef HAVE_QRUPDATE

  F77_RET_T
  F77_FUNC (cch1up, CCH1UP) (const octave_idx_type&, FloatComplex*,
                             const octave_idx_type&, FloatComplex*, float*);

  F77_RET_T
  F77_FUNC (cch1dn, CCH1DN) (const octave_idx_type&, FloatComplex*,
                             const octave_idx_type&, FloatComplex*,
                             float*, octave_idx_type&);

  F77_RET_T
  F77_FUNC (cchinx, CCHINX) (const octave_idx_type&, FloatComplex*,
                             const octave_idx_type&, const octave_idx_type&,
                             FloatComplex*, float*, octave_idx_type&);

  F77_RET_T
  F77_FUNC (cchdex, CCHDEX) (const octave_idx_type&, FloatComplex*,
                             const octave_idx_type&, const octave_idx_type&,
                             float*);

  F77_RET_T
  F77_FUNC (cchshx, CCHSHX) (const octave_idx_type&, FloatComplex*,
                             const octave_idx_type&, const octave_idx_type&,
                             const octave_idx_type&, FloatComplex*, float*);
#endif
}

static Matrix
chol2inv_internal (const Matrix& r, bool is_upper = true)
{
  Matrix retval;

  octave_idx_type r_nr = r.rows ();
  octave_idx_type r_nc = r.cols ();

  if (r_nr != r_nc)
    (*current_liboctave_error_handler) ("chol2inv requires square matrix");

  octave_idx_type n = r_nc;
  octave_idx_type info = 0;

  Matrix tmp = r;
  double *v = tmp.fortran_vec ();

  if (info == 0)
    {
      if (is_upper)
        F77_XFCN (dpotri, DPOTRI, (F77_CONST_CHAR_ARG2 ("U", 1), n,
                                   v, n, info
                                   F77_CHAR_ARG_LEN (1)));
      else
        F77_XFCN (dpotri, DPOTRI, (F77_CONST_CHAR_ARG2 ("L", 1), n,
                                   v, n, info
                                   F77_CHAR_ARG_LEN (1)));

      // If someone thinks of a more graceful way of doing this (or
      // faster for that matter :-)), please let me know!

      if (n > 1)
        {
          if (is_upper)
            for (octave_idx_type j = 0; j < r_nc; j++)
              for (octave_idx_type i = j+1; i < r_nr; i++)
                tmp.xelem (i, j) = tmp.xelem (j, i);
          else
            for (octave_idx_type j = 0; j < r_nc; j++)
              for (octave_idx_type i = j+1; i < r_nr; i++)
                tmp.xelem (j, i) = tmp.xelem (i, j);
        }

      retval = tmp;
    }

  return retval;
}

static FloatMatrix
chol2inv_internal (const FloatMatrix& r, bool is_upper = true)
{
  FloatMatrix retval;

  octave_idx_type r_nr = r.rows ();
  octave_idx_type r_nc = r.cols ();

  if (r_nr != r_nc)
    (*current_liboctave_error_handler) ("chol2inv requires square matrix");

  octave_idx_type n = r_nc;
  octave_idx_type info = 0;

  FloatMatrix tmp = r;
  float *v = tmp.fortran_vec ();

  if (info == 0)
    {
      if (is_upper)
        F77_XFCN (spotri, SPOTRI, (F77_CONST_CHAR_ARG2 ("U", 1), n,
                                   v, n, info
                                   F77_CHAR_ARG_LEN (1)));
      else
        F77_XFCN (spotri, SPOTRI, (F77_CONST_CHAR_ARG2 ("L", 1), n,
                                   v, n, info
                                   F77_CHAR_ARG_LEN (1)));

      // If someone thinks of a more graceful way of doing this (or
      // faster for that matter :-)), please let me know!

      if (n > 1)
        {
          if (is_upper)
            for (octave_idx_type j = 0; j < r_nc; j++)
              for (octave_idx_type i = j+1; i < r_nr; i++)
                tmp.xelem (i, j) = tmp.xelem (j, i);
          else
            for (octave_idx_type j = 0; j < r_nc; j++)
              for (octave_idx_type i = j+1; i < r_nr; i++)
                tmp.xelem (j, i) = tmp.xelem (i, j);
        }

      retval = tmp;
    }

  return retval;
}

static ComplexMatrix
chol2inv_internal (const ComplexMatrix& r, bool is_upper = true)
{
  ComplexMatrix retval;

  octave_idx_type r_nr = r.rows ();
  octave_idx_type r_nc = r.cols ();

  if (r_nr != r_nc)
    (*current_liboctave_error_handler) ("chol2inv requires square matrix");

  octave_idx_type n = r_nc;
  octave_idx_type info;

  ComplexMatrix tmp = r;

  if (is_upper)
    F77_XFCN (zpotri, ZPOTRI, (F77_CONST_CHAR_ARG2 ("U", 1), n,
                               tmp.fortran_vec (), n, info
                               F77_CHAR_ARG_LEN (1)));
  else
    F77_XFCN (zpotri, ZPOTRI, (F77_CONST_CHAR_ARG2 ("L", 1), n,
                               tmp.fortran_vec (), n, info
                               F77_CHAR_ARG_LEN (1)));

  // If someone thinks of a more graceful way of doing this (or
  // faster for that matter :-)), please let me know!

  if (n > 1)
    {
      if (is_upper)
        for (octave_idx_type j = 0; j < r_nc; j++)
          for (octave_idx_type i = j+1; i < r_nr; i++)
            tmp.xelem (i, j) = tmp.xelem (j, i);
      else
        for (octave_idx_type j = 0; j < r_nc; j++)
          for (octave_idx_type i = j+1; i < r_nr; i++)
            tmp.xelem (j, i) = tmp.xelem (i, j);
    }

  retval = tmp;

  return retval;
}

static FloatComplexMatrix
chol2inv_internal (const FloatComplexMatrix& r, bool is_upper = true)
{
  FloatComplexMatrix retval;

  octave_idx_type r_nr = r.rows ();
  octave_idx_type r_nc = r.cols ();

  if (r_nr != r_nc)
    (*current_liboctave_error_handler) ("chol2inv requires square matrix");

  octave_idx_type n = r_nc;
  octave_idx_type info;

  FloatComplexMatrix tmp = r;

  if (is_upper)
    F77_XFCN (cpotri, CPOTRI, (F77_CONST_CHAR_ARG2 ("U", 1), n,
                               tmp.fortran_vec (), n, info
                               F77_CHAR_ARG_LEN (1)));
  else
    F77_XFCN (cpotri, CPOTRI, (F77_CONST_CHAR_ARG2 ("L", 1), n,
                               tmp.fortran_vec (), n, info
                               F77_CHAR_ARG_LEN (1)));

  // If someone thinks of a more graceful way of doing this (or
  // faster for that matter :-)), please let me know!

  if (n > 1)
    {
      if (is_upper)
        for (octave_idx_type j = 0; j < r_nc; j++)
          for (octave_idx_type i = j+1; i < r_nr; i++)
            tmp.xelem (i, j) = tmp.xelem (j, i);
      else
        for (octave_idx_type j = 0; j < r_nc; j++)
          for (octave_idx_type i = j+1; i < r_nr; i++)
            tmp.xelem (j, i) = tmp.xelem (i, j);
    }

  retval = tmp;

  return retval;
}

template <typename T>
T
chol2inv (const T& r)
{
  return chol2inv_internal (r);
}

// Compute the inverse of a matrix using the Cholesky factorization.
template <typename T>
T
chol<T>::inverse (void) const
{
  return chol2inv_internal (chol_mat, is_upper);
}

template <typename T>
void
chol<T>::set (const T& R)
{
  if (! R.is_square ())
    (*current_liboctave_error_handler) ("chol: requires square matrix");

  chol_mat = R;
}

#if ! defined (HAVE_QRUPDATE)

template <typename T>
void
chol<T>::update (const VT& u)
{
  warn_qrupdate_once ();

  octave_idx_type n = chol_mat.rows ();

  if (u.numel () != n)
    (*current_liboctave_error_handler) ("cholupdate: dimension mismatch");

  init (chol_mat.hermitian () * chol_mat + T (u) * T (u).hermitian (),
        true, false);
}

template <typename T>
static bool
singular (const T& a)
{
  static typename T::element_type zero (0);
  for (octave_idx_type i = 0; i < a.rows (); i++)
    if (a(i,i) == zero) return true;
  return false;
}

template <typename T>
octave_idx_type
chol<T>::downdate (const VT& u)
{
  warn_qrupdate_once ();

  octave_idx_type info = -1;

  octave_idx_type n = chol_mat.rows ();

  if (u.numel () != n)
    (*current_liboctave_error_handler) ("cholupdate: dimension mismatch");

  if (singular (chol_mat))
    info = 2;
  else
    {
      info = init (chol_mat.hermitian () * chol_mat
                   - T (u) * T (u).hermitian (), true, false);
      if (info) info = 1;
    }

  return info;
}

template <typename T>
octave_idx_type
chol<T>::insert_sym (const VT& u, octave_idx_type j)
{
  static typename T::element_type zero (0);

  warn_qrupdate_once ();

  octave_idx_type info = -1;

  octave_idx_type n = chol_mat.rows ();

  if (u.numel () != n + 1)
    (*current_liboctave_error_handler) ("cholinsert: dimension mismatch");
  if (j < 0 || j > n)
    (*current_liboctave_error_handler) ("cholinsert: index out of range");

  if (singular (chol_mat))
    info = 2;
  else if (imag (u(j)) != zero)
    info = 3;
  else
    {
      T a = chol_mat.hermitian () * chol_mat;
      T a1 (n+1, n+1);
      for (octave_idx_type k = 0; k < n+1; k++)
        for (octave_idx_type l = 0; l < n+1; l++)
          {
            if (l == j)
              a1(k, l) = u(k);
            else if (k == j)
              a1(k, l) = conj (u(l));
            else
              a1(k, l) = a(k < j ? k : k-1, l < j ? l : l-1);
          }
      info = init (a1, true, false);
      if (info) info = 1;
    }

  return info;
}

template <typename T>
void
chol<T>::delete_sym (octave_idx_type j)
{
  warn_qrupdate_once ();

  octave_idx_type n = chol_mat.rows ();

  if (j < 0 || j > n-1)
    (*current_liboctave_error_handler) ("choldelete: index out of range");

  T a = chol_mat.hermitian () * chol_mat;
  a.delete_elements (1, idx_vector (j));
  a.delete_elements (0, idx_vector (j));
  init (a, true, false);
}

template <typename T>
void
chol<T>::shift_sym (octave_idx_type i, octave_idx_type j)
{
  warn_qrupdate_once ();

  octave_idx_type n = chol_mat.rows ();

  if (i < 0 || i > n-1 || j < 0 || j > n-1)
    (*current_liboctave_error_handler) ("cholshift: index out of range");

  T a = chol_mat.hermitian () * chol_mat;
  Array<octave_idx_type> p (dim_vector (n, 1));
  for (octave_idx_type k = 0; k < n; k++) p(k) = k;
  if (i < j)
    {
      for (octave_idx_type k = i; k < j; k++) p(k) = k+1;
      p(j) = i;
    }
  else if (j < i)
    {
      p(j) = i;
      for (octave_idx_type k = j+1; k < i+1; k++) p(k) = k-1;
    }

  init (a.index (idx_vector (p), idx_vector (p)), true, false);
}

#endif

// Specializations.

template <>
octave_idx_type
chol<Matrix>::init (const Matrix& a, bool upper, bool calc_cond)
{
  octave_idx_type a_nr = a.rows ();
  octave_idx_type a_nc = a.cols ();

  if (a_nr != a_nc)
    (*current_liboctave_error_handler) ("chol: requires square matrix");

  octave_idx_type n = a_nc;
  octave_idx_type info;

  is_upper = upper;

  chol_mat.clear (n, n);
  if (is_upper)
    for (octave_idx_type j = 0; j < n; j++)
      {
        for (octave_idx_type i = 0; i <= j; i++)
          chol_mat.xelem (i, j) = a(i, j);
        for (octave_idx_type i = j+1; i < n; i++)
          chol_mat.xelem (i, j) = 0.0;
      }
  else
    for (octave_idx_type j = 0; j < n; j++)
      {
        for (octave_idx_type i = 0; i < j; i++)
          chol_mat.xelem (i, j) = 0.0;
        for (octave_idx_type i = j; i < n; i++)
          chol_mat.xelem (i, j) = a(i, j);
      }
  double *h = chol_mat.fortran_vec ();

  // Calculate the norm of the matrix, for later use.
  double anorm = 0;
  if (calc_cond)
    anorm = xnorm (a, 1);

  if (is_upper)
    F77_XFCN (dpotrf, DPOTRF, (F77_CONST_CHAR_ARG2 ("U", 1), n, h, n, info
                               F77_CHAR_ARG_LEN (1)));
  else
    F77_XFCN (dpotrf, DPOTRF, (F77_CONST_CHAR_ARG2 ("L", 1), n, h, n, info
                               F77_CHAR_ARG_LEN (1)));

  xrcond = 0.0;
  if (info > 0)
    chol_mat.resize (info - 1, info - 1);
  else if (calc_cond)
    {
      octave_idx_type dpocon_info = 0;

      // Now calculate the condition number for non-singular matrix.
      Array<double> z (dim_vector (3*n, 1));
      double *pz = z.fortran_vec ();
      Array<octave_idx_type> iz (dim_vector (n, 1));
      octave_idx_type *piz = iz.fortran_vec ();
      if (is_upper)
        F77_XFCN (dpocon, DPOCON, (F77_CONST_CHAR_ARG2 ("U", 1), n, h,
                                   n, anorm, xrcond, pz, piz, dpocon_info
                                   F77_CHAR_ARG_LEN (1)));
      else
        F77_XFCN (dpocon, DPOCON, (F77_CONST_CHAR_ARG2 ("L", 1), n, h,
                                   n, anorm, xrcond, pz, piz, dpocon_info
                                   F77_CHAR_ARG_LEN (1)));

      if (dpocon_info != 0)
        info = -1;
    }

  return info;
}

#if defined (HAVE_QRUPDATE)

template <>
void
chol<Matrix>::update (const ColumnVector& u)
{
  octave_idx_type n = chol_mat.rows ();

  if (u.numel () != n)
    (*current_liboctave_error_handler) ("cholupdate: dimension mismatch");

  ColumnVector utmp = u;

  OCTAVE_LOCAL_BUFFER (double, w, n);

  F77_XFCN (dch1up, DCH1UP, (n, chol_mat.fortran_vec (), chol_mat.rows (),
                             utmp.fortran_vec (), w));
}

template <>
octave_idx_type
chol<Matrix>::downdate (const ColumnVector& u)
{
  octave_idx_type info = -1;

  octave_idx_type n = chol_mat.rows ();

  if (u.numel () != n)
    (*current_liboctave_error_handler) ("cholupdate: dimension mismatch");

  ColumnVector utmp = u;

  OCTAVE_LOCAL_BUFFER (double, w, n);

  F77_XFCN (dch1dn, DCH1DN, (n, chol_mat.fortran_vec (), chol_mat.rows (),
                             utmp.fortran_vec (), w, info));

  return info;
}

template <>
octave_idx_type
chol<Matrix>::insert_sym (const ColumnVector& u, octave_idx_type j)
{
  octave_idx_type info = -1;

  octave_idx_type n = chol_mat.rows ();

  if (u.numel () != n + 1)
    (*current_liboctave_error_handler) ("cholinsert: dimension mismatch");
  if (j < 0 || j > n)
    (*current_liboctave_error_handler) ("cholinsert: index out of range");

  ColumnVector utmp = u;

  OCTAVE_LOCAL_BUFFER (double, w, n);

  chol_mat.resize (n+1, n+1);

  F77_XFCN (dchinx, DCHINX, (n, chol_mat.fortran_vec (), chol_mat.rows (),
                             j + 1, utmp.fortran_vec (), w, info));

  return info;
}

template <>
void
chol<Matrix>::delete_sym (octave_idx_type j)
{
  octave_idx_type n = chol_mat.rows ();

  if (j < 0 || j > n-1)
    (*current_liboctave_error_handler) ("choldelete: index out of range");

  OCTAVE_LOCAL_BUFFER (double, w, n);

  F77_XFCN (dchdex, DCHDEX, (n, chol_mat.fortran_vec (), chol_mat.rows (),
                             j + 1, w));

  chol_mat.resize (n-1, n-1);
}

template <>
void
chol<Matrix>::shift_sym (octave_idx_type i, octave_idx_type j)
{
  octave_idx_type n = chol_mat.rows ();

  if (i < 0 || i > n-1 || j < 0 || j > n-1)
    (*current_liboctave_error_handler) ("cholshift: index out of range");

  OCTAVE_LOCAL_BUFFER (double, w, 2*n);

  F77_XFCN (dchshx, DCHSHX, (n, chol_mat.fortran_vec (), chol_mat.rows (),
                             i + 1, j + 1, w));
}

#endif

template <>
octave_idx_type
chol<FloatMatrix>::init (const FloatMatrix& a, bool upper, bool calc_cond)
{
  octave_idx_type a_nr = a.rows ();
  octave_idx_type a_nc = a.cols ();

  if (a_nr != a_nc)
    (*current_liboctave_error_handler) ("chol: requires square matrix");

  octave_idx_type n = a_nc;
  octave_idx_type info;

  is_upper = upper;

  chol_mat.clear (n, n);
  if (is_upper)
    for (octave_idx_type j = 0; j < n; j++)
      {
        for (octave_idx_type i = 0; i <= j; i++)
          chol_mat.xelem (i, j) = a(i, j);
        for (octave_idx_type i = j+1; i < n; i++)
          chol_mat.xelem (i, j) = 0.0f;
      }
  else
    for (octave_idx_type j = 0; j < n; j++)
      {
        for (octave_idx_type i = 0; i < j; i++)
          chol_mat.xelem (i, j) = 0.0f;
        for (octave_idx_type i = j; i < n; i++)
          chol_mat.xelem (i, j) = a(i, j);
      }
  float *h = chol_mat.fortran_vec ();

  // Calculate the norm of the matrix, for later use.
  float anorm = 0;
  if (calc_cond)
    anorm = xnorm (a, 1);

  if (is_upper)
    F77_XFCN (spotrf, SPOTRF, (F77_CONST_CHAR_ARG2 ("U", 1), n, h, n, info
                               F77_CHAR_ARG_LEN (1)));
  else
    F77_XFCN (spotrf, SPOTRF, (F77_CONST_CHAR_ARG2 ("L", 1), n, h, n, info
                               F77_CHAR_ARG_LEN (1)));

  xrcond = 0.0;
  if (info > 0)
    chol_mat.resize (info - 1, info - 1);
  else if (calc_cond)
    {
      octave_idx_type spocon_info = 0;

      // Now calculate the condition number for non-singular matrix.
      Array<float> z (dim_vector (3*n, 1));
      float *pz = z.fortran_vec ();
      Array<octave_idx_type> iz (dim_vector (n, 1));
      octave_idx_type *piz = iz.fortran_vec ();
      if (is_upper)
        F77_XFCN (spocon, SPOCON, (F77_CONST_CHAR_ARG2 ("U", 1), n, h,
                                   n, anorm, xrcond, pz, piz, spocon_info
                                   F77_CHAR_ARG_LEN (1)));
      else
        F77_XFCN (spocon, SPOCON, (F77_CONST_CHAR_ARG2 ("L", 1), n, h,
                                   n, anorm, xrcond, pz, piz, spocon_info
                                   F77_CHAR_ARG_LEN (1)));

      if (spocon_info != 0)
        info = -1;
    }

  return info;
}

#ifdef HAVE_QRUPDATE

template <>
void
chol<FloatMatrix>::update (const FloatColumnVector& u)
{
  octave_idx_type n = chol_mat.rows ();

  if (u.numel () != n)
    (*current_liboctave_error_handler) ("cholupdate: dimension mismatch");

  FloatColumnVector utmp = u;

  OCTAVE_LOCAL_BUFFER (float, w, n);

  F77_XFCN (sch1up, SCH1UP, (n, chol_mat.fortran_vec (), chol_mat.rows (),
                             utmp.fortran_vec (), w));
}

template <>
octave_idx_type
chol<FloatMatrix>::downdate (const FloatColumnVector& u)
{
  octave_idx_type info = -1;

  octave_idx_type n = chol_mat.rows ();

  if (u.numel () != n)
    (*current_liboctave_error_handler) ("cholupdate: dimension mismatch");

  FloatColumnVector utmp = u;

  OCTAVE_LOCAL_BUFFER (float, w, n);

  F77_XFCN (sch1dn, SCH1DN, (n, chol_mat.fortran_vec (), chol_mat.rows (),
                             utmp.fortran_vec (), w, info));

  return info;
}

template <>
octave_idx_type
chol<FloatMatrix>::insert_sym (const FloatColumnVector& u, octave_idx_type j)
{
  octave_idx_type info = -1;

  octave_idx_type n = chol_mat.rows ();

  if (u.numel () != n + 1)
    (*current_liboctave_error_handler) ("cholinsert: dimension mismatch");
  if (j < 0 || j > n)
    (*current_liboctave_error_handler) ("cholinsert: index out of range");

  FloatColumnVector utmp = u;

  OCTAVE_LOCAL_BUFFER (float, w, n);

  chol_mat.resize (n+1, n+1);

  F77_XFCN (schinx, SCHINX, (n, chol_mat.fortran_vec (), chol_mat.rows (),
                             j + 1, utmp.fortran_vec (), w, info));

  return info;
}

template <>
void
chol<FloatMatrix>::delete_sym (octave_idx_type j)
{
  octave_idx_type n = chol_mat.rows ();

  if (j < 0 || j > n-1)
    (*current_liboctave_error_handler) ("choldelete: index out of range");

  OCTAVE_LOCAL_BUFFER (float, w, n);

  F77_XFCN (schdex, SCHDEX, (n, chol_mat.fortran_vec (), chol_mat.rows (),
                             j + 1, w));

  chol_mat.resize (n-1, n-1);
}

template <>
void
chol<FloatMatrix>::shift_sym (octave_idx_type i, octave_idx_type j)
{
  octave_idx_type n = chol_mat.rows ();

  if (i < 0 || i > n-1 || j < 0 || j > n-1)
    (*current_liboctave_error_handler) ("cholshift: index out of range");

  OCTAVE_LOCAL_BUFFER (float, w, 2*n);

  F77_XFCN (schshx, SCHSHX, (n, chol_mat.fortran_vec (), chol_mat.rows (),
                             i + 1, j + 1, w));
}

#endif

template <>
octave_idx_type
chol<ComplexMatrix>::init (const ComplexMatrix& a, bool upper, bool calc_cond)
{
  octave_idx_type a_nr = a.rows ();
  octave_idx_type a_nc = a.cols ();

  if (a_nr != a_nc)
    (*current_liboctave_error_handler) ("chol: requires square matrix");

  octave_idx_type n = a_nc;
  octave_idx_type info;

  is_upper = upper;

  chol_mat.clear (n, n);
  if (is_upper)
    for (octave_idx_type j = 0; j < n; j++)
      {
        for (octave_idx_type i = 0; i <= j; i++)
          chol_mat.xelem (i, j) = a(i, j);
        for (octave_idx_type i = j+1; i < n; i++)
          chol_mat.xelem (i, j) = 0.0;
      }
  else
    for (octave_idx_type j = 0; j < n; j++)
      {
        for (octave_idx_type i = 0; i < j; i++)
          chol_mat.xelem (i, j) = 0.0;
        for (octave_idx_type i = j; i < n; i++)
          chol_mat.xelem (i, j) = a(i, j);
      }
  Complex *h = chol_mat.fortran_vec ();

  // Calculate the norm of the matrix, for later use.
  double anorm = 0;
  if (calc_cond)
    anorm = xnorm (a, 1);

  if (is_upper)
    F77_XFCN (zpotrf, ZPOTRF, (F77_CONST_CHAR_ARG2 ("U", 1), n, h, n, info
                               F77_CHAR_ARG_LEN (1)));
  else
    F77_XFCN (zpotrf, ZPOTRF, (F77_CONST_CHAR_ARG2 ("L", 1), n, h, n, info
                               F77_CHAR_ARG_LEN (1)));

  xrcond = 0.0;
  if (info > 0)
    chol_mat.resize (info - 1, info - 1);
  else if (calc_cond)
    {
      octave_idx_type zpocon_info = 0;

      // Now calculate the condition number for non-singular matrix.
      Array<Complex> z (dim_vector (2*n, 1));
      Complex *pz = z.fortran_vec ();
      Array<double> rz (dim_vector (n, 1));
      double *prz = rz.fortran_vec ();
      F77_XFCN (zpocon, ZPOCON, (F77_CONST_CHAR_ARG2 ("U", 1), n, h,
                                 n, anorm, xrcond, pz, prz, zpocon_info
                                 F77_CHAR_ARG_LEN (1)));

      if (zpocon_info != 0)
        info = -1;
    }

  return info;
}

#ifdef HAVE_QRUPDATE

template <>
void
chol<ComplexMatrix>::update (const ComplexColumnVector& u)
{
  octave_idx_type n = chol_mat.rows ();

  if (u.numel () != n)
    (*current_liboctave_error_handler) ("cholupdate: dimension mismatch");

  ComplexColumnVector utmp = u;

  OCTAVE_LOCAL_BUFFER (double, rw, n);

  F77_XFCN (zch1up, ZCH1UP, (n, chol_mat.fortran_vec (), chol_mat.rows (),
                             utmp.fortran_vec (), rw));
}

template <>
octave_idx_type
chol<ComplexMatrix>::downdate (const ComplexColumnVector& u)
{
  octave_idx_type info = -1;

  octave_idx_type n = chol_mat.rows ();

  if (u.numel () != n)
    (*current_liboctave_error_handler) ("cholupdate: dimension mismatch");

  ComplexColumnVector utmp = u;

  OCTAVE_LOCAL_BUFFER (double, rw, n);

  F77_XFCN (zch1dn, ZCH1DN, (n, chol_mat.fortran_vec (), chol_mat.rows (),
                             utmp.fortran_vec (), rw, info));

  return info;
}

template <>
octave_idx_type
chol<ComplexMatrix>::insert_sym (const ComplexColumnVector& u,
                                 octave_idx_type j)
{
  octave_idx_type info = -1;

  octave_idx_type n = chol_mat.rows ();

  if (u.numel () != n + 1)
    (*current_liboctave_error_handler) ("cholinsert: dimension mismatch");
  if (j < 0 || j > n)
    (*current_liboctave_error_handler) ("cholinsert: index out of range");

  ComplexColumnVector utmp = u;

  OCTAVE_LOCAL_BUFFER (double, rw, n);

  chol_mat.resize (n+1, n+1);

  F77_XFCN (zchinx, ZCHINX, (n, chol_mat.fortran_vec (), chol_mat.rows (),
                             j + 1, utmp.fortran_vec (), rw, info));

  return info;
}

template <>
void
chol<ComplexMatrix>::delete_sym (octave_idx_type j)
{
  octave_idx_type n = chol_mat.rows ();

  if (j < 0 || j > n-1)
    (*current_liboctave_error_handler) ("choldelete: index out of range");

  OCTAVE_LOCAL_BUFFER (double, rw, n);

  F77_XFCN (zchdex, ZCHDEX, (n, chol_mat.fortran_vec (), chol_mat.rows (),
                             j + 1, rw));

  chol_mat.resize (n-1, n-1);
}

template <>
void
chol<ComplexMatrix>::shift_sym (octave_idx_type i, octave_idx_type j)
{
  octave_idx_type n = chol_mat.rows ();

  if (i < 0 || i > n-1 || j < 0 || j > n-1)
    (*current_liboctave_error_handler) ("cholshift: index out of range");

  OCTAVE_LOCAL_BUFFER (Complex, w, n);
  OCTAVE_LOCAL_BUFFER (double, rw, n);

  F77_XFCN (zchshx, ZCHSHX, (n, chol_mat.fortran_vec (), chol_mat.rows (),
                             i + 1, j + 1, w, rw));
}

#endif

template <>
octave_idx_type
chol<FloatComplexMatrix>::init (const FloatComplexMatrix& a, bool upper,
                                bool calc_cond)
{
  octave_idx_type a_nr = a.rows ();
  octave_idx_type a_nc = a.cols ();

  if (a_nr != a_nc)
    (*current_liboctave_error_handler) ("chol: requires square matrix");

  octave_idx_type n = a_nc;
  octave_idx_type info;

  is_upper = upper;

  chol_mat.clear (n, n);
  if (is_upper)
    for (octave_idx_type j = 0; j < n; j++)
      {
        for (octave_idx_type i = 0; i <= j; i++)
          chol_mat.xelem (i, j) = a(i, j);
        for (octave_idx_type i = j+1; i < n; i++)
          chol_mat.xelem (i, j) = 0.0f;
      }
  else
    for (octave_idx_type j = 0; j < n; j++)
      {
        for (octave_idx_type i = 0; i < j; i++)
          chol_mat.xelem (i, j) = 0.0f;
        for (octave_idx_type i = j; i < n; i++)
          chol_mat.xelem (i, j) = a(i, j);
      }
  FloatComplex *h = chol_mat.fortran_vec ();

  // Calculate the norm of the matrix, for later use.
  float anorm = 0;
  if (calc_cond)
    anorm = xnorm (a, 1);

  if (is_upper)
    F77_XFCN (cpotrf, CPOTRF, (F77_CONST_CHAR_ARG2 ("U", 1), n, h, n, info
                               F77_CHAR_ARG_LEN (1)));
  else
    F77_XFCN (cpotrf, CPOTRF, (F77_CONST_CHAR_ARG2 ("L", 1), n, h, n, info
                               F77_CHAR_ARG_LEN (1)));

  xrcond = 0.0;
  if (info > 0)
    chol_mat.resize (info - 1, info - 1);
  else if (calc_cond)
    {
      octave_idx_type cpocon_info = 0;

      // Now calculate the condition number for non-singular matrix.
      Array<FloatComplex> z (dim_vector (2*n, 1));
      FloatComplex *pz = z.fortran_vec ();
      Array<float> rz (dim_vector (n, 1));
      float *prz = rz.fortran_vec ();
      F77_XFCN (cpocon, CPOCON, (F77_CONST_CHAR_ARG2 ("U", 1), n, h,
                                 n, anorm, xrcond, pz, prz, cpocon_info
                                 F77_CHAR_ARG_LEN (1)));

      if (cpocon_info != 0)
        info = -1;
    }

  return info;
}

#ifdef HAVE_QRUPDATE

template <>
void
chol<FloatComplexMatrix>::update (const FloatComplexColumnVector& u)
{
  octave_idx_type n = chol_mat.rows ();

  if (u.numel () != n)
    (*current_liboctave_error_handler) ("cholupdate: dimension mismatch");

  FloatComplexColumnVector utmp = u;

  OCTAVE_LOCAL_BUFFER (float, rw, n);

  F77_XFCN (cch1up, CCH1UP, (n, chol_mat.fortran_vec (), chol_mat.rows (),
                             utmp.fortran_vec (), rw));
}

template <>
octave_idx_type
chol<FloatComplexMatrix>::downdate (const FloatComplexColumnVector& u)
{
  octave_idx_type info = -1;

  octave_idx_type n = chol_mat.rows ();

  if (u.numel () != n)
    (*current_liboctave_error_handler) ("cholupdate: dimension mismatch");

  FloatComplexColumnVector utmp = u;

  OCTAVE_LOCAL_BUFFER (float, rw, n);

  F77_XFCN (cch1dn, CCH1DN, (n, chol_mat.fortran_vec (), chol_mat.rows (),
                             utmp.fortran_vec (), rw, info));

  return info;
}

template <>
octave_idx_type
chol<FloatComplexMatrix>::insert_sym (const FloatComplexColumnVector& u,
                                      octave_idx_type j)
{
  octave_idx_type info = -1;

  octave_idx_type n = chol_mat.rows ();

  if (u.numel () != n + 1)
    (*current_liboctave_error_handler) ("cholinsert: dimension mismatch");
  if (j < 0 || j > n)
    (*current_liboctave_error_handler) ("cholinsert: index out of range");

  FloatComplexColumnVector utmp = u;

  OCTAVE_LOCAL_BUFFER (float, rw, n);

  chol_mat.resize (n+1, n+1);

  F77_XFCN (cchinx, CCHINX, (n, chol_mat.fortran_vec (), chol_mat.rows (),
                             j + 1, utmp.fortran_vec (), rw, info));

  return info;
}

template <>
void
chol<FloatComplexMatrix>::delete_sym (octave_idx_type j)
{
  octave_idx_type n = chol_mat.rows ();

  if (j < 0 || j > n-1)
    (*current_liboctave_error_handler) ("choldelete: index out of range");

  OCTAVE_LOCAL_BUFFER (float, rw, n);

  F77_XFCN (cchdex, CCHDEX, (n, chol_mat.fortran_vec (), chol_mat.rows (),
                             j + 1, rw));

  chol_mat.resize (n-1, n-1);
}

template <>
void
chol<FloatComplexMatrix>::shift_sym (octave_idx_type i, octave_idx_type j)
{
  octave_idx_type n = chol_mat.rows ();

  if (i < 0 || i > n-1 || j < 0 || j > n-1)
    (*current_liboctave_error_handler) ("cholshift: index out of range");

  OCTAVE_LOCAL_BUFFER (FloatComplex, w, n);
  OCTAVE_LOCAL_BUFFER (float, rw, n);

  F77_XFCN (cchshx, CCHSHX, (n, chol_mat.fortran_vec (), chol_mat.rows (),
                             i + 1, j + 1, w, rw));
}

#endif

// Instantiations we need.

template class chol<Matrix>;

template class chol<FloatMatrix>;

template class chol<ComplexMatrix>;

template class chol<FloatComplexMatrix>;

template Matrix
chol2inv<Matrix> (const Matrix& r);

template ComplexMatrix
chol2inv<ComplexMatrix> (const ComplexMatrix& r);

template FloatMatrix
chol2inv<FloatMatrix> (const FloatMatrix& r);

template FloatComplexMatrix
chol2inv<FloatComplexMatrix> (const FloatComplexMatrix& r);