changeset 6390:334499d75c5c

[project @ 2007-03-07 18:11:28 by jwe]
author jwe
date Wed, 07 Mar 2007 18:11:29 +0000
parents f427b33aeb4c
children 3f3e86e9fb57
files liboctave/CMatrix.cc liboctave/ChangeLog liboctave/dMatrix.cc scripts/ChangeLog scripts/plot/__uiobject_adopt__.m scripts/plot/__uiobject_draw_axes__.m scripts/plot/__uiobject_make_handle__.in scripts/plot/mesh.m src/DLD-FUNCTIONS/urlwrite.cc
diffstat 9 files changed, 65 insertions(+), 38 deletions(-) [+]
line wrap: on
line diff
--- a/liboctave/CMatrix.cc	Wed Mar 07 09:13:50 2007 +0000
+++ b/liboctave/CMatrix.cc	Wed Mar 07 18:11:29 2007 +0000
@@ -1683,7 +1683,7 @@
 	  if (typ == MatrixType::Permuted_Upper)
 	    {
 	      (*current_liboctave_error_handler)
-		("Permuted triangular matrix not implemented");
+		("permuted triangular matrix not implemented");
 	    }
 	  else
 	    {
@@ -1790,7 +1790,7 @@
 	  if (typ == MatrixType::Permuted_Lower)
 	    {
 	      (*current_liboctave_error_handler)
-		("Permuted triangular matrix not implemented");
+		("permuted triangular matrix not implemented");
 	    }
 	  else
 	    {
@@ -3768,22 +3768,30 @@
 	      if (nr == 1)
 		F77_FUNC (xzdotu, XZDOTU) (nc, m.data (), 1, a.data (), 1, *c);
 	      else
-		F77_XFCN (zgemv, ZGEMV, (F77_CONST_CHAR_ARG2 ("N", 1),
-					 nr, nc, 1.0,  m.data (), ld,
-					 a.data (), 1, 0.0, c, 1
-					 F77_CHAR_ARG_LEN (1)));
+		{
+		  F77_XFCN (zgemv, ZGEMV, (F77_CONST_CHAR_ARG2 ("N", 1),
+					   nr, nc, 1.0,  m.data (), ld,
+					   a.data (), 1, 0.0, c, 1
+					   F77_CHAR_ARG_LEN (1)));
+
+		  if (f77_exception_encountered)
+		    (*current_liboctave_error_handler)
+		      ("unrecoverable error in zgemv");
+		}
 	    }
 	  else
-	    F77_XFCN (zgemm, ZGEMM, (F77_CONST_CHAR_ARG2 ("N", 1),
-				     F77_CONST_CHAR_ARG2 ("N", 1),
-				     nr, a_nc, nc, 1.0, m.data (),
-				     ld, a.data (), lda, 0.0, c, nr
-				     F77_CHAR_ARG_LEN (1)
-				     F77_CHAR_ARG_LEN (1)));
-
-	  if (f77_exception_encountered)
-	    (*current_liboctave_error_handler)
-	      ("unrecoverable error in zgemm");
+	    {
+	      F77_XFCN (zgemm, ZGEMM, (F77_CONST_CHAR_ARG2 ("N", 1),
+				       F77_CONST_CHAR_ARG2 ("N", 1),
+				       nr, a_nc, nc, 1.0, m.data (),
+				       ld, a.data (), lda, 0.0, c, nr
+				       F77_CHAR_ARG_LEN (1)
+				       F77_CHAR_ARG_LEN (1)));
+
+	      if (f77_exception_encountered)
+		(*current_liboctave_error_handler)
+		  ("unrecoverable error in zgemm");
+	    }
 	}
     }
 
--- a/liboctave/ChangeLog	Wed Mar 07 09:13:50 2007 +0000
+++ b/liboctave/ChangeLog	Wed Mar 07 18:11:29 2007 +0000
@@ -1,5 +1,8 @@
 2007-03-07  John W. Eaton  <jwe@octave.org>
 
+	* dMatrix.cc, CMatrix.cc (operator *): Only check
+	f77_exception_encountered immediately after calls that use F77_XFCN.
+
 	* Array.cc (assign1 (Array<LT>&, const Array<RT>&, const LT&)):
 	Only allow resizing empty LHS if it is 0x0.
 
--- a/liboctave/dMatrix.cc	Wed Mar 07 09:13:50 2007 +0000
+++ b/liboctave/dMatrix.cc	Wed Mar 07 18:11:29 2007 +0000
@@ -1345,7 +1345,7 @@
 	  if (typ == MatrixType::Permuted_Upper)
 	    {
 	      (*current_liboctave_error_handler)
-		("Permuted triangular matrix not implemented");
+		("permuted triangular matrix not implemented");
 	    }
 	  else
 	    {
@@ -1451,7 +1451,7 @@
 	  if (typ == MatrixType::Permuted_Lower)
 	    {
 	      (*current_liboctave_error_handler)
-		("Permuted triangular matrix not implemented");
+		("permuted triangular matrix not implemented");
 	    }
 	  else
 	    {
@@ -3145,22 +3145,30 @@
 	      if (nr == 1)
 		F77_FUNC (xddot, XDDOT) (nc, m.data (), 1, a.data (), 1, *c);
 	      else
-		F77_XFCN (dgemv, DGEMV, (F77_CONST_CHAR_ARG2 ("N", 1),
-					 nr, nc, 1.0,  m.data (), ld,
-					 a.data (), 1, 0.0, c, 1
-					 F77_CHAR_ARG_LEN (1)));
+		{
+		  F77_XFCN (dgemv, DGEMV, (F77_CONST_CHAR_ARG2 ("N", 1),
+					   nr, nc, 1.0,  m.data (), ld,
+					   a.data (), 1, 0.0, c, 1
+					   F77_CHAR_ARG_LEN (1)));
+
+		  if (f77_exception_encountered)
+		    (*current_liboctave_error_handler)
+		      ("unrecoverable error in dgemv");
+		}
             }
 	  else
-	    F77_XFCN (dgemm, DGEMM, (F77_CONST_CHAR_ARG2 ("N", 1),
-				     F77_CONST_CHAR_ARG2 ("N", 1),
-				     nr, a_nc, nc, 1.0, m.data (),
-				     ld, a.data (), lda, 0.0, c, nr
-				     F77_CHAR_ARG_LEN (1)
-				     F77_CHAR_ARG_LEN (1)));
-
-	  if (f77_exception_encountered)
-	    (*current_liboctave_error_handler)
-	      ("unrecoverable error in dgemm");
+	    {
+	      F77_XFCN (dgemm, DGEMM, (F77_CONST_CHAR_ARG2 ("N", 1),
+				       F77_CONST_CHAR_ARG2 ("N", 1),
+				       nr, a_nc, nc, 1.0, m.data (),
+				       ld, a.data (), lda, 0.0, c, nr
+				       F77_CHAR_ARG_LEN (1)
+				       F77_CHAR_ARG_LEN (1)));
+
+	      if (f77_exception_encountered)
+		(*current_liboctave_error_handler)
+		  ("unrecoverable error in dgemm");
+	    }
 	}
     }
 
--- a/scripts/ChangeLog	Wed Mar 07 09:13:50 2007 +0000
+++ b/scripts/ChangeLog	Wed Mar 07 18:11:29 2007 +0000
@@ -1,5 +1,11 @@
 2007-03-07  John W. Eaton  <jwe@octave.org>
 
+	* plot/mesh.m: Call newplot before doing anything.
+
+	* plot/__uiobject_draw_axes__.m: Send "e\n" at end of data, not
+	just "e".  Only flush plot stream once.
+	From Daniel J Sebald <daniel.sebald@ieee.org>.
+
 	* strings/blanks.m: Omit first index in assignment.
 
 2007-03-07  Paul Kienzle  <pkienzle@users.sf.net>
--- a/scripts/plot/__uiobject_adopt__.m	Wed Mar 07 09:13:50 2007 +0000
+++ b/scripts/plot/__uiobject_adopt__.m	Wed Mar 07 18:11:29 2007 +0000
@@ -33,7 +33,7 @@
       ## Put this child at the end of the list.  If it is already in
       ## the list, move it.
       kids(kids == child) = [];
-      kids = [kids, child];
+      kids = [kids, child]
       set (parent, "children", kids);
     else
       error ("__uiobject_adopt__: expecting parent to be a handle");
--- a/scripts/plot/__uiobject_draw_axes__.m	Wed Mar 07 09:13:50 2007 +0000
+++ b/scripts/plot/__uiobject_draw_axes__.m	Wed Mar 07 18:11:29 2007 +0000
@@ -642,7 +642,6 @@
     endif
 
     fputs (plot_stream, "set style data lines;\n");
-    fflush (plot_stream);
 
     if (! use_gnuplot_for_images)
       for i = 1:ximg_data_idx
@@ -688,13 +687,11 @@
 	      endfor
 	    endif
 	  endif
-	  fputs (plot_stream, "e");
-	  fflush (plot_stream);
+	  fputs (plot_stream, "e\n");
 	endif
       endfor
     endif
 
-    fputs (plot_stream, "\n");
     fflush (plot_stream);
 
   else
--- a/scripts/plot/__uiobject_make_handle__.in	Wed Mar 07 09:13:50 2007 +0000
+++ b/scripts/plot/__uiobject_make_handle__.in	Wed Mar 07 18:11:29 2007 +0000
@@ -30,7 +30,7 @@
 
   if (nargin == 1)
     idx = __uiobject_alloc__ ();
-    h = __uiobject_get_handle__ ();
+    h = __uiobject_get_handle__ ()
     __uiobject_list__(idx).handle = h;
     __uiobject_list__(idx).object = s;
     __request_drawnow__ ();
--- a/scripts/plot/mesh.m	Wed Mar 07 09:13:50 2007 +0000
+++ b/scripts/plot/mesh.m	Wed Mar 07 18:11:29 2007 +0000
@@ -32,6 +32,8 @@
 
 function h = mesh (x, y, z)
 
+  newplot ();
+
   if (nargin == 1)
     z = x;
     if (ismatrix (z))
--- a/src/DLD-FUNCTIONS/urlwrite.cc	Wed Mar 07 09:13:50 2007 +0000
+++ b/src/DLD-FUNCTIONS/urlwrite.cc	Wed Mar 07 18:11:29 2007 +0000
@@ -153,6 +153,9 @@
   // Set a pointer to our struct to pass to the callback.
   curl_easy_setopt (curl, CURLOPT_WRITEDATA, static_cast<void*> (&stream));
 
+  // Follow redirects.
+  curl_easy_setopt (curl, CURLOPT_FOLLOWLOCATION, 1);
+
   curl_easy_setopt (curl, CURLOPT_NOPROGRESS, false);
   curl_easy_setopt (curl, CURLOPT_PROGRESSFUNCTION, progress_func);
   curl_easy_setopt (curl, CURLOPT_PROGRESSDATA, url.c_str ());