diff --git a/dune/istl/io.hh b/dune/istl/io.hh index 63dc737aff8808041c5a6e082a3ed531195c3f09..61bcf8d58e2deecdcf20c8e270aa3465dd3178b5 100644 --- a/dune/istl/io.hh +++ b/dune/istl/io.hh @@ -15,6 +15,7 @@ #include "istlexception.hh" #include <dune/common/fvector.hh> #include <dune/common/fmatrix.hh> +#include <dune/common/dynmatrix.hh> #include <dune/common/diagonalmatrix.hh> #include <dune/common/unused.hh> @@ -473,6 +474,27 @@ namespace Dune { } } + /** + * \brief Helper method for the writeMatrixToMatlab routine. + * + * \code + * #include <dune/istl/io.hh> + * \endcode + * + * This specialization for DynamicMatrices ends the recursion + */ + template <class FieldType> + void writeMatrixToMatlabHelper(const DynamicMatrix<FieldType>& matrix, int rowOffset, + int colOffset, std::ostream& s) + { + for (int i=0; i<matrix.N(); i++) + for (int j=0; j<matrix.M(); j++) { + //+1 for Matlab numbering + s << rowOffset + i + 1 << " " << colOffset + j + 1 << " "; + MatlabPODWriter<FieldType>::write(matrix[i][j], s)<< std::endl; + } + } + /** * \brief Helper method for the writeMatrixToMatlab routine. * diff --git a/dune/istl/matrixutils.hh b/dune/istl/matrixutils.hh index fe888eae3aebd1e7ed609d2c506505ed23b763fc..22f808b88dd25d01e0e0885a19c9df675f5dfd29 100644 --- a/dune/istl/matrixutils.hh +++ b/dune/istl/matrixutils.hh @@ -8,6 +8,7 @@ #include <limits> #include <dune/common/typetraits.hh> #include <dune/common/fmatrix.hh> +#include <dune/common/dynmatrix.hh> #include <dune/common/diagonalmatrix.hh> #include <dune/common/unused.hh> #include <dune/istl/scaledidmatrix.hh> @@ -361,6 +362,33 @@ namespace Dune } }; + template <class T> + struct MatrixDimension<Dune::DynamicMatrix<T> > + { + typedef Dune::DynamicMatrix<T> MatrixType; + typedef typename MatrixType::size_type size_type; + + static size_type rowdim(const MatrixType& /*A*/, size_type /*r*/) + { + return 1; + } + + static size_type coldim(const MatrixType& /*A*/, size_type /*r*/) + { + return 1; + } + + static size_type rowdim(const MatrixType& A) + { + return A.N(); + } + + static size_type coldim(const MatrixType& A) + { + return A.M(); + } + }; + template<typename K, int n, int m, typename TA> struct MatrixDimension<Matrix<FieldMatrix<K,n,m>, TA> > {