view liboctave/numeric/schur.cc @ 31607:aac27ad79be6 stable

maint: Re-indent code after switch to using namespace macros. * build-env.h, build-env.in.cc, Cell.h, __betainc__.cc, __eigs__.cc, __ftp__.cc, __ichol__.cc, __ilu__.cc, __isprimelarge__.cc, __magick_read__.cc, __pchip_deriv__.cc, amd.cc, base-text-renderer.cc, base-text-renderer.h, besselj.cc, bitfcns.cc, bsxfun.cc, c-file-ptr-stream.h, call-stack.cc, call-stack.h, ccolamd.cc, cellfun.cc, chol.cc, colamd.cc, dasrt.cc, data.cc, debug.cc, defaults.cc, defaults.h, det.cc, display.cc, display.h, dlmread.cc, dynamic-ld.cc, dynamic-ld.h, 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, file-io.cc, filter.cc, find.cc, ft-text-renderer.cc, ft-text-renderer.h, gcd.cc, gl-render.cc, gl-render.h, gl2ps-print.cc, gl2ps-print.h, graphics-toolkit.cc, graphics-toolkit.h, graphics.cc, gsvd.cc, gtk-manager.cc, gtk-manager.h, help.cc, help.h, 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, latex-text-renderer.cc, latex-text-renderer.h, load-path.cc, load-path.h, load-save.cc, load-save.h, lookup.cc, ls-hdf5.cc, ls-mat4.cc, ls-mat5.cc, lsode.cc, lu.cc, mappers.cc, matrix_type.cc, max.cc, mex.cc, mexproto.h, mxarray.h, mxtypes.in.h, oct-errno.in.cc, oct-hdf5-types.cc, oct-hist.cc, oct-hist.h, oct-map.cc, oct-map.h, oct-opengl.h, oct-prcstrm.h, oct-process.cc, oct-process.h, oct-stdstrm.h, oct-stream.cc, oct-stream.h, oct-strstrm.h, octave-default-image.h, ordqz.cc, ordschur.cc, pager.cc, pager.h, pinv.cc, pow2.cc, pr-output.cc, psi.cc, qr.cc, quadcc.cc, rand.cc, regexp.cc, settings.cc, settings.h, sighandlers.cc, sighandlers.h, sparse-xpow.cc, sqrtm.cc, stack-frame.cc, stack-frame.h, stream-euler.cc, strfns.cc, svd.cc, syminfo.cc, syminfo.h, symrcm.cc, symrec.cc, symrec.h, symscope.cc, symscope.h, symtab.cc, symtab.h, sysdep.cc, sysdep.h, text-engine.cc, text-engine.h, text-renderer.cc, text-renderer.h, time.cc, toplev.cc, typecast.cc, url-handle-manager.cc, url-handle-manager.h, urlwrite.cc, utils.cc, utils.h, variables.cc, variables.h, xdiv.cc, __delaunayn__.cc, __init_fltk__.cc, __init_gnuplot__.cc, __ode15__.cc, __voronoi__.cc, audioread.cc, convhulln.cc, gzip.cc, 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-diag.cc, ov-base-int.cc, ov-base-mat.cc, ov-base-mat.h, ov-base-scalar.cc, ov-base.cc, ov-base.h, ov-bool-mat.cc, ov-bool-mat.h, ov-bool-sparse.cc, ov-bool.cc, ov-builtin.h, ov-cell.cc, ov-ch-mat.cc, ov-class.cc, ov-class.h, ov-classdef.cc, ov-classdef.h, ov-complex.cc, ov-cx-diag.cc, ov-cx-mat.cc, ov-cx-sparse.cc, ov-dld-fcn.cc, ov-dld-fcn.h, ov-fcn-handle.cc, ov-fcn-handle.h, ov-fcn.h, 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-flt-re-mat.h, ov-intx.h, ov-java.cc, ov-lazy-idx.cc, ov-legacy-range.cc, ov-magic-int.cc, ov-mex-fcn.cc, ov-mex-fcn.h, ov-null-mat.cc, ov-perm.cc, ov-range.cc, ov-re-diag.cc, ov-re-mat.cc, ov-re-mat.h, ov-re-sparse.cc, ov-scalar.cc, ov-str-mat.cc, ov-struct.cc, ov-typeinfo.cc, ov-typeinfo.h, ov-usr-fcn.cc, ov-usr-fcn.h, ov.cc, ov.h, ovl.h, octave.cc, octave.h, op-b-sbm.cc, op-bm-sbm.cc, op-cs-scm.cc, op-fm-fcm.cc, op-fs-fcm.cc, op-s-scm.cc, op-scm-cs.cc, op-scm-s.cc, op-sm-cs.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, oct-lvalue.cc, oct-lvalue.h, 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-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: Re-indent code after switch to using namespace macros.
author Rik <rik@octave.org>
date Thu, 01 Dec 2022 18:02:15 -0800
parents e88a07dec498
children 597f3ee61a48
line wrap: on
line source

////////////////////////////////////////////////////////////////////////
//
// Copyright (C) 1994-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 "Array.h"
#include "CMatrix.h"
#include "dMatrix.h"
#include "fCMatrix.h"
#include "fMatrix.h"
#include "lo-error.h"
#include "lo-lapack-proto.h"
#include "oct-locbuf.h"
#include "schur.h"

OCTAVE_BEGIN_NAMESPACE(octave)

OCTAVE_BEGIN_NAMESPACE(math)

// For real types.

static F77_INT
select_ana (const double& a, const double&)
{
  return (a < 0.0);
}

static F77_INT
select_dig (const double& a, const double& b)
{
  return (hypot (a, b) < 1.0);
}

static F77_INT
select_ana (const float& a, const float&)
{
  return (a < 0.0);
}

static F77_INT
select_dig (const float& a, const float& b)
{
  return (hypot (a, b) < 1.0);
}

// For complex types.

static F77_INT
select_ana (const F77_DBLE_CMPLX& a_arg)
{
  const Complex a = reinterpret_cast<const Complex&> (a_arg);
  return a.real () < 0.0;
}

static F77_INT
select_dig (const F77_DBLE_CMPLX& a_arg)
{
  const Complex& a = reinterpret_cast<const Complex&> (a_arg);
  return (abs (a) < 1.0);
}

static F77_INT
select_ana (const F77_CMPLX& a_arg)
{
  const FloatComplex& a = reinterpret_cast<const FloatComplex&> (a_arg);
  return a.real () < 0.0;
}

static F77_INT
select_dig (const F77_CMPLX& a_arg)
{
  const FloatComplex& a = reinterpret_cast<const FloatComplex&> (a_arg);
  return (abs (a) < 1.0);
}

template <>
OCTAVE_API F77_INT
schur<Matrix>::init (const Matrix& a, const std::string& ord,
                     bool calc_unitary)
{
  F77_INT a_nr = to_f77_int (a.rows ());
  F77_INT a_nc = to_f77_int (a.cols ());

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

  if (a_nr == 0)
    {
      m_schur_mat.clear ();
      m_unitary_schur_mat.clear ();
      return 0;
    }

  // Workspace requirements may need to be fixed if any of the
  // following change.

  char jobvs;
  char sense = 'N';
  char sort = 'N';

  if (calc_unitary)
    jobvs = 'V';
  else
    jobvs = 'N';

  char ord_char = (ord.empty () ? 'U' : ord[0]);

  if (ord_char == 'A' || ord_char == 'D'
      || ord_char == 'a' || ord_char == 'd')
    sort = 'S';

  volatile double_selector selector = nullptr;
  if (ord_char == 'A' || ord_char == 'a')
    selector = select_ana;
  else if (ord_char == 'D' || ord_char == 'd')
    selector = select_dig;

  F77_INT n = a_nc;
  F77_INT lwork = 8 * n;
  F77_INT liwork = 1;
  F77_INT info;
  F77_INT sdim;
  double rconde;
  double rcondv;

  m_schur_mat = a;

  if (calc_unitary)
    m_unitary_schur_mat.clear (n, n);

  double *s = m_schur_mat.fortran_vec ();
  double *q = m_unitary_schur_mat.fortran_vec ();

  Array<double> wr (dim_vector (n, 1));
  double *pwr = wr.fortran_vec ();

  Array<double> wi (dim_vector (n, 1));
  double *pwi = wi.fortran_vec ();

  Array<double> work (dim_vector (lwork, 1));
  double *pwork = work.fortran_vec ();

  // BWORK is not referenced for the non-ordered Schur routine.
  F77_INT ntmp = (ord_char == 'N' || ord_char == 'n') ? 0 : n;
  Array<F77_INT> bwork (dim_vector (ntmp, 1));
  F77_INT *pbwork = bwork.fortran_vec ();

  Array<F77_INT> iwork (dim_vector (liwork, 1));
  F77_INT *piwork = iwork.fortran_vec ();

  F77_XFCN (dgeesx, DGEESX, (F77_CONST_CHAR_ARG2 (&jobvs, 1),
                             F77_CONST_CHAR_ARG2 (&sort, 1),
                             selector,
                             F77_CONST_CHAR_ARG2 (&sense, 1),
                             n, s, n, sdim, pwr, pwi, q, n, rconde, rcondv,
                             pwork, lwork, piwork, liwork, pbwork, info
                             F77_CHAR_ARG_LEN (1)
                             F77_CHAR_ARG_LEN (1)
                             F77_CHAR_ARG_LEN (1)));

  return info;
}

template <>
OCTAVE_API F77_INT
schur<FloatMatrix>::init (const FloatMatrix& a, const std::string& ord,
                          bool calc_unitary)
{
  F77_INT a_nr = to_f77_int (a.rows ());
  F77_INT a_nc = to_f77_int (a.cols ());

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

  if (a_nr == 0)
    {
      m_schur_mat.clear ();
      m_unitary_schur_mat.clear ();
      return 0;
    }

  // Workspace requirements may need to be fixed if any of the
  // following change.

  char jobvs;
  char sense = 'N';
  char sort = 'N';

  if (calc_unitary)
    jobvs = 'V';
  else
    jobvs = 'N';

  char ord_char = (ord.empty () ? 'U' : ord[0]);

  if (ord_char == 'A' || ord_char == 'D'
      || ord_char == 'a' || ord_char == 'd')
    sort = 'S';

  volatile float_selector selector = nullptr;
  if (ord_char == 'A' || ord_char == 'a')
    selector = select_ana;
  else if (ord_char == 'D' || ord_char == 'd')
    selector = select_dig;

  F77_INT n = a_nc;
  F77_INT lwork = 8 * n;
  F77_INT liwork = 1;
  F77_INT info;
  F77_INT sdim;
  float rconde;
  float rcondv;

  m_schur_mat = a;

  if (calc_unitary)
    m_unitary_schur_mat.clear (n, n);

  float *s = m_schur_mat.fortran_vec ();
  float *q = m_unitary_schur_mat.fortran_vec ();

  Array<float> wr (dim_vector (n, 1));
  float *pwr = wr.fortran_vec ();

  Array<float> wi (dim_vector (n, 1));
  float *pwi = wi.fortran_vec ();

  Array<float> work (dim_vector (lwork, 1));
  float *pwork = work.fortran_vec ();

  // BWORK is not referenced for the non-ordered Schur routine.
  F77_INT ntmp = (ord_char == 'N' || ord_char == 'n') ? 0 : n;
  Array<F77_INT> bwork (dim_vector (ntmp, 1));
  F77_INT *pbwork = bwork.fortran_vec ();

  Array<F77_INT> iwork (dim_vector (liwork, 1));
  F77_INT *piwork = iwork.fortran_vec ();

  F77_XFCN (sgeesx, SGEESX, (F77_CONST_CHAR_ARG2 (&jobvs, 1),
                             F77_CONST_CHAR_ARG2 (&sort, 1),
                             selector,
                             F77_CONST_CHAR_ARG2 (&sense, 1),
                             n, s, n, sdim, pwr, pwi, q, n, rconde, rcondv,
                             pwork, lwork, piwork, liwork, pbwork, info
                             F77_CHAR_ARG_LEN (1)
                             F77_CHAR_ARG_LEN (1)
                             F77_CHAR_ARG_LEN (1)));

  return info;
}

template <>
OCTAVE_API F77_INT
schur<ComplexMatrix>::init (const ComplexMatrix& a, const std::string& ord,
                            bool calc_unitary)
{
  F77_INT a_nr = to_f77_int (a.rows ());
  F77_INT a_nc = to_f77_int (a.cols ());

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

  if (a_nr == 0)
    {
      m_schur_mat.clear ();
      m_unitary_schur_mat.clear ();
      return 0;
    }

  // Workspace requirements may need to be fixed if any of the
  // following change.

  char jobvs;
  char sense = 'N';
  char sort = 'N';

  if (calc_unitary)
    jobvs = 'V';
  else
    jobvs = 'N';

  char ord_char = (ord.empty () ? 'U' : ord[0]);

  if (ord_char == 'A' || ord_char == 'D'
      || ord_char == 'a' || ord_char == 'd')
    sort = 'S';

  volatile complex_selector selector = nullptr;
  if (ord_char == 'A' || ord_char == 'a')
    selector = select_ana;
  else if (ord_char == 'D' || ord_char == 'd')
    selector = select_dig;

  F77_INT n = a_nc;
  F77_INT lwork = 8 * n;
  F77_INT info;
  F77_INT sdim;
  double rconde;
  double rcondv;

  m_schur_mat = a;
  if (calc_unitary)
    m_unitary_schur_mat.clear (n, n);

  Complex *s = m_schur_mat.fortran_vec ();
  Complex *q = m_unitary_schur_mat.fortran_vec ();

  Array<double> rwork (dim_vector (n, 1));
  double *prwork = rwork.fortran_vec ();

  Array<Complex> w (dim_vector (n, 1));
  Complex *pw = w.fortran_vec ();

  Array<Complex> work (dim_vector (lwork, 1));
  Complex *pwork = work.fortran_vec ();

  // BWORK is not referenced for non-ordered Schur.
  F77_INT ntmp = (ord_char == 'N' || ord_char == 'n') ? 0 : n;
  Array<F77_INT> bwork (dim_vector (ntmp, 1));
  F77_INT *pbwork = bwork.fortran_vec ();

  F77_XFCN (zgeesx, ZGEESX,
            (F77_CONST_CHAR_ARG2 (&jobvs, 1),
             F77_CONST_CHAR_ARG2 (&sort, 1),
             selector,
             F77_CONST_CHAR_ARG2 (&sense, 1),
             n, F77_DBLE_CMPLX_ARG (s), n, sdim, F77_DBLE_CMPLX_ARG (pw),
             F77_DBLE_CMPLX_ARG (q), n, rconde, rcondv,
             F77_DBLE_CMPLX_ARG (pwork), lwork, prwork, pbwork, info
             F77_CHAR_ARG_LEN (1)
             F77_CHAR_ARG_LEN (1)
             F77_CHAR_ARG_LEN (1)));

  return info;
}

template <>
OCTAVE_API schur<ComplexMatrix>
rsf2csf<ComplexMatrix, Matrix> (const Matrix& s_arg, const Matrix& u_arg)
{
  ComplexMatrix s (s_arg);
  ComplexMatrix u (u_arg);

  F77_INT n = to_f77_int (s.rows ());

  if (s.columns () != n || u.rows () != n || u.columns () != n)
    (*current_liboctave_error_handler)
      ("rsf2csf: inconsistent matrix dimensions");

  if (n > 0)
    {
      OCTAVE_LOCAL_BUFFER (double, c, n-1);
      OCTAVE_LOCAL_BUFFER (double, sx, n-1);

      F77_XFCN (zrsf2csf, ZRSF2CSF,
                (n, F77_DBLE_CMPLX_ARG (s.fortran_vec ()),
                 F77_DBLE_CMPLX_ARG (u.fortran_vec ()), c, sx));
    }

  return schur<ComplexMatrix> (s, u);
}

template <>
OCTAVE_API F77_INT
schur<FloatComplexMatrix>::init (const FloatComplexMatrix& a,
                                 const std::string& ord, bool calc_unitary)
{
  F77_INT a_nr = to_f77_int (a.rows ());
  F77_INT a_nc = to_f77_int (a.cols ());

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

  if (a_nr == 0)
    {
      m_schur_mat.clear ();
      m_unitary_schur_mat.clear ();
      return 0;
    }

  // Workspace requirements may need to be fixed if any of the
  // following change.

  char jobvs;
  char sense = 'N';
  char sort = 'N';

  if (calc_unitary)
    jobvs = 'V';
  else
    jobvs = 'N';

  char ord_char = (ord.empty () ? 'U' : ord[0]);

  if (ord_char == 'A' || ord_char == 'D'
      || ord_char == 'a' || ord_char == 'd')
    sort = 'S';

  volatile float_complex_selector selector = nullptr;
  if (ord_char == 'A' || ord_char == 'a')
    selector = select_ana;
  else if (ord_char == 'D' || ord_char == 'd')
    selector = select_dig;

  F77_INT n = a_nc;
  F77_INT lwork = 8 * n;
  F77_INT info;
  F77_INT sdim;
  float rconde;
  float rcondv;

  m_schur_mat = a;
  if (calc_unitary)
    m_unitary_schur_mat.clear (n, n);

  FloatComplex *s = m_schur_mat.fortran_vec ();
  FloatComplex *q = m_unitary_schur_mat.fortran_vec ();

  Array<float> rwork (dim_vector (n, 1));
  float *prwork = rwork.fortran_vec ();

  Array<FloatComplex> w (dim_vector (n, 1));
  FloatComplex *pw = w.fortran_vec ();

  Array<FloatComplex> work (dim_vector (lwork, 1));
  FloatComplex *pwork = work.fortran_vec ();

  // BWORK is not referenced for non-ordered Schur.
  F77_INT ntmp = (ord_char == 'N' || ord_char == 'n') ? 0 : n;
  Array<F77_INT> bwork (dim_vector (ntmp, 1));
  F77_INT *pbwork = bwork.fortran_vec ();

  F77_XFCN (cgeesx, CGEESX,
            (F77_CONST_CHAR_ARG2 (&jobvs, 1),
             F77_CONST_CHAR_ARG2 (&sort, 1),
             selector,
             F77_CONST_CHAR_ARG2 (&sense, 1),
             n, F77_CMPLX_ARG (s), n, sdim, F77_CMPLX_ARG (pw),
             F77_CMPLX_ARG (q), n,
             rconde, rcondv,
             F77_CMPLX_ARG (pwork), lwork, prwork, pbwork, info
             F77_CHAR_ARG_LEN (1)
             F77_CHAR_ARG_LEN (1)
             F77_CHAR_ARG_LEN (1)));

  return info;
}

template <>
OCTAVE_API schur<FloatComplexMatrix>
rsf2csf<FloatComplexMatrix, FloatMatrix> (const FloatMatrix& s_arg,
    const FloatMatrix& u_arg)
{
  FloatComplexMatrix s (s_arg);
  FloatComplexMatrix u (u_arg);

  F77_INT n = to_f77_int (s.rows ());

  if (s.columns () != n || u.rows () != n || u.columns () != n)
    (*current_liboctave_error_handler)
      ("rsf2csf: inconsistent matrix dimensions");

  if (n > 0)
    {
      OCTAVE_LOCAL_BUFFER (float, c, n-1);
      OCTAVE_LOCAL_BUFFER (float, sx, n-1);

      F77_XFCN (crsf2csf, CRSF2CSF,
                (n, F77_CMPLX_ARG (s.fortran_vec ()),
                 F77_CMPLX_ARG (u.fortran_vec ()), c, sx));
    }

  return schur<FloatComplexMatrix> (s, u);
}

// Instantiations we need.

template class schur<ComplexMatrix>;

template class schur<FloatComplexMatrix>;

template class schur<FloatMatrix>;

template class schur<Matrix>;

OCTAVE_END_NAMESPACE(math)
OCTAVE_END_NAMESPACE(octave)