Mercurial > hg > octave-lyh
changeset 9469:c6edba80dfae
sanity checks for loading sparse matrices
author | John W. Eaton <jwe@octave.org> |
---|---|
date | Wed, 29 Jul 2009 12:15:27 -0400 |
parents | 5af462716bff |
children | bcdf878e2686 |
files | liboctave/CMatrix.cc liboctave/CNDArray.cc liboctave/CSparse.cc liboctave/ChangeLog liboctave/Sparse.cc liboctave/Sparse.h liboctave/boolSparse.cc liboctave/dMatrix.cc liboctave/dNDArray.cc liboctave/dSparse.cc liboctave/fCMatrix.cc liboctave/fCNDArray.cc liboctave/fMatrix.cc liboctave/fNDArray.cc liboctave/lo-utils.cc liboctave/lo-utils.h liboctave/sparse-util.cc liboctave/sparse-util.h src/ChangeLog src/ls-mat-ascii.cc src/ov-bool-sparse.cc src/ov-bool.cc src/ov-complex.cc src/ov-cx-sparse.cc src/ov-float.cc src/ov-flt-complex.cc src/ov-re-sparse.cc src/ov-scalar.cc |
diffstat | 28 files changed, 365 insertions(+), 240 deletions(-) [+] |
line wrap: on
line diff
--- a/liboctave/CMatrix.cc +++ b/liboctave/CMatrix.cc @@ -3646,7 +3646,7 @@ for (octave_idx_type i = 0; i < nr; i++) for (octave_idx_type j = 0; j < nc; j++) { - tmp = octave_read_complex (is); + tmp = octave_read_value<Complex> (is); if (is) a.elem (i, j) = tmp; else
--- a/liboctave/CNDArray.cc +++ b/liboctave/CNDArray.cc @@ -921,7 +921,7 @@ Complex tmp; for (octave_idx_type i = 0; i < nel; i++) { - tmp = octave_read_complex (is); + tmp = octave_read_value<Complex> (is); if (is) a.elem (i) = tmp; else
--- a/liboctave/CSparse.cc +++ b/liboctave/CSparse.cc @@ -7466,48 +7466,9 @@ std::istream& operator >> (std::istream& is, SparseComplexMatrix& a) { - octave_idx_type nr = a.rows (); - octave_idx_type nc = a.cols (); - octave_idx_type nz = a.nzmax (); - - if (nr > 0 && nc > 0) - { - octave_idx_type itmp, jtmp, jold = 0; - Complex tmp; - octave_idx_type ii = 0; - - a.cidx (0) = 0; - for (octave_idx_type i = 0; i < nz; i++) - { - is >> itmp; - itmp--; - is >> jtmp; - jtmp--; - tmp = octave_read_complex (is); - - if (is) - { - if (jold != jtmp) - { - for (octave_idx_type j = jold; j < jtmp; j++) - a.cidx(j+1) = ii; - - jold = jtmp; - } - a.data (ii) = tmp; - a.ridx (ii++) = itmp; - } - else - goto done; - } - - for (octave_idx_type j = jold; j < nc; j++) - a.cidx(j+1) = ii; - } - - done: - - return is; + typedef SparseComplexMatrix::element_type elt_type; + + return read_sparse_matrix<elt_type> (is, a, octave_read_value<Complex>); } SparseComplexMatrix
--- a/liboctave/ChangeLog +++ b/liboctave/ChangeLog @@ -1,3 +1,27 @@ +2009-07-29 John W. Eaton <jwe@octave.org> + + * fMatrix.cc (operator >>): Use template function to read value. + * fNDArray.cc (operator >>): Likeise. + * fCMatrix.cc (operator >>): Use template function to read value. + * fCNDArray.cc (operator >>): Likeise. + * dMatrix.cc (operator >>): Use template function to read value. + * dNDArray.cc (operator >>): Likeise. + * CMatrix.cc (operator >>): Use template function to read value. + * CNDArray.cc (operator >>): Likeise. + + * lo-utils.cc, lo-utils.h (octave_read_value): New template + (octave_read_value<double>, octave_read_value<Complex>): + Provide specializations. + (octave_read_double, octave_read_complex, octave_read_float, + octave_rread_float_complex): Define in terms of template functions. + * Sparse.h (read_sparse_matrix): New template function. + * dSparse.cc (operator >>): Call read_sparse_matrix. + * CSparse.cc (operator >>): Likewise. + * boolSparse.cc (operator >>): Likewise. + * sparse-util.cc, sparse-util.h (sparse_indices_ok): New function. + * Sparse.cc (Sparse<T>::indices_ok, Sparse<T>::SparseRep::indices_ok): + New member functions. + 2009-07-20 John W. Eaton <jwe@octave.org> * lo-ieee.cc (octave_ieee_init) [__NetBSD__]: Call nan to
--- a/liboctave/Sparse.cc +++ b/liboctave/Sparse.cc @@ -43,6 +43,7 @@ #include "Sparse.h" #include "sparse-sort.h" +#include "sparse-util.h" #include "oct-spparms.h" template <class T> @@ -199,6 +200,13 @@ } template <class T> +bool +Sparse<T>::SparseRep::indices_ok (void) const +{ + return sparse_indices_ok (r, c, nrows, ncols, nnz ()); +} + +template <class T> template <class U> Sparse<T>::Sparse (const Sparse<U>& a) : dimensions (a.dimensions), idx (0), idx_count (0)
--- a/liboctave/Sparse.h +++ b/liboctave/Sparse.h @@ -33,6 +33,7 @@ #include "Array.h" #include "Array2.h" #include "dim-vector.h" +#include "lo-error.h" #include "lo-utils.h" #include "oct-sort.h" @@ -131,6 +132,8 @@ void change_length (octave_idx_type nz); + bool indices_ok (void) const; + private: // No assignment! @@ -565,6 +568,8 @@ return result; } + + bool indices_ok (void) const { return rep->indices_ok (); } }; // NOTE: these functions should be friends of the Sparse<T> class and @@ -579,6 +584,98 @@ /* friend */ int assign1 (Sparse<LT>& lhs, const Sparse<RT>& rhs); +template<typename T> +std::istream& +read_sparse_matrix (std::istream& is, Sparse<T>& a, + T (*read_fcn) (std::istream&)) +{ + octave_idx_type nr = a.rows (); + octave_idx_type nc = a.cols (); + octave_idx_type nz = a.nzmax (); + + if (nr > 0 && nc > 0) + { + octave_idx_type itmp; + octave_idx_type jtmp; + octave_idx_type iold = 0; + octave_idx_type jold = 0; + octave_idx_type ii = 0; + T tmp; + + a.cidx (0) = 0; + for (octave_idx_type i = 0; i < nz; i++) + { + is >> itmp; + itmp--; + + if (itmp < iold) + { + (*current_liboctave_error_handler) + ("invalid sparse matrix: row indices must appear in ascending order in each column"); + is.setstate (std::ios::failbit); + goto done; + } + + if (itmp < 0 || itmp >= nr) + { + (*current_liboctave_error_handler) + ("invalid sparse matrix: row index = %d out of range", + itmp + 1); + is.setstate (std::ios::failbit); + goto done; + } + + + iold = itmp; + + is >> jtmp; + jtmp--; + + if (jtmp < jold) + { + (*current_liboctave_error_handler) + ("invalid sparse matrix: column indices must appear in ascending order"); + is.setstate (std::ios::failbit); + goto done; + } + + if (jtmp < 0 || jtmp >= nc) + { + (*current_liboctave_error_handler) + ("invalid sparse matrix: column index = %d out of range", + jtmp + 1); + is.setstate (std::ios::failbit); + goto done; + } + + tmp = read_fcn (is); + + if (is) + { + if (jold != jtmp) + { + for (octave_idx_type j = jold; j < jtmp; j++) + a.cidx(j+1) = ii; + + jold = jtmp; + iold = 0; + } + a.data (ii) = tmp; + a.ridx (ii++) = itmp; + } + else + goto done; + } + + for (octave_idx_type j = jold; j < nc; j++) + a.cidx(j+1) = ii; + } + + done: + + return is; +} + #define INSTANTIATE_SPARSE_ASSIGN(LT, RT, API) \ template API int assign (Sparse<LT>&, const Sparse<RT>&); \ template API int assign1 (Sparse<LT>&, const Sparse<RT>&);
--- a/liboctave/boolSparse.cc +++ b/liboctave/boolSparse.cc @@ -180,47 +180,9 @@ std::istream& operator >> (std::istream& is, SparseBoolMatrix& a) { - octave_idx_type nr = a.rows (); - octave_idx_type nc = a.cols (); - octave_idx_type nz = a.nzmax (); + typedef SparseBoolMatrix::element_type elt_type; - if (nr > 0 && nc > 0) - { - octave_idx_type itmp, jtmp, jold = 0; - bool tmp; - octave_idx_type ii = 0; - - a.cidx (0) = 0; - for (octave_idx_type i = 0; i < nz; i++) - { - is >> itmp; - itmp--; - is >> jtmp; - jtmp--; - is >> tmp; - if (is) - { - if (jold != jtmp) - { - for (octave_idx_type j = jold; j < jtmp; j++) - a.cidx(j+1) = ii; - - jold = jtmp; - } - a.data (ii) = tmp; - a.ridx (ii++) = itmp; - } - else - goto done; - } - - for (octave_idx_type j = jold; j < nc; j++) - a.cidx(j+1) = ii; - } - - done: - - return is; + return read_sparse_matrix<elt_type> (is, a, octave_read_value<bool>); } SparseBoolMatrix
--- a/liboctave/dMatrix.cc +++ b/liboctave/dMatrix.cc @@ -3064,7 +3064,7 @@ for (octave_idx_type i = 0; i < nr; i++) for (octave_idx_type j = 0; j < nc; j++) { - tmp = octave_read_double (is); + tmp = octave_read_value<double> (is); if (is) a.elem (i, j) = tmp; else
--- a/liboctave/dNDArray.cc +++ b/liboctave/dNDArray.cc @@ -970,7 +970,7 @@ double tmp; for (octave_idx_type i = 0; i < nel; i++) { - tmp = octave_read_double (is); + tmp = octave_read_value<double> (is); if (is) a.elem (i) = tmp; else
--- a/liboctave/dSparse.cc +++ b/liboctave/dSparse.cc @@ -7584,48 +7584,9 @@ std::istream& operator >> (std::istream& is, SparseMatrix& a) { - octave_idx_type nr = a.rows (); - octave_idx_type nc = a.cols (); - octave_idx_type nz = a.nzmax (); - - if (nr > 0 && nc > 0) - { - octave_idx_type itmp, jtmp, jold = 0; - double tmp; - octave_idx_type ii = 0; - - a.cidx (0) = 0; - for (octave_idx_type i = 0; i < nz; i++) - { - is >> itmp; - itmp--; - is >> jtmp; - jtmp--; - tmp = octave_read_double (is); - - if (is) - { - if (jold != jtmp) - { - for (octave_idx_type j = jold; j < jtmp; j++) - a.cidx(j+1) = ii; - - jold = jtmp; - } - a.data (ii) = tmp; - a.ridx (ii++) = itmp; - } - else - goto done; - } - - for (octave_idx_type j = jold; j < nc; j++) - a.cidx(j+1) = ii; - } - - done: - - return is; + typedef SparseMatrix::element_type elt_type; + + return read_sparse_matrix<elt_type> (is, a, octave_read_value<double>); } SparseMatrix
--- a/liboctave/fCMatrix.cc +++ b/liboctave/fCMatrix.cc @@ -3639,7 +3639,7 @@ for (octave_idx_type i = 0; i < nr; i++) for (octave_idx_type j = 0; j < nc; j++) { - tmp = octave_read_complex (is); + tmp = octave_read_value<FloatComplex> (is); if (is) a.elem (i, j) = tmp; else
--- a/liboctave/fCNDArray.cc +++ b/liboctave/fCNDArray.cc @@ -916,7 +916,7 @@ FloatComplex tmp; for (octave_idx_type i = 0; i < nel; i++) { - tmp = octave_read_complex (is); + tmp = octave_read_value<FloatComplex> (is); if (is) a.elem (i) = tmp; else
--- a/liboctave/fMatrix.cc +++ b/liboctave/fMatrix.cc @@ -3063,7 +3063,7 @@ for (octave_idx_type i = 0; i < nr; i++) for (octave_idx_type j = 0; j < nc; j++) { - tmp = octave_read_float (is); + tmp = octave_read_value<float> (is); if (is) a.elem (i, j) = tmp; else
--- a/liboctave/fNDArray.cc +++ b/liboctave/fNDArray.cc @@ -925,7 +925,7 @@ float tmp; for (octave_idx_type i = 0; i < nel; i++) { - tmp = octave_read_float (is); + tmp = octave_read_value<float> (is); if (is) a.elem (i) = tmp; else
--- a/liboctave/lo-utils.cc +++ b/liboctave/lo-utils.cc @@ -290,8 +290,9 @@ return d; } +template <> double -octave_read_double (std::istream& is) +octave_read_value (std::istream& is) { double d = 0.0; @@ -345,8 +346,9 @@ return d; } +template <> Complex -octave_read_complex (std::istream& is) +octave_read_value (std::istream& is) { double re = 0.0, im = 0.0; @@ -359,12 +361,12 @@ if (ch == '(') { - re = octave_read_double (is); + re = octave_read_value<double> (is); ch = is.get (); if (ch == ',') { - im = octave_read_double (is); + im = octave_read_value<double> (is); ch = is.get (); if (ch == ')') @@ -380,57 +382,13 @@ else { is.putback (ch); - cx = octave_read_double (is); + cx = octave_read_value<double> (is); } return cx; } -void -octave_write_double (std::ostream& os, double d) -{ - if (lo_ieee_is_NA (d)) - os << "NA"; - else if (lo_ieee_isnan (d)) - os << "NaN"; - else if (lo_ieee_isinf (d)) - os << (d < 0 ? "-Inf" : "Inf"); - else - os << d; -} - -void -octave_write_complex (std::ostream& os, const Complex& c) -{ - os << "("; - octave_write_double (os, real (c)); - os << ","; - octave_write_double (os, imag (c)); - os << ")"; -} - - - - - - - - - - - - - - - - - - - - - - static inline float read_float_inf_nan_na (std::istream& is, char c, char sign = '+') { @@ -480,8 +438,9 @@ return d; } +template <> float -octave_read_float (std::istream& is) +octave_read_value (std::istream& is) { float d = 0.0; @@ -535,8 +494,9 @@ return d; } +template <> FloatComplex -octave_read_float_complex (std::istream& is) +octave_read_value (std::istream& is) { float re = 0.0, im = 0.0; @@ -549,12 +509,12 @@ if (ch == '(') { - re = octave_read_float (is); + re = octave_read_value<float> (is); ch = is.get (); if (ch == ',') { - im = octave_read_float (is); + im = octave_read_value<float> (is); ch = is.get (); if (ch == ')') @@ -570,7 +530,7 @@ else { is.putback (ch); - cx = octave_read_float (is); + cx = octave_read_value<float> (is); } return cx; @@ -578,6 +538,29 @@ } void +octave_write_double (std::ostream& os, double d) +{ + if (lo_ieee_is_NA (d)) + os << "NA"; + else if (lo_ieee_isnan (d)) + os << "NaN"; + else if (lo_ieee_isinf (d)) + os << (d < 0 ? "-Inf" : "Inf"); + else + os << d; +} + +void +octave_write_complex (std::ostream& os, const Complex& c) +{ + os << "("; + octave_write_double (os, real (c)); + os << ","; + octave_write_double (os, imag (c)); + os << ")"; +} + +void octave_write_float (std::ostream& os, float d) { if (lo_ieee_is_NA (d))
--- a/liboctave/lo-utils.h +++ b/liboctave/lo-utils.h @@ -26,7 +26,7 @@ #include <cstdio> -#include <iosfwd> +#include <iostream> #include <string> #include "oct-cmplx.h" @@ -62,15 +62,48 @@ extern "C" OCTINTERP_API int octave_strncasecmp (const char *s1, const char *s2, size_t n); -extern OCTAVE_API double octave_read_double (std::istream& is); -extern OCTAVE_API Complex octave_read_complex (std::istream& is); +template <typename T> +T +octave_read_value (std::istream& is) +{ + T retval; + is >> retval; + return retval; +} + +template <> OCTAVE_API double octave_read_value (std::istream& is); +template <> OCTAVE_API Complex octave_read_value (std::istream& is); +template <> OCTAVE_API float octave_read_value (std::istream& is); +template <> OCTAVE_API FloatComplex octave_read_value (std::istream& is); + +// The next four functions are provided for backward compatibility. +inline double +octave_read_double (std::istream& is) +{ + return octave_read_value<double> (is); +} + +inline Complex +octave_read_complex (std::istream& is) +{ + return octave_read_value<Complex> (is); +} + +inline float +octave_read_float (std::istream& is) +{ + return octave_read_value<float> (is); +} + +inline FloatComplex +octave_read_float_complex (std::istream& is) +{ + return octave_read_value<FloatComplex> (is); +} extern OCTAVE_API void octave_write_double (std::ostream& os, double dval); extern OCTAVE_API void octave_write_complex (std::ostream& os, const Complex& cval); -extern OCTAVE_API float octave_read_float (std::istream& is); -extern OCTAVE_API FloatComplex octave_read_float_complex (std::istream& is); - extern OCTAVE_API void octave_write_float (std::ostream& os, float dval); extern OCTAVE_API void octave_write_float_complex (std::ostream& os, const FloatComplex& cval);
--- a/liboctave/sparse-util.cc +++ b/liboctave/sparse-util.cc @@ -1,6 +1,6 @@ /* -Copyright (C) 2005, 2007, 2008 David Bateman +Copyright (C) 2005, 2007, 2008, 2009 David Bateman Copyright (C) 1998, 1999, 2000, 2001, 2002, 2003, 2004, 2005 Andy Adler This file is part of Octave. @@ -58,6 +58,66 @@ } +bool +sparse_indices_ok (octave_idx_type *r, octave_idx_type *c, + octave_idx_type nrows, octave_idx_type ncols, + octave_idx_type nnz) +{ + if (nnz > 0) + { + if (c[0] != 0) + { + (*current_liboctave_error_handler) + ("invalid sparse matrix: cidx[0] must be zero"); + return false; + } + + octave_idx_type jold = 0; + + for (octave_idx_type j = 1; j < ncols+1; j++) + { + if (c[j] < c[j-1]) + { + (*current_liboctave_error_handler) + ("invalid sparse matrix: cidx elements must appear in ascending order"); + return false; + } + + if (c[j] > nnz) + { + (*current_liboctave_error_handler) + ("invalid sparse matrix: cidx[%d] = %d exceeds number of nonzero elements", j, c[j]+1); + return false; + } + + if (c[j] != jold) + { + for (octave_idx_type i = jold+1; i < c[j]; i++) + { + if (r[i] < r[i-1]) + { + (*current_liboctave_error_handler) + ("invalid sparse matrix: ridx elements must appear in ascending order for each column"); + return false; + } + + if (r[i] >= nrows) + { + (*current_liboctave_error_handler) + ("invalid sparse matrix: ridx[%d] = %d out of range", + i, r[i]+1); + return false; + } + } + + jold = c[j]; + } + } + } + + return true; +} + /* ;;; Local Variables: *** ;;; mode: C++ ***
--- a/liboctave/sparse-util.h +++ b/liboctave/sparse-util.h @@ -1,6 +1,6 @@ /* -Copyright (C) 2005, 2007, 2008 David Bateman +Copyright (C) 2005, 2007, 2008, 2009 David Bateman Copyright (C) 1998, 1999, 2000, 2001, 2002, 2003, 2004, 2005 Andy Adler This file is part of Octave. @@ -31,6 +31,11 @@ int line, const char *message); extern OCTAVE_API int SparseCholPrint (const char *fmt, ...); +extern OCTAVE_API bool +sparse_indices_ok (octave_idx_type *r, octave_idx_type *c, + octave_idx_type nrows, octave_idx_type ncols, + octave_idx_type nnz); + #endif /*
--- a/src/ChangeLog +++ b/src/ChangeLog @@ -1,3 +1,20 @@ +2009-07-29 John W. Eaton <jwe@octave.org> + + * ov-bool.cc (octave_bool::load_ascii): Call template function to + read value. + * ov-scalar.cc (octave_scalar::load_ascii): Likewise. + * ov-complex.cc (octave_complex::load_ascii): Likewise. + * ov-float.cc (octave_float_scalar::load_ascii): Likewise. + * ov-flt-complex.cc (octave_float_complex::load_ascii): Likewise. + * ls-mat-ascii.cc (read_mat_ascii_data): Likewise. + + * ov-re-sparse.cc (octave_sparse_matrix::load_binary, + octave_sparse_matrix::load_hdf5): Perform sanity check on indices. + * ov-cx-sparse.cc (octave_sparse_complex_matrix::load_binary, + octave_sparse_complex_matrix::load_hdf5): Likewise. + * ov-bool-sparse.cc (octave_sparse_bool_matrix::load_binary, + octave_sparse_bool_matrix::load_hdf5): Likewise. + 2009-07-29 Jaroslav Hajek <highegg@gmail.com> * ov-fcn-handle.cc (octave_fcn_handle::do_multi_index_op):
--- a/src/ls-mat-ascii.cc +++ b/src/ls-mat-ascii.cc @@ -256,7 +256,7 @@ { OCTAVE_QUIT; - d = octave_read_double (tmp_stream); + d = octave_read_value<double> (tmp_stream); if (tmp_stream || tmp_stream.eof ()) {
--- a/src/ov-bool-sparse.cc +++ b/src/ov-bool-sparse.cc @@ -83,9 +83,9 @@ retval = new octave_bool (tmp (0)); } - else if (matrix.cols () > 0 && matrix.rows () > 0 && - double (matrix.byte_size ()) > double (matrix.rows ()) * - double (matrix.cols ()) * sizeof (bool)) + else if (matrix.cols () > 0 && matrix.rows () > 0 + && (double (matrix.byte_size ()) > double (matrix.rows ()) + * double (matrix.cols ()) * sizeof (bool))) retval = new octave_bool_matrix (matrix.matrix_value ()); } @@ -271,7 +271,7 @@ swap_bytes<4> (&tmp); if (tmp != -2) { - error("load: only 2D sparse matrices are supported"); + error ("load: only 2D sparse matrices are supported"); return false; } @@ -324,6 +324,9 @@ for (int i = 0; i < nz; i++) m.data(i) = (htmp[i] ? 1 : 0); + if (! m.indices_ok ()) + return false; + matrix = m; return true; @@ -605,8 +608,8 @@ H5Sget_simple_extent_dims (space_hid, hdims, maxdims); - if (static_cast<int> (hdims[0]) != nc + 1 || - static_cast<int> (hdims[1]) != 1) + if (static_cast<int> (hdims[0]) != nc + 1 + || static_cast<int> (hdims[1]) != 1) { H5Sclose (space_hid); H5Dclose (data_hid); @@ -640,8 +643,8 @@ H5Sget_simple_extent_dims (space_hid, hdims, maxdims); - if (static_cast<int> (hdims[0]) != nz || - static_cast<int> (hdims[1]) != 1) + if (static_cast<int> (hdims[0]) != nz + || static_cast<int> (hdims[1]) != 1) { H5Sclose (space_hid); H5Dclose (data_hid); @@ -650,7 +653,8 @@ } itmp = m.xridx (); - if (H5Dread (data_hid, H5T_NATIVE_IDX, H5S_ALL, H5S_ALL, H5P_DEFAULT, itmp) < 0) + if (H5Dread (data_hid, H5T_NATIVE_IDX, H5S_ALL, H5S_ALL, + H5P_DEFAULT, itmp) < 0) { H5Sclose (space_hid); H5Dclose (data_hid); @@ -675,8 +679,8 @@ H5Sget_simple_extent_dims (space_hid, hdims, maxdims); - if (static_cast<int> (hdims[0]) != nz || - static_cast<int> (hdims[1]) != 1) + if (static_cast<int> (hdims[0]) != nz + || static_cast<int> (hdims[1]) != 1) { H5Sclose (space_hid); H5Dclose (data_hid); @@ -686,7 +690,9 @@ OCTAVE_LOCAL_BUFFER (hbool_t, htmp, nz); bool retval = false; - if (H5Dread (data_hid, H5T_NATIVE_HBOOL, H5S_ALL, H5S_ALL, H5P_DEFAULT, htmp) >= 0) + if (H5Dread (data_hid, H5T_NATIVE_HBOOL, H5S_ALL, H5S_ALL, + H5P_DEFAULT, htmp) >= 0 + && m.indices_ok ()) { retval = true;
--- a/src/ov-bool.cc +++ b/src/ov-bool.cc @@ -127,7 +127,7 @@ bool octave_bool::load_ascii (std::istream& is) { - scalar = (octave_read_double (is) != 0.); + scalar = (octave_read_value<double> (is) != 0.); if (!is) {
--- a/src/ov-complex.cc +++ b/src/ov-complex.cc @@ -259,7 +259,7 @@ bool octave_complex::load_ascii (std::istream& is) { - scalar = octave_read_complex (is); + scalar = octave_read_value<Complex> (is); if (!is) {
--- a/src/ov-cx-sparse.cc +++ b/src/ov-cx-sparse.cc @@ -81,15 +81,15 @@ else if (nr == 0 || nc == 0) retval = new octave_matrix (Matrix (nr, nc)); else if (matrix.all_elements_are_real ()) - if (matrix.cols () > 0 && matrix.rows () > 0 && - double (matrix.byte_size ()) > double (matrix.rows ()) * - double (matrix.cols ()) * sizeof (double)) + if (matrix.cols () > 0 && matrix.rows () > 0 + && (double (matrix.byte_size ()) > double (matrix.rows ()) + * double (matrix.cols ()) * sizeof (double))) retval = new octave_matrix (::real (matrix.matrix_value ())); else retval = new octave_sparse_matrix (::real (matrix)); - else if (matrix.cols () > 0 && matrix.rows () > 0 && - double (matrix.byte_size ()) > double (matrix.rows ()) * - double (matrix.cols ()) * sizeof (Complex)) + else if (matrix.cols () > 0 && matrix.rows () > 0 + && (double (matrix.byte_size ()) > double (matrix.rows ()) + * double (matrix.cols ()) * sizeof (Complex))) retval = new octave_complex_matrix (matrix.matrix_value ()); } else @@ -312,7 +312,7 @@ swap_bytes<4> (&tmp); if (tmp != -2) { - error("load: only 2D sparse matrices are supported"); + error ("load: only 2D sparse matrices are supported"); return false; } @@ -362,6 +362,10 @@ if (error_state || ! is) return false; + + if (! m.indices_ok ()) + return false; + matrix = m; return true; @@ -681,8 +685,8 @@ H5Sget_simple_extent_dims (space_hid, hdims, maxdims); - if (static_cast<int> (hdims[0]) != nc + 1 || - static_cast<int> (hdims[1]) != 1) + if (static_cast<int> (hdims[0]) != nc + 1 + || static_cast<int> (hdims[1]) != 1) { H5Sclose (space_hid); H5Dclose (data_hid); @@ -716,8 +720,8 @@ H5Sget_simple_extent_dims (space_hid, hdims, maxdims); - if (static_cast<int> (hdims[0]) != nz || - static_cast<int> (hdims[1]) != 1) + if (static_cast<int> (hdims[0]) != nz + || static_cast<int> (hdims[1]) != 1) { H5Sclose (space_hid); H5Dclose (data_hid); @@ -763,8 +767,8 @@ H5Sget_simple_extent_dims (space_hid, hdims, maxdims); - if (static_cast<int> (hdims[0]) != nz || - static_cast<int> (hdims[1]) != 1) + if (static_cast<int> (hdims[0]) != nz + || static_cast<int> (hdims[1]) != 1) { H5Sclose (space_hid); H5Dclose (data_hid); @@ -774,7 +778,9 @@ Complex *ctmp = m.xdata (); bool retval = false; - if (H5Dread (data_hid, complex_type, H5S_ALL, H5S_ALL, H5P_DEFAULT, ctmp) >= 0) + if (H5Dread (data_hid, complex_type, H5S_ALL, H5S_ALL, + H5P_DEFAULT, ctmp) >= 0 + && m.indices_ok ()) { retval = true; matrix = m;
--- a/src/ov-float.cc +++ b/src/ov-float.cc @@ -139,7 +139,7 @@ bool octave_float_scalar::load_ascii (std::istream& is) { - scalar = octave_read_float (is); + scalar = octave_read_value<float> (is); if (!is) { error ("load: failed to load scalar constant");
--- a/src/ov-flt-complex.cc +++ b/src/ov-flt-complex.cc @@ -244,7 +244,7 @@ bool octave_float_complex::load_ascii (std::istream& is) { - scalar = octave_read_float_complex (is); + scalar = octave_read_value<FloatComplex> (is); if (!is) {
--- a/src/ov-re-sparse.cc +++ b/src/ov-re-sparse.cc @@ -83,9 +83,9 @@ retval = new octave_scalar (tmp (0)); } - else if (matrix.cols () > 0 && matrix.rows () > 0 && - double (matrix.byte_size ()) > double (matrix.rows ()) * - double (matrix.cols ()) * sizeof (double)) + else if (matrix.cols () > 0 && matrix.rows () > 0 + && (double (matrix.byte_size ()) > double (matrix.rows ()) + * double (matrix.cols ()) * sizeof (double))) retval = new octave_matrix (matrix.matrix_value ()); } @@ -326,7 +326,7 @@ swap_bytes<4> (&tmp); if (tmp != -2) { - error("load: only 2D sparse matrices are supported"); + error ("load: only 2D sparse matrices are supported"); return false; } @@ -375,6 +375,10 @@ if (error_state || ! is) return false; + + if (! m.indices_ok ()) + return false; + matrix = m; return true; @@ -682,8 +686,8 @@ H5Sget_simple_extent_dims (space_hid, hdims, maxdims); - if (static_cast<int> (hdims[0]) != nc + 1 || - static_cast<int> (hdims[1]) != 1) + if (static_cast<int> (hdims[0]) != nc + 1 + || static_cast<int> (hdims[1]) != 1) { H5Sclose (space_hid); H5Dclose (data_hid); @@ -762,22 +766,20 @@ } double *dtmp = m.xdata (); + bool retval = false; if (H5Dread (data_hid, H5T_NATIVE_DOUBLE, H5S_ALL, H5S_ALL, - H5P_DEFAULT, dtmp) < 0) + H5P_DEFAULT, dtmp) >= 0 + && m.indices_ok ()) { - H5Sclose (space_hid); - H5Dclose (data_hid); - H5Gclose (group_hid); - return false; + retval = true; + matrix = m; } H5Sclose (space_hid); H5Dclose (data_hid); H5Gclose (group_hid); - matrix = m; - - return true; + return retval; } #endif