From 142b88f265a3aa1b2302d234bed63615f6e3b50f Mon Sep 17 00:00:00 2001 From: Oliver Sander <oliver.sander@tu-dresden.de> Date: Sat, 29 Dec 2018 09:39:15 +0100 Subject: [PATCH] Implement the 'printmatrix' method for scalar matrix entries Simplify the implementation while we're there. The IsNumber trait allows much simpler recursion termination. --- dune/istl/io.hh | 82 +++++++++++++--------------------------- dune/istl/test/iotest.cc | 17 +++++++++ 2 files changed, 44 insertions(+), 55 deletions(-) diff --git a/dune/istl/io.hh b/dune/istl/io.hh index ccb64a556..8c5c5b4c2 100644 --- a/dune/istl/io.hh +++ b/dune/istl/io.hh @@ -151,6 +151,31 @@ namespace Dune { } } + /** + * \brief Print one row of a matrix, specialization for number types + * + * \code + * #include <dune/istl/io.hh> + * \endcode + */ + template<class K> + void print_row (std::ostream& s, const K& value, + typename FieldMatrix<K,1,1>::size_type I, + typename FieldMatrix<K,1,1>::size_type J, + typename FieldMatrix<K,1,1>::size_type therow, + int width, int precision, + typename std::enable_if_t<Dune::IsNumber<K>::value>* sfinae = nullptr) + { + DUNE_UNUSED_PARAMETER(I); + DUNE_UNUSED_PARAMETER(J); + DUNE_UNUSED_PARAMETER(therow); + DUNE_UNUSED_PARAMETER(precision); + + s << " "; // space in front of each entry + s.width(width); // set width for each entry anew + s << value; + } + /** * \brief Print one row of a matrix * @@ -161,7 +186,8 @@ namespace Dune { template<class M> void print_row (std::ostream& s, const M& A, typename M::size_type I, typename M::size_type J, typename M::size_type therow, - int width, int precision) + int width, int precision, + typename std::enable_if_t<!Dune::IsNumber<M>::value>* sfinae = nullptr) { typename M::size_type i0=I; for (typename M::size_type i=0; i<A.N(); i++) @@ -190,60 +216,6 @@ namespace Dune { } } - /** - * \brief Print one row of a matrix, specialization for FieldMatrix - * - * \code - * #include <dune/istl/io.hh> - * \endcode - */ - template<class K, int n, int m> - void print_row (std::ostream& s, const FieldMatrix<K,n,m>& A, - typename FieldMatrix<K,n,m>::size_type I, - typename FieldMatrix<K,n,m>::size_type J, - typename FieldMatrix<K,n,m>::size_type therow, int width, - int precision) - { - DUNE_UNUSED_PARAMETER(J); - DUNE_UNUSED_PARAMETER(precision); - - typedef typename FieldMatrix<K,n,m>::size_type size_type; - - for (size_type i=0; i<n; i++) - if (I+i==therow) - for (int j=0; j<m; j++) - { - s << " "; // space in front of each entry - s.width(width); // set width for each entry anew - s << A[i][j]; // yeah, the number ! - } - } - - /** - * \brief Print one row of a matrix, specialization for FieldMatrix<K,1,1> - * - * \code - * #include <dune/istl/io.hh> - * \endcode - */ - template<class K> - void print_row (std::ostream& s, const FieldMatrix<K,1,1>& A, - typename FieldMatrix<K,1,1>::size_type I, - typename FieldMatrix<K,1,1>::size_type J, - typename FieldMatrix<K,1,1>::size_type therow, - int width, int precision) - { - DUNE_UNUSED_PARAMETER(J); - DUNE_UNUSED_PARAMETER(precision); - - if (I==therow) - { - s << " "; // space in front of each entry - s.width(width); // set width for each entry anew - s << static_cast<K>(A); // yeah, the number ! - } - } - /** * \brief Print a generic block matrix * diff --git a/dune/istl/test/iotest.cc b/dune/istl/test/iotest.cc index aa471835c..f82666f48 100644 --- a/dune/istl/test/iotest.cc +++ b/dune/istl/test/iotest.cc @@ -83,4 +83,21 @@ int main(int argc, char** argv) v3.push_back(1); v3.push_back(2); Dune::writeVectorToMatlabHelper(v3, std::cout); + + // Test the printmatrix method + { + Dune::BCRSMatrix<double> matrix; + setupLaplacian(matrix, 3); + Dune::printmatrix(std::cout, matrix, "matrix", "--"); + } + { + Dune::BCRSMatrix<Dune::FieldMatrix<double,1,1> > matrix; + setupLaplacian(matrix, 3); + Dune::printmatrix(std::cout, matrix, "matrix", "--"); + } + { + Dune::BCRSMatrix<Dune::FieldMatrix<double,2,3> > matrix; + setupLaplacian(matrix, 3); + Dune::printmatrix(std::cout, matrix, "matrix", "--"); + } } -- GitLab