changeset 10607:f7501986e42d

make schur more Matlab compatible
author Jaroslav Hajek <highegg@gmail.com>
date Thu, 06 May 2010 09:49:36 +0200
parents ec34c7acd057
children f9860b622680
files liboctave/ChangeLog liboctave/CmplxSCHUR.cc liboctave/dbleSCHUR.cc liboctave/fCmplxSCHUR.cc liboctave/floatSCHUR.cc src/ChangeLog src/DLD-FUNCTIONS/schur.cc
diffstat 7 files changed, 68 insertions(+), 26 deletions(-) [+]
line wrap: on
line diff
--- a/liboctave/ChangeLog	Wed May 05 21:14:18 2010 -0700
+++ b/liboctave/ChangeLog	Thu May 06 09:49:36 2010 +0200
@@ -1,3 +1,11 @@
+2010-05-06  Jaroslav Hajek  <highegg@gmail.com>
+
+	* dbleSCHUR.cc (SCHUR::init): Handle empty matrix case.
+	Use clear rather than resize to realloc matrix.
+	* floatSCHUR.cc (FloatSCHUR::init): Ditto.
+	* CmplxSCHUR.cc (ComplexSCHUR::init): Ditto.
+	* fCmplxSCHUR.cc (FloatComplexSCHUR::init): Ditto.
+
 2010-05-04  John W. Eaton  <jwe@octave.org>
 
 	* CollocWt.cc (diff, jcobi, dfopr): New functions, based on
--- a/liboctave/CmplxSCHUR.cc	Wed May 05 21:14:18 2010 -0700
+++ b/liboctave/CmplxSCHUR.cc	Thu May 06 09:49:36 2010 +0200
@@ -70,6 +70,12 @@
         ("ComplexSCHUR requires square matrix");
       return -1;
     }
+  else if (a_nr == 0)
+    {
+      schur_mat.clear ();
+      unitary_mat.clear ();
+      return 0;
+    }
 
   // Workspace requirements may need to be fixed if any of the
   // following change.
@@ -104,7 +110,7 @@
 
   schur_mat = a;
   if (calc_unitary)
-    unitary_mat.resize (n, n);
+    unitary_mat.clear (n, n);
 
   Complex *s = schur_mat.fortran_vec ();
   Complex *q = unitary_mat.fortran_vec ();
--- a/liboctave/dbleSCHUR.cc	Wed May 05 21:14:18 2010 -0700
+++ b/liboctave/dbleSCHUR.cc	Thu May 06 09:49:36 2010 +0200
@@ -70,6 +70,12 @@
       (*current_liboctave_error_handler) ("SCHUR requires square matrix");
       return -1;
     }
+  else if (a_nr == 0)
+    {
+      schur_mat.clear ();
+      unitary_mat.clear ();
+      return 0;
+    }
 
   // Workspace requirements may need to be fixed if any of the
   // following change.
@@ -106,7 +112,7 @@
   schur_mat = a;
 
   if (calc_unitary)
-    unitary_mat.resize (n, n);
+    unitary_mat.clear (n, n);
 
   double *s = schur_mat.fortran_vec ();
   double *q = unitary_mat.fortran_vec ();
--- a/liboctave/fCmplxSCHUR.cc	Wed May 05 21:14:18 2010 -0700
+++ b/liboctave/fCmplxSCHUR.cc	Thu May 06 09:49:36 2010 +0200
@@ -70,6 +70,12 @@
         ("FloatComplexSCHUR requires square matrix");
       return -1;
     }
+  else if (a_nr == 0)
+    {
+      schur_mat.clear ();
+      unitary_mat.clear ();
+      return 0;
+    }
 
   // Workspace requirements may need to be fixed if any of the
   // following change.
@@ -104,7 +110,7 @@
 
   schur_mat = a;
   if (calc_unitary)
-    unitary_mat.resize (n, n);
+    unitary_mat.clear (n, n);
 
   FloatComplex *s = schur_mat.fortran_vec ();
   FloatComplex *q = unitary_mat.fortran_vec ();
--- a/liboctave/floatSCHUR.cc	Wed May 05 21:14:18 2010 -0700
+++ b/liboctave/floatSCHUR.cc	Thu May 06 09:49:36 2010 +0200
@@ -70,6 +70,12 @@
       (*current_liboctave_error_handler) ("FloatSCHUR requires square matrix");
       return -1;
     }
+  else if (a_nr == 0)
+    {
+      schur_mat.clear ();
+      unitary_mat.clear ();
+      return 0;
+    }
 
   // Workspace requirements may need to be fixed if any of the
   // following change.
@@ -106,7 +112,7 @@
   schur_mat = a;
 
   if (calc_unitary)
-    unitary_mat.resize (n, n);
+    unitary_mat.clear (n, n);
 
   float *s = schur_mat.fortran_vec ();
   float *q = unitary_mat.fortran_vec ();
--- a/src/ChangeLog	Wed May 05 21:14:18 2010 -0700
+++ b/src/ChangeLog	Thu May 06 09:49:36 2010 +0200
@@ -1,3 +1,9 @@
+2010-05-06  Jaroslav Hajek  <highegg@gmail.com>
+
+	* DLD-FUNCTIONS/schur.cc (Fschur): Recognize "complex" option for
+	Matlab compatibility. Simplify argument handling and improve error
+	messages.
+
 2010-05-05  John W. Eaton  <jwe@octave.org>
 
 	* utils.cc (reset_warning_state): New function.
--- a/src/DLD-FUNCTIONS/schur.cc	Wed May 05 21:14:18 2010 -0700
+++ b/src/DLD-FUNCTIONS/schur.cc	Thu May 06 09:49:36 2010 +0200
@@ -41,6 +41,7 @@
 DEFUN_DLD (schur, args, nargout,
   "-*- texinfo -*-\n\
 @deftypefn {Loadable Function} {@var{s} =} schur (@var{a})\n\
+@deftypefnx {Loadable Function} {@var{s} =} schur (@var{a}, \"complex\")\n\
 @deftypefnx {Loadable Function} {[@var{u}, @var{s}] =} schur (@var{a}, @var{opt})\n\
 @cindex Schur decomposition\n\
 The Schur decomposition is used to compute eigenvalues of a\n\
@@ -148,6 +149,8 @@
 @code{s}.\n\
 @end ifnottex\n\
 \n\
+A complex decomposition may be forced by passing \"complex\" as @var{opt}.\n\
+\n\
 The eigenvalues are optionally ordered along the diagonal according to\n\
 the value of @code{opt}.  @code{opt = \"a\"} indicates that all\n\
 eigenvalues with negative real parts should be moved to the leading\n\
@@ -229,35 +232,40 @@
         }
     }
 
-  char ord_char = ord.empty () ? 'U' : ord[0];
+  bool force_complex = false;
 
-  if (ord_char != 'U' && ord_char != 'A' && ord_char != 'D'
-      && ord_char != 'u' && ord_char != 'a' && ord_char != 'd')
+  if (ord == "complex")
+    {
+      force_complex = true;
+      ord = std::string ();
+    }
+  else
     {
-      warning ("schur: incorrect ordered schur argument `%c'",
-               ord.c_str ());
-      return retval;
+      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 `%c'",
+                   ord.c_str ());
+          return retval;
+        }
     }
 
   octave_idx_type nr = arg.rows ();
   octave_idx_type nc = arg.columns ();
 
-  int arg_is_empty = empty_arg ("schur", nr, nc);
-
-  if (arg_is_empty < 0)
-    return retval;
-  else if (arg_is_empty > 0)
-    return octave_value_list (2, Matrix ());
-
   if (nr != nc)
     {
       gripe_square_matrix_required ("schur");
       return retval;
     }
 
-  if (arg.is_single_type ())
+  if (! arg.is_numeric_type ())
+    gripe_wrong_type_arg ("schur", arg);
+  else if (arg.is_single_type ())
     {
-      if (arg.is_real_type ())
+      if (! force_complex && arg.is_real_type ())
         {
           FloatMatrix tmp = arg.float_matrix_value ();
 
@@ -276,7 +284,7 @@
                 }
             }
         }
-      else if (arg.is_complex_type ())
+      else
         {
           FloatComplexMatrix ctmp = arg.float_complex_matrix_value ();
 
@@ -299,7 +307,7 @@
     }
   else
     {
-      if (arg.is_real_type ())
+      if (! force_complex && arg.is_real_type ())
         {
           Matrix tmp = arg.matrix_value ();
 
@@ -318,7 +326,7 @@
                 }
             }
         }
-      else if (arg.is_complex_type ())
+      else
         {
           ComplexMatrix ctmp = arg.complex_matrix_value ();
 
@@ -338,10 +346,6 @@
                 }
             }
         }
-      else
-        {
-          gripe_wrong_type_arg ("schur", arg);
-        }
     }
  
   return retval;