view libinterp/corefcn/eig.cc @ 17787:175b392e91fe

Use GNU style coding conventions for code in libinterp/ * libinterp/corefcn/Cell.h, libinterp/corefcn/__contourc__.cc, libinterp/corefcn/__dispatch__.cc, libinterp/corefcn/__lin_interpn__.cc, libinterp/corefcn/__pchip_deriv__.cc, libinterp/corefcn/__qp__.cc, libinterp/corefcn/balance.cc, libinterp/corefcn/besselj.cc, libinterp/corefcn/betainc.cc, libinterp/corefcn/bitfcns.cc, libinterp/corefcn/bsxfun.cc, libinterp/corefcn/c-file-ptr-stream.cc, libinterp/corefcn/c-file-ptr-stream.h, libinterp/corefcn/cellfun.cc, libinterp/corefcn/colloc.cc, libinterp/corefcn/comment-list.h, libinterp/corefcn/conv2.cc, libinterp/corefcn/daspk.cc, libinterp/corefcn/dasrt.cc, libinterp/corefcn/dassl.cc, libinterp/corefcn/data.cc, libinterp/corefcn/debug.cc, libinterp/corefcn/defaults.cc, libinterp/corefcn/defaults.in.h, libinterp/corefcn/defun-int.h, libinterp/corefcn/defun.cc, libinterp/corefcn/det.cc, libinterp/corefcn/dirfns.cc, libinterp/corefcn/display.cc, libinterp/corefcn/dlmread.cc, libinterp/corefcn/dot.cc, libinterp/corefcn/dynamic-ld.cc, libinterp/corefcn/dynamic-ld.h, libinterp/corefcn/eig.cc, libinterp/corefcn/ellipj.cc, libinterp/corefcn/error.cc, libinterp/corefcn/error.h, libinterp/corefcn/event-queue.h, libinterp/corefcn/fft.cc, libinterp/corefcn/fft2.cc, libinterp/corefcn/fftn.cc, libinterp/corefcn/file-io.cc, libinterp/corefcn/filter.cc, libinterp/corefcn/find.cc, libinterp/corefcn/gammainc.cc, libinterp/corefcn/gcd.cc, libinterp/corefcn/getgrent.cc, libinterp/corefcn/getpwent.cc, libinterp/corefcn/getrusage.cc, libinterp/corefcn/givens.cc, libinterp/corefcn/gl-render.cc, libinterp/corefcn/gl2ps-renderer.cc, libinterp/corefcn/gl2ps-renderer.h, libinterp/corefcn/graphics.cc, libinterp/corefcn/graphics.in.h, libinterp/corefcn/gripes.cc, libinterp/corefcn/gripes.h, libinterp/corefcn/help.cc, libinterp/corefcn/hess.cc, libinterp/corefcn/hex2num.cc, libinterp/corefcn/input.cc, libinterp/corefcn/input.h, libinterp/corefcn/inv.cc, libinterp/corefcn/jit-ir.h, libinterp/corefcn/jit-typeinfo.cc, libinterp/corefcn/jit-typeinfo.h, libinterp/corefcn/jit-util.h, libinterp/corefcn/kron.cc, libinterp/corefcn/load-path.cc, libinterp/corefcn/load-path.h, libinterp/corefcn/load-save.cc, libinterp/corefcn/load-save.h, libinterp/corefcn/lookup.cc, libinterp/corefcn/ls-ascii-helper.cc, libinterp/corefcn/ls-hdf5.cc, libinterp/corefcn/ls-hdf5.h, libinterp/corefcn/ls-mat-ascii.cc, libinterp/corefcn/ls-mat-ascii.h, libinterp/corefcn/ls-mat4.cc, libinterp/corefcn/ls-mat5.cc, libinterp/corefcn/ls-mat5.h, libinterp/corefcn/ls-oct-ascii.cc, libinterp/corefcn/lsode.cc, libinterp/corefcn/lu.cc, libinterp/corefcn/luinc.cc, libinterp/corefcn/mappers.cc, libinterp/corefcn/matrix_type.cc, libinterp/corefcn/max.cc, libinterp/corefcn/md5sum.cc, libinterp/corefcn/mex.cc, libinterp/corefcn/mexproto.h, libinterp/corefcn/mgorth.cc, libinterp/corefcn/mxarray.in.h, libinterp/corefcn/nproc.cc, libinterp/corefcn/oct-hist.cc, libinterp/corefcn/oct-lvalue.h, libinterp/corefcn/oct-map.cc, libinterp/corefcn/oct-map.h, libinterp/corefcn/oct-obj.h, libinterp/corefcn/oct-prcstrm.h, libinterp/corefcn/oct-stdstrm.h, libinterp/corefcn/oct-stream.cc, libinterp/corefcn/oct-stream.h, libinterp/corefcn/octave-link.cc, libinterp/corefcn/octave-link.h, libinterp/corefcn/pager.cc, libinterp/corefcn/pinv.cc, libinterp/corefcn/pr-output.cc, libinterp/corefcn/procstream.h, libinterp/corefcn/profiler.cc, libinterp/corefcn/pt-jit.cc, libinterp/corefcn/pt-jit.h, libinterp/corefcn/quad.cc, libinterp/corefcn/quadcc.cc, libinterp/corefcn/qz.cc, libinterp/corefcn/rand.cc, libinterp/corefcn/rcond.cc, libinterp/corefcn/regexp.cc, libinterp/corefcn/schur.cc, libinterp/corefcn/sighandlers.cc, libinterp/corefcn/sighandlers.h, libinterp/corefcn/sparse-xdiv.cc, libinterp/corefcn/sparse-xdiv.h, libinterp/corefcn/sparse-xpow.cc, libinterp/corefcn/sparse.cc, libinterp/corefcn/spparms.cc, libinterp/corefcn/sqrtm.cc, libinterp/corefcn/str2double.cc, libinterp/corefcn/strfind.cc, libinterp/corefcn/strfns.cc, libinterp/corefcn/sub2ind.cc, libinterp/corefcn/svd.cc, libinterp/corefcn/syl.cc, libinterp/corefcn/symtab.cc, libinterp/corefcn/symtab.h, libinterp/corefcn/syscalls.cc, libinterp/corefcn/sysdep.cc, libinterp/corefcn/sysdep.h, libinterp/corefcn/time.cc, libinterp/corefcn/toplev.cc, libinterp/corefcn/toplev.h, libinterp/corefcn/tril.cc, libinterp/corefcn/txt-eng-ft.cc, libinterp/corefcn/txt-eng-ft.h, libinterp/corefcn/txt-eng.h, libinterp/corefcn/typecast.cc, libinterp/corefcn/urlwrite.cc, libinterp/corefcn/utils.cc, libinterp/corefcn/variables.cc, libinterp/corefcn/variables.h, libinterp/corefcn/xdiv.cc, libinterp/corefcn/xdiv.h, libinterp/corefcn/xnorm.h, libinterp/corefcn/xpow.cc, libinterp/corefcn/xpow.h, libinterp/corefcn/zfstream.cc, libinterp/corefcn/zfstream.h, libinterp/dldfcn/__delaunayn__.cc, libinterp/dldfcn/__dsearchn__.cc, libinterp/dldfcn/__eigs__.cc, libinterp/dldfcn/__fltk_uigetfile__.cc, libinterp/dldfcn/__glpk__.cc, libinterp/dldfcn/__init_fltk__.cc, libinterp/dldfcn/__init_gnuplot__.cc, libinterp/dldfcn/__magick_read__.cc, libinterp/dldfcn/__voronoi__.cc, libinterp/dldfcn/amd.cc, libinterp/dldfcn/ccolamd.cc, libinterp/dldfcn/chol.cc, libinterp/dldfcn/colamd.cc, libinterp/dldfcn/convhulln.cc, libinterp/dldfcn/dmperm.cc, libinterp/dldfcn/fftw.cc, libinterp/dldfcn/qr.cc, libinterp/dldfcn/symbfact.cc, libinterp/dldfcn/symrcm.cc, libinterp/dldfcn/tsearch.cc, libinterp/octave-value/ov-base-diag.cc, libinterp/octave-value/ov-base-diag.h, libinterp/octave-value/ov-base-int.cc, libinterp/octave-value/ov-base-int.h, libinterp/octave-value/ov-base-mat.h, libinterp/octave-value/ov-base-scalar.cc, libinterp/octave-value/ov-base-scalar.h, libinterp/octave-value/ov-base-sparse.cc, libinterp/octave-value/ov-base-sparse.h, libinterp/octave-value/ov-base.cc, libinterp/octave-value/ov-base.h, libinterp/octave-value/ov-bool-mat.cc, libinterp/octave-value/ov-bool-mat.h, libinterp/octave-value/ov-bool-sparse.cc, libinterp/octave-value/ov-bool-sparse.h, libinterp/octave-value/ov-bool.cc, libinterp/octave-value/ov-bool.h, libinterp/octave-value/ov-builtin.cc, libinterp/octave-value/ov-builtin.h, libinterp/octave-value/ov-cell.cc, libinterp/octave-value/ov-cell.h, libinterp/octave-value/ov-ch-mat.cc, libinterp/octave-value/ov-ch-mat.h, libinterp/octave-value/ov-class.cc, libinterp/octave-value/ov-class.h, libinterp/octave-value/ov-colon.h, libinterp/octave-value/ov-complex.cc, libinterp/octave-value/ov-complex.h, libinterp/octave-value/ov-cx-diag.cc, libinterp/octave-value/ov-cx-diag.h, libinterp/octave-value/ov-cx-mat.cc, libinterp/octave-value/ov-cx-mat.h, libinterp/octave-value/ov-cx-sparse.cc, libinterp/octave-value/ov-cx-sparse.h, libinterp/octave-value/ov-dld-fcn.h, libinterp/octave-value/ov-fcn-handle.cc, libinterp/octave-value/ov-fcn-handle.h, libinterp/octave-value/ov-fcn-inline.cc, libinterp/octave-value/ov-fcn-inline.h, libinterp/octave-value/ov-fcn.h, libinterp/octave-value/ov-float.cc, libinterp/octave-value/ov-float.h, libinterp/octave-value/ov-flt-complex.cc, libinterp/octave-value/ov-flt-complex.h, libinterp/octave-value/ov-flt-cx-diag.cc, libinterp/octave-value/ov-flt-cx-diag.h, libinterp/octave-value/ov-flt-cx-mat.cc, libinterp/octave-value/ov-flt-cx-mat.h, libinterp/octave-value/ov-flt-re-diag.cc, libinterp/octave-value/ov-flt-re-diag.h, libinterp/octave-value/ov-flt-re-mat.cc, libinterp/octave-value/ov-flt-re-mat.h, libinterp/octave-value/ov-int16.cc, libinterp/octave-value/ov-int32.cc, libinterp/octave-value/ov-int64.cc, libinterp/octave-value/ov-int8.cc, libinterp/octave-value/ov-intx.h, libinterp/octave-value/ov-java.cc, libinterp/octave-value/ov-lazy-idx.h, libinterp/octave-value/ov-mex-fcn.cc, libinterp/octave-value/ov-mex-fcn.h, libinterp/octave-value/ov-null-mat.cc, libinterp/octave-value/ov-null-mat.h, libinterp/octave-value/ov-oncleanup.cc, libinterp/octave-value/ov-perm.cc, libinterp/octave-value/ov-perm.h, libinterp/octave-value/ov-range.cc, libinterp/octave-value/ov-range.h, libinterp/octave-value/ov-re-diag.cc, libinterp/octave-value/ov-re-diag.h, libinterp/octave-value/ov-re-mat.cc, libinterp/octave-value/ov-re-mat.h, libinterp/octave-value/ov-re-sparse.cc, libinterp/octave-value/ov-re-sparse.h, libinterp/octave-value/ov-scalar.cc, libinterp/octave-value/ov-scalar.h, libinterp/octave-value/ov-str-mat.cc, libinterp/octave-value/ov-str-mat.h, libinterp/octave-value/ov-struct.cc, libinterp/octave-value/ov-struct.h, libinterp/octave-value/ov-type-conv.h, libinterp/octave-value/ov-typeinfo.cc, libinterp/octave-value/ov-typeinfo.h, libinterp/octave-value/ov-uint16.cc, libinterp/octave-value/ov-uint32.cc, libinterp/octave-value/ov-uint64.cc, libinterp/octave-value/ov-uint8.cc, libinterp/octave-value/ov-usr-fcn.cc, libinterp/octave-value/ov-usr-fcn.h, libinterp/octave-value/ov.cc, libinterp/octave-value/ov.h, libinterp/octave.cc, libinterp/operators/op-b-bm.cc, libinterp/operators/op-b-sbm.cc, libinterp/operators/op-bm-b.cc, libinterp/operators/op-bm-bm.cc, libinterp/operators/op-cdm-cdm.cc, libinterp/operators/op-chm.cc, libinterp/operators/op-class.cc, libinterp/operators/op-cm-cm.cc, libinterp/operators/op-cm-cs.cc, libinterp/operators/op-cm-s.cc, libinterp/operators/op-cm-scm.cc, libinterp/operators/op-cm-sm.cc, libinterp/operators/op-cs-cm.cc, libinterp/operators/op-cs-cs.cc, libinterp/operators/op-cs-scm.cc, libinterp/operators/op-cs-sm.cc, libinterp/operators/op-dm-dm.cc, libinterp/operators/op-dm-scm.cc, libinterp/operators/op-double-conv.cc, libinterp/operators/op-fcdm-fcdm.cc, libinterp/operators/op-fcm-fcm.cc, libinterp/operators/op-fcm-fcs.cc, libinterp/operators/op-fcm-fm.cc, libinterp/operators/op-fcm-fs.cc, libinterp/operators/op-fcs-fcm.cc, libinterp/operators/op-fcs-fcs.cc, libinterp/operators/op-fcs-fm.cc, libinterp/operators/op-fcs-fs.cc, libinterp/operators/op-fdm-fdm.cc, libinterp/operators/op-float-conv.cc, libinterp/operators/op-fm-fcm.cc, libinterp/operators/op-fm-fcs.cc, libinterp/operators/op-fm-fm.cc, libinterp/operators/op-fm-fs.cc, libinterp/operators/op-fs-fcm.cc, libinterp/operators/op-fs-fcs.cc, libinterp/operators/op-fs-fm.cc, libinterp/operators/op-fs-fs.cc, libinterp/operators/op-m-cm.cc, libinterp/operators/op-m-cs.cc, libinterp/operators/op-m-m.cc, libinterp/operators/op-m-s.cc, libinterp/operators/op-m-scm.cc, libinterp/operators/op-m-sm.cc, libinterp/operators/op-pm-scm.cc, libinterp/operators/op-range.cc, libinterp/operators/op-s-cm.cc, libinterp/operators/op-s-cs.cc, libinterp/operators/op-s-scm.cc, libinterp/operators/op-sbm-b.cc, libinterp/operators/op-sbm-bm.cc, libinterp/operators/op-sbm-sbm.cc, libinterp/operators/op-scm-cm.cc, libinterp/operators/op-scm-cs.cc, libinterp/operators/op-scm-m.cc, libinterp/operators/op-scm-s.cc, libinterp/operators/op-scm-scm.cc, libinterp/operators/op-scm-sm.cc, libinterp/operators/op-sm-cm.cc, libinterp/operators/op-sm-m.cc, libinterp/operators/op-sm-s.cc, libinterp/operators/op-sm-scm.cc, libinterp/operators/op-sm-sm.cc, libinterp/operators/op-str-m.cc, libinterp/operators/op-str-s.cc, libinterp/operators/op-str-str.cc, libinterp/operators/ops.h, libinterp/parse-tree/lex.h, libinterp/parse-tree/parse.h, libinterp/parse-tree/pt-arg-list.cc, libinterp/parse-tree/pt-arg-list.h, libinterp/parse-tree/pt-assign.cc, libinterp/parse-tree/pt-assign.h, libinterp/parse-tree/pt-binop.cc, libinterp/parse-tree/pt-binop.h, libinterp/parse-tree/pt-bp.h, libinterp/parse-tree/pt-cbinop.cc, libinterp/parse-tree/pt-check.cc, libinterp/parse-tree/pt-colon.cc, libinterp/parse-tree/pt-colon.h, libinterp/parse-tree/pt-const.cc, libinterp/parse-tree/pt-decl.cc, libinterp/parse-tree/pt-decl.h, libinterp/parse-tree/pt-eval.cc, libinterp/parse-tree/pt-except.h, libinterp/parse-tree/pt-exp.h, libinterp/parse-tree/pt-fcn-handle.cc, libinterp/parse-tree/pt-id.cc, libinterp/parse-tree/pt-id.h, libinterp/parse-tree/pt-idx.cc, libinterp/parse-tree/pt-idx.h, libinterp/parse-tree/pt-loop.h, libinterp/parse-tree/pt-mat.cc, libinterp/parse-tree/pt-misc.cc, libinterp/parse-tree/pt-misc.h, libinterp/parse-tree/pt-pr-code.cc, libinterp/parse-tree/pt-select.h, libinterp/parse-tree/pt-stmt.h, libinterp/parse-tree/token.h, libinterp/version.cc: Use GNU style coding conventions for code in libinterp/
author Rik <rik@octave.org>
date Mon, 28 Oct 2013 19:51:46 -0700
parents d63878346099
children 6a71e5030df5
line wrap: on
line source

/*

Copyright (C) 1996-2013 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 "EIG.h"
#include "fEIG.h"

#include "defun.h"
#include "error.h"
#include "gripes.h"
#include "oct-obj.h"
#include "utils.h"

DEFUN (eig, args, nargout,
       "-*- texinfo -*-\n\
@deftypefn  {Built-in Function} {@var{lambda} =} eig (@var{A})\n\
@deftypefnx {Built-in Function} {@var{lambda} =} eig (@var{A}, @var{B})\n\
@deftypefnx {Built-in Function} {[@var{V}, @var{lambda}] =} eig (@var{A})\n\
@deftypefnx {Built-in Function} {[@var{V}, @var{lambda}] =} eig (@var{A}, @var{B})\n\
Compute the eigenvalues (and optionally the eigenvectors) of a matrix\n\
or a pair of matrices\n\
\n\
The algorithm used depends on whether there are one or two input\n\
matrices, if they are real or complex and if they are symmetric\n\
(Hermitian if complex) or non-symmetric.\n\
\n\
The eigenvalues returned by @code{eig} are not ordered.\n\
@seealso{eigs, svd}\n\
@end deftypefn")
{
  octave_value_list retval;

  int nargin = args.length ();

  if (nargin > 2 || nargin == 0 || nargout > 2)
    {
      print_usage ();
      return retval;
    }

  octave_value arg_a, arg_b;

  octave_idx_type nr_a = 0, nr_b = 0;
  octave_idx_type nc_a = 0, nc_b = 0;

  arg_a = args(0);
  nr_a = arg_a.rows ();
  nc_a = arg_a.columns ();

  int arg_is_empty = empty_arg ("eig", nr_a, nc_a);
  if (arg_is_empty < 0)
    return retval;
  else if (arg_is_empty > 0)
    return octave_value_list (2, Matrix ());

  if (!(arg_a.is_single_type () || arg_a.is_double_type ()))
    {
      gripe_wrong_type_arg ("eig", arg_a);
      return retval;
    }

  if (nargin == 2)
    {
      arg_b = args(1);
      nr_b = arg_b.rows ();
      nc_b = arg_b.columns ();

      arg_is_empty = empty_arg ("eig", nr_b, nc_b);
      if (arg_is_empty < 0)
        return retval;
      else if (arg_is_empty > 0)
        return octave_value_list (2, Matrix ());

      if (!(arg_b.is_single_type () || arg_b.is_double_type ()))
        {
          gripe_wrong_type_arg ("eig", arg_b);
          return retval;
        }
    }

  if (nr_a != nc_a)
    {
      gripe_square_matrix_required ("eig");
      return retval;
    }

  if (nargin == 2 && nr_b != nc_b)
    {
      gripe_square_matrix_required ("eig");
      return retval;
    }

  Matrix tmp_a, tmp_b;
  ComplexMatrix ctmp_a, ctmp_b;
  FloatMatrix ftmp_a, ftmp_b;
  FloatComplexMatrix fctmp_a, fctmp_b;

  if (arg_a.is_single_type ())
    {
      FloatEIG result;

      if (nargin == 1)
        {
          if (arg_a.is_real_type ())
            {
              ftmp_a = arg_a.float_matrix_value ();

              if (error_state)
                return retval;
              else
                result = FloatEIG (ftmp_a, nargout > 1);
            }
          else
            {
              fctmp_a = arg_a.float_complex_matrix_value ();

              if (error_state)
                return retval;
              else
                result = FloatEIG (fctmp_a, nargout > 1);
            }
        }
      else if (nargin == 2)
        {
          if (arg_a.is_real_type () && arg_b.is_real_type ())
            {
              ftmp_a = arg_a.float_matrix_value ();
              ftmp_b = arg_b.float_matrix_value ();

              if (error_state)
                return retval;
              else
                result = FloatEIG (ftmp_a, ftmp_b, nargout > 1);
            }
          else
            {
              fctmp_a = arg_a.float_complex_matrix_value ();
              fctmp_b = arg_b.float_complex_matrix_value ();

              if (error_state)
                return retval;
              else
                result = FloatEIG (fctmp_a, fctmp_b, nargout > 1);
            }
        }

      if (! error_state)
        {
          if (nargout == 0 || nargout == 1)
            {
              retval(0) = result.eigenvalues ();
            }
          else
            {
              // Blame it on Matlab.

              FloatComplexDiagMatrix d (result.eigenvalues ());

              retval(1) = d;
              retval(0) = result.eigenvectors ();
            }
        }
    }
  else
    {
      EIG result;

      if (nargin == 1)
        {
          if (arg_a.is_real_type ())
            {
              tmp_a = arg_a.matrix_value ();

              if (error_state)
                return retval;
              else
                result = EIG (tmp_a, nargout > 1);
            }
          else
            {
              ctmp_a = arg_a.complex_matrix_value ();

              if (error_state)
                return retval;
              else
                result = EIG (ctmp_a, nargout > 1);
            }
        }
      else if (nargin == 2)
        {
          if (arg_a.is_real_type () && arg_b.is_real_type ())
            {
              tmp_a = arg_a.matrix_value ();
              tmp_b = arg_b.matrix_value ();

              if (error_state)
                return retval;
              else
                result = EIG (tmp_a, tmp_b, nargout > 1);
            }
          else
            {
              ctmp_a = arg_a.complex_matrix_value ();
              ctmp_b = arg_b.complex_matrix_value ();

              if (error_state)
                return retval;
              else
                result = EIG (ctmp_a, ctmp_b, nargout > 1);
            }
        }

      if (! error_state)
        {
          if (nargout == 0 || nargout == 1)
            {
              retval(0) = result.eigenvalues ();
            }
          else
            {
              // Blame it on Matlab.

              ComplexDiagMatrix d (result.eigenvalues ());

              retval(1) = d;
              retval(0) = result.eigenvectors ();
            }
        }
    }

  return retval;
}

/*
%!assert (eig ([1, 2; 2, 1]), [-1; 3], sqrt (eps))

%!test
%! [v, d] = eig ([1, 2; 2, 1]);
%! x = 1 / sqrt (2);
%! assert (d, [-1, 0; 0, 3], sqrt (eps));
%! assert (v, [-x, x; x, x], sqrt (eps));

%!assert (eig (single ([1, 2; 2, 1])), single ([-1; 3]), sqrt (eps ("single")))

%!test
%! [v, d] = eig (single ([1, 2; 2, 1]));
%! x = single (1 / sqrt (2));
%! assert (d, single ([-1, 0; 0, 3]), sqrt (eps ("single")));
%! assert (v, [-x, x; x, x], sqrt (eps ("single")));

%!test
%! A = [1, 2; -1, 1];  B = [3, 3; 1, 2];
%! [v, d] = eig (A, B);
%! assert (A * v(:, 1), d(1, 1) * B * v(:, 1), sqrt (eps));
%! assert (A * v(:, 2), d(2, 2) * B * v(:, 2), sqrt (eps));

%!test
%! A = single ([1, 2; -1, 1]);  B = single ([3, 3; 1, 2]);
%! [v, d] = eig (A, B);
%! assert (A * v(:, 1), d(1, 1) * B * v(:, 1), sqrt (eps ("single")));
%! assert (A * v(:, 2), d(2, 2) * B * v(:, 2), sqrt (eps ("single")));

%!test
%! A = [1, 2; 2, 1];  B = [3, -2; -2, 3];
%! [v, d] = eig (A, B);
%! assert (A * v(:, 1), d(1, 1) * B * v(:, 1), sqrt (eps));
%! assert (A * v(:, 2), d(2, 2) * B * v(:, 2), sqrt (eps));

%!test
%! A = single ([1, 2; 2, 1]);  B = single ([3, -2; -2, 3]);
%! [v, d] = eig (A, B);
%! assert (A * v(:, 1), d(1, 1) * B * v(:, 1), sqrt (eps ("single")));
%! assert (A * v(:, 2), d(2, 2) * B * v(:, 2), sqrt (eps ("single")));

%!test
%! A = [1+3i, 2+i; 2-i, 1+3i];  B = [5+9i, 2+i; 2-i, 5+9i];
%! [v, d] = eig (A, B);
%! assert (A * v(:, 1), d(1, 1) * B * v(:, 1), sqrt (eps));
%! assert (A * v(:, 2), d(2, 2) * B * v(:, 2), sqrt (eps));

%!test
%! A = single ([1+3i, 2+i; 2-i, 1+3i]);  B = single ([5+9i, 2+i; 2-i, 5+9i]);
%! [v, d] = eig (A, B);
%! assert (A * v(:, 1), d(1, 1) * B * v(:, 1), sqrt (eps ("single")));
%! assert (A * v(:, 2), d(2, 2) * B * v(:, 2), sqrt (eps ("single")));

%!test
%! A = [1+3i, 2+3i; 3-8i, 8+3i];  B = [8+i, 3+i; 4-9i, 3+i];
%! [v, d] = eig (A, B);
%! assert (A * v(:, 1), d(1, 1) * B * v(:, 1), sqrt (eps));
%! assert (A * v(:, 2), d(2, 2) * B * v(:, 2), sqrt (eps));

%!test
%! A = single ([1+3i, 2+3i; 3-8i, 8+3i]);  B = single ([8+i, 3+i; 4-9i, 3+i]);
%! [v, d] = eig (A, B);
%! assert (A * v(:, 1), d(1, 1) * B * v(:, 1), sqrt (eps ("single")));
%! assert (A * v(:, 2), d(2, 2) * B * v(:, 2), sqrt (eps ("single")));

%!test
%! A = [1, 2; 3, 8];  B = [8, 3; 4, 3];
%! [v, d] = eig (A, B);
%! assert (A * v(:, 1), d(1, 1) * B * v(:, 1), sqrt (eps));
%! assert (A * v(:, 2), d(2, 2) * B * v(:, 2), sqrt (eps));

%!error eig ()
%!error eig ([1, 2; 3, 4], [4, 3; 2, 1], 1)
%!error <EIG requires same size matrices> eig ([1, 2; 3, 4], 2)
%!error <argument must be a square matrix> eig ([1, 2; 3, 4; 5, 6])
%!error <wrong type argument> eig ("abcd")
%!error <wrong type argument> eig ([1 2 ; 2 3], "abcd")
%!error <wrong type argument> eig (false, [1 2 ; 2 3])
*/