view liboctave/numeric/lu.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 7e67c7f82fc1
children aba2e6293dd8
line wrap: on
line source

/*

Copyright (C) 1996-2015 John W. Eaton
Copyright (C) 2009 VZLU Prague

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 "CColVector.h"
#include "CMatrix.h"
#include "dColVector.h"
#include "dMatrix.h"
#include "f77-fcn.h"
#include "fCColVector.h"
#include "fCMatrix.h"
#include "fColVector.h"
#include "fMatrix.h"
#include "lo-error.h"
#include "lu.h"
#include "oct-locbuf.h"

extern "C"
{
  F77_RET_T
  F77_FUNC (dgetrf, DGETRF) (const octave_idx_type&, const octave_idx_type&,
                             double*, const octave_idx_type&,
                             octave_idx_type*, octave_idx_type&);

#ifdef HAVE_QRUPDATE_LUU
  F77_RET_T
  F77_FUNC (dlu1up, DLU1UP) (const octave_idx_type&, const octave_idx_type&,
                             double *, const octave_idx_type&,
                             double *, const octave_idx_type&,
                             double *, double *);

  F77_RET_T
  F77_FUNC (dlup1up, DLUP1UP) (const octave_idx_type&, const octave_idx_type&,
                               double *, const octave_idx_type&,
                               double *, const octave_idx_type&,
                               octave_idx_type *, const double *,
                               const double *, double *);
#endif

  F77_RET_T
  F77_FUNC (sgetrf, SGETRF) (const octave_idx_type&, const octave_idx_type&,
                             float*, const octave_idx_type&, octave_idx_type*,
                             octave_idx_type&);

#ifdef HAVE_QRUPDATE_LUU
  F77_RET_T
  F77_FUNC (slu1up, SLU1UP) (const octave_idx_type&, const octave_idx_type&,
                             float *, const octave_idx_type&,
                             float *, const octave_idx_type&,
                             float *, float *);

  F77_RET_T
  F77_FUNC (slup1up, SLUP1UP) (const octave_idx_type&, const octave_idx_type&,
                               float *, const octave_idx_type&,
                               float *, const octave_idx_type&,
                               octave_idx_type *, const float *,
                               const float *, float *);
#endif

  F77_RET_T
  F77_FUNC (zgetrf, ZGETRF) (const octave_idx_type&, const octave_idx_type&,
                             Complex*, const octave_idx_type&,
                             octave_idx_type*, octave_idx_type&);

#ifdef HAVE_QRUPDATE_LUU
  F77_RET_T
  F77_FUNC (zlu1up, ZLU1UP) (const octave_idx_type&, const octave_idx_type&,
                             Complex *, const octave_idx_type&,
                             Complex *, const octave_idx_type&,
                             Complex *, Complex *);

  F77_RET_T
  F77_FUNC (zlup1up, ZLUP1UP) (const octave_idx_type&, const octave_idx_type&,
                               Complex *, const octave_idx_type&,
                               Complex *, const octave_idx_type&,
                               octave_idx_type *, const Complex *,
                               const Complex *, Complex *);
#endif

  F77_RET_T
  F77_FUNC (cgetrf, CGETRF) (const octave_idx_type&, const octave_idx_type&,
                             FloatComplex*, const octave_idx_type&,
                             octave_idx_type*, octave_idx_type&);

#ifdef HAVE_QRUPDATE_LUU
  F77_RET_T
  F77_FUNC (clu1up, CLU1UP) (const octave_idx_type&, const octave_idx_type&,
                             FloatComplex *, const octave_idx_type&,
                             FloatComplex *, const octave_idx_type&,
                             FloatComplex *, FloatComplex *);

  F77_RET_T
  F77_FUNC (clup1up, CLUP1UP) (const octave_idx_type&, const octave_idx_type&,
                               FloatComplex *, const octave_idx_type&,
                               FloatComplex *, const octave_idx_type&,
                               octave_idx_type *, const FloatComplex *,
                               const FloatComplex *, FloatComplex *);
#endif
}

template <typename T>
lu<T>::lu (const T& l, const T& u,
                           const PermMatrix& p)
  : a_fact (u), l_fact (l), ipvt (p.transpose ().col_perm_vec ())
{
  if (l.columns () != u.rows ())
    (*current_liboctave_error_handler) ("lu: dimension mismatch");
}

template <typename T>
bool
lu<T>::packed (void) const
{
  return l_fact.dims () == dim_vector ();
}

template <typename T>
void
lu<T>::unpack (void)
{
  if (packed ())
    {
      l_fact = L ();
      a_fact = U (); // FIXME: sub-optimal
      ipvt = getp ();
    }
}

template <typename T>
T
lu<T>::L (void) const
{
  if (packed ())
    {
      octave_idx_type a_nr = a_fact.rows ();
      octave_idx_type a_nc = a_fact.cols ();
      octave_idx_type mn = (a_nr < a_nc ? a_nr : a_nc);

      T l (a_nr, mn, ELT_T (0.0));

      for (octave_idx_type i = 0; i < a_nr; i++)
        {
          if (i < a_nc)
            l.xelem (i, i) = 1.0;

          for (octave_idx_type j = 0; j < (i < a_nc ? i : a_nc); j++)
            l.xelem (i, j) = a_fact.xelem (i, j);
        }

      return l;
    }
  else
    return l_fact;
}

template <typename T>
T
lu<T>::U (void) const
{
  if (packed ())
    {
      octave_idx_type a_nr = a_fact.rows ();
      octave_idx_type a_nc = a_fact.cols ();
      octave_idx_type mn = (a_nr < a_nc ? a_nr : a_nc);

      T u (mn, a_nc, ELT_T (0.0));

      for (octave_idx_type i = 0; i < mn; i++)
        {
          for (octave_idx_type j = i; j < a_nc; j++)
            u.xelem (i, j) = a_fact.xelem (i, j);
        }

      return u;
    }
  else
    return a_fact;
}

template <typename T>
T
lu<T>::Y (void) const
{
  if (! packed ())
    (*current_liboctave_error_handler)
      ("lu: Y () not implemented for unpacked form");

  return a_fact;
}

template <typename T>
Array<octave_idx_type>
lu<T>::getp (void) const
{
  if (packed ())
    {
      octave_idx_type a_nr = a_fact.rows ();

      Array<octave_idx_type> pvt (dim_vector (a_nr, 1));

      for (octave_idx_type i = 0; i < a_nr; i++)
        pvt.xelem (i) = i;

      for (octave_idx_type i = 0; i < ipvt.numel (); i++)
        {
          octave_idx_type k = ipvt.xelem (i);

          if (k != i)
            {
              octave_idx_type tmp = pvt.xelem (k);
              pvt.xelem (k) = pvt.xelem (i);
              pvt.xelem (i) = tmp;
            }
        }

      return pvt;
    }
  else
    return ipvt;
}

template <typename T>
PermMatrix
lu<T>::P (void) const
{
  return PermMatrix (getp (), false);
}

template <typename T>
ColumnVector
lu<T>::P_vec (void) const
{
  octave_idx_type a_nr = a_fact.rows ();

  ColumnVector p (a_nr);

  Array<octave_idx_type> pvt = getp ();

  for (octave_idx_type i = 0; i < a_nr; i++)
    p.xelem (i) = static_cast<double> (pvt.xelem (i) + 1);

  return p;
}

template <typename T>
bool
lu<T>::regular (void) const
{
  bool retval = true;

  octave_idx_type k = std::min (a_fact.rows (), a_fact.columns ());

  for (octave_idx_type i = 0; i < k; i++)
    {
      if (a_fact(i, i) == ELT_T ())
        {
          retval = false;
          break;
        }
    }

  return retval;
}

#if ! defined (HAVE_QRUPDATE_LUU)

template <typename T>
void
lu<T>::update (const VT&, const VT&)
{
  (*current_liboctave_error_handler)
    ("luupdate: support for qrupdate with LU updates "
     "was unavailable or disabled when liboctave was built");
}

template <typename T>
void
lu<T>::update (const T&, const T&)
{
  (*current_liboctave_error_handler)
    ("luupdate: support for qrupdate with LU updates "
     "was unavailable or disabled when liboctave was built");
}

template <typename T>
void
lu<T>::update_piv (const VT&, const VT&)
{
  (*current_liboctave_error_handler)
    ("luupdate: support for qrupdate with LU updates "
     "was unavailable or disabled when liboctave was built");
}

template <typename T>
void
lu<T>::update_piv (const T&, const T&)
{
  (*current_liboctave_error_handler)
    ("luupdate: support for qrupdate with LU updates "
     "was unavailable or disabled when liboctave was built");
}

#endif

// Specializations.

template <>
lu<Matrix>::lu (const Matrix& a)
{
  octave_idx_type a_nr = a.rows ();
  octave_idx_type a_nc = a.cols ();
  octave_idx_type mn = (a_nr < a_nc ? a_nr : a_nc);

  ipvt.resize (dim_vector (mn, 1));
  octave_idx_type *pipvt = ipvt.fortran_vec ();

  a_fact = a;
  double *tmp_data = a_fact.fortran_vec ();

  octave_idx_type info = 0;

  F77_XFCN (dgetrf, DGETRF, (a_nr, a_nc, tmp_data, a_nr, pipvt, info));

  for (octave_idx_type i = 0; i < mn; i++)
    pipvt[i] -= 1;
}

#ifdef HAVE_QRUPDATE_LUU

template <>
void
lu<Matrix>::update (const ColumnVector& u, const ColumnVector& v)
{
  if (packed ())
    unpack ();

  Matrix& l = l_fact;
  Matrix& r = a_fact;

  octave_idx_type m = l.rows ();
  octave_idx_type n = r.columns ();
  octave_idx_type k = l.columns ();

  if (u.numel () != m || v.numel () != n)
    (*current_liboctave_error_handler) ("luupdate: dimensions mismatch");

  ColumnVector utmp = u;
  ColumnVector vtmp = v;
  F77_XFCN (dlu1up, DLU1UP, (m, n, l.fortran_vec (), m, r.fortran_vec (), k,
                             utmp.fortran_vec (), vtmp.fortran_vec ()));
}

template <>
void
lu<Matrix>::update (const Matrix& u, const Matrix& v)
{
  if (packed ())
    unpack ();

  Matrix& l = l_fact;
  Matrix& r = a_fact;

  octave_idx_type m = l.rows ();
  octave_idx_type n = r.columns ();
  octave_idx_type k = l.columns ();

  if (u.rows () != m || v.rows () != n || u.cols () != v.cols ())
    (*current_liboctave_error_handler) ("luupdate: dimensions mismatch");

  for (volatile octave_idx_type i = 0; i < u.cols (); i++)
    {
      ColumnVector utmp = u.column (i);
      ColumnVector vtmp = v.column (i);
      F77_XFCN (dlu1up, DLU1UP, (m, n, l.fortran_vec (),
                                 m, r.fortran_vec (), k,
                                 utmp.fortran_vec (), vtmp.fortran_vec ()));
    }
}

template <>
void
lu<Matrix>::update_piv (const ColumnVector& u, const ColumnVector& v)
{
  if (packed ())
    unpack ();

  Matrix& l = l_fact;
  Matrix& r = a_fact;

  octave_idx_type m = l.rows ();
  octave_idx_type n = r.columns ();
  octave_idx_type k = l.columns ();

  if (u.numel () != m || v.numel () != n)
    (*current_liboctave_error_handler) ("luupdate: dimensions mismatch");

  ColumnVector utmp = u;
  ColumnVector vtmp = v;
  OCTAVE_LOCAL_BUFFER (double, w, m);
  for (octave_idx_type i = 0; i < m; i++) ipvt(i) += 1; // increment
  F77_XFCN (dlup1up, DLUP1UP, (m, n, l.fortran_vec (),
                               m, r.fortran_vec (), k,
                               ipvt.fortran_vec (),
                               utmp.data (), vtmp.data (), w));
  for (octave_idx_type i = 0; i < m; i++) ipvt(i) -= 1; // decrement
}

template <>
void
lu<Matrix>::update_piv (const Matrix& u, const Matrix& v)
{
  if (packed ())
    unpack ();

  Matrix& l = l_fact;
  Matrix& r = a_fact;

  octave_idx_type m = l.rows ();
  octave_idx_type n = r.columns ();
  octave_idx_type k = l.columns ();

  if (u.rows () != m || v.rows () != n || u.cols () != v.cols ())
    (*current_liboctave_error_handler) ("luupdate: dimensions mismatch");

  OCTAVE_LOCAL_BUFFER (double, w, m);
  for (octave_idx_type i = 0; i < m; i++) ipvt(i) += 1; // increment
  for (volatile octave_idx_type i = 0; i < u.cols (); i++)
    {
      ColumnVector utmp = u.column (i);
      ColumnVector vtmp = v.column (i);
      F77_XFCN (dlup1up, DLUP1UP, (m, n, l.fortran_vec (),
                                   m, r.fortran_vec (), k,
                                   ipvt.fortran_vec (),
                                   utmp.data (), vtmp.data (), w));
    }
  for (octave_idx_type i = 0; i < m; i++) ipvt(i) -= 1; // decrement
}

#endif

template <>
lu<FloatMatrix>::lu (const FloatMatrix& a)
{
  octave_idx_type a_nr = a.rows ();
  octave_idx_type a_nc = a.cols ();
  octave_idx_type mn = (a_nr < a_nc ? a_nr : a_nc);

  ipvt.resize (dim_vector (mn, 1));
  octave_idx_type *pipvt = ipvt.fortran_vec ();

  a_fact = a;
  float *tmp_data = a_fact.fortran_vec ();

  octave_idx_type info = 0;

  F77_XFCN (sgetrf, SGETRF, (a_nr, a_nc, tmp_data, a_nr, pipvt, info));

  for (octave_idx_type i = 0; i < mn; i++)
    pipvt[i] -= 1;
}

#ifdef HAVE_QRUPDATE_LUU

template <>
void
lu<FloatMatrix>::update (const FloatColumnVector& u, const FloatColumnVector& v)
{
  if (packed ())
    unpack ();

  FloatMatrix& l = l_fact;
  FloatMatrix& r = a_fact;

  octave_idx_type m = l.rows ();
  octave_idx_type n = r.columns ();
  octave_idx_type k = l.columns ();

  if (u.numel () != m || v.numel () != n)
    (*current_liboctave_error_handler) ("luupdate: dimensions mismatch");

  FloatColumnVector utmp = u;
  FloatColumnVector vtmp = v;
  F77_XFCN (slu1up, SLU1UP, (m, n, l.fortran_vec (),
                             m, r.fortran_vec (), k,
                             utmp.fortran_vec (), vtmp.fortran_vec ()));
}

template <>
void
lu<FloatMatrix>::update (const FloatMatrix& u, const FloatMatrix& v)
{
  if (packed ())
    unpack ();

  FloatMatrix& l = l_fact;
  FloatMatrix& r = a_fact;

  octave_idx_type m = l.rows ();
  octave_idx_type n = r.columns ();
  octave_idx_type k = l.columns ();

  if (u.rows () != m || v.rows () != n || u.cols () != v.cols ())
    (*current_liboctave_error_handler) ("luupdate: dimensions mismatch");

  for (volatile octave_idx_type i = 0; i < u.cols (); i++)
    {
      FloatColumnVector utmp = u.column (i);
      FloatColumnVector vtmp = v.column (i);
      F77_XFCN (slu1up, SLU1UP, (m, n, l.fortran_vec (),
                                 m, r.fortran_vec (), k,
                                 utmp.fortran_vec (), vtmp.fortran_vec ()));
    }
}

template <>
void
lu<FloatMatrix>::update_piv (const FloatColumnVector& u,
                          const FloatColumnVector& v)
{
  if (packed ())
    unpack ();

  FloatMatrix& l = l_fact;
  FloatMatrix& r = a_fact;

  octave_idx_type m = l.rows ();
  octave_idx_type n = r.columns ();
  octave_idx_type k = l.columns ();

  if (u.numel () != m || v.numel () != n)
    (*current_liboctave_error_handler) ("luupdate: dimensions mismatch");

  FloatColumnVector utmp = u;
  FloatColumnVector vtmp = v;
  OCTAVE_LOCAL_BUFFER (float, w, m);
  for (octave_idx_type i = 0; i < m; i++) ipvt(i) += 1; // increment
  F77_XFCN (slup1up, SLUP1UP, (m, n, l.fortran_vec (),
                               m, r.fortran_vec (), k,
                               ipvt.fortran_vec (),
                               utmp.data (), vtmp.data (), w));
  for (octave_idx_type i = 0; i < m; i++) ipvt(i) -= 1; // decrement
}

template <>
void
lu<FloatMatrix>::update_piv (const FloatMatrix& u, const FloatMatrix& v)
{
  if (packed ())
    unpack ();

  FloatMatrix& l = l_fact;
  FloatMatrix& r = a_fact;

  octave_idx_type m = l.rows ();
  octave_idx_type n = r.columns ();
  octave_idx_type k = l.columns ();

  if (u.rows () != m || v.rows () != n || u.cols () != v.cols ())
    (*current_liboctave_error_handler) ("luupdate: dimensions mismatch");

  OCTAVE_LOCAL_BUFFER (float, w, m);
  for (octave_idx_type i = 0; i < m; i++) ipvt(i) += 1; // increment
  for (volatile octave_idx_type i = 0; i < u.cols (); i++)
    {
      FloatColumnVector utmp = u.column (i);
      FloatColumnVector vtmp = v.column (i);
      F77_XFCN (slup1up, SLUP1UP, (m, n, l.fortran_vec (),
                                   m, r.fortran_vec (), k,
                                   ipvt.fortran_vec (),
                                   utmp.data (), vtmp.data (), w));
    }
  for (octave_idx_type i = 0; i < m; i++) ipvt(i) -= 1; // decrement
}

#endif

template <>
lu<ComplexMatrix>::lu (const ComplexMatrix& a)
{
  octave_idx_type a_nr = a.rows ();
  octave_idx_type a_nc = a.cols ();
  octave_idx_type mn = (a_nr < a_nc ? a_nr : a_nc);

  ipvt.resize (dim_vector (mn, 1));
  octave_idx_type *pipvt = ipvt.fortran_vec ();

  a_fact = a;
  Complex *tmp_data = a_fact.fortran_vec ();

  octave_idx_type info = 0;

  F77_XFCN (zgetrf, ZGETRF, (a_nr, a_nc, tmp_data, a_nr, pipvt, info));

  for (octave_idx_type i = 0; i < mn; i++)
    pipvt[i] -= 1;
}

#ifdef HAVE_QRUPDATE_LUU

template <>
void
lu<ComplexMatrix>::update (const ComplexColumnVector& u,
                           const ComplexColumnVector& v)
{
  if (packed ())
    unpack ();

  ComplexMatrix& l = l_fact;
  ComplexMatrix& r = a_fact;

  octave_idx_type m = l.rows ();
  octave_idx_type n = r.columns ();
  octave_idx_type k = l.columns ();

  if (u.numel () != m || v.numel () != n)
    (*current_liboctave_error_handler) ("luupdate: dimensions mismatch");

  ComplexColumnVector utmp = u;
  ComplexColumnVector vtmp = v;
  F77_XFCN (zlu1up, ZLU1UP, (m, n, l.fortran_vec (), m, r.fortran_vec (), k,
                             utmp.fortran_vec (), vtmp.fortran_vec ()));
}

template <>
void
lu<ComplexMatrix>::update (const ComplexMatrix& u, const ComplexMatrix& v)
{
  if (packed ())
    unpack ();

  ComplexMatrix& l = l_fact;
  ComplexMatrix& r = a_fact;

  octave_idx_type m = l.rows ();
  octave_idx_type n = r.columns ();
  octave_idx_type k = l.columns ();

  if (u.rows () != m || v.rows () != n || u.cols () != v.cols ())
    (*current_liboctave_error_handler) ("luupdate: dimensions mismatch");

  for (volatile octave_idx_type i = 0; i < u.cols (); i++)
    {
      ComplexColumnVector utmp = u.column (i);
      ComplexColumnVector vtmp = v.column (i);
      F77_XFCN (zlu1up, ZLU1UP, (m, n, l.fortran_vec (),
                                 m, r.fortran_vec (), k,
                                 utmp.fortran_vec (), vtmp.fortran_vec ()));
    }
}

template <>
void
lu<ComplexMatrix>::update_piv (const ComplexColumnVector& u,
                               const ComplexColumnVector& v)
{
  if (packed ())
    unpack ();

  ComplexMatrix& l = l_fact;
  ComplexMatrix& r = a_fact;

  octave_idx_type m = l.rows ();
  octave_idx_type n = r.columns ();
  octave_idx_type k = l.columns ();

  if (u.numel () != m || v.numel () != n)
    (*current_liboctave_error_handler) ("luupdate: dimensions mismatch");

  ComplexColumnVector utmp = u;
  ComplexColumnVector vtmp = v;
  OCTAVE_LOCAL_BUFFER (Complex, w, m);
  for (octave_idx_type i = 0; i < m; i++) ipvt(i) += 1; // increment
  F77_XFCN (zlup1up, ZLUP1UP, (m, n, l.fortran_vec (),
                               m, r.fortran_vec (), k,
                               ipvt.fortran_vec (),
                               utmp.data (), vtmp.data (), w));
  for (octave_idx_type i = 0; i < m; i++) ipvt(i) -= 1; // decrement
}

template <>
void
lu<ComplexMatrix>::update_piv (const ComplexMatrix& u, const ComplexMatrix& v)
{
  if (packed ())
    unpack ();

  ComplexMatrix& l = l_fact;
  ComplexMatrix& r = a_fact;

  octave_idx_type m = l.rows ();
  octave_idx_type n = r.columns ();
  octave_idx_type k = l.columns ();

  if (u.rows () != m || v.rows () != n || u.cols () != v.cols ())
    (*current_liboctave_error_handler) ("luupdate: dimensions mismatch");

  OCTAVE_LOCAL_BUFFER (Complex, w, m);
  for (octave_idx_type i = 0; i < m; i++) ipvt(i) += 1; // increment
  for (volatile octave_idx_type i = 0; i < u.cols (); i++)
    {
      ComplexColumnVector utmp = u.column (i);
      ComplexColumnVector vtmp = v.column (i);
      F77_XFCN (zlup1up, ZLUP1UP, (m, n, l.fortran_vec (),
                                   m, r.fortran_vec (), k,
                                   ipvt.fortran_vec (),
                                   utmp.data (), vtmp.data (), w));
    }
  for (octave_idx_type i = 0; i < m; i++) ipvt(i) -= 1; // decrement
}

#endif

template <>
lu<FloatComplexMatrix>::lu (const FloatComplexMatrix& a)
{
  octave_idx_type a_nr = a.rows ();
  octave_idx_type a_nc = a.cols ();
  octave_idx_type mn = (a_nr < a_nc ? a_nr : a_nc);

  ipvt.resize (dim_vector (mn, 1));
  octave_idx_type *pipvt = ipvt.fortran_vec ();

  a_fact = a;
  FloatComplex *tmp_data = a_fact.fortran_vec ();

  octave_idx_type info = 0;

  F77_XFCN (cgetrf, CGETRF, (a_nr, a_nc, tmp_data, a_nr, pipvt, info));

  for (octave_idx_type i = 0; i < mn; i++)
    pipvt[i] -= 1;
}

#ifdef HAVE_QRUPDATE_LUU

template <>
void
lu<FloatComplexMatrix>::update (const FloatComplexColumnVector& u,
                             const FloatComplexColumnVector& v)
{
  if (packed ())
    unpack ();

  FloatComplexMatrix& l = l_fact;
  FloatComplexMatrix& r = a_fact;

  octave_idx_type m = l.rows ();
  octave_idx_type n = r.columns ();
  octave_idx_type k = l.columns ();

  if (u.numel () == m && v.numel () == n)
    {
      FloatComplexColumnVector utmp = u;
      FloatComplexColumnVector vtmp = v;
      F77_XFCN (clu1up, CLU1UP, (m, n, l.fortran_vec (), m, r.fortran_vec (), k,
                                 utmp.fortran_vec (), vtmp.fortran_vec ()));
    }
  else
    (*current_liboctave_error_handler) ("luupdate: dimensions mismatch");
}

template <>
void
lu<FloatComplexMatrix>::update (const FloatComplexMatrix& u,
                             const FloatComplexMatrix& v)
{
  if (packed ())
    unpack ();

  FloatComplexMatrix& l = l_fact;
  FloatComplexMatrix& r = a_fact;

  octave_idx_type m = l.rows ();
  octave_idx_type n = r.columns ();
  octave_idx_type k = l.columns ();

  if (u.rows () != m || v.rows () != n || u.cols () != v.cols ())
    (*current_liboctave_error_handler) ("luupdate: dimensions mismatch");

  for (volatile octave_idx_type i = 0; i < u.cols (); i++)
    {
      FloatComplexColumnVector utmp = u.column (i);
      FloatComplexColumnVector vtmp = v.column (i);
      F77_XFCN (clu1up, CLU1UP, (m, n, l.fortran_vec (),
                                 m, r.fortran_vec (), k,
                                 utmp.fortran_vec (), vtmp.fortran_vec ()));
    }
}

template <>
void
lu<FloatComplexMatrix>::update_piv (const FloatComplexColumnVector& u,
                                 const FloatComplexColumnVector& v)
{
  if (packed ())
    unpack ();

  FloatComplexMatrix& l = l_fact;
  FloatComplexMatrix& r = a_fact;

  octave_idx_type m = l.rows ();
  octave_idx_type n = r.columns ();
  octave_idx_type k = l.columns ();

  if (u.numel () != m || v.numel () != n)
    (*current_liboctave_error_handler) ("luupdate: dimensions mismatch");

  FloatComplexColumnVector utmp = u;
  FloatComplexColumnVector vtmp = v;
  OCTAVE_LOCAL_BUFFER (FloatComplex, w, m);
  for (octave_idx_type i = 0; i < m; i++) ipvt(i) += 1; // increment
  F77_XFCN (clup1up, CLUP1UP, (m, n, l.fortran_vec (),
                               m, r.fortran_vec (), k,
                               ipvt.fortran_vec (),
                               utmp.data (), vtmp.data (), w));
  for (octave_idx_type i = 0; i < m; i++) ipvt(i) -= 1; // decrement
}

template <>
void
lu<FloatComplexMatrix>::update_piv (const FloatComplexMatrix& u,
                                 const FloatComplexMatrix& v)
{
  if (packed ())
    unpack ();

  FloatComplexMatrix& l = l_fact;
  FloatComplexMatrix& r = a_fact;

  octave_idx_type m = l.rows ();
  octave_idx_type n = r.columns ();
  octave_idx_type k = l.columns ();

  if (u.rows () != m || v.rows () != n || u.cols () != v.cols ())
    (*current_liboctave_error_handler) ("luupdate: dimensions mismatch");

  OCTAVE_LOCAL_BUFFER (FloatComplex, w, m);
  for (octave_idx_type i = 0; i < m; i++) ipvt(i) += 1; // increment
  for (volatile octave_idx_type i = 0; i < u.cols (); i++)
    {
      FloatComplexColumnVector utmp = u.column (i);
      FloatComplexColumnVector vtmp = v.column (i);
      F77_XFCN (clup1up, CLUP1UP, (m, n, l.fortran_vec (),
                                   m, r.fortran_vec (), k,
                                   ipvt.fortran_vec (),
                                   utmp.data (), vtmp.data (), w));
    }
  for (octave_idx_type i = 0; i < m; i++) ipvt(i) -= 1; // decrement
}

#endif

// Instantiations we need.

template class lu<Matrix>;

template class lu<FloatMatrix>;

template class lu<ComplexMatrix>;

template class lu<FloatComplexMatrix>;