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++)
152 [[maybe_unused]]
typename FieldMatrix<K,1,1>::size_type I,
153 [[maybe_unused]]
typename FieldMatrix<K,1,1>::size_type J,
154 [[maybe_unused]]
typename FieldMatrix<K,1,1>::size_type therow,
156 [[maybe_unused]]
int precision,
172 void print_row (std::ostream& s,
const M& A,
typename M::size_type I,
173 typename M::size_type J,
typename M::size_type therow,
174 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);
270 template<
class B,
int n,
int m,
class A>
273 std::string title, std::string rowtext,
274 int width=3,
int precision=2)
278 std::ios_base::fmtflags oldflags = s.flags();
280 s.setf(std::ios_base::scientific, std::ios_base::floatfield);
281 int oldprec = s.precision();
282 s.precision(precision);
287 <<
",rowdim=" << MatrixDimension<Matrix>::rowdim(mat)
288 <<
",coldim=" << MatrixDimension<Matrix>::coldim(mat)
293 for(Row row=mat.
begin(); row != mat.
end(); ++row) {
295 bool reachedEnd=
false;
298 for(
int innerrow=0; innerrow<n; ++innerrow) {
301 Col col=row->begin();
302 for(; col != row->end(); ++col,++count) {
305 if(count>=skipcols+width)
308 if(count==skipcols) {
312 s << row.index()<<
": ";
315 s<<col.index()<<
": |";
317 if(count==skipcols) {
318 for(
typename std::string::size_type i=0; i < rowtext.length(); i++)
324 for(
int innercol=0; innercol < m; ++innercol) {
326 s<<(*col)[innerrow][innercol]<<
" ";
331 if(innerrow==n-1 && col==row->end())
344 s.precision(oldprec);
350 struct MatlabPODWriter
352 static std::ostream& write(
const T& t, std::ostream& s)
359 struct MatlabPODWriter<
std::complex<T> >
361 static std::ostream& write(
const std::complex<T>& t, std::ostream& s)
363 s << t.real() <<
" " << t.imag();
378 template <
class FieldType>
380 int rowOffset,
int colOffset,
385 s << rowOffset + 1 <<
" " << colOffset + 1 <<
" ";
386 MatlabPODWriter<FieldType>::write(value, s)<< std::endl;
396 template <
class MatrixType>
398 int externalRowOffset,
int externalColOffset,
403 std::vector<typename MatrixType::size_type> colOffset(matrix.M());
404 if (colOffset.size() > 0)
407 for (
typename MatrixType::size_type i=0; i<matrix.M()-1; i++)
408 colOffset[i+1] = colOffset[i] +
409 MatrixDimension<MatrixType>::coldim(matrix,i);
411 typename MatrixType::size_type rowOffset = 0;
414 for (
typename MatrixType::size_type rowIdx=0; rowIdx<matrix.N(); rowIdx++)
416 auto cIt = matrix[rowIdx].begin();
417 auto cEndIt = matrix[rowIdx].end();
420 for (; cIt!=cEndIt; ++cIt)
422 externalRowOffset+rowOffset,
423 externalColOffset + colOffset[cIt.index()],
426 rowOffset += MatrixDimension<MatrixType>::rowdim(matrix, rowIdx);
450 template <
class MatrixType>
452 const std::string& filename,
int outputPrecision = 18)
454 std::ofstream outStream(filename.c_str());
455 int oldPrecision = outStream.precision();
456 outStream.precision(outputPrecision);
459 outStream.precision(oldPrecision);
464 void writeVectorToMatlabHelper (
const V& v, std::ostream& stream)
466 if constexpr (IsNumber<V>()) {
467 stream << v << std::endl;
469 for (
const auto& entry : v)
470 writeVectorToMatlabHelper(entry, stream);
491 template <
class VectorType>
493 const std::string& filename,
int outputPrecision = 18)
495 std::ofstream outStream(filename.c_str());
496 int oldPrecision = outStream.precision();
497 outStream.precision(outputPrecision);
499 writeVectorToMatlabHelper(vector, outStream);
500 outStream.precision(oldPrecision);
508 friend NullStream &operator<<(NullStream &dev0, Any &&) {
515 template <
class Stream,
class SVGMatrixOptions>
516 void writeSVGMatrixHeader(Stream &out,
const SVGMatrixOptions &opts,
517 std::pair<std::size_t, size_t> offsets) {
518 auto [col_offset, row_offset] = offsets;
519 double width = opts.width;
520 double height = opts.height;
522 if (opts.width == 0 and opts.height == 0)
523 width = height = 500;
525 width = opts.height * (double(col_offset) / row_offset);
526 if (opts.height == 0)
527 height = opts.width * (double(row_offset) / col_offset);
530 double scale_width = width / col_offset;
531 double scale_height = height / row_offset;
534 out <<
"<svg xmlns='http://www.w3.org/2000/svg' width='" << std::ceil(width)
535 <<
"' height='" << std::ceil(height) <<
"' version='1.1'>\n"
537 << opts.style <<
"</style>\n"
538 <<
"<g transform='scale(" << scale_width <<
" " << scale_height
543 template <
class Mat,
class Stream,
class SVGMatrixOptions,
544 class RowPrefix,
class ColPrefix>
545 std::pair<std::size_t, size_t>
546 writeSVGMatrix(
const Mat &mat, Stream &out, SVGMatrixOptions opts,
547 RowPrefix row_prefix, ColPrefix col_prefix) {
549 const auto& block_size = opts.block_size;
550 const auto& interspace = opts.interspace;
552 const std::size_t rows = mat.N();
553 const std::size_t cols = mat.M();
556 const bool write_header = opts.write_header;
557 opts.write_header =
false;
560 std::size_t row_offset = interspace;
561 std::size_t col_offset = interspace;
564 auto for_each_entry = [&mat](
const auto &call_back) {
565 for (
auto row_it = mat.begin(); row_it != mat.end(); ++row_it) {
566 for (
auto col_it = row_it->begin(); col_it != row_it->end(); ++col_it) {
567 call_back(row_it.index(), col_it.index(), *col_it);
573 std::stringstream ss;
576 row_prefix.push_back(0);
577 col_prefix.push_back(0);
580 if constexpr (Dune::blockLevel<typename Mat::block_type>() == 0) {
582 for_each_entry([&](
const auto &row,
const auto &col,
const auto &val) {
583 std::size_t x_off = interspace + col * (interspace + block_size);
584 std::size_t y_off = interspace + row * (interspace + block_size);
585 row_prefix.back() = row;
586 col_prefix.back() = col;
587 opts.writeSVGBlock(ss, row_prefix, col_prefix, val,
588 {x_off, y_off, block_size, block_size});
590 col_offset += cols * (block_size + interspace);
591 row_offset += rows * (block_size + interspace);
596 std::vector<std::size_t> col_offsets(cols + 1, null_offset);
597 std::vector<std::size_t> row_offsets(rows + 1, null_offset);
598 for_each_entry([&](
const auto &row,
const auto &col,
const auto &val) {
605 if (col_offsets[col + 1] == null_offset)
606 col_offsets[col + 1] = sub_size.first;
609 if (row_offsets[row + 1] == null_offset)
610 row_offsets[row + 1] = sub_size.second;
614 auto min_row_offset = *std::min_element(begin(row_offsets), end(row_offsets));
615 std::replace(begin(row_offsets), end(row_offsets), null_offset, min_row_offset);
616 auto min_col_offset = *std::min_element(begin(col_offsets), end(col_offsets));
617 std::replace(begin(col_offsets), end(col_offsets), null_offset, min_col_offset);
620 col_offsets[0] = interspace;
621 row_offsets[0] = interspace;
622 for (std::size_t i = 1; i < col_offsets.size(); i++)
623 col_offsets[i] += col_offsets[i - 1] + interspace;
624 for (std::size_t i = 1; i < row_offsets.size(); i++)
625 row_offsets[i] += row_offsets[i - 1] + interspace;
627 for_each_entry([&](
const auto &row,
const auto &col,
const auto &val) {
630 col_offsets[col + 1] - col_offsets[col] - interspace;
632 row_offsets[row + 1] - row_offsets[row] - interspace;
633 row_prefix.back() = row;
634 col_prefix.back() = col;
636 ss <<
"<svg x='" << col_offsets[col] <<
"' y='" << row_offsets[row]
637 <<
"' width='" << width <<
"' height='" << height <<
"'>\n";
642 col_offset = col_offsets.back();
643 row_offset = row_offsets.back();
649 writeSVGMatrixHeader(out, opts, {col_offset, row_offset});
651 col_prefix.pop_back();
652 row_prefix.pop_back();
654 opts.writeSVGBlock(out, row_prefix, col_prefix, mat,
655 {0, 0, col_offset, row_offset});
660 out <<
"</g>\n</svg>\n";
663 return {col_offset, row_offset};
686 std::string
style =
" .matrix-block {\n"
687 " fill: cornflowerblue;\n"
688 " fill-opacity: 0.4;\n"
689 " stroke-width: 2;\n"
691 " stroke-opacity: 0.5;\n"
693 " .matrix-block:hover {\n"
694 " fill: lightcoral;\n"
695 " fill-opacity: 0.4;\n"
696 " stroke-opacity: 1;\n"
717 template <
class RowPrefix,
class ColPrefix>
719 const ColPrefix &col_prefix)
const {
721 return "matrix-block";
732 template <
class Stream,
class RowPrefix,
class ColPrefix,
class Block>
734 const ColPrefix &col_prefix,
735 const Block &block)
const {
736 if (this->write_block_title) {
738 assert(row_prefix.size() == col_prefix.size());
739 for (std::size_t i = 0; i < row_prefix.size(); ++i)
740 out <<
"[" << row_prefix[i] <<
", "<< col_prefix[i] <<
"]";
741 if constexpr (Dune::blockLevel<Block>() == 0)
742 out <<
": " << block;
768 template <
class Stream,
class RowPrefix,
class ColPrefix,
class Block>
770 const RowPrefix &row_prefix,
771 const ColPrefix &col_prefix,
const Block block,
772 const std::array<std::size_t, 4> &svg_box)
const {
776 std::string block_class = this->
blockStyleClass(row_prefix, col_prefix);
778 out <<
"<rect class='" << block_class <<
"' x='" << x_off <<
"' y='"
779 << y_off <<
"' width='" <<
width <<
"' height='" <<
height <<
"'";
780 if constexpr (Dune::blockLevel<Block>() == 0 and std::is_convertible<Block,double>{})
782 out <<
" style='fill-opacity: 1;fill:" <<
color_fill(
double(block)) <<
"'";
806 template <
class Mat,
class SVGOptions = DefaultSVGMatrixOptions>
808 SVGOptions opts = {}) {
Implementation of the BCRSMatrix class.
Helper functions for determining the vector/matrix block level.
A sparse block matrix with compressed row storage.
Definition: bcrsmatrix.hh:466
A dense n x m matrix.
Definition: fmatrix.hh:117
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:620
RowIterator begin()
Get iterator to first row.
Definition: matrix.hh:614
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:700
size_type N() const
Return the number of rows.
Definition: matrix.hh:695
A Vector class with statically reserved memory.
Definition: reservedvector.hh:47
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.
void writeSVGMatrix(const Mat &mat, std::ostream &out, SVGOptions opts={})
Writes the visualization of matrix in the SVG format.
Definition: io.hh:807
void writeMatrixToMatlab(const MatrixType &matrix, const std::string &filename, int outputPrecision=18)
Writes sparse matrix in a Matlab-readable format.
Definition: io.hh:451
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)
Print one row of a matrix, specialization for number types.
Definition: io.hh:151
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 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 writeMatrixToMatlabHelper(const FieldType &value, int rowOffset, int colOffset, std::ostream &s, typename std::enable_if_t< Dune::IsNumber< FieldType >::value > *sfinae=nullptr)
Helper method for the writeMatrixToMatlab routine.
Definition: io.hh:379
void writeVectorToMatlab(const VectorType &vector, const std::string &filename, int outputPrecision=18)
Writes vectors in a Matlab-readable format.
Definition: io.hh:492
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 printSparseMatrix(std::ostream &s, const BCRSMatrix< FieldMatrix< B, n, m >, A > &mat, std::string title, std::string rowtext, int width=3, int precision=2)
Prints a BCRSMatrix with fixed sized blocks.
Definition: io.hh:271
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
auto max(ADLTag< 0 >, const V &v1, const V &v2)
implements binary Simd::max()
Definition: defaults.hh:81
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:674
std::string style
CSS style block to write in header.
Definition: io.hh:686
std::size_t width
Final width size (pixels) of the SVG header. If 0, size is automatic.
Definition: io.hh:680
std::function< std::string(const double &)> color_fill
Color fill for default options.
Definition: io.hh:710
std::size_t block_size
size (pixels) of the deepst block/value of the matrix
Definition: io.hh:676
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:769
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:733
std::size_t interspace
size (pixels) of the interspace between blocks
Definition: io.hh:678
bool write_block_title
(Helper) Whether to write a title on the rectangle value
Definition: io.hh:725
std::size_t height
Final height size (pixels) of the SVG header. If 0, size is automatic.
Definition: io.hh:682
bool write_header
Whether to write the SVG header.
Definition: io.hh:684
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:718
Whether this type acts as a scalar in the context of (hierarchically blocked) containers.
Definition: typetraits.hh:194