From 671852e5acb0049cd2a6c70ac8cc3bac4fb41d47 Mon Sep 17 00:00:00 2001 From: Santiago Ospina De Los Rios <sospinar@gmail.com> Date: Thu, 22 Dec 2022 14:32:18 +0530 Subject: [PATCH] Swap first two arguments in the writeSVGMatrix function This is to be consistent with other writers --- dune/istl/io.hh | 33 ++++++++++++++++++++++++++------- dune/istl/test/iotest.cc | 2 +- 2 files changed, 27 insertions(+), 8 deletions(-) diff --git a/dune/istl/io.hh b/dune/istl/io.hh index f0f22c034..c6c5ea7e9 100644 --- a/dune/istl/io.hh +++ b/dune/istl/io.hh @@ -540,10 +540,10 @@ namespace Dune { } //! Writes (recursively) the svg content for an sparse matrix - template <class Mat, class Stream, class SVGMatrixOptions, + template <class Stream, class Mat, class SVGMatrixOptions, class RowPrefix, class ColPrefix> std::pair<std::size_t, size_t> - writeSVGMatrix(const Mat &mat, Stream &out, SVGMatrixOptions opts, + writeSVGMatrix(Stream &out, const Mat &mat, SVGMatrixOptions opts, RowPrefix row_prefix, ColPrefix col_prefix) { // get values to fill the offsets const auto& block_size = opts.block_size; @@ -599,7 +599,7 @@ namespace Dune { NullStream dev0; // get size of sub-block auto sub_size = - writeSVGMatrix(val, dev0, opts, row_prefix, col_prefix); + writeSVGMatrix(dev0, val, opts, row_prefix, col_prefix); // if we didn't see col size before if (col_offsets[col + 1] == null_offset) // write it in the offset vector @@ -636,7 +636,7 @@ namespace Dune { ss << "<svg x='" << col_offsets[col] << "' y='" << row_offsets[row] << "' width='" << width << "' height='" << height << "'>\n"; // write a nested svg with the contents of the sub-block - writeSVGMatrix(val, ss, opts, row_prefix, col_prefix); + writeSVGMatrix(ss, val, opts, row_prefix, col_prefix); ss << "</svg>\n"; }); col_offset = col_offsets.back(); @@ -804,12 +804,31 @@ namespace Dune { * @param opts SVG Options object */ template <class Mat, class SVGOptions = DefaultSVGMatrixOptions> - void writeSVGMatrix(const Mat &mat, std::ostream &out, - SVGOptions opts = {}) { + void writeSVGMatrix(std::ostream &out, const Mat &mat, SVGOptions opts = {}) { // We need a vector that can fit all the multi-indices for rows and columns using IndexPrefix = Dune::ReservedVector<std::size_t, blockLevel<Mat>()>; // Call overload for Mat type - Impl::writeSVGMatrix(mat, out, opts, IndexPrefix{}, IndexPrefix{}); + Impl::writeSVGMatrix(out, mat, opts, IndexPrefix{}, IndexPrefix{}); + } + + /** + * @brief Writes the visualization of matrix in the SVG format + * @details The default visualization writes a rectangle for every block + * bounding box. This is enough to visualize patterns. If you need a + * more advance SVG object in each matrix block (e.g. color scale, or + * write the value in text), just provide a custom SVGOptions that + * fulfills the DefaultSVGMatrixOptions interface. + * + * @tparam Mat Matrix type to write + * @tparam SVGOptions Options object type (see DefaultSVGMatrixOptions) + * @param mat The matrix to write + * @param out A output stream to write SVG to + * @param opts SVG Options object + */ + template <class Mat, class SVGOptions = DefaultSVGMatrixOptions> + [[deprecated("Use signature where std::stream is the first argument.")]] + void writeSVGMatrix(const Mat &mat, std::ostream &out, SVGOptions opts = {}) { + writeSVGMatrix(out, mat, opts); } /** @} end documentation */ diff --git a/dune/istl/test/iotest.cc b/dune/istl/test/iotest.cc index 8f2b9745d..8704bb85a 100644 --- a/dune/istl/test/iotest.cc +++ b/dune/istl/test/iotest.cc @@ -25,7 +25,7 @@ void testWriteMatrix(BlockType b=BlockType(0.0)) A[0][0] += b; writeMatrixToMatlabHelper(A, 0, 0, std::cout); - writeSVGMatrix(A, std::cout); + writeSVGMatrix(std::cout, A); } /* uses the writeVectorToMatlab method, filled with dummy data */ -- GitLab