changeset 14394:ed8c4921bf61 stable

correctly fill result for M * v for Nx0 * 0x1 operations * dColVector.cc (operator * (const Matrix&, const ColumnVector&)): Fill result if NC is 0. * CColVector.cc (operator * (const ComplexMatrix&, const ComplexColumnVector&)): Likewise. * fCColVector.cc (const FloatComplexMatrix&, const FloatComplexColumnVector&)): Likewise. * fColVector.cc (const FloatMatrix&, const FloatColumnVector&)): Likewise.
author John W. Eaton <jwe@octave.org>
date Thu, 23 Feb 2012 13:50:26 -0500
parents ba4d6343524b
children 08e48e7a4c8a
files liboctave/CColVector.cc liboctave/dColVector.cc liboctave/fCColVector.cc liboctave/fColVector.cc
diffstat 4 files changed, 41 insertions(+), 20 deletions(-) [+]
line wrap: on
line diff
--- a/liboctave/CColVector.cc	Wed Feb 22 14:39:41 2012 -0500
+++ b/liboctave/CColVector.cc	Thu Feb 23 13:50:26 2012 -0500
@@ -346,13 +346,19 @@
 
       if (nr != 0)
         {
-          Complex *y = retval.fortran_vec ();
+          if (nc == 0)
+            retval.fill (0.0);
+          else
+            {
+              Complex *y = retval.fortran_vec ();
 
-          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)));
+              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/dColVector.cc	Wed Feb 22 14:39:41 2012 -0500
+++ b/liboctave/dColVector.cc	Thu Feb 23 13:50:26 2012 -0500
@@ -212,12 +212,17 @@
 
       if (nr != 0)
         {
-          double *y = retval.fortran_vec ();
+          if (nc == 0)
+            retval.fill (0.0);
+          else
+            {
+              double *y = retval.fortran_vec ();
 
-          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)));
+              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)));
+            }
         }
     }
 
--- a/liboctave/fCColVector.cc	Wed Feb 22 14:39:41 2012 -0500
+++ b/liboctave/fCColVector.cc	Thu Feb 23 13:50:26 2012 -0500
@@ -346,12 +346,17 @@
 
       if (nr != 0)
         {
-          FloatComplex *y = retval.fortran_vec ();
+          if (nc == 0)
+            retval.fill (0.0);
+          else
+            {
+              FloatComplex *y = retval.fortran_vec ();
 
-          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)));
+              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)));
+            }
         }
     }
 
--- a/liboctave/fColVector.cc	Wed Feb 22 14:39:41 2012 -0500
+++ b/liboctave/fColVector.cc	Thu Feb 23 13:50:26 2012 -0500
@@ -211,12 +211,17 @@
 
       if (nr != 0)
         {
-          float *y = retval.fortran_vec ();
+          if (nc == 0)
+            retval.fill (0.0);
+          else
+            {
+              float *y = retval.fortran_vec ();
 
-          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)));
+              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)));
+            }
         }
     }