changeset 718:e81d3a66725e

[project @ 1994-09-21 14:58:18 by jwe]
author jwe
date Wed, 21 Sep 1994 14:58:18 +0000
parents b14a2dda50e7
children 02814aea28c8
files src/balance.cc src/chol.cc src/dassl.cc src/det.cc src/eig.cc src/expm.cc src/fsolve.cc src/hess.cc src/inv.cc src/log.cc src/lsode.cc src/lu.cc src/npsol.cc src/pinv.cc src/qr.cc src/quad.cc src/qzval.cc src/schur.cc src/svd.cc src/syl.cc
diffstat 20 files changed, 101 insertions(+), 45 deletions(-) [+]
line wrap: on
line diff
--- a/src/balance.cc
+++ b/src/balance.cc
@@ -96,8 +96,12 @@
 
 // Check argument 1 dimensions.
 
-  if (empty_arg ("balance", a_nr, a_nc) < 0)
+  int arg_is_empty = empty_arg ("balance", a_nr, a_nc);
+
+  if (arg_is_empty < 0)
     return retval;
+  if (arg_is_empty > 0)
+    return Octave_object (2, Matrix ());
 
   if (a_nr != a_nc)
     {
--- a/src/chol.cc
+++ b/src/chol.cc
@@ -54,8 +54,12 @@
   int nr = arg.rows ();
   int nc = arg.columns ();
 
-  if (empty_arg ("chol", nr, nc) < 0)
+  int arg_is_empty = empty_arg ("chol", nr, nc);
+
+  if (arg_is_empty < 0)
     return retval;
+  if (arg_is_empty > 0)
+    return Matrix ();
 
   if (arg.is_real_type ())
     {
--- a/src/dassl.cc
+++ b/src/dassl.cc
@@ -103,7 +103,7 @@
   return retval;
 }
 
-DEFUN_DLD_BUILTIN ("dassl", Fdassl, Sdassl, 5, 1,
+DEFUN_DLD_BUILTIN ("dassl", Fdassl, Sdassl, 5, 2,
   "dassl (\"function_name\", x_0, xdot_0, t_out)\n\
 dassl (F, X_0, XDOT_0, T_OUT, T_CRIT)\n\
 \n\
@@ -125,7 +125,7 @@
     }
 
   dassl_fcn = is_valid_function (args(0), "dassl", 1);
-  if (! dassl_fcn || takes_correct_nargs (dassl_fcn, 4, "dassl", 1) != 1)
+  if (! dassl_fcn || takes_correct_nargs (dassl_fcn, 3, "dassl", 1) != 1)
     return retval;
 
   ColumnVector state = args(1).vector_value ();
--- a/src/det.cc
+++ b/src/det.cc
@@ -60,8 +60,11 @@
       return retval;
     }
 
-  if (empty_arg ("det", nr, nc) < 0)
+  int arg_is_empty = empty_arg ("det", nr, nc);
+  if (arg_is_empty < 0)
     return retval;
+  if (arg_is_empty > 0)
+    return Matrix (1, 1, 1.0);
 
   if (nr != nc)
     {
--- a/src/eig.cc
+++ b/src/eig.cc
@@ -53,8 +53,11 @@
   int nr = arg.rows ();
   int nc = arg.columns ();
 
-  if (empty_arg ("eig", nr, nc) < 0)
+  int arg_is_empty = empty_arg ("eig", nr, nc);
+  if (arg_is_empty < 0)
     return retval;
+  else if (arg_is_empty > 0)
+    return Octave_object (2, Matrix ());
 
   if (nr != nc)
     {
--- a/src/expm.cc
+++ b/src/expm.cc
@@ -85,8 +85,12 @@
   int nr = arg.rows ();
   int nc = arg.columns ();
 
-  if (empty_arg ("expm", nr, nc) < 0)
+  int arg_is_empty = empty_arg ("expm", nr, nc);
+
+  if (arg_is_empty < 0)
     return retval;
+  if (arg_is_empty > 0)
+    return Matrix ();
 
   if (nr != nc)
     {
--- a/src/fsolve.cc
+++ b/src/fsolve.cc
@@ -46,6 +46,8 @@
 int
 hybrd_info_to_fsolve_info (int info)
 {
+  cerr << "info: " << info << "\n";
+
   switch (info)
     {
     case -1:
@@ -119,7 +121,7 @@
   return retval;
 }
 
-DEFUN_DLD_BUILTIN ("fsolve", Ffsolve, Sfsolve, 5, 1,
+DEFUN_DLD_BUILTIN ("fsolve", Ffsolve, Sfsolve, 2, 1,
   "Solve nonlinear equations using Minpack.  Usage:\n\
 \n\
   [X, INFO] = fsolve (F, X0)\n\
@@ -135,14 +137,14 @@
 
   int nargin = args.length ();
 
-  if (nargin < 2 || nargin > 6 || nargout > 3)
+  if (nargin != 2 || nargout > 3)
     {
       print_usage ("fsolve");
       return retval;
     }
 
   fsolve_fcn = is_valid_function (args(0), "fsolve", 1);
-  if (! fsolve_fcn || takes_correct_nargs (fsolve_fcn, 2, "fsolve", 1) != 1)
+  if (! fsolve_fcn || takes_correct_nargs (fsolve_fcn, 1, "fsolve", 1) != 1)
     return retval;
 
   ColumnVector x = args(1).vector_value ();
--- a/src/hess.cc
+++ b/src/hess.cc
@@ -54,8 +54,12 @@
   int nr = arg.rows ();
   int nc = arg.columns ();
 
-  if (empty_arg ("hess", nr, nc) < 0)
+  int arg_is_empty = empty_arg ("hess", nr, nc);
+
+  if (arg_is_empty < 0)
     return retval;
+  else if (arg_is_empty > 0)
+    return Octave_object (2, Matrix ());
 
   if (nr != nc)
     {
--- a/src/inv.cc
+++ b/src/inv.cc
@@ -54,8 +54,12 @@
   int nr = arg.rows ();
   int nc = arg.columns ();
 
-  if (empty_arg ("inverse", nr, nc) < 0)
+  int arg_is_empty = empty_arg ("inverse", nr, nc);
+
+  if (arg_is_empty < 0)
     return retval;
+  else if (arg_is_empty > 0)
+    return Matrix ();
 
   if (nr != nc)
     {
--- a/src/log.cc
+++ b/src/log.cc
@@ -51,9 +51,13 @@
     }
 
   tree_constant arg = args(0);
-    
-  if (empty_arg ("logm", arg.rows (), arg.columns ()) < 0)
+
+  int arg_is_empty = empty_arg ("logm", arg.rows (), arg.columns ());
+
+  if (arg_is_empty < 0)
     return retval;
+  else if (arg_is_empty > 0)
+    return Matrix ();
 
   if (arg.is_real_scalar ())
     {
@@ -159,9 +163,13 @@
     }
 
   tree_constant arg = args(0);
-    
-  if (empty_arg ("sqrtm", arg.rows (), arg.columns ()))
+
+  int arg_is_empty = empty_arg ("sqrtm", arg.rows (), arg.columns ());
+
+  if (arg_is_empty < 0)
     return retval;
+  else if (arg_is_empty > 0)
+    return Matrix ();
 
   if (arg.is_real_scalar ())
     {
--- a/src/lsode.cc
+++ b/src/lsode.cc
@@ -92,7 +92,7 @@
   return retval;
 }
 
-DEFUN_DLD_BUILTIN ("lsode", Flsode, Slsode, 6, 1,
+DEFUN_DLD_BUILTIN ("lsode", Flsode, Slsode, 4, 1,
   "lsode (F, X0, T_OUT, T_CRIT)\n\
 \n\
 The first argument is the name of the function to call to\n\
@@ -113,7 +113,7 @@
     }
 
   lsode_fcn = is_valid_function (args(0), "lsode", 1);
-  if (! lsode_fcn || takes_correct_nargs (lsode_fcn, 3, "lsode", 1) != 1)
+  if (! lsode_fcn || takes_correct_nargs (lsode_fcn, 2, "lsode", 1) != 1)
     return retval;
 
   ColumnVector state = args(1).vector_value ();
--- a/src/lu.cc
+++ b/src/lu.cc
@@ -53,8 +53,12 @@
   int nr = arg.rows ();
   int nc = arg.columns ();
 
-  if (empty_arg ("lu", nr, nc) < 0)
+  int arg_is_empty = empty_arg ("lu", nr, nc);
+
+  if (arg_is_empty < 0)
     return retval;
+  else if (arg_is_empty > 0)
+    return Octave_object (3, Matrix ());
 
   if (nr != nc)
     {
--- a/src/npsol.cc
+++ b/src/npsol.cc
@@ -65,10 +65,8 @@
       decision_vars = d;
     }
 
-//  tree_constant name = npsol_objective->name ();
   Octave_object args;
-  args(1) = decision_vars;
-//  args(0) = name;
+  args(0) = decision_vars;
 
   static double retval;
   retval = 0.0;
@@ -140,10 +138,8 @@
       decision_vars = d;
     }
 
-//  tree_constant name = npsol_constraints->name ();
   Octave_object args;
-  args(1) = decision_vars;
-//  args(0) = name;
+  args(0) = decision_vars;
 
   if (npsol_constraints)
     {
@@ -232,12 +228,12 @@
 #endif
 
 #if defined (NPSOL_MISSING)
-DEFUN_DLD_BUILTIN ("npsol", Fnpsol, Snpsol, 11, 3,
+DEFUN_DLD_BUILTIN ("npsol", Fnpsol, Snpsol, 10, 4,
   "This function requires NPSOL, which is not freely\n\
 redistributable.  For more information, read the file\n\
 libcruft/npsol/README.MISSING in the source distribution.")
 #else
-DEFUN_DLD_BUILTIN ("npsol", Fnpsol, Snpsol, 11, 3,
+DEFUN_DLD_BUILTIN ("npsol", Fnpsol, Snpsol, 10, 4,
   "[X, OBJ, INFO, LAMBDA] = npsol (X, PHI [, LB, UB] [, LB, A, UB] [, LB, G, UB])\n\
 \n\
 Groups of arguments surrounded in `[]' are optional, but\n\
@@ -306,7 +302,7 @@
 
   npsol_objective = is_valid_function (args(1), "npsol", 1);
   if (! npsol_objective
-      || takes_correct_nargs (npsol_objective, 2, "npsol", 1) != 1)
+      || takes_correct_nargs (npsol_objective, 1, "npsol", 1) != 1)
     return retval;
 
   Objective func (npsol_objective_function);
@@ -410,7 +406,7 @@
 	}
       else
 	{
-	  if (takes_correct_nargs (npsol_constraints, 2, "npsol", 1))
+	  if (takes_correct_nargs (npsol_constraints, 1, "npsol", 1))
 	    {
 	      ColumnVector nlub = args(nargin).vector_value ();
 	      ColumnVector nllb = args(nargin-2).vector_value ();
@@ -453,7 +449,7 @@
 	}
       else
 	{
-	  if (takes_correct_nargs (npsol_constraints, 2, "npsol", 1))
+	  if (takes_correct_nargs (npsol_constraints, 1, "npsol", 1))
 	    {
 	      ColumnVector nlub = args(nargin).vector_value ();
 	      ColumnVector nllb = args(nargin-2).vector_value ();
--- a/src/pinv.cc
+++ b/src/pinv.cc
@@ -69,7 +69,7 @@
 
   if (arg_is_empty < 0)
     return retval;
-  else if (arg_is_empty)
+  else if (arg_is_empty > 0)
     return Matrix ();
 
   if (arg.is_real_type ())
--- a/src/qr.cc
+++ b/src/qr.cc
@@ -67,9 +67,13 @@
     }
 
   tree_constant arg = args(0);
-    
-  if (empty_arg ("qr", arg.rows (), arg.columns ()) < 0)
+
+  int arg_is_empty = empty_arg ("qr", arg.rows (), arg.columns ());
+
+  if (arg_is_empty < 0)
     return retval;
+  else if (arg_is_empty > 0)
+    return Octave_object (3, Matrix ());
 
   QR::type type = nargout == 1 ? QR::raw
     : (nargin == 2 ? QR::economy : QR::std);
--- a/src/quad.cc
+++ b/src/quad.cc
@@ -83,7 +83,7 @@
   return retval;
 }
 
-DEFUN_DLD_BUILTIN ("quad", Fquad, Squad, 6, 3,
+DEFUN_DLD_BUILTIN ("quad", Fquad, Squad, 5, 3,
   "[V, IER, NFUN] = quad (F, A, B [, TOL] [, SING])\n\
 \n\
 Where the first argument is the name of the  function to call to\n\
@@ -108,8 +108,8 @@
       return retval;
     }
 
-  quad_fcn = is_valid_function (args(0), "fsolve", 1);
-  if (! quad_fcn || takes_correct_nargs (quad_fcn, 2, "fsolve", 1) != 1)
+  quad_fcn = is_valid_function (args(0), "quad", 1);
+  if (! quad_fcn || takes_correct_nargs (quad_fcn, 1, "quad", 1) != 1)
     return retval;
 
   double a = args(1).double_value ();
--- a/src/qzval.cc
+++ b/src/qzval.cc
@@ -78,9 +78,13 @@
 
   int b_nr = arg_b.rows();
   int b_nc = arg_b.columns();
-  
-  if (empty_arg ("qzvalue", a_nr, a_nc) < 0
-      || empty_arg ("qzvalue", b_nr, b_nc) < 0)
+
+  int arg_a_is_empty = empty_arg ("qzvalue", a_nr, a_nc);
+  int arg_b_is_empty = empty_arg ("qzvalue", b_nr, b_nc);
+
+  if (arg_a_is_empty > 0 && arg_b_is_empty > 0)
+    return Matrix ();
+  else if (arg_a_is_empty || arg_b_is_empty)
     return retval;
 
 // Arguments are not empty, so check for correct dimensions.
--- a/src/schur.cc
+++ b/src/schur.cc
@@ -83,8 +83,12 @@
   int nr = arg.rows ();
   int nc = arg.columns ();
 
-  if (empty_arg ("schur", nr, nc) < 0)
+  int arg_is_empty = empty_arg ("schur", nr, nc);
+
+  if (arg_is_empty < 0)
     return retval;
+  else if (arg_is_empty > 0)
+    return Octave_object (2, Matrix ());
 
   if (nr != nc)
     {
--- a/src/svd.cc
+++ b/src/svd.cc
@@ -55,8 +55,12 @@
 
   tree_constant arg = args(0);
 
-  if (empty_arg ("svd", arg.rows (), arg.columns ()) < 0)
+  int arg_is_empty = empty_arg ("svd", arg.rows (), arg.columns ());
+
+  if (arg_is_empty < 0)
     return retval;
+  else if (arg_is_empty > 0)
+    return Octave_object (3, Matrix ());
 
   SVD::type type = (nargin == 2) ? SVD::economy : SVD::std;
 
--- a/src/syl.cc
+++ b/src/syl.cc
@@ -81,10 +81,14 @@
 
   int c_nr = arg_c.rows ();
   int c_nc = arg_c.columns ();
-  
-  if (empty_arg ("syl", a_nr, a_nc) < 0
-      || empty_arg ("syl", b_nr, b_nc) < 0
-      || empty_arg ("syl", c_nr, c_nc) < 0)
+
+  int arg_a_is_empty = empty_arg ("syl", a_nr, a_nc);
+  int arg_b_is_empty = empty_arg ("syl", b_nr, b_nc);
+  int arg_c_is_empty = empty_arg ("syl", c_nr, c_nc);
+
+  if (arg_a_is_empty > 0 && arg_b_is_empty > 0 && arg_c_is_empty > 0)
+    return Matrix ();
+  else if (arg_a_is_empty || arg_b_is_empty || arg_c_is_empty)
     return retval;
 
 // Arguments are not empty, so check for correct dimensions.