changeset 13884:22e16fd68b8a gui

Merge default onto gui
author Jordi Gutiérrez Hermoso <jordigh@octave.org>
date Thu, 17 Nov 2011 21:58:56 -0500
parents b09321ab8464 (current diff) ecf0c6bca0c9 (diff)
children 36f90899e058
files
diffstat 16 files changed, 264 insertions(+), 212 deletions(-) [+]
line wrap: on
line diff
--- /dev/null	Thu Jan 01 00:00:00 1970 +0000
+++ b/.hgsub	Thu Nov 17 21:58:56 2011 -0500
@@ -0,0 +1,1 @@
+gnulib = [git]git://git.sv.gnu.org/gnulib
--- /dev/null	Thu Jan 01 00:00:00 1970 +0000
+++ b/.hgsubstate	Thu Nov 17 21:58:56 2011 -0500
@@ -0,0 +1,1 @@
+84c3f9cf54eaa688c5a1a26925886535079a91de gnulib
--- a/build-aux/bootstrap	Thu Nov 17 23:04:02 2011 +0000
+++ b/build-aux/bootstrap	Thu Nov 17 21:58:56 2011 -0500
@@ -169,8 +169,10 @@
 # default.
 bootstrap_sync=false
 
-# Use git to update gnulib sources
-use_git=true
+# Don't use git to update gnulib sources. We keep gnulib under a
+# Mercurial subrepository instead
+use_git=false
+GNULIB_SRCDIR=gnulib/
 
 # find_tool ENVVAR NAMES...
 # -------------------------
--- a/scripts/geometry/voronoi.m	Thu Nov 17 23:04:02 2011 +0000
+++ b/scripts/geometry/voronoi.m	Thu Nov 17 21:58:56 2011 -0500
@@ -70,8 +70,6 @@
 
   if (nargin < 1)
     print_usage ();
-  elseif (nargout > 2)
-    error ("voronoi: No more than two output arguments supported");
   endif
 
   narg = 1;
@@ -131,7 +129,9 @@
   ybox = [xmin - scale * xdelta; xmax + scale * xdelta; ...
           xmax + scale * xdelta; xmin - scale * xdelta];
 
-  [p, c, infi] = __voronoi__ ([[x(:) ; xbox(:)], [y(:); ybox(:)]], opts{:});
+  [p, c, infi] = __voronoi__ ("voronoi",
+                              [[x(:) ; xbox(:)], [y(:); ybox(:)]],
+                              opts{:});
 
   idx = find (! infi);
   ll = length (idx);
--- a/scripts/geometry/voronoin.m	Thu Nov 17 23:04:02 2011 +0000
+++ b/scripts/geometry/voronoin.m	Thu Nov 17 21:58:56 2011 -0500
@@ -48,14 +48,14 @@
 
   if (np <= dim)
     error ("voronoin: number of points must be greater than their dimension");
-  elseif (nargin == 2 && ! (ischar (options) || iscellstr (options)))
-    error ("voronoin: OPTIONS argument must be a string or cell array of strings");
   endif
 
+  caller = "voronoin";
+
   if (nargin == 1)
-    [C, F] = __voronoi__ (pts);
+    [C, F] = __voronoi__ (caller, pts);
   else
-    [C, F] = __voronoi__ (pts, options);
+    [C, F] = __voronoi__ (caller, pts, options);
   endif
 
 endfunction
--- a/scripts/special-matrix/hilb.m	Thu Nov 17 23:04:02 2011 +0000
+++ b/scripts/special-matrix/hilb.m	Thu Nov 17 21:58:56 2011 -0500
@@ -54,28 +54,26 @@
 
 function retval = hilb (n)
 
-
   if (nargin != 1)
     print_usage ();
+  elseif (! isscalar (n))
+    error ("hilb: N must be a scalar integer");
   endif
 
-  nmax = length (n);
-  if (nmax == 1)
-    retval = zeros (n);
-    tmp = 1:n;
-    for i = 1:n
-      retval (i, :) = 1.0 ./ (tmp + (i - 1));
-    endfor
-  else
-    error ("hilb: expecting scalar argument, found something else");
-  endif
+  retval = zeros (n);
+  tmp = 1:n;
+  for i = 1:n
+    retval(i, :) = 1.0 ./ tmp;
+    tmp++;
+  endfor
 
 endfunction
 
-%!assert((hilb (2) == [1, 1/2; 1/2, 1/3]
-%! && hilb (3) == [1, 1/2, 1/3; 1/2, 1/3, 1/4; 1/3, 1/4, 1/5]));
+
+%!assert (hilb (2), [1, 1/2; 1/2, 1/3])
+%!assert (hilb (3), [1, 1/2, 1/3; 1/2, 1/3, 1/4; 1/3, 1/4, 1/5])
 
-%!error hilb ();
+%!error hilb ()
+%!error hilb (1, 2)
+%!error <N must be a scalar integer> hilb (ones(2))
 
-%!error hilb (1, 2);
-
--- a/scripts/special-matrix/vander.m	Thu Nov 17 23:04:02 2011 +0000
+++ b/scripts/special-matrix/vander.m	Thu Nov 17 21:58:56 2011 -0500
@@ -60,37 +60,36 @@
     print_usage ();
   endif
 
-  if (isvector (c))
-    retval = zeros (length (c), n, class (c));
-    ## avoiding many ^s appears to be faster for n >= 100.
-    d = 1;
-    c = c(:);
-    for i = n:-1:1
-      retval(:,i) = d;
-      d = c .* d;
-    endfor
-  else
-    error ("vander: argument must be a vector");
+  if (! isvector (c))
+    error ("vander: polynomial C must be a vector");
   endif
 
+  ## avoiding many ^s appears to be faster for n >= 100.
+  retval = zeros (length (c), n, class (c));
+  d = 1;
+  c = c(:);
+  for i = n:-1:1
+    retval(:,i) = d;
+    d .*= c;
+  endfor
+
 endfunction
 
+
 %!test
 %! c = [0,1,2,3];
 %! expect = [0,0,0,1; 1,1,1,1; 8,4,2,1; 27,9,3,1];
-%! result = vander(c);
-%! assert(expect, result);
+%! assert(vander (c), expect);
 
-%!assert((vander (1) == 1 && vander ([1, 2, 3]) == vander ([1; 2; 3])
-%! && vander ([1, 2, 3]) == [1, 1, 1; 4, 2, 1; 9, 3, 1]
-%! && vander ([1, 2, 3]*i) == [-1, i, 1; -4, 2i, 1; -9, 3i, 1]));
+%!assert (vander (1), 1)
+%!assert (vander ([1, 2, 3]), vander ([1; 2; 3]))
+%!assert (vander ([1, 2, 3]), [1, 1, 1; 4, 2, 1; 9, 3, 1])
+%!assert (vander ([1, 2, 3]*i), [-1, i, 1; -4, 2i, 1; -9, 3i, 1])
 
 %!assert(vander (2, 3), [4, 2, 1])
 %!assert(vander ([2, 3], 3), [4, 2, 1; 9, 3, 1])
 
-%!error vander ([1, 2; 3, 4]);
+%!error vander ();
+%!error vander (1, 2, 3);
+%!error <polynomial C must be a vector> vander ([1, 2; 3, 4]);
 
-%!error vander ();
-
-%!error vander (1, 2, 3);
-
--- a/src/DLD-FUNCTIONS/__voronoi__.cc	Thu Nov 17 23:04:02 2011 +0000
+++ b/src/DLD-FUNCTIONS/__voronoi__.cc	Thu Nov 17 21:58:56 2011 -0500
@@ -57,8 +57,8 @@
 
 DEFUN_DLD (__voronoi__, args, ,
         "-*- texinfo -*-\n\
-@deftypefn  {Loadable Function} {@var{C}, @var{F} =} __voronoi__ (@var{pts})\n\
-@deftypefnx {Loadable Function} {@var{C}, @var{F} =} __voronoi__ (@var{pts}, @var{options})\n\
+@deftypefn  {Loadable Function} {@var{C}, @var{F} =} __voronoi__ (@var{caller}, @var{pts})\n\
+@deftypefnx {Loadable Function} {@var{C}, @var{F} =} __voronoi__ (@var{caller}, @var{pts}, @var{options})\n\
 @deftypefnx {Loadable Function} {@var{C}, @var{F}, @var{Inf_Pts} =} __voronoi__ (@dots{})\n\
 Undocumented internal function.\n\
 @end deftypefn")
@@ -70,107 +70,125 @@
   retval(0) = 0.0;
 
   int nargin = args.length ();
-  if (nargin < 1 || nargin > 2)
+  if (nargin < 2 || nargin > 3)
     {
       print_usage ();
       return retval;
     }
 
-  std::string options = "";
+  std::string caller = args(0).string_value ();
+
+  Matrix points = args(1).matrix_value ();
+  const octave_idx_type dim = points.columns ();
+  const octave_idx_type num_points = points.rows ();
+
+  points = points.transpose ();
 
-  if (nargin == 2)
+  std::string options;
+
+  if (dim <= 4)
+    options = " Qbb";
+  else
+    options = " Qbb Qx";
+
+  if (nargin == 3)
     {
-      if (args(1).is_string ())
-        options = args(1).string_value ();
-      else if (args(1).is_empty ())
-        ;  // Use default options
-      else if (args(1).is_cellstr ())
+      octave_value opt_arg = args(2);
+
+      if (opt_arg.is_string ())
+        options = " " + opt_arg.string_value ();
+      else if (opt_arg.is_empty ())
+        ; // Use default options.
+      else if (opt_arg.is_cellstr ())
         {
           options = "";
-          Array<std::string> tmp = args(1).cellstr_value ();
+
+          Array<std::string> tmp = opt_arg.cellstr_value ();
 
           for (octave_idx_type i = 0; i < tmp.numel (); i++)
-            options += tmp(i) + " ";
+            options += " " + tmp(i);
         }
       else
         {
-          error ("__voronoi__: OPTIONS argument must be a string, cell array of strings, or empty");
+          error ("%s: OPTIONS must be a string, cell array of strings, or empty",
+                 caller.c_str ());
           return retval;
         }
     }
 
-  Matrix p (args(0).matrix_value ());
-  const octave_idx_type dim = p.columns ();
-  const octave_idx_type np = p.rows ();
-
-  p = p.transpose ();
-  double *pt_array = p.fortran_vec ();
-
   boolT ismalloc = false;
 
   // Replace the 0 pointer with stdout for debugging information
   FILE *outfile = 0;
   FILE *errfile = stderr;
 
-  // Qhull flags argument is not const char*
-  OCTAVE_LOCAL_BUFFER (char, flags, 12 + options.length ());
+  // qh_new_qhull command and points arguments are not const...
+
+  std::string cmd = "qhull v" + options;
 
-  sprintf (flags, "qhull v Fv %s", options.c_str ());
+  OCTAVE_LOCAL_BUFFER (char, cmd_str, cmd.length () + 1);
 
-  int exitcode = qh_new_qhull (dim, np, pt_array, 
-                               ismalloc, flags, outfile, errfile);
+  strcpy (cmd_str, cmd.c_str ());
+
+  int exitcode = qh_new_qhull (dim, num_points, points.fortran_vec (),
+                               ismalloc, cmd_str, outfile, errfile);
   if (! exitcode) 
     {
-      facetT *facet;
-      vertexT *vertex;
+      // Calling findgood_all provides the number of Voronoi vertices
+      // (sets qh num_good).
 
-      octave_idx_type i = 0, n = 1, k = 0, m = 0, j = 0, r = 0;
+      qh_findgood_all (qh facet_list);
 
-      // Count number of vertices for size of NI.  FIXME -- does Qhull
-      // have a way to query this value directly?
-      octave_idx_type nv = 0;
-      FORALLvertices
-        {
-          nv++;
-        }
+      octave_idx_type num_voronoi_regions
+        = qh num_vertices - qh_setsize (qh del_vertices);
 
-      OCTAVE_LOCAL_BUFFER (octave_idx_type, ni, nv);
+      octave_idx_type num_voronoi_vertices = qh num_good;
 
-      for (i = 0; i < nv; i++)
-        ni[i] = 0;
+      // Find the voronoi centers for all facets.
 
       qh_setvoronoi_all ();
 
-      bool infinity_seen = false;
-      facetT *neighbor, **neighborp;
-      coordT *voronoi_vertex;
+      facetT *facet;
+      vertexT *vertex;
+      octave_idx_type k;
+
+      // Find the number of Voronoi vertices for each Voronoi cell and
+      // store them in NI so we can use them later to set the dimensions
+      // of the RowVector objects used to collect them.
 
       FORALLfacets
         {
           facet->seen = false;
         }
+      
+      OCTAVE_LOCAL_BUFFER (octave_idx_type, ni, num_voronoi_regions);
+      for (octave_idx_type i = 0; i < num_voronoi_regions; i++)
+        ni[i] = 0;
+
+      k = 0;
 
       FORALLvertices
         {
           if (qh hull_dim == 3)
             qh_order_vertexneighbors (vertex);
           
-          infinity_seen = false;
+          bool infinity_seen = false;
+
+          facetT *neighbor, **neighborp;
 
           FOREACHneighbor_ (vertex)
             {
-              if (! neighbor->upperdelaunay)
+              if (neighbor->upperdelaunay)
                 {
-                  if (! neighbor->seen)
+                  if (! infinity_seen)
                     {
-                      n++;
-                      neighbor->seen = true;
+                      infinity_seen = true;
+                      ni[k]++;
                     }
-                  ni[k]++;
                 }
-              else if (! infinity_seen)
+              else
                 {
-                  infinity_seen = true;
+                  neighbor->seen = true;
                   ni[k]++;
                 }
             }
@@ -178,38 +196,65 @@
           k++;
         }
 
-      Matrix v (n, dim);
-      for (octave_idx_type d = 0; d < dim; d++)
-        v(0,d) = octave_Inf;
+      // If Qhull finds fewer regions than points, we will pad the end
+      // of the at_inf and C arrays so that they always contain at least
+      // as many elements as the given points array.
+
+      // FIXME -- is it possible (or does it make sense) for
+      // num_voronoi_regions to ever be larger than num_points?
+
+      octave_idx_type nr = (num_points > num_voronoi_regions
+                            ? num_points : num_voronoi_regions);
+
+      boolMatrix at_inf (nr, 1, false);
 
-      boolMatrix AtInf (nv, 1, false);
-      std::list<RowVector> F;
-      k = 0;
-      i = 0;
+      // The list of Voronoi vertices.  The first element is always
+      // Inf.
+      Matrix F (num_voronoi_vertices+1, dim);
+
+      for (octave_idx_type d = 0; d < dim; d++)
+        F(0,d) = octave_Inf;
+
+      // The cell array of vectors of indices into F that represent the
+      // vertices of the Voronoi regions (cells).
+
+      Cell C (nr, 1);
+
+      // Now loop through the list of vertices again and store the
+      // coordinates of the Voronoi vertices and the lists of indices
+      // for the cells.
 
       FORALLfacets
         {
           facet->seen = false;
         }
 
+      octave_idx_type i = 0;
+      k = 0;
+
       FORALLvertices
         {
           if (qh hull_dim == 3)
             qh_order_vertexneighbors (vertex);
 
-          infinity_seen = false;
+          bool infinity_seen = false;
+
+          octave_idx_type idx = qh_pointid (vertex->point);
 
-          octave_idx_type n_vertices = ni[k++];
+          octave_idx_type num_vertices = ni[k++];
 
-          // Qhull seems to sometimes produce "facets" with a single
-          // vertex.  Is that a bug?  How can a facet have just one
+          // Qhull seems to sometimes produces regions with a single
+          // vertex.  Is that a bug?  How can a region have just one
           // vertex?  Let's skip it.
 
-          if (n_vertices == 1)
+          if (num_vertices == 1)
             continue;
 
-          RowVector facet_list (n_vertices);
-          m = 0;
+          RowVector facet_list (num_vertices);
+
+          octave_idx_type m = 0;
+
+          facetT *neighbor, **neighborp;
 
           FOREACHneighbor_(vertex)
             {
@@ -219,62 +264,47 @@
                     {
                       infinity_seen = true;
                       facet_list(m++) = 1;
-                      AtInf(j) = true;
+                      at_inf(idx) = true;
                     }
                 }
               else
                 {
                   if (! neighbor->seen)
                     {
-                      voronoi_vertex = neighbor->center;
                       i++;
                       for (octave_idx_type d = 0; d < dim; d++)
-                        {
-                          v(i,d) = *(voronoi_vertex++);
-                        }
+                        F(i,d) = neighbor->center[d];
+
                       neighbor->seen = true;
                       neighbor->visitid = i;
                     }
+
                   facet_list(m++) = neighbor->visitid + 1;
                 }
             }
-          F.push_back (facet_list);
-          j++;
+
+          C(idx) = facet_list;
         }
 
-      // For compatibility with Matlab, pad the cell array of vertex
-      // lists with empty matrices if there are fewer facets than
-      // points.
-      octave_idx_type f_len = F.size ();
-      Cell C(np > f_len ? np : f_len, 1);
-
-      i = 0;
-      for (std::list<RowVector>::const_iterator it = F.begin ();
-           it != F.end (); it++)
-        C.elem (i++) = *it;
-
-      F.clear ();
-
-      AtInf.resize (f_len, 1);
-      retval(2) = AtInf;
+      retval(2) = at_inf;
       retval(1) = C;
-      retval(0) = v;
+      retval(0) = F;
     }
   else
-    error ("__voronoi__: qhull failed");
+    error ("%s: qhull failed", caller.c_str ());
 
-  // free memory from Qhull
+  // Free memory from Qhull
   qh_freeqhull (! qh_ALL);
 
   int curlong, totlong;
   qh_memfreeshort (&curlong, &totlong);
 
   if (curlong || totlong)
-    warning ("__voronoi__: did not free %d bytes of long memory (%d pieces)",
-             totlong, curlong);
+    warning ("%s: qhull did not free %d bytes of long memory (%d pieces)",
+             caller.c_str (), totlong, curlong);
 
 #else
-  error ("__voronoi__: not available in this version of Octave");
+  error ("%s: not available in this version of Octave", caller.c_str ());
 #endif
 
   return retval;
--- a/src/DLD-FUNCTIONS/convhulln.cc	Thu Nov 17 23:04:02 2011 +0000
+++ b/src/DLD-FUNCTIONS/convhulln.cc	Thu Nov 17 21:58:56 2011 -0500
@@ -94,30 +94,33 @@
       return retval;
     }
 
-  Matrix p (args(0).matrix_value ());
-  const octave_idx_type dim = p.columns ();
-  const octave_idx_type n = p.rows ();
+  Matrix points (args(0).matrix_value ());
+  const octave_idx_type dim = points.columns ();
+  const octave_idx_type num_points = points.rows ();
+
+  points = points.transpose ();
 
-  // Default options
   std::string options;
+
   if (dim <= 4)
-    options = "Qt";
+    options = " Qt";
   else
-    options = "Qt Qx";
+    options = " Qt Qx";
 
   if (nargin == 2)
     {
       if (args(1).is_string ())
-        options = args(1).string_value ();
+        options = " " + args(1).string_value ();
       else if (args(1).is_empty ())
-        ; // Use default options
+        ; // Use default options.
       else if (args(1).is_cellstr ())
         {
           options = "";
+
           Array<std::string> tmp = args(1).cellstr_value ();
 
           for (octave_idx_type i = 0; i < tmp.numel (); i++)
-            options += tmp(i) + " ";
+            options += " " + tmp(i);
         }
       else
         {
@@ -126,53 +129,56 @@
         }
      }
 
-
-  p = p.transpose ();
-  double *pt_array = p.fortran_vec ();
   boolT ismalloc = false;
 
-  // FIXME: we can't just pass options.c_str () to qh_new_qhull
-  // because the argument is not declared const.  Ugh.  Unless qh_new_qhull
-  // really needs to modify this argument, someone should fix QHULL.
-  OCTAVE_LOCAL_BUFFER (char, flags, 7 + options.length ());
-
-  sprintf (flags, "qhull %s", options.c_str ());
-
-  // Replace the 0 pointer with stdout for debugging information
+  // Replace the 0 pointer with stdout for debugging information.
   FILE *outfile = 0;
   FILE *errfile = stderr;
       
-  int exitcode = qh_new_qhull (dim, n, pt_array, 
-                               ismalloc, flags, outfile, errfile);
+  // qh_new_qhull command and points arguments are not const...
+
+  std::string cmd = "qhull" + options;
+
+  OCTAVE_LOCAL_BUFFER (char, cmd_str, cmd.length () + 1);
+
+  strcpy (cmd_str, cmd.c_str ());
+
+  int exitcode = qh_new_qhull (dim, num_points, points.fortran_vec (),
+                               ismalloc, cmd_str, outfile, errfile);
   if (! exitcode)
     {
-      vertexT *vertex, **vertexp;
-      facetT *facet;
-      setT *vertices;
       bool nonsimp_seen = false;
+
       octave_idx_type nf = qh num_facets;
 
       Matrix idx (nf, dim + 1);
 
-      octave_idx_type i = 0, j;
+      facetT *facet;
+
+      octave_idx_type i = 0;
+
       FORALLfacets
         {
-          j = 0;
+          octave_idx_type j = 0;
 
           if (! nonsimp_seen && ! facet->simplicial)
             {
               nonsimp_seen = true;
 
-              if (options.find ("QJ") != std::string::npos)
+              if (cmd.find ("QJ") != std::string::npos)
                 {
-                  // should never happen with QJ
-                  error ("convhulln: qhull failed.  Option 'QJ' returned non-simplicial facet");
-                  break;
+                  // Should never happen with QJ.
+                  error ("convhulln: qhull failed: option 'QJ' returned non-simplicial facet");
+                  return retval;
                 }
             }
+
           if (dim == 3)
             {
-              vertices = qh_facet3vertex (facet);
+              setT *vertices = qh_facet3vertex (facet);
+
+              vertexT *vertex, **vertexp;
+
               FOREACHvertex_ (vertices)
                 idx(i, j++) = 1 + qh_pointid(vertex->point);
 
@@ -182,29 +188,34 @@
             {
               if (facet->toporient ^ qh_ORIENTclock)
                 {
+                  vertexT *vertex, **vertexp;
+
                   FOREACHvertex_ (facet->vertices)
                     idx(i, j++) = 1 + qh_pointid(vertex->point);
                 }
               else
                 {
+                  vertexT *vertex, **vertexp;
+
                   FOREACHvertexreverse12_ (facet->vertices)
                     idx(i, j++) = 1 + qh_pointid(vertex->point);
                 }
             }
           if (j < dim)
-            warning ("facet %d only has %d vertices", i, j);
+            warning ("convhulln: facet %d only has %d vertices", i, j);
 
           i++;
         }
 
-      // Remove extra dimension if all facets were simplicial
+      // Remove extra dimension if all facets were simplicial.
+
       if (! nonsimp_seen)
         idx.resize (nf, dim, 0.0);
 
       if (nargout == 2)
-        // calculate volume of convex hull
-        // taken from qhull src/geom2.c
         {
+          // Calculate volume of convex hull, taken from qhull src/geom2.c.
+
           realT area;
           realT dist;
 
@@ -240,7 +251,7 @@
   else
     error ("convhulln: qhull failed");
 
-  // free memory from Qhull
+  // Free memory from Qhull
   qh_freeqhull (! qh_ALL);
 
   int curlong, totlong;
--- a/src/graphics.cc	Thu Nov 17 23:04:02 2011 +0000
+++ b/src/graphics.cc	Thu Nov 17 21:58:56 2011 -0500
@@ -564,7 +564,7 @@
   graphics_object go = gh_manager::get_object (props.get___myhandle__ ());
   graphics_object ax = go.get_ancestor ("axes");
 
-  Matrix retval (1, pos.numel (), 0);
+  Matrix retval;
 
   if (ax.valid_object ())
     {
@@ -583,6 +583,8 @@
                            v2 = ax_xform.transform (pos(0) + pos(2),
                                                     pos(1) + pos(3), 0);
 
+              retval.resize (1, 4);
+
               retval(0) = v1(0) - ax_bbox(0) + 1;
               retval(1) = ax_bbox(1) + ax_bbox(3) - v1(1) + 1;
               retval(2) = v2(0) - v1(0);
@@ -592,6 +594,8 @@
             {
               ColumnVector v = ax_xform.transform (pos(0), pos(1), pos(2));
 
+              retval.resize (1, 3);
+
               retval(0) = v(0) - ax_bbox(0) + 1;
               retval(1) = ax_bbox(1) + ax_bbox(3) - v(1) + 1;
               retval(2) = 0;
@@ -611,6 +615,8 @@
                                v2 = ax_xform.untransform (retval(0) + retval(2) + ax_bbox(0) - 1,
                                                           ax_bbox(1) + ax_bbox(3)  - (retval(1) + retval(3)) + 1);
 
+                  retval.resize (1, 4);
+
                   retval(0) = v1(0);
                   retval(1) = v1(1);
                   retval(2) = v2(0) - v1(0);
@@ -621,6 +627,8 @@
                   ColumnVector v = ax_xform.untransform (retval(0) + ax_bbox(0) - 1,
                                                          ax_bbox(1) + ax_bbox(3)  - retval(1) + 1);
 
+                  retval.resize (1, 3);
+
                   retval(0) = v(0);
                   retval(1) = v(1);
                   retval(2) = v(2);
--- a/src/ls-mat5.cc	Thu Nov 17 23:04:02 2011 +0000
+++ b/src/ls-mat5.cc	Thu Nov 17 21:58:56 2011 -0500
@@ -1224,7 +1224,9 @@
               }
             else
               {
-                octave_class* cls = new octave_class (m, classname);
+                octave_class* cls
+                  = new octave_class (m, classname,
+                                      std::list<std::string> ());
 
                 if (cls->reconstruct_exemplar ())
                   {
--- a/src/ov-class.cc	Thu Nov 17 23:04:02 2011 +0000
+++ b/src/ov-class.cc	Thu Nov 17 21:58:56 2011 -0500
@@ -515,7 +515,7 @@
               }
             else
               retval(0) = octave_value (map.index (idx.front ()),
-                                        class_name ());
+                                        c_name, parent_list);
           }
           break;
 
@@ -591,8 +591,8 @@
       else
         {
           if (type.length () == 1 && type[0] == '(')
-            retval(0) = octave_value (map.index (idx.front ()), class_name (),
-                                      parent_class_name_list ());
+            retval(0) = octave_value (map.index (idx.front ()), c_name,
+                                      parent_list);
           else
             gripe_invalid_index1 ();
         }
@@ -993,7 +993,7 @@
   if (meth.is_defined ())
     {
       octave_value_list args;
-      args(0) = octave_value (new octave_class (map, c_name));
+      args(0) = octave_value (new octave_class (map, c_name, parent_list));
 
       octave_value_list tmp = feval (meth.function_value (), args, 1);
 
@@ -1123,7 +1123,7 @@
   if (meth.is_defined ())
     {
       octave_value_list args;
-      args(0) = octave_value (new octave_class (map, c_name));
+      args(0) = octave_value (new octave_class (map, c_name, parent_list));
 
       octave_value_list tmp = feval (meth.function_value (), args, 1);
 
@@ -1954,7 +1954,9 @@
                   if (! error_state)
                     {
                       if (nargin == 2)
-                        retval = octave_value (new octave_class (m, id));
+                        retval
+                          = octave_value (new octave_class
+                                          (m, id, std::list<std::string> ()));
                       else
                         {
                           octave_value_list parents = args.slice (2, nargin-2);
--- a/src/ov-class.h	Thu Nov 17 23:04:02 2011 +0000
+++ b/src/ov-class.h	Thu Nov 17 21:58:56 2011 -0500
@@ -55,19 +55,18 @@
     { }
 
   octave_class (const octave_map& m, const std::string& id,
-                const std::list<std::string>& plist
-                  = std::list<std::string> ())
+                const std::list<std::string>& plist)
     : octave_base_value (), map (m), c_name (id),
       parent_list (plist), obsolete_copies (0)
     { }
 
+  octave_class (const octave_map& m, const std::string& id,
+                const octave_value_list& parents);
+
   octave_class (const octave_class& s)
     : octave_base_value (s), map (s.map), c_name (s.c_name),
       parent_list (s.parent_list), obsolete_copies (0)  { }
 
-  octave_class (const octave_map& m, const std::string& id,
-                const octave_value_list& parents);
-
   ~octave_class (void) { }
 
   octave_base_value *clone (void) const { return new octave_class (*this); }
@@ -76,7 +75,7 @@
 
   octave_base_value *empty_clone (void) const
   {
-    return new octave_class (octave_map (map.keys ()), class_name ());
+    return new octave_class (octave_map (map.keys ()), c_name, parent_list);
   }
 
   Cell dotref (const octave_value_list& idx);
--- a/src/ov.h	Thu Nov 17 23:04:02 2011 +0000
+++ b/src/ov.h	Thu Nov 17 21:58:56 2011 -0500
@@ -280,8 +280,7 @@
   octave_value (const octave_scalar_map& m);
   octave_value (const Octave_map& m);
   octave_value (const Octave_map& m, const std::string& id,
-                const std::list<std::string>& plist
-                  = std::list<std::string> ());
+                const std::list<std::string>& plist);
   octave_value (const octave_value_list& m, bool = false);
   octave_value (octave_value::magic_colon);
 
--- a/src/sparse-xpow.cc	Thu Nov 17 23:04:02 2011 +0000
+++ b/src/sparse-xpow.cc	Thu Nov 17 21:58:56 2011 -0500
@@ -63,7 +63,7 @@
   octave_idx_type nc = a.cols ();
 
   if (nr == 0 || nc == 0 || nr != nc)
-    error ("for A^b, A must be square");
+    error ("for A^b, A must be a square matrix");
   else
     {
       if (static_cast<int> (b) == b)
@@ -136,7 +136,7 @@
   octave_idx_type nc = a.cols ();
 
   if (nr == 0 || nc == 0 || nr != nc)
-    error ("for A^b, A must be square");
+    error ("for A^b, A must be a square matrix");
   else
     {
       if (static_cast<int> (b) == b)
--- a/src/xpow.cc	Thu Nov 17 23:04:02 2011 +0000
+++ b/src/xpow.cc	Thu Nov 17 21:58:56 2011 -0500
@@ -104,7 +104,7 @@
   octave_idx_type nc = b.cols ();
 
   if (nr == 0 || nc == 0 || nr != nc)
-    error ("for x^A, A must be square");
+    error ("for x^A, A must be a square matrix");
   else
     {
       EIG b_eig (b);
@@ -155,7 +155,7 @@
   octave_idx_type nc = b.cols ();
 
   if (nr == 0 || nc == 0 || nr != nc)
-    error ("for x^A, A must be square");
+    error ("for x^A, A must be a square matrix");
   else
     {
       EIG b_eig (b);
@@ -194,7 +194,7 @@
   octave_idx_type nc = a.cols ();
 
   if (nr == 0 || nc == 0 || nr != nc)
-    error ("for A^b, A must be square");
+    error ("for A^b, A must be a square matrix");
   else
     {
       if (static_cast<int> (b) == b)
@@ -280,7 +280,7 @@
   octave_idx_type nc = a.cols ();
 
   if (nr == 0 || nc == 0 || nr != nc)
-    error ("for A^b, A must be square");
+    error ("for A^b, A must be a square matrix");
   else
     {
       if (static_cast<int> (b) == b)
@@ -324,7 +324,7 @@
   octave_idx_type nc = a.cols ();
 
   if (nr == 0 || nc == 0 || nr != nc)
-    error ("for A^b, A must be square");
+    error ("for A^b, A must be a square matrix");
   else
     {
       EIG a_eig (a);
@@ -372,7 +372,7 @@
   octave_idx_type nc = b.cols ();
 
   if (nr == 0 || nc == 0 || nr != nc)
-    error ("for x^A, A must be square");
+    error ("for x^A, A must be a square matrix");
   else
     {
       EIG b_eig (b);
@@ -420,7 +420,7 @@
   octave_idx_type nc = b.cols ();
 
   if (nr == 0 || nc == 0 || nr != nc)
-    error ("for x^A, A must be square");
+    error ("for x^A, A must be a square matrix");
   else
     {
       EIG b_eig (b);
@@ -459,7 +459,7 @@
   octave_idx_type nc = a.cols ();
 
   if (nr == 0 || nc == 0 || nr != nc)
-    error ("for A^b, A must be square");
+    error ("for A^b, A must be a square matrix");
   else
     {
       if (static_cast<int> (b) == b)
@@ -545,7 +545,7 @@
   octave_idx_type nc = a.cols ();
 
   if (nr == 0 || nc == 0 || nr != nc)
-    error ("for A^b, A must be square");
+    error ("for A^b, A must be a square matrix");
   else
     {
       EIG a_eig (a);
@@ -579,7 +579,7 @@
   octave_idx_type nc = a.cols ();
 
   if (nr == 0 || nc == 0 || nr != nc)
-    error ("for A^b, A must be square");
+    error ("for A^b, A must be a square matrix");
   else
     {
       ComplexDiagMatrix r (nr, nc);
@@ -1553,7 +1553,7 @@
   octave_idx_type nc = b.cols ();
 
   if (nr == 0 || nc == 0 || nr != nc)
-    error ("for x^A, A must be square");
+    error ("for x^A, A must be a square matrix");
   else
     {
       FloatEIG b_eig (b);
@@ -1605,7 +1605,7 @@
   octave_idx_type nc = b.cols ();
 
   if (nr == 0 || nc == 0 || nr != nc)
-    error ("for x^A, A must be square");
+    error ("for x^A, A must be a square matrix");
   else
     {
       FloatEIG b_eig (b);
@@ -1644,7 +1644,7 @@
   octave_idx_type nc = a.cols ();
 
   if (nr == 0 || nc == 0 || nr != nc)
-    error ("for A^b, A must be square");
+    error ("for A^b, A must be a square matrix");
   else
     {
       if (static_cast<int> (b) == b)
@@ -1730,7 +1730,7 @@
   octave_idx_type nc = a.cols ();
 
   if (nr == 0 || nc == 0 || nr != nc)
-    error ("for A^b, A must be square");
+    error ("for A^b, A must be a square matrix");
   else
     {
       if (static_cast<int> (b) == b)
@@ -1762,7 +1762,7 @@
   octave_idx_type nc = a.cols ();
 
   if (nr == 0 || nc == 0 || nr != nc)
-    error ("for A^b, A must be square");
+    error ("for A^b, A must be a square matrix");
   else
     {
       FloatEIG a_eig (a);
@@ -1810,7 +1810,7 @@
   octave_idx_type nc = b.cols ();
 
   if (nr == 0 || nc == 0 || nr != nc)
-    error ("for x^A, A must be square");
+    error ("for x^A, A must be a square matrix");
   else
     {
       FloatEIG b_eig (b);
@@ -1858,7 +1858,7 @@
   octave_idx_type nc = b.cols ();
 
   if (nr == 0 || nc == 0 || nr != nc)
-    error ("for x^A, A must be square");
+    error ("for x^A, A must be a square matrix");
   else
     {
       FloatEIG b_eig (b);
@@ -1897,7 +1897,7 @@
   octave_idx_type nc = a.cols ();
 
   if (nr == 0 || nc == 0 || nr != nc)
-    error ("for A^b, A must be square");
+    error ("for A^b, A must be a square matrix");
   else
     {
       if (static_cast<int> (b) == b)
@@ -1983,7 +1983,7 @@
   octave_idx_type nc = a.cols ();
 
   if (nr == 0 || nc == 0 || nr != nc)
-    error ("for A^b, A must be square");
+    error ("for A^b, A must be a square matrix");
   else
     {
       FloatEIG a_eig (a);
@@ -2017,7 +2017,7 @@
   octave_idx_type nc = a.cols ();
 
   if (nr == 0 || nc == 0 || nr != nc)
-    error ("for A^b, A must be square");
+    error ("for A^b, A must be a square matrix");
   else
     {
       FloatComplexDiagMatrix r (nr, nc);