17#include "istlexception.hh"
20#include <dune/common/hybridutilities.hh>
53 int& counter,
int columns,
int width)
58 if (counter%columns==0)
69 if (counter%columns==0)
75 for (
const auto& entry : v)
89 void printvector (std::ostream& s,
const V& v, std::string title,
90 std::string rowtext,
int columns=1,
int width=10,
97 std::ios_base::fmtflags oldflags = s.flags();
100 s.setf(std::ios_base::scientific, std::ios_base::floatfield);
101 int oldprec = s.precision();
102 s.precision(precision);
105 s << title <<
" [blocks=" << v.N() <<
",dimension=" << v.dim() <<
"]"
112 if (counter%columns!=0)
117 s.precision(oldprec);
133 inline void fill_row (std::ostream& s,
int m,
int width, [[maybe_unused]]
int precision)
135 for (
int j=0; j<m; j++)
151 std::enable_if_t<Dune::IsNumber<K>::value,
int> = 0>
153 [[maybe_unused]]
typename FieldMatrix<K,1,1>::size_type I,
154 [[maybe_unused]]
typename FieldMatrix<K,1,1>::size_type J,
155 [[maybe_unused]]
typename FieldMatrix<K,1,1>::size_type therow,
157 [[maybe_unused]]
int precision)
172 std::enable_if_t<not Dune::IsNumber<M>::value,
int> = 0>
173 void print_row (std::ostream& s,
const M& A,
typename M::size_type I,
174 typename M::size_type J,
typename M::size_type therow,
175 int width,
int precision)
177 typename M::size_type i0=I;
178 for (
typename M::size_type i=0; i<A.
N(); i++)
180 if (therow>=i0 && therow<i0+MatrixDimension<M>::rowdim(A,i))
183 typename M::size_type j0=J;
184 for (
typename M::size_type j=0; j<A.
M(); j++)
187 typename M::ConstColIterator it = A[i].find(j);
191 print_row(s,*it,i0,j0,therow,width,precision);
193 fill_row(s,MatrixDimension<M>::coldim(A,j),width,precision);
196 j0 += MatrixDimension<M>::coldim(A,j);
200 i0 += MatrixDimension<M>::rowdim(A,i);
214 std::string rowtext,
int width=10,
int precision=2)
218 std::ios_base::fmtflags oldflags = s.flags();
221 s.setf(std::ios_base::scientific, std::ios_base::floatfield);
222 int oldprec = s.precision();
223 s.precision(precision);
229 <<
",rowdim=" << MatrixDimension<M>::rowdim(A)
230 <<
",coldim=" << MatrixDimension<M>::coldim(A)
234 for (
typename M::size_type i=0; i<MatrixDimension<M>::rowdim(A); i++)
246 s.precision(oldprec);
252 void printInnerMatrixElement(std::ostream& s,
253 const B& innerMatrixElement,
254 int innerrow,
int innercol)
256 s<<innerMatrixElement<<
" ";
259 template<
class B,
int n>
260 void printInnerMatrixElement(std::ostream& s,
261 const ScaledIdentityMatrix<B,n> innerMatrixElement,
262 int innerrow,
int innercol)
264 if (innerrow == innercol)
265 s<<innerMatrixElement.scalar()<<
" ";
270 template<
class B,
int n,
int m>
271 void printInnerMatrixElement(std::ostream& s,
272 const FieldMatrix<B,n,m> innerMatrixElement,
273 int innerrow,
int innercol)
275 s<<innerMatrixElement[innerrow][innercol]<<
" ";
300 template<
class A,
class InnerMatrixType>
303 std::string title, std::string rowtext,
304 int width=3,
int precision=2)
308 std::ios_base::fmtflags oldflags = s.flags();
310 s.setf(std::ios_base::scientific, std::ios_base::floatfield);
311 int oldprec = s.precision();
312 s.precision(precision);
317 <<
",rowdim=" << MatrixDimension<Matrix>::rowdim(mat)
318 <<
",coldim=" << MatrixDimension<Matrix>::coldim(mat)
323 constexpr int n = std::decay_t<decltype(Impl::asMatrix(std::declval<InnerMatrixType>()))>::rows;
324 constexpr int m = std::decay_t<decltype(Impl::asMatrix(std::declval<InnerMatrixType>()))>::cols;
325 for(Row row=mat.
begin(); row != mat.
end(); ++row) {
327 bool reachedEnd=
false;
330 for(
int innerrow=0; innerrow<n; ++innerrow) {
333 Col col=row->begin();
334 for(; col != row->end(); ++col,++count) {
337 if(count>=skipcols+width)
340 if(count==skipcols) {
344 s << row.index()<<
": ";
347 s<<col.index()<<
": |";
349 if(count==skipcols) {
350 for(
typename std::string::size_type i=0; i < rowtext.length(); i++)
356 for(
int innercol=0; innercol < m; ++innercol) {
358 Impl::printInnerMatrixElement(s,*col,innerrow,innercol);
363 if(innerrow==n-1 && col==row->end())
376 s.precision(oldprec);
382 struct MatlabPODWriter
384 static std::ostream& write(
const T& t, std::ostream& s)
391 struct MatlabPODWriter<
std::complex<T> >
393 static std::ostream& write(
const std::complex<T>& t, std::ostream& s)
395 s << t.real() <<
" " << t.imag();
410 template <
class FieldType,
411 std::enable_if_t<Dune::IsNumber<FieldType>::value,
int> = 0>
413 int rowOffset,
int colOffset,
417 s << rowOffset + 1 <<
" " << colOffset + 1 <<
" ";
418 MatlabPODWriter<FieldType>::write(value, s)<< std::endl;
428 template <
class MatrixType,
429 std::enable_if_t<not Dune::IsNumber<MatrixType>::value,
int> = 0>
431 int externalRowOffset,
int externalColOffset,
435 std::vector<typename MatrixType::size_type> colOffset(matrix.M());
436 if (colOffset.size() > 0)
439 for (
typename MatrixType::size_type i=0; i<matrix.M()-1; i++)
440 colOffset[i+1] = colOffset[i] +
441 MatrixDimension<MatrixType>::coldim(matrix,i);
443 typename MatrixType::size_type rowOffset = 0;
446 for (
typename MatrixType::size_type rowIdx=0; rowIdx<matrix.N(); rowIdx++)
448 auto cIt = matrix[rowIdx].begin();
449 auto cEndIt = matrix[rowIdx].end();
452 for (; cIt!=cEndIt; ++cIt)
454 externalRowOffset+rowOffset,
455 externalColOffset + colOffset[cIt.index()],
458 rowOffset += MatrixDimension<MatrixType>::rowdim(matrix, rowIdx);
482 template <
class MatrixType>
484 const std::string& filename,
int outputPrecision = 18)
486 std::ofstream outStream(filename.c_str());
487 int oldPrecision = outStream.precision();
488 outStream.precision(outputPrecision);
491 outStream.precision(oldPrecision);
496 void writeVectorToMatlabHelper (
const V& v, std::ostream& stream)
498 if constexpr (IsNumber<V>()) {
499 stream << v << std::endl;
501 for (
const auto& entry : v)
502 writeVectorToMatlabHelper(entry, stream);
523 template <
class VectorType>
525 const std::string& filename,
int outputPrecision = 18)
527 std::ofstream outStream(filename.c_str());
528 int oldPrecision = outStream.precision();
529 outStream.precision(outputPrecision);
531 writeVectorToMatlabHelper(vector, outStream);
532 outStream.precision(oldPrecision);
540 friend NullStream &operator<<(NullStream &dev0, Any &&) {
547 template <
class Stream,
class SVGMatrixOptions>
548 void writeSVGMatrixHeader(Stream &out,
const SVGMatrixOptions &opts,
549 std::pair<std::size_t, size_t> offsets) {
550 auto [col_offset, row_offset] = offsets;
551 double width = opts.width;
552 double height = opts.height;
554 if (opts.width == 0 and opts.height == 0)
555 width = height = 500;
557 width = opts.height * (double(col_offset) / row_offset);
558 if (opts.height == 0)
559 height = opts.width * (double(row_offset) / col_offset);
562 double scale_width = width / col_offset;
563 double scale_height = height / row_offset;
566 out <<
"<svg xmlns='http://www.w3.org/2000/svg' width='" << std::ceil(width)
567 <<
"' height='" << std::ceil(height) <<
"' version='1.1'>\n"
569 << opts.style <<
"</style>\n"
570 <<
"<g transform='scale(" << scale_width <<
" " << scale_height
575 template <
class Stream,
class Mat,
class SVGMatrixOptions,
576 class RowPrefix,
class ColPrefix>
577 std::pair<std::size_t, size_t>
578 writeSVGMatrix(Stream &out,
const Mat &mat, SVGMatrixOptions opts,
579 RowPrefix row_prefix, ColPrefix col_prefix) {
581 const auto& block_size = opts.block_size;
582 const auto& interspace = opts.interspace;
584 const std::size_t rows = mat.N();
585 const std::size_t cols = mat.M();
588 const bool write_header = opts.write_header;
589 opts.write_header =
false;
592 std::size_t row_offset = interspace;
593 std::size_t col_offset = interspace;
596 auto for_each_entry = [&mat](
const auto &call_back) {
597 for (
auto row_it = mat.begin(); row_it != mat.end(); ++row_it) {
598 for (
auto col_it = row_it->begin(); col_it != row_it->end(); ++col_it) {
599 call_back(row_it.index(), col_it.index(), *col_it);
605 std::stringstream ss;
608 row_prefix.push_back(0);
609 col_prefix.push_back(0);
612 if constexpr (Dune::blockLevel<typename Mat::block_type>() == 0) {
614 for_each_entry([&](
const auto &row,
const auto &col,
const auto &val) {
615 std::size_t x_off = interspace + col * (interspace + block_size);
616 std::size_t y_off = interspace + row * (interspace + block_size);
617 row_prefix.back() = row;
618 col_prefix.back() = col;
619 opts.writeSVGBlock(ss, row_prefix, col_prefix, val,
620 {x_off, y_off, block_size, block_size});
622 col_offset += cols * (block_size + interspace);
623 row_offset += rows * (block_size + interspace);
628 std::vector<std::size_t> col_offsets(cols + 1, null_offset);
629 std::vector<std::size_t> row_offsets(rows + 1, null_offset);
630 for_each_entry([&](
const auto &row,
const auto &col,
const auto &val) {
637 if (col_offsets[col + 1] == null_offset)
638 col_offsets[col + 1] = sub_size.first;
641 if (row_offsets[row + 1] == null_offset)
642 row_offsets[row + 1] = sub_size.second;
646 auto min_row_offset = *std::min_element(begin(row_offsets), end(row_offsets));
647 std::replace(begin(row_offsets), end(row_offsets), null_offset, min_row_offset);
648 auto min_col_offset = *std::min_element(begin(col_offsets), end(col_offsets));
649 std::replace(begin(col_offsets), end(col_offsets), null_offset, min_col_offset);
652 col_offsets[0] = interspace;
653 row_offsets[0] = interspace;
654 for (std::size_t i = 1; i < col_offsets.size(); i++)
655 col_offsets[i] += col_offsets[i - 1] + interspace;
656 for (std::size_t i = 1; i < row_offsets.size(); i++)
657 row_offsets[i] += row_offsets[i - 1] + interspace;
659 for_each_entry([&](
const auto &row,
const auto &col,
const auto &val) {
662 col_offsets[col + 1] - col_offsets[col] - interspace;
664 row_offsets[row + 1] - row_offsets[row] - interspace;
665 row_prefix.back() = row;
666 col_prefix.back() = col;
668 ss <<
"<svg x='" << col_offsets[col] <<
"' y='" << row_offsets[row]
669 <<
"' width='" << width <<
"' height='" << height <<
"'>\n";
674 col_offset = col_offsets.back();
675 row_offset = row_offsets.back();
681 writeSVGMatrixHeader(out, opts, {col_offset, row_offset});
683 col_prefix.pop_back();
684 row_prefix.pop_back();
686 opts.writeSVGBlock(out, row_prefix, col_prefix, mat,
687 {0, 0, col_offset, row_offset});
692 out <<
"</g>\n</svg>\n";
695 return {col_offset, row_offset};
718 std::string
style =
" .matrix-block {\n"
719 " fill: cornflowerblue;\n"
720 " fill-opacity: 0.4;\n"
721 " stroke-width: 2;\n"
723 " stroke-opacity: 0.5;\n"
725 " .matrix-block:hover {\n"
726 " fill: lightcoral;\n"
727 " fill-opacity: 0.4;\n"
728 " stroke-opacity: 1;\n"
749 template <
class RowPrefix,
class ColPrefix>
751 const ColPrefix &col_prefix)
const {
753 return "matrix-block";
764 template <
class Stream,
class RowPrefix,
class ColPrefix,
class Block>
766 const ColPrefix &col_prefix,
767 const Block &block)
const {
768 if (this->write_block_title) {
770 assert(row_prefix.size() == col_prefix.size());
771 for (std::size_t i = 0; i < row_prefix.size(); ++i)
772 out <<
"[" << row_prefix[i] <<
", "<< col_prefix[i] <<
"]";
773 if constexpr (Dune::blockLevel<Block>() == 0)
774 out <<
": " << block;
800 template <
class Stream,
class RowPrefix,
class ColPrefix,
class Block>
802 const RowPrefix &row_prefix,
803 const ColPrefix &col_prefix,
const Block block,
804 const std::array<std::size_t, 4> &svg_box)
const {
808 std::string block_class = this->
blockStyleClass(row_prefix, col_prefix);
810 out <<
"<rect class='" << block_class <<
"' x='" << x_off <<
"' y='"
811 << y_off <<
"' width='" <<
width <<
"' height='" <<
height <<
"'";
812 if constexpr (Dune::blockLevel<Block>() == 0 and std::is_convertible<Block,double>{})
814 out <<
" style='fill-opacity: 1;fill:" <<
color_fill(
double(block)) <<
"'";
838 template <
class Mat,
class SVGOptions = DefaultSVGMatrixOptions>
863 template <
class Mat,
class SVGOptions = DefaultSVGMatrixOptions>
864 [[deprecated(
"Use signature where std::stream is the first argument. This will be removed after Dune 2.10.")]]
Helper functions for determining the vector/matrix block level.
A sparse block matrix with compressed row storage.
Definition: bcrsmatrix.hh:466
ConstIterator class for sequential access.
Definition: matrix.hh:404
A generic dynamic dense matrix.
Definition: matrix.hh:561
RowIterator end()
Get iterator to one beyond last row.
Definition: matrix.hh:616
RowIterator begin()
Get iterator to first row.
Definition: matrix.hh:610
row_type::const_iterator ConstColIterator
Const iterator for the entries of each row.
Definition: matrix.hh:589
size_type M() const
Return the number of columns.
Definition: matrix.hh:696
size_type N() const
Return the number of rows.
Definition: matrix.hh:691
A Vector class with statically reserved memory.
Definition: reservedvector.hh:47
Implementation of the BCRSMatrix class.
Implements a matrix constructed from a given type representing a field and compile-time given number ...
Implements a vector constructed from a given type representing a field and a compile-time given size.
constexpr auto max
Function object that returns the greater of the given values.
Definition: hybridutilities.hh:484
void writeSVGMatrix(const Mat &mat, std::ostream &out, SVGOptions opts={})
Writes the visualization of matrix in the SVG format.
Definition: io.hh:865
void writeMatrixToMatlab(const MatrixType &matrix, const std::string &filename, int outputPrecision=18)
Writes sparse matrix in a Matlab-readable format.
Definition: io.hh:483
void writeMatrixToMatlabHelper(const FieldType &value, int rowOffset, int colOffset, std::ostream &s)
Helper method for the writeMatrixToMatlab routine.
Definition: io.hh:412
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)
Print one row of a matrix, specialization for number types.
Definition: io.hh:152
void printmatrix(std::ostream &s, const M &A, std::string title, std::string rowtext, int width=10, int precision=2)
Print a generic block matrix.
Definition: io.hh:213
void printSparseMatrix(std::ostream &s, const BCRSMatrix< InnerMatrixType, A > &mat, std::string title, std::string rowtext, int width=3, int precision=2)
Prints a BCRSMatrix with fixed sized blocks.
Definition: io.hh:301
void printvector(std::ostream &s, const V &v, std::string title, std::string rowtext, int columns=1, int width=10, int precision=2)
Print an ISTL vector.
Definition: io.hh:89
void writeSVGMatrix(std::ostream &out, const Mat &mat, SVGOptions opts={})
Writes the visualization of matrix in the SVG format.
Definition: io.hh:839
void writeVectorToMatlab(const VectorType &vector, const std::string &filename, int outputPrecision=18)
Writes vectors in a Matlab-readable format.
Definition: io.hh:524
void recursive_printvector(std::ostream &s, const V &v, std::string rowtext, int &counter, int columns, int width)
Recursively print a vector.
Definition: io.hh:52
void fill_row(std::ostream &s, int m, int width, int precision)
Print a row of zeros for a non-existing block.
Definition: io.hh:133
Some handy generic functions for ISTL matrices.
Dune namespace.
Definition: alignedallocator.hh:13
An stl-compliant random-access container which stores everything on the stack.
Default options class to write SVG matrices.
Definition: io.hh:706
std::string style
CSS style block to write in header.
Definition: io.hh:718
std::size_t width
Final width size (pixels) of the SVG header. If 0, size is automatic.
Definition: io.hh:712
std::function< std::string(const double &)> color_fill
Color fill for default options.
Definition: io.hh:742
std::size_t block_size
size (pixels) of the deepst block/value of the matrix
Definition: io.hh:708
void writeSVGBlock(Stream &out, const RowPrefix &row_prefix, const ColPrefix &col_prefix, const Block block, const std::array< std::size_t, 4 > &svg_box) const
Write an SVG object for a given block/value in the matrix.
Definition: io.hh:801
void writeBlockTitle(Stream &out, const RowPrefix &row_prefix, const ColPrefix &col_prefix, const Block &block) const
Helper function writes a title for a given block and prefix.
Definition: io.hh:765
std::size_t interspace
size (pixels) of the interspace between blocks
Definition: io.hh:710
bool write_block_title
(Helper) Whether to write a title on the rectangle value
Definition: io.hh:757
std::size_t height
Final height size (pixels) of the SVG header. If 0, size is automatic.
Definition: io.hh:714
bool write_header
Whether to write the SVG header.
Definition: io.hh:716
std::string blockStyleClass(const RowPrefix &row_prefix, const ColPrefix &col_prefix) const
Helper function that returns an style class for a given prefix.
Definition: io.hh:750
Whether this type acts as a scalar in the context of (hierarchically blocked) containers.
Definition: typetraits.hh:194