view libinterp/corefcn/schur.cc @ 21200:fcac5dbbf9ed

maint: Indent #ifdef blocks in libinterp. * builtins.h, 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, c-file-ptr-stream.h, cellfun.cc, colloc.cc, comment-list.cc, conv2.cc, daspk.cc, dasrt.cc, dassl.cc, data.cc, debug.cc, defaults.cc, defaults.in.h, defun-dld.h, defun.cc, defun.h, 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, gammainc.cc, gcd.cc, getgrent.cc, getpwent.cc, getrusage.cc, givens.cc, gl-render.cc, gl2ps-print.cc, graphics.cc, graphics.in.h, gripes.cc, hash.cc, help.cc, hess.cc, hex2num.cc, input.cc, inv.cc, jit-ir.cc, jit-typeinfo.cc, jit-util.cc, jit-util.h, 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-oct-text.h, ls-utils.cc, ls-utils.h, lsode.cc, lu.cc, luinc.cc, mappers.cc, matrix_type.cc, max.cc, mex.h, mexproto.h, mgorth.cc, nproc.cc, oct-errno.in.cc, oct-fstrm.cc, oct-hdf5-types.cc, oct-hdf5.h, oct-hist.cc, oct-iostrm.cc, oct-lvalue.cc, oct-map.cc, oct-prcstrm.cc, oct-procbuf.cc, oct-stream.cc, oct-strstrm.cc, 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, 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, sysdep.h, time.cc, toplev.cc, tril.cc, tsearch.cc, txt-eng-ft.cc, txt-eng.cc, typecast.cc, urlwrite.cc, utils.cc, variables.cc, xdiv.cc, 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, oct-qhull.h, qr.cc, symbfact.cc, symrcm.cc, oct-conf.in.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-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, 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-colon.h, 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: Indent #ifdef blocks in libinterp.
author Rik <rik@octave.org>
date Fri, 05 Feb 2016 16:29:08 -0800
parents 538b57866b90
children e69eaee28737
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 <typename 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", "A");

  if (! arg.is_numeric_type ())
    err_wrong_type_arg ("schur", arg);

  octave_value_list retval;

  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 <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 ();

  if (! args(0).is_numeric_type ())
    err_wrong_type_arg ("rsf2csf", args(0));
  if (! args(1).is_numeric_type ())
    err_wrong_type_arg ("rsf2csf", args(1));
  if (args(0).is_complex_type () || args(1).is_complex_type ())
    error ("rsf2csf: UR and TR must be real matrices");

  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));

      return 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));

      return ovl (cs.unitary_matrix (), cs.schur_matrix ());
    }
}

/*
%!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);
*/