Dune Core Modules (unstable)

io.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright © DUNE Project contributors, see file LICENSE.md in module root
2// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
3// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
4// vi: set et ts=4 sw=2 sts=2:
5#ifndef DUNE_ISTL_IO_HH
6#define DUNE_ISTL_IO_HH
7
8#include <cmath>
9#include <complex>
10#include <limits>
11#include <ios>
12#include <iomanip>
13#include <fstream>
14#include <string>
15
16#include "matrixutils.hh"
17#include "istlexception.hh"
20#include <dune/common/hybridutilities.hh>
22
25
26namespace Dune {
27
40 //
41 // pretty printing of vectors
42 //
43
51 template<class V>
52 void recursive_printvector (std::ostream& s, const V& v, std::string rowtext,
53 int& counter, int columns, int width)
54 {
55 if constexpr (IsNumber<V>())
56 {
57 // Print one number
58 if (counter%columns==0)
59 {
60 s << rowtext; // start a new row
61 s << " "; // space in front of each entry
62 s.width(4); // set width for counter
63 s << counter; // number of first entry in a line
64 }
65 s << " "; // space in front of each entry
66 s.width(width); // set width for each entry anew
67 s << v; // yeah, the number !
68 counter++; // increment the counter
69 if (counter%columns==0)
70 s << std::endl; // start a new line
71 }
72 else
73 {
74 // Recursively print a vector
75 for (const auto& entry : v)
76 recursive_printvector(s,entry,rowtext,counter,columns,width);
77 }
78 }
79
80
88 template<class V>
89 void printvector (std::ostream& s, const V& v, std::string title,
90 std::string rowtext, int columns=1, int width=10,
91 int precision=2)
92 {
93 // count the numbers printed to make columns
94 int counter=0;
95
96 // remember old flags
97 std::ios_base::fmtflags oldflags = s.flags();
98
99 // set the output format
100 s.setf(std::ios_base::scientific, std::ios_base::floatfield);
101 int oldprec = s.precision();
102 s.precision(precision);
103
104 // print title
105 s << title << " [blocks=" << v.N() << ",dimension=" << v.dim() << "]"
106 << std::endl;
107
108 // print data from all blocks
109 recursive_printvector(s,v,rowtext,counter,columns,width);
110
111 // check if new line is required
112 if (counter%columns!=0)
113 s << std::endl;
114
115 // reset the output format
116 s.flags(oldflags);
117 s.precision(oldprec);
118 }
119
120
122 //
123 // pretty printing of matrices
124 //
125
133 inline void fill_row (std::ostream& s, int m, int width, [[maybe_unused]] int precision)
134 {
135 for (int j=0; j<m; j++)
136 {
137 s << " "; // space in front of each entry
138 s.width(width); // set width for each entry anew
139 s << "."; // yeah, the number !
140 }
141 }
142
150 template<class K,
151 std::enable_if_t<Dune::IsNumber<K>::value, int> = 0>
152 void print_row (std::ostream& s, const K& value,
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,
156 int width,
157 [[maybe_unused]] int precision)
158 {
159 s << " "; // space in front of each entry
160 s.width(width); // set width for each entry anew
161 s << value;
162 }
163
171 template<class M,
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)
176 {
177 typename M::size_type i0=I;
178 for (typename M::size_type i=0; i<A.N(); i++)
179 {
180 if (therow>=i0 && therow<i0+MatrixDimension<M>::rowdim(A,i))
181 {
182 // the row is in this block row !
183 typename M::size_type j0=J;
184 for (typename M::size_type j=0; j<A.M(); j++)
185 {
186 // find this block
187 typename M::ConstColIterator it = A[i].find(j);
188
189 // print row or filler
190 if (it!=A[i].end())
191 print_row(s,*it,i0,j0,therow,width,precision);
192 else
193 fill_row(s,MatrixDimension<M>::coldim(A,j),width,precision);
194
195 // advance columns
196 j0 += MatrixDimension<M>::coldim(A,j);
197 }
198 }
199 // advance rows
200 i0 += MatrixDimension<M>::rowdim(A,i);
201 }
202 }
203
212 template<class M>
213 void printmatrix (std::ostream& s, const M& A, std::string title,
214 std::string rowtext, int width=10, int precision=2)
215 {
216
217 // remember old flags
218 std::ios_base::fmtflags oldflags = s.flags();
219
220 // set the output format
221 s.setf(std::ios_base::scientific, std::ios_base::floatfield);
222 int oldprec = s.precision();
223 s.precision(precision);
224
225 // print title
226 s << title
227 << " [n=" << A.N()
228 << ",m=" << A.M()
229 << ",rowdim=" << MatrixDimension<M>::rowdim(A)
230 << ",coldim=" << MatrixDimension<M>::coldim(A)
231 << "]" << std::endl;
232
233 // print all rows
234 for (typename M::size_type i=0; i<MatrixDimension<M>::rowdim(A); i++)
235 {
236 s << rowtext; // start a new row
237 s << " "; // space in front of each entry
238 s.width(4); // set width for counter
239 s << i; // number of first entry in a line
240 print_row(s,A,0,0,i,width,precision); // generic print
241 s << std::endl; // start a new line
242 }
243
244 // reset the output format
245 s.flags(oldflags);
246 s.precision(oldprec);
247 }
248
249 namespace Impl
250 {
251 template<class B>
252 void printInnerMatrixElement(std::ostream& s,
253 const B& innerMatrixElement,
254 int innerrow, int innercol)
255 {
256 s<<innerMatrixElement<<" ";
257 }
258
259 template<class B, int n>
260 void printInnerMatrixElement(std::ostream& s,
261 const ScaledIdentityMatrix<B,n> innerMatrixElement,
262 int innerrow, int innercol)
263 {
264 if (innerrow == innercol)
265 s<<innerMatrixElement.scalar()<<" ";
266 else
267 s<<"-";
268 }
269
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)
274 {
275 s<<innerMatrixElement[innerrow][innercol]<<" ";
276 }
277 }
278
300 template<class A, class InnerMatrixType>
301 void printSparseMatrix(std::ostream& s,
303 std::string title, std::string rowtext,
304 int width=3, int precision=2)
305 {
307 // remember old flags
308 std::ios_base::fmtflags oldflags = s.flags();
309 // set the output format
310 s.setf(std::ios_base::scientific, std::ios_base::floatfield);
311 int oldprec = s.precision();
312 s.precision(precision);
313 // print title
314 s << title
315 << " [n=" << mat.N()
316 << ",m=" << mat.M()
317 << ",rowdim=" << MatrixDimension<Matrix>::rowdim(mat)
318 << ",coldim=" << MatrixDimension<Matrix>::coldim(mat)
319 << "]" << std::endl;
320
321 typedef typename Matrix::ConstRowIterator Row;
322
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) {
326 int skipcols=0;
327 bool reachedEnd=false;
328
329 while(!reachedEnd) {
330 for(int innerrow=0; innerrow<n; ++innerrow) {
331 int count=0;
332 typedef typename Matrix::ConstColIterator Col;
333 Col col=row->begin();
334 for(; col != row->end(); ++col,++count) {
335 if(count<skipcols)
336 continue;
337 if(count>=skipcols+width)
338 break;
339 if(innerrow==0) {
340 if(count==skipcols) {
341 s << rowtext; // start a new row
342 s << " "; // space in front of each entry
343 s.width(4); // set width for counter
344 s << row.index()<<": "; // number of first entry in a line
345 }
346 s.width(4);
347 s<<col.index()<<": |";
348 } else {
349 if(count==skipcols) {
350 for(typename std::string::size_type i=0; i < rowtext.length(); i++)
351 s<<" ";
352 s<<" ";
353 }
354 s<<" |";
355 }
356 for(int innercol=0; innercol < m; ++innercol) {
357 s.width(9);
358 Impl::printInnerMatrixElement(s,*col,innerrow,innercol);
359 }
360
361 s<<"|";
362 }
363 if(innerrow==n-1 && col==row->end())
364 reachedEnd = true;
365 else
366 s << std::endl;
367 }
368 skipcols += width;
369 s << std::endl;
370 }
371 s << std::endl;
372 }
373
374 // reset the output format
375 s.flags(oldflags);
376 s.precision(oldprec);
377 }
378
379 namespace
380 {
381 template<typename T>
382 struct MatlabPODWriter
383 {
384 static std::ostream& write(const T& t, std::ostream& s)
385 {
386 s << t;
387 return s;
388 }
389 };
390 template<typename T>
391 struct MatlabPODWriter<std::complex<T> >
392 {
393 static std::ostream& write(const std::complex<T>& t, std::ostream& s)
394 {
395 s << t.real() << " " << t.imag();
396 return s;
397 }
398 };
399 } // anonymous namespace
400
410 template <class FieldType,
411 std::enable_if_t<Dune::IsNumber<FieldType>::value, int> = 0>
412 void writeMatrixToMatlabHelper(const FieldType& value,
413 int rowOffset, int colOffset,
414 std::ostream& s)
415 {
416 //+1 for Matlab numbering
417 s << rowOffset + 1 << " " << colOffset + 1 << " ";
418 MatlabPODWriter<FieldType>::write(value, s)<< std::endl;
419 }
420
428 template <class MatrixType,
429 std::enable_if_t<not Dune::IsNumber<MatrixType>::value, int> = 0>
430 void writeMatrixToMatlabHelper(const MatrixType& matrix,
431 int externalRowOffset, int externalColOffset,
432 std::ostream& s)
433 {
434 // Precompute the accumulated sizes of the columns
435 std::vector<typename MatrixType::size_type> colOffset(matrix.M());
436 if (colOffset.size() > 0)
437 colOffset[0] = 0;
438
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);
442
443 typename MatrixType::size_type rowOffset = 0;
444
445 // Loop over all matrix rows
446 for (typename MatrixType::size_type rowIdx=0; rowIdx<matrix.N(); rowIdx++)
447 {
448 auto cIt = matrix[rowIdx].begin();
449 auto cEndIt = matrix[rowIdx].end();
450
451 // Loop over all columns in this row
452 for (; cIt!=cEndIt; ++cIt)
454 externalRowOffset+rowOffset,
455 externalColOffset + colOffset[cIt.index()],
456 s);
457
458 rowOffset += MatrixDimension<MatrixType>::rowdim(matrix, rowIdx);
459 }
460
461 }
462
482 template <class MatrixType>
483 void writeMatrixToMatlab(const MatrixType& matrix,
484 const std::string& filename, int outputPrecision = 18)
485 {
486 std::ofstream outStream(filename.c_str());
487 int oldPrecision = outStream.precision();
488 outStream.precision(outputPrecision);
489
490 writeMatrixToMatlabHelper(matrix, 0, 0, outStream);
491 outStream.precision(oldPrecision);
492 }
493
494 // Recursively write vector entries to a stream
495 template<class V>
496 void writeVectorToMatlabHelper (const V& v, std::ostream& stream)
497 {
498 if constexpr (IsNumber<V>()) {
499 stream << v << std::endl;
500 } else {
501 for (const auto& entry : v)
502 writeVectorToMatlabHelper(entry, stream);
503 }
504 }
505
523 template <class VectorType>
524 void writeVectorToMatlab(const VectorType& vector,
525 const std::string& filename, int outputPrecision = 18)
526 {
527 std::ofstream outStream(filename.c_str());
528 int oldPrecision = outStream.precision();
529 outStream.precision(outputPrecision);
530
531 writeVectorToMatlabHelper(vector, outStream);
532 outStream.precision(oldPrecision);
533 }
534
535 namespace Impl {
536
538 struct NullStream {
539 template <class Any>
540 friend NullStream &operator<<(NullStream &dev0, Any &&) {
541 return dev0;
542 }
543 };
544
546 // svg shall be closed with a group and an svg. i.e. "</g></svg>"
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;
553 // if empty, we try to figure out a sensible value of width and height
554 if (opts.width == 0 and opts.height == 0)
555 width = height = 500;
556 if (opts.width == 0)
557 width = opts.height * (double(col_offset) / row_offset);
558 if (opts.height == 0)
559 height = opts.width * (double(row_offset) / col_offset);
560
561 // scale group w.r.t final offsets
562 double scale_width = width / col_offset;
563 double scale_height = height / row_offset;
564
565 // write the header text
566 out << "<svg xmlns='http://www.w3.org/2000/svg' width='" << std::ceil(width)
567 << "' height='" << std::ceil(height) << "' version='1.1'>\n"
568 << "<style>\n"
569 << opts.style << "</style>\n"
570 << "<g transform='scale(" << scale_width << " " << scale_height
571 << ")'>\n";
572 }
573
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) {
580 // get values to fill the offsets
581 const auto& block_size = opts.block_size;
582 const auto& interspace = opts.interspace;
583
584 const std::size_t rows = mat.N();
585 const std::size_t cols = mat.M();
586
587 // disable header write for recursive calls
588 const bool write_header = opts.write_header;
589 opts.write_header = false;
590
591 // counter of offsets for every block
592 std::size_t row_offset = interspace;
593 std::size_t col_offset = interspace;
594
595 // lambda helper: for-each value
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);
600 }
601 }
602 };
603
604 // accumulate content in another stream so that we write in correct order
605 std::stringstream ss;
606
607 // we need to append current row and col values to the prefixes
608 row_prefix.push_back(0);
609 col_prefix.push_back(0);
610
611 // do we need to write nested matrix blocks?
612 if constexpr (Dune::blockLevel<typename Mat::block_type>() == 0) {
613 // simple case: write svg block content to stream for each value
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});
621 });
622 col_offset += cols * (block_size + interspace);
623 row_offset += rows * (block_size + interspace);
624 } else {
625 // before we write anything, we need to calculate the
626 // offset for every {row,col} index
627 const auto null_offset = std::numeric_limits<std::size_t>::max();
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) {
631 NullStream dev0;
632 // get size of sub-block
633 auto sub_size =
634 writeSVGMatrix(dev0, val, opts, row_prefix, col_prefix);
635
636 // if we didn't see col size before
637 if (col_offsets[col + 1] == null_offset) // write it in the offset vector
638 col_offsets[col + 1] = sub_size.first;
639
640 // repeat process for row sizes
641 if (row_offsets[row + 1] == null_offset)
642 row_offsets[row + 1] = sub_size.second;
643 });
644
645 // if some rows/cols were not visited, make an educated guess with the minimum offset
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);
650
651 // we have sizes for every block: to get offsets we make a partial sum
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;
658
659 for_each_entry([&](const auto &row, const auto &col, const auto &val) {
660 // calculate svg view from offsets
661 std::size_t width =
662 col_offsets[col + 1] - col_offsets[col] - interspace;
663 std::size_t height =
664 row_offsets[row + 1] - row_offsets[row] - interspace;
665 row_prefix.back() = row;
666 col_prefix.back() = col;
667 // content of the sub-block has origin at {0,0}: shift it to the correct place
668 ss << "<svg x='" << col_offsets[col] << "' y='" << row_offsets[row]
669 << "' width='" << width << "' height='" << height << "'>\n";
670 // write a nested svg with the contents of the sub-block
671 writeSVGMatrix(ss, val, opts, row_prefix, col_prefix);
672 ss << "</svg>\n";
673 });
674 col_offset = col_offsets.back();
675 row_offset = row_offsets.back();
676 }
677
678 // write content in order!
679 // (i) if required, first header
680 if (write_header)
681 writeSVGMatrixHeader(out, opts, {col_offset, row_offset});
682
683 col_prefix.pop_back();
684 row_prefix.pop_back();
685 // (ii) an svg block for this level
686 opts.writeSVGBlock(out, row_prefix, col_prefix, mat,
687 {0, 0, col_offset, row_offset});
688 // (iii) the content of the matrix
689 out << ss.str();
690 // (iv) if required, close the header
691 if (write_header)
692 out << "</g>\n</svg>\n";
693
694 // return the total required for this block
695 return {col_offset, row_offset};
696 }
697 } // namespace Impl
698
699
708 std::size_t block_size = 10;
710 std::size_t interspace = 5;
712 std::size_t width = 500;
714 std::size_t height = 0;
716 bool write_header = true;
718 std::string style = " .matrix-block {\n"
719 " fill: cornflowerblue;\n"
720 " fill-opacity: 0.4;\n"
721 " stroke-width: 2;\n"
722 " stroke: black;\n"
723 " stroke-opacity: 0.5;\n"
724 " }\n"
725 " .matrix-block:hover {\n"
726 " fill: lightcoral;\n"
727 " fill-opacity: 0.4;\n"
728 " stroke-opacity: 1;\n"
729 " }\n";
730
742 std::function<std::string(const double&)> color_fill;
743
749 template <class RowPrefix, class ColPrefix>
750 std::string blockStyleClass(const RowPrefix &row_prefix,
751 const ColPrefix &col_prefix) const {
752 // here, you can potentially give a different style to each block
753 return "matrix-block";
754 }
755
757 bool write_block_title = true;
758
764 template <class Stream, class RowPrefix, class ColPrefix, class Block>
765 void writeBlockTitle(Stream& out, const RowPrefix &row_prefix,
766 const ColPrefix &col_prefix,
767 const Block &block) const {
768 if (this->write_block_title) {
769 out << "<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;
775 out << "</title>\n";
776 }
777 }
778
800 template <class Stream, class RowPrefix, class ColPrefix, class Block>
801 void writeSVGBlock(Stream &out,
802 const RowPrefix &row_prefix,
803 const ColPrefix &col_prefix, const Block block,
804 const std::array<std::size_t, 4> &svg_box) const {
805 // get bounding box values
806 auto &[x_off, y_off, width, height] = svg_box;
807 // get style class
808 std::string block_class = this->blockStyleClass(row_prefix, col_prefix);
809 // write a rectangle on the bounding box
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>{})
813 if (color_fill)
814 out << " style='fill-opacity: 1;fill:" << color_fill(double(block)) << "'";
815
816 out << ">\n";
817 // give the rectangle a title (in html this shows info about the block)
818 this->writeBlockTitle(out,row_prefix, col_prefix, block);
819 // close rectangle
820 out << "</rect>\n";
821 }
822 };
823
838 template <class Mat, class SVGOptions = DefaultSVGMatrixOptions>
839 void writeSVGMatrix(std::ostream &out, const Mat &mat, SVGOptions opts = {}) {
840 // We need a vector that can fit all the multi-indices for rows and columns
842 // Call overload for Mat type
843 Impl::writeSVGMatrix(out, mat, opts, IndexPrefix{}, IndexPrefix{});
844 }
845
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.")]]
865 void writeSVGMatrix(const Mat &mat, std::ostream &out, SVGOptions opts = {}) {
866 writeSVGMatrix(out, mat, opts);
867 }
868
871} // namespace Dune
872
873#endif
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
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
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
STL namespace.
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
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 21, 23:30, 2024)