changeset 9625:cbabf50315ca

optimize Matrix*ColumnVector
author Jaroslav Hajek <highegg@gmail.com>
date Sun, 06 Sep 2009 11:00:19 +0200
parents 3fc7272937ce
children bccba774af8b
files liboctave/CColVector.cc liboctave/ChangeLog liboctave/dColVector.cc liboctave/fCColVector.cc liboctave/fColVector.cc
diffstat 5 files changed, 51 insertions(+), 52 deletions(-) [+]
line wrap: on
line diff
--- a/liboctave/CColVector.cc	Sun Sep 06 10:54:36 2009 +0200
+++ b/liboctave/CColVector.cc	Sun Sep 06 11:00:19 2009 +0200
@@ -340,20 +340,17 @@
     gripe_nonconformant ("operator *", nr, nc, a_len, 1);
   else
     {
-      if (nc == 0 || nr == 0)
-	retval.resize (nr, 0.0);
-      else
-	{
-	  octave_idx_type ld = nr;
+      retval.clear (nr);
+
+      if (nr != 0)
+        {
+          Complex *y = retval.fortran_vec ();
 
-	  retval.resize (nr);
-	  Complex *y = retval.fortran_vec ();
-
-	  F77_XFCN (zgemv, ZGEMV, (F77_CONST_CHAR_ARG2 ("N", 1),
-				   nr, nc, 1.0, m.data (), ld,
-				   a.data (), 1, 0.0, y, 1
-				   F77_CHAR_ARG_LEN (1)));
-	}
+          F77_XFCN (zgemv, ZGEMV, (F77_CONST_CHAR_ARG2 ("N", 1),
+                                   nr, nc, 1.0, m.data (), nr,
+                                   a.data (), 1, 0.0, y, 1
+                                   F77_CHAR_ARG_LEN (1)));
+        }
     }
 
   return retval;
--- a/liboctave/ChangeLog	Sun Sep 06 10:54:36 2009 +0200
+++ b/liboctave/ChangeLog	Sun Sep 06 11:00:19 2009 +0200
@@ -1,3 +1,14 @@
+2009-09-06  Jaroslav Hajek  <highegg@gmail.com>
+
+	* dColVector.h (operator *(const Matrix&, const ColumnVector)):
+	Optimize.
+	* fColVector.h (operator *(const FloatMatrix&, const
+	FloatColumnVector)): Optimize.
+	* CColVector.h (operator *(const ComplexMatrix&, const
+	ComplexColumnVector)): Optimize.
+	* fCColVector.h (operator *(const FloatComplexMatrix&, const
+	FloatComplexColumnVector)): Optimize.
+
 2009-09-04  Jaroslav Hajek  <highegg@gmail.com>
 
 	* Array.cc (Array<T>::clear (const dim_vector&)): new method.
--- a/liboctave/dColVector.cc	Sun Sep 06 10:54:36 2009 +0200
+++ b/liboctave/dColVector.cc	Sun Sep 06 11:00:19 2009 +0200
@@ -209,20 +209,17 @@
     gripe_nonconformant ("operator *", nr, nc, a_len, 1);
   else
     {
-      if (nr == 0 || nc == 0)
-	retval.resize (nr, 0.0);
-      else
-	{
-	  octave_idx_type ld = nr;
+      retval.clear (nr);
+
+      if (nr != 0)
+        {
+          double *y = retval.fortran_vec ();
 
-	  retval.resize (nr);
-	  double *y = retval.fortran_vec ();
-
-	  F77_XFCN (dgemv, DGEMV, (F77_CONST_CHAR_ARG2 ("N", 1),
-				   nr, nc, 1.0, m.data (), ld,
-				   a.data (), 1, 0.0, y, 1
-				   F77_CHAR_ARG_LEN (1)));
-	}
+          F77_XFCN (dgemv, DGEMV, (F77_CONST_CHAR_ARG2 ("N", 1),
+                                   nr, nc, 1.0, m.data (), nr,
+                                   a.data (), 1, 0.0, y, 1
+                                   F77_CHAR_ARG_LEN (1)));
+        }
     }
 
   return retval;
--- a/liboctave/fCColVector.cc	Sun Sep 06 10:54:36 2009 +0200
+++ b/liboctave/fCColVector.cc	Sun Sep 06 11:00:19 2009 +0200
@@ -340,20 +340,17 @@
     gripe_nonconformant ("operator *", nr, nc, a_len, 1);
   else
     {
-      if (nc == 0 || nr == 0)
-	retval.resize (nr, 0.0);
-      else
-	{
-	  octave_idx_type ld = nr;
+      retval.clear (nr);
+
+      if (nr != 0)
+        {
+          FloatComplex *y = retval.fortran_vec ();
 
-	  retval.resize (nr);
-	  FloatComplex *y = retval.fortran_vec ();
-
-	  F77_XFCN (cgemv, CGEMV, (F77_CONST_CHAR_ARG2 ("N", 1),
-				   nr, nc, 1.0, m.data (), ld,
-				   a.data (), 1, 0.0, y, 1
-				   F77_CHAR_ARG_LEN (1)));
-	}
+          F77_XFCN (cgemv, CGEMV, (F77_CONST_CHAR_ARG2 ("N", 1),
+                                   nr, nc, 1.0f, m.data (), nr,
+                                   a.data (), 1, 0.0f, y, 1
+                                   F77_CHAR_ARG_LEN (1)));
+        }
     }
 
   return retval;
--- a/liboctave/fColVector.cc	Sun Sep 06 10:54:36 2009 +0200
+++ b/liboctave/fColVector.cc	Sun Sep 06 11:00:19 2009 +0200
@@ -209,20 +209,17 @@
     gripe_nonconformant ("operator *", nr, nc, a_len, 1);
   else
     {
-      if (nr == 0 || nc == 0)
-	retval.resize (nr, 0.0);
-      else
-	{
-	  octave_idx_type ld = nr;
+      retval.clear (nr);
+
+      if (nr != 0)
+        {
+          float *y = retval.fortran_vec ();
 
-	  retval.resize (nr);
-	  float *y = retval.fortran_vec ();
-
-	  F77_XFCN (sgemv, SGEMV, (F77_CONST_CHAR_ARG2 ("N", 1),
-				   nr, nc, 1.0, m.data (), ld,
-				   a.data (), 1, 0.0, y, 1
-				   F77_CHAR_ARG_LEN (1)));
-	}
+          F77_XFCN (sgemv, SGEMV, (F77_CONST_CHAR_ARG2 ("N", 1),
+                                   nr, nc, 1.0f, m.data (), nr,
+                                   a.data (), 1, 0.0f, y, 1
+                                   F77_CHAR_ARG_LEN (1)));
+        }
     }
 
   return retval;