view libinterp/corefcn/schur.cc @ 21100:e39e05d90788

Switch gripe_XXX to either err_XXX or warn_XXX naming scheme. * libinterp/corefcn/errwarn.h, libinterp/corefcn/errwarn.cc: New header and .cc file with common errors and warnings for libinterp. * libinterp/corefcn/module.mk: Add errwarn.h, errwarn.cc to build system. * liboctave/util/lo-array-errwarn.h, liboctave/util/lo-array-errwarn.cc: New header and .cc file with common errors and warnings for liboctave. * liboctave/util/module.mk: Add lo-array-errwarn.h, lo-array-errwarn.cc to build system. * lo-array-gripes.h: #include "lo-array-errwarn.h" for access to class index_exception. Remove const char *error_id_XXX prototypes. * lo-array-gripes.cc: Remove const char *error_id_XXX initializations. Remove index_exception method definitions. * Cell.cc, __pchip_deriv__.cc, __qp__.cc, balance.cc, betainc.cc, cellfun.cc, daspk.cc, dasrt.cc, dassl.cc, data.cc, debug.cc, defaults.cc, det.cc, dirfns.cc, eig.cc, fft.cc, fft2.cc, fftn.cc, find.cc, gammainc.cc, gcd.cc, getgrent.cc, getpwent.cc, graphics.in.h, help.cc, hess.cc, hex2num.cc, input.cc, inv.cc, jit-typeinfo.cc, load-save.cc, lookup.cc, ls-hdf5.cc, ls-mat-ascii.cc, ls-mat4.cc, ls-mat5.cc, ls-oct-binary.cc, ls-oct-text.cc, lsode.cc, lu.cc, luinc.cc, max.cc, mgorth.cc, oct-hist.cc, oct-procbuf.cc, oct-stream.cc, oct.h, pager.cc, pinv.cc, pr-output.cc, quad.cc, qz.cc, rand.cc, rcond.cc, regexp.cc, schur.cc, sparse-xdiv.cc, sparse-xpow.cc, sparse.cc, spparms.cc, sqrtm.cc, str2double.cc, strfind.cc, strfns.cc, sub2ind.cc, svd.cc, sylvester.cc, syscalls.cc, typecast.cc, utils.cc, variables.cc, xdiv.cc, xnorm.cc, xpow.cc, __eigs__.cc, __glpk__.cc, __magick_read__.cc, __osmesa_print__.cc, audiodevinfo.cc, audioread.cc, chol.cc, dmperm.cc, fftw.cc, qr.cc, symbfact.cc, symrcm.cc, ov-base-diag.cc, ov-base-int.cc, ov-base-mat.cc, ov-base-scalar.cc, ov-base-sparse.cc, ov-base.cc, ov-bool-mat.cc, ov-bool-sparse.cc, ov-bool.cc, ov-builtin.cc, ov-cell.cc, ov-ch-mat.cc, ov-class.cc, ov-complex.cc, ov-complex.h, ov-cs-list.cc, ov-cx-diag.cc, ov-cx-mat.cc, ov-cx-sparse.cc, ov-fcn-handle.cc, ov-fcn-inline.cc, ov-float.cc, ov-float.h, ov-flt-complex.cc, ov-flt-complex.h, ov-flt-cx-diag.cc, ov-flt-cx-mat.cc, ov-flt-re-mat.cc, ov-int16.cc, ov-int32.cc, ov-int64.cc, ov-int8.cc, ov-intx.h, ov-mex-fcn.cc, ov-perm.cc, ov-range.cc, ov-re-mat.cc, ov-re-sparse.cc, ov-scalar.cc, ov-scalar.h, ov-str-mat.cc, ov-struct.cc, ov-type-conv.h, ov-uint16.cc, ov-uint32.cc, ov-uint64.cc, ov-uint8.cc, ov-usr-fcn.cc, ov.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-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-int.h, 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-range.cc, op-s-cm.cc, op-s-cs.cc, op-s-m.cc, op-s-s.cc, op-s-scm.cc, op-s-sm.cc, op-sbm-b.cc, op-sbm-bm.cc, op-sbm-sbm.cc, op-scm-cm.cc, op-scm-cs.cc, op-scm-m.cc, op-scm-s.cc, op-scm-scm.cc, op-scm-sm.cc, op-sm-cm.cc, op-sm-cs.cc, op-sm-m.cc, op-sm-s.cc, op-sm-scm.cc, op-sm-sm.cc, op-str-m.cc, op-str-s.cc, op-str-str.cc, op-struct.cc, op-ui16-ui16.cc, op-ui32-ui32.cc, op-ui64-ui64.cc, op-ui8-ui8.cc, ops.h, lex.ll, pt-assign.cc, pt-eval.cc, pt-idx.cc, pt-loop.cc, pt-mat.cc, pt-stmt.cc, Array-util.cc, Array-util.h, Array.cc, CColVector.cc, CDiagMatrix.cc, CMatrix.cc, CNDArray.cc, CRowVector.cc, CSparse.cc, DiagArray2.cc, MDiagArray2.cc, MSparse.cc, PermMatrix.cc, Range.cc, Sparse.cc, dColVector.cc, dDiagMatrix.cc, dMatrix.cc, dNDArray.cc, dRowVector.cc, dSparse.cc, fCColVector.cc, fCDiagMatrix.cc, fCMatrix.cc, fCNDArray.cc, fCRowVector.cc, fColVector.cc, fDiagMatrix.cc, fMatrix.cc, fNDArray.cc, fRowVector.cc, idx-vector.cc, CmplxGEPBAL.cc, dbleGEPBAL.cc, fCmplxGEPBAL.cc, floatGEPBAL.cc, Sparse-diag-op-defs.h, Sparse-op-defs.h, Sparse-perm-op-defs.h, mx-inlines.cc, mx-op-defs.h, oct-binmap.h: Replace 'include "gripes.h"' with 'include "errwarn.h". Change all gripe_XXX to err_XXX or warn_XXX or errwarn_XXX.
author Rik <rik@octave.org>
date Mon, 18 Jan 2016 18:28:06 -0800
parents 5e00ed38a58b
children 3d0d84305600
line wrap: on
line source

/*

Copyright (C) 1996-2015 John W. Eaton

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 <string>

#include "CmplxSCHUR.h"
#include "dbleSCHUR.h"
#include "fCmplxSCHUR.h"
#include "floatSCHUR.h"

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

template <class Matrix>
static octave_value
mark_upper_triangular (const Matrix& a)
{
  octave_value retval = a;

  octave_idx_type n = a.rows ();
  assert (a.columns () == n);

  const typename Matrix::element_type zero = typename Matrix::element_type ();

  for (octave_idx_type i = 0; i < n; i++)
    if (a(i,i) == zero)
      return retval;

  retval.matrix_type (MatrixType::Upper);

  return retval;
}

DEFUN (schur, args, nargout,
       "-*- texinfo -*-\n\
@deftypefn  {} {@var{S} =} schur (@var{A})\n\
@deftypefnx {} {@var{S} =} schur (@var{A}, \"real\")\n\
@deftypefnx {} {@var{S} =} schur (@var{A}, \"complex\")\n\
@deftypefnx {} {@var{S} =} schur (@var{A}, @var{opt})\n\
@deftypefnx {} {[@var{U}, @var{S}] =} schur (@dots{})\n\
@cindex Schur decomposition\n\
Compute the Schur@tie{}decomposition of @var{A}.\n\
\n\
The Schur@tie{}decomposition is defined as\n\
@tex\n\
$$\n\
 S = U^T A U\n\
$$\n\
@end tex\n\
@ifnottex\n\
\n\
@example\n\
@code{@var{S} = @var{U}' * @var{A} * @var{U}}\n\
@end example\n\
\n\
@end ifnottex\n\
where @var{U} is a unitary matrix\n\
@tex\n\
($U^T U$ is identity)\n\
@end tex\n\
@ifnottex\n\
(@code{@var{U}'* @var{U}} is identity)\n\
@end ifnottex\n\
and @var{S} is upper triangular.  The eigenvalues of @var{A} (and @var{S})\n\
are the diagonal elements of @var{S}.  If the matrix @var{A} is real, then\n\
the real Schur@tie{}decomposition is computed, in which the matrix @var{U}\n\
is orthogonal and @var{S} is block upper triangular with blocks of size at\n\
most\n\
@tex\n\
$2 \\times 2$\n\
@end tex\n\
@ifnottex\n\
@code{2 x 2}\n\
@end ifnottex\n\
along the diagonal.  The diagonal elements of @var{S}\n\
(or the eigenvalues of the\n\
@tex\n\
$2 \\times 2$\n\
@end tex\n\
@ifnottex\n\
@code{2 x 2}\n\
@end ifnottex\n\
blocks, when appropriate) are the eigenvalues of @var{A} and @var{S}.\n\
\n\
The default for real matrices is a real Schur@tie{}decomposition.\n\
A complex decomposition may be forced by passing the flag\n\
@qcode{\"complex\"}.\n\
\n\
The eigenvalues are optionally ordered along the diagonal according to the\n\
value of @var{opt}.  @code{@var{opt} = \"a\"} indicates that all eigenvalues\n\
with negative real parts should be moved to the leading block of @var{S}\n\
(used in @code{are}), @code{@var{opt} = \"d\"} indicates that all\n\
eigenvalues with magnitude less than one should be moved to the leading\n\
block of @var{S} (used in @code{dare}), and @code{@var{opt} = \"u\"}, the\n\
default, indicates that no ordering of eigenvalues should occur.  The\n\
leading @var{k} columns of @var{U} always span the @var{A}-invariant\n\
subspace corresponding to the @var{k} leading eigenvalues of @var{S}.\n\
\n\
The Schur@tie{}decomposition is used to compute eigenvalues of a square\n\
matrix, and has applications in the solution of algebraic Riccati equations\n\
in control (see @code{are} and @code{dare}).\n\
@seealso{rsf2csf, ordschur, lu, chol, hess, qr, qz, svd}\n\
@end deftypefn")
{
  int nargin = args.length ();

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

  octave_value arg = args(0);

  std::string ord;
  if (nargin == 2)
    ord = args(1).xstring_value ("schur: second argument must be a string");

  bool force_complex = false;

  if (ord == "real")
    {
      ord = "";
    }
  else if (ord == "complex")
    {
      force_complex = true;
      ord = "";
    }
  else
    {
      char ord_char = ord.empty () ? 'U' : ord[0];

      if (ord_char != 'U' && ord_char != 'A' && ord_char != 'D'
          && ord_char != 'u' && ord_char != 'a' && ord_char != 'd')
        {
          warning ("schur: incorrect ordered schur argument '%s'",
                   ord.c_str ());
          return ovl ();
        }
    }

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

  if (nr != nc)
    err_square_matrix_required ("schur");

  octave_value_list retval;

  if (! arg.is_numeric_type ())
    err_wrong_type_arg ("schur", arg);
  else if (arg.is_single_type ())
    {
      if (! force_complex && arg.is_real_type ())
        {
          FloatMatrix tmp = arg.float_matrix_value ();

          if (nargout <= 1)
            {
              FloatSCHUR result (tmp, ord, false);
              retval = ovl (result.schur_matrix ());
            }
          else
            {
              FloatSCHUR result (tmp, ord, true);
              retval = ovl (result.unitary_matrix (),
                            result.schur_matrix ());
            }
        }
      else
        {
          FloatComplexMatrix ctmp = arg.float_complex_matrix_value ();

          if (nargout <= 1)
            {
              FloatComplexSCHUR result (ctmp, ord, false);
              retval = ovl (mark_upper_triangular (result.schur_matrix ()));
            }
          else
            {
              FloatComplexSCHUR result (ctmp, ord, true);
              retval = ovl (result.unitary_matrix (),
                            mark_upper_triangular (result.schur_matrix ()));
            }
        }
    }
  else
    {
      if (! force_complex && arg.is_real_type ())
        {
          Matrix tmp = arg.matrix_value ();

          if (nargout <= 1)
            {
              SCHUR result (tmp, ord, false);
              retval = ovl (result.schur_matrix ());
            }
          else
            {
              SCHUR result (tmp, ord, true);
              retval = ovl (result.unitary_matrix (),
                            result.schur_matrix ());
            }
        }
      else
        {
          ComplexMatrix ctmp = arg.complex_matrix_value ();

          if (nargout <= 1)
            {
              ComplexSCHUR result (ctmp, ord, false);
              retval = ovl (mark_upper_triangular (result.schur_matrix ()));
            }
          else
            {
              ComplexSCHUR result (ctmp, ord, true);
              retval = ovl (result.unitary_matrix (),
                            mark_upper_triangular (result.schur_matrix ()));
            }
        }
    }

  return retval;
}

/*
%!test
%! a = [1, 2, 3; 4, 5, 9; 7, 8, 6];
%! [u, s] = schur (a);
%! assert (u' * a * u, s, sqrt (eps));

%!test
%! a = single ([1, 2, 3; 4, 5, 9; 7, 8, 6]);
%! [u, s] = schur (a);
%! assert (u' * a * u, s, sqrt (eps ("single")));

%!error schur ()
%!error schur (1,2,3)
%!error [a,b,c] = schur (1)
%!error <argument must be a square matrix> schur ([1, 2, 3; 4, 5, 6])
%!error <wrong type argument 'cell'> schur ({1})
%!warning <incorrect ordered schur argument> schur ([1, 2; 3, 4], "bad_opt");

*/

DEFUN (rsf2csf, args, nargout,
       "-*- texinfo -*-\n\
@deftypefn {} {[@var{U}, @var{T}] =} rsf2csf (@var{UR}, @var{TR})\n\
Convert a real, upper quasi-triangular Schur@tie{}form @var{TR} to a complex,\n\
upper triangular Schur@tie{}form @var{T}.\n\
\n\
Note that the following relations hold:\n\
\n\
@tex\n\
$UR \\cdot TR \\cdot {UR}^T = U T U^{\\dagger}$ and\n\
$U^{\\dagger} U$ is the identity matrix I.\n\
@end tex\n\
@ifnottex\n\
@tcode{@var{UR} * @var{TR} * @var{UR}' = @var{U} * @var{T} * @var{U}'} and\n\
@code{@var{U}' * @var{U}} is the identity matrix I.\n\
@end ifnottex\n\
\n\
Note also that @var{U} and @var{T} are not unique.\n\
@seealso{schur}\n\
@end deftypefn")
{
  if (args.length () != 2 || nargout > 2)
    print_usage ();

  octave_value_list retval;

  if (! args(0).is_numeric_type ())
    err_wrong_type_arg ("rsf2csf", args(0));
  else if (! args(1).is_numeric_type ())
    err_wrong_type_arg ("rsf2csf", args(1));
  else if (args(0).is_complex_type () || args(1).is_complex_type ())
    error ("rsf2csf: UR and TR must be real matrices");
  else
    {
      if (args(0).is_single_type () || args(1).is_single_type ())
        {
          FloatMatrix u = args(0).float_matrix_value ();
          FloatMatrix t = args(1).float_matrix_value ();

          FloatComplexSCHUR cs (FloatSCHUR (t, u));

          retval = ovl (cs.unitary_matrix (), cs.schur_matrix ());
        }
      else
        {
          Matrix u = args(0).matrix_value ();
          Matrix t = args(1).matrix_value ();

          ComplexSCHUR cs (SCHUR (t, u));

          retval = ovl (cs.unitary_matrix (), cs.schur_matrix ());
        }
    }

  return retval;
}

/*
%!test
%! A = [1, 1, 1, 2; 1, 2, 1, 1; 1, 1, 3, 1; -2, 1, 1, 1];
%! [u, t] = schur (A);
%! [U, T] = rsf2csf (u, t);
%! assert (norm (u * t * u' - U * T * U'), 0, 1e-12);
%! assert (norm (A - U * T * U'), 0, 1e-12);

%!test
%! A = rand (10);
%! [u, t] = schur (A);
%! [U, T] = rsf2csf (u, t);
%! assert (norm (tril (T, -1)), 0);
%! assert (norm (U * U'), 1, 1e-14);

%!test
%! A = [0, 1;-1, 0];
%! [u, t] = schur (A);
%! [U, T] = rsf2csf (u,t);
%! assert (U * T * U', A, 1e-14);
*/