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*/,
255 int /*innercol*/)
256 {
257 s<<innerMatrixElement<<" ";
258 }
259
260 template<class B, int n>
261 void printInnerMatrixElement(std::ostream& s,
262 const ScaledIdentityMatrix<B,n> innerMatrixElement,
263 int innerrow, int innercol)
264 {
265 if (innerrow == innercol)
266 s<<innerMatrixElement.scalar()<<" ";
267 else
268 s<<"-";
269 }
270
271 template<class B, int n, int m>
272 void printInnerMatrixElement(std::ostream& s,
273 const FieldMatrix<B,n,m> innerMatrixElement,
274 int innerrow, int innercol)
275 {
276 s<<innerMatrixElement[innerrow][innercol]<<" ";
277 }
278 }
279
301 template<class A, class InnerMatrixType>
302 void printSparseMatrix(std::ostream& s,
304 std::string title, std::string rowtext,
305 int width=3, int precision=2)
306 {
308 // remember old flags
309 std::ios_base::fmtflags oldflags = s.flags();
310 // set the output format
311 s.setf(std::ios_base::scientific, std::ios_base::floatfield);
312 int oldprec = s.precision();
313 s.precision(precision);
314 // print title
315 s << title
316 << " [n=" << mat.N()
317 << ",m=" << mat.M()
318 << ",rowdim=" << MatrixDimension<Matrix>::rowdim(mat)
319 << ",coldim=" << MatrixDimension<Matrix>::coldim(mat)
320 << "]" << std::endl;
321
322 typedef typename Matrix::ConstRowIterator Row;
323
324 constexpr int n = std::decay_t<decltype(Impl::asMatrix(std::declval<InnerMatrixType>()))>::rows;
325 constexpr int m = std::decay_t<decltype(Impl::asMatrix(std::declval<InnerMatrixType>()))>::cols;
326 for(Row row=mat.begin(); row != mat.end(); ++row) {
327 int skipcols=0;
328 bool reachedEnd=false;
329
330 while(!reachedEnd) {
331 for(int innerrow=0; innerrow<n; ++innerrow) {
332 int count=0;
333 typedef typename Matrix::ConstColIterator Col;
334 Col col=row->begin();
335 for(; col != row->end(); ++col,++count) {
336 if(count<skipcols)
337 continue;
338 if(count>=skipcols+width)
339 break;
340 if(innerrow==0) {
341 if(count==skipcols) {
342 s << rowtext; // start a new row
343 s << " "; // space in front of each entry
344 s.width(4); // set width for counter
345 s << row.index()<<": "; // number of first entry in a line
346 }
347 s.width(4);
348 s<<col.index()<<": |";
349 } else {
350 if(count==skipcols) {
351 for(typename std::string::size_type i=0; i < rowtext.length(); i++)
352 s<<" ";
353 s<<" ";
354 }
355 s<<" |";
356 }
357 for(int innercol=0; innercol < m; ++innercol) {
358 s.width(9);
359 Impl::printInnerMatrixElement(s,*col,innerrow,innercol);
360 }
361
362 s<<"|";
363 }
364 if(innerrow==n-1 && col==row->end())
365 reachedEnd = true;
366 else
367 s << std::endl;
368 }
369 skipcols += width;
370 s << std::endl;
371 }
372 s << std::endl;
373 }
374
375 // reset the output format
376 s.flags(oldflags);
377 s.precision(oldprec);
378 }
379
380 namespace
381 {
382 template<typename T>
383 struct MatlabPODWriter
384 {
385 static std::ostream& write(const T& t, std::ostream& s)
386 {
387 s << t;
388 return s;
389 }
390 };
391 template<typename T>
392 struct MatlabPODWriter<std::complex<T> >
393 {
394 static std::ostream& write(const std::complex<T>& t, std::ostream& s)
395 {
396 s << t.real() << " " << t.imag();
397 return s;
398 }
399 };
400 } // anonymous namespace
401
411 template <class FieldType,
412 std::enable_if_t<Dune::IsNumber<FieldType>::value, int> = 0>
413 void writeMatrixToMatlabHelper(const FieldType& value,
414 int rowOffset, int colOffset,
415 std::ostream& s)
416 {
417 //+1 for Matlab numbering
418 s << rowOffset + 1 << " " << colOffset + 1 << " ";
419 MatlabPODWriter<FieldType>::write(value, s)<< std::endl;
420 }
421
429 template <class MatrixType,
430 std::enable_if_t<not Dune::IsNumber<MatrixType>::value, int> = 0>
431 void writeMatrixToMatlabHelper(const MatrixType& matrix,
432 int externalRowOffset, int externalColOffset,
433 std::ostream& s)
434 {
435 // Precompute the accumulated sizes of the columns
436 std::vector<typename MatrixType::size_type> colOffset(matrix.M());
437 if (colOffset.size() > 0)
438 colOffset[0] = 0;
439
440 for (typename MatrixType::size_type i=0; i<matrix.M()-1; i++)
441 colOffset[i+1] = colOffset[i] +
442 MatrixDimension<MatrixType>::coldim(matrix,i);
443
444 typename MatrixType::size_type rowOffset = 0;
445
446 // Loop over all matrix rows
447 for (typename MatrixType::size_type rowIdx=0; rowIdx<matrix.N(); rowIdx++)
448 {
449 auto cIt = matrix[rowIdx].begin();
450 auto cEndIt = matrix[rowIdx].end();
451
452 // Loop over all columns in this row
453 for (; cIt!=cEndIt; ++cIt)
455 externalRowOffset+rowOffset,
456 externalColOffset + colOffset[cIt.index()],
457 s);
458
459 rowOffset += MatrixDimension<MatrixType>::rowdim(matrix, rowIdx);
460 }
461
462 }
463
483 template <class MatrixType>
484 void writeMatrixToMatlab(const MatrixType& matrix,
485 const std::string& filename, int outputPrecision = 18)
486 {
487 std::ofstream outStream(filename.c_str());
488 int oldPrecision = outStream.precision();
489 outStream.precision(outputPrecision);
490
491 writeMatrixToMatlabHelper(matrix, 0, 0, outStream);
492 outStream.precision(oldPrecision);
493 }
494
495 // Recursively write vector entries to a stream
496 template<class V>
497 void writeVectorToMatlabHelper (const V& v, std::ostream& stream)
498 {
499 if constexpr (IsNumber<V>()) {
500 stream << v << std::endl;
501 } else {
502 for (const auto& entry : v)
503 writeVectorToMatlabHelper(entry, stream);
504 }
505 }
506
524 template <class VectorType>
525 void writeVectorToMatlab(const VectorType& vector,
526 const std::string& filename, int outputPrecision = 18)
527 {
528 std::ofstream outStream(filename.c_str());
529 int oldPrecision = outStream.precision();
530 outStream.precision(outputPrecision);
531
532 writeVectorToMatlabHelper(vector, outStream);
533 outStream.precision(oldPrecision);
534 }
535
536 namespace Impl {
537
539 struct NullStream {
540 template <class Any>
541 friend NullStream &operator<<(NullStream &dev0, Any &&) {
542 return dev0;
543 }
544 };
545
547 // svg shall be closed with a group and an svg. i.e. "</g></svg>"
548 template <class Stream, class SVGMatrixOptions>
549 void writeSVGMatrixHeader(Stream &out, const SVGMatrixOptions &opts,
550 std::pair<std::size_t, size_t> offsets) {
551 auto [col_offset, row_offset] = offsets;
552 double width = opts.width;
553 double height = opts.height;
554 // if empty, we try to figure out a sensible value of width and height
555 if (opts.width == 0 and opts.height == 0)
556 width = height = 500;
557 if (opts.width == 0)
558 width = opts.height * (double(col_offset) / row_offset);
559 if (opts.height == 0)
560 height = opts.width * (double(row_offset) / col_offset);
561
562 // scale group w.r.t final offsets
563 double scale_width = width / col_offset;
564 double scale_height = height / row_offset;
565
566 // write the header text
567 out << "<svg xmlns='http://www.w3.org/2000/svg' width='" << std::ceil(width)
568 << "' height='" << std::ceil(height) << "' version='1.1'>\n"
569 << "<style>\n"
570 << opts.style << "</style>\n"
571 << "<g transform='scale(" << scale_width << " " << scale_height
572 << ")'>\n";
573 }
574
576 template <class Stream, class Mat, class SVGMatrixOptions,
577 class RowPrefix, class ColPrefix>
578 std::pair<std::size_t, size_t>
579 writeSVGMatrix(Stream &out, const Mat &mat, SVGMatrixOptions opts,
580 RowPrefix row_prefix, ColPrefix col_prefix) {
581 // get values to fill the offsets
582 const auto& block_size = opts.block_size;
583 const auto& interspace = opts.interspace;
584
585 const std::size_t rows = mat.N();
586 const std::size_t cols = mat.M();
587
588 // disable header write for recursive calls
589 const bool write_header = opts.write_header;
590 opts.write_header = false;
591
592 // counter of offsets for every block
593 std::size_t row_offset = interspace;
594 std::size_t col_offset = interspace;
595
596 // lambda helper: for-each value
597 auto for_each_entry = [&mat](const auto &call_back) {
598 for (auto row_it = mat.begin(); row_it != mat.end(); ++row_it) {
599 for (auto col_it = row_it->begin(); col_it != row_it->end(); ++col_it) {
600 call_back(row_it.index(), col_it.index(), *col_it);
601 }
602 }
603 };
604
605 // accumulate content in another stream so that we write in correct order
606 std::stringstream ss;
607
608 // we need to append current row and col values to the prefixes
609 row_prefix.push_back(0);
610 col_prefix.push_back(0);
611
612 // do we need to write nested matrix blocks?
613 if constexpr (Dune::blockLevel<typename Mat::block_type>() == 0) {
614 // simple case: write svg block content to stream for each value
615 for_each_entry([&](const auto &row, const auto &col, const auto &val) {
616 std::size_t x_off = interspace + col * (interspace + block_size);
617 std::size_t y_off = interspace + row * (interspace + block_size);
618 row_prefix.back() = row;
619 col_prefix.back() = col;
620 opts.writeSVGBlock(ss, row_prefix, col_prefix, val,
621 {x_off, y_off, block_size, block_size});
622 });
623 col_offset += cols * (block_size + interspace);
624 row_offset += rows * (block_size + interspace);
625 } else {
626 // before we write anything, we need to calculate the
627 // offset for every {row,col} index
628 const auto null_offset = std::numeric_limits<std::size_t>::max();
629 std::vector<std::size_t> col_offsets(cols + 1, null_offset);
630 std::vector<std::size_t> row_offsets(rows + 1, null_offset);
631 for_each_entry([&](const auto &row, const auto &col, const auto &val) {
632 NullStream dev0;
633 // get size of sub-block
634 auto sub_size =
635 writeSVGMatrix(dev0, val, opts, row_prefix, col_prefix);
636
637 // if we didn't see col size before
638 if (col_offsets[col + 1] == null_offset) // write it in the offset vector
639 col_offsets[col + 1] = sub_size.first;
640
641 // repeat process for row sizes
642 if (row_offsets[row + 1] == null_offset)
643 row_offsets[row + 1] = sub_size.second;
644 });
645
646 // if some rows/cols were not visited, make an educated guess with the minimum offset
647 auto min_row_offset = *std::min_element(begin(row_offsets), end(row_offsets));
648 std::replace(begin(row_offsets), end(row_offsets), null_offset, min_row_offset);
649 auto min_col_offset = *std::min_element(begin(col_offsets), end(col_offsets));
650 std::replace(begin(col_offsets), end(col_offsets), null_offset, min_col_offset);
651
652 // we have sizes for every block: to get offsets we make a partial sum
653 col_offsets[0] = interspace;
654 row_offsets[0] = interspace;
655 for (std::size_t i = 1; i < col_offsets.size(); i++)
656 col_offsets[i] += col_offsets[i - 1] + interspace;
657 for (std::size_t i = 1; i < row_offsets.size(); i++)
658 row_offsets[i] += row_offsets[i - 1] + interspace;
659
660 for_each_entry([&](const auto &row, const auto &col, const auto &val) {
661 // calculate svg view from offsets
662 std::size_t width =
663 col_offsets[col + 1] - col_offsets[col] - interspace;
664 std::size_t height =
665 row_offsets[row + 1] - row_offsets[row] - interspace;
666 row_prefix.back() = row;
667 col_prefix.back() = col;
668 // content of the sub-block has origin at {0,0}: shift it to the correct place
669 ss << "<svg x='" << col_offsets[col] << "' y='" << row_offsets[row]
670 << "' width='" << width << "' height='" << height << "'>\n";
671 // write a nested svg with the contents of the sub-block
672 writeSVGMatrix(ss, val, opts, row_prefix, col_prefix);
673 ss << "</svg>\n";
674 });
675 col_offset = col_offsets.back();
676 row_offset = row_offsets.back();
677 }
678
679 // write content in order!
680 // (i) if required, first header
681 if (write_header)
682 writeSVGMatrixHeader(out, opts, {col_offset, row_offset});
683
684 col_prefix.pop_back();
685 row_prefix.pop_back();
686 // (ii) an svg block for this level
687 opts.writeSVGBlock(out, row_prefix, col_prefix, mat,
688 {0, 0, col_offset, row_offset});
689 // (iii) the content of the matrix
690 out << ss.str();
691 // (iv) if required, close the header
692 if (write_header)
693 out << "</g>\n</svg>\n";
694
695 // return the total required for this block
696 return {col_offset, row_offset};
697 }
698 } // namespace Impl
699
700
709 std::size_t block_size = 10;
711 std::size_t interspace = 5;
713 std::size_t width = 500;
715 std::size_t height = 0;
717 bool write_header = true;
719 std::string style = " .matrix-block {\n"
720 " fill: cornflowerblue;\n"
721 " fill-opacity: 0.4;\n"
722 " stroke-width: 2;\n"
723 " stroke: black;\n"
724 " stroke-opacity: 0.5;\n"
725 " }\n"
726 " .matrix-block:hover {\n"
727 " fill: lightcoral;\n"
728 " fill-opacity: 0.4;\n"
729 " stroke-opacity: 1;\n"
730 " }\n";
731
743 std::function<std::string(const double&)> color_fill;
744
750 template <class RowPrefix, class ColPrefix>
751 std::string blockStyleClass(const RowPrefix & /*row_prefix*/,
752 const ColPrefix & /*col_prefix*/) const {
753 // here, you can potentially give a different style to each block
754 return "matrix-block";
755 }
756
758 bool write_block_title = true;
759
765 template <class Stream, class RowPrefix, class ColPrefix, class Block>
766 void writeBlockTitle(Stream& out, const RowPrefix &row_prefix,
767 const ColPrefix &col_prefix,
768 const Block &block) const {
769 if (this->write_block_title) {
770 out << "<title>";
771 assert(row_prefix.size() == col_prefix.size());
772 for (std::size_t i = 0; i < row_prefix.size(); ++i)
773 out << "[" << row_prefix[i] << ", "<< col_prefix[i] << "]";
774 if constexpr (Dune::blockLevel<Block>() == 0)
775 out << ": " << block;
776 out << "</title>\n";
777 }
778 }
779
801 template <class Stream, class RowPrefix, class ColPrefix, class Block>
802 void writeSVGBlock(Stream &out,
803 const RowPrefix &row_prefix,
804 const ColPrefix &col_prefix, const Block block,
805 const std::array<std::size_t, 4> &svg_box) const {
806 // get bounding box values
807 auto &[x_off, y_off, block_width, block_height] = svg_box;
808 // get style class
809 std::string block_class = this->blockStyleClass(row_prefix, col_prefix);
810 // write a rectangle on the bounding box
811 out << "<rect class='" << block_class << "' x='" << x_off << "' y='"
812 << y_off << "' width='" << block_width << "' height='" << block_height << "'";
813 if constexpr (Dune::blockLevel<Block>() == 0 and std::is_convertible<Block,double>{})
814 if (color_fill)
815 out << " style='fill-opacity: 1;fill:" << color_fill(double(block)) << "'";
816
817 out << ">\n";
818 // give the rectangle a title (in html this shows info about the block)
819 this->writeBlockTitle(out,row_prefix, col_prefix, block);
820 // close rectangle
821 out << "</rect>\n";
822 }
823 };
824
839 template <class Mat, class SVGOptions = DefaultSVGMatrixOptions>
840 void writeSVGMatrix(std::ostream &out, const Mat &mat, SVGOptions opts = {}) {
841 // We need a vector that can fit all the multi-indices for rows and columns
843 // Call overload for Mat type
844 Impl::writeSVGMatrix(out, mat, opts, IndexPrefix{}, IndexPrefix{});
845 }
846
849} // namespace Dune
850
851#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:467
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:49
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:489
void writeMatrixToMatlab(const MatrixType &matrix, const std::string &filename, int outputPrecision=18)
Writes sparse matrix in a Matlab-readable format.
Definition: io.hh:484
void writeMatrixToMatlabHelper(const FieldType &value, int rowOffset, int colOffset, std::ostream &s)
Helper method for the writeMatrixToMatlab routine.
Definition: io.hh:413
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:302
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:840
void writeVectorToMatlab(const VectorType &vector, const std::string &filename, int outputPrecision=18)
Writes vectors in a Matlab-readable format.
Definition: io.hh:525
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:707
std::string style
CSS style block to write in header.
Definition: io.hh:719
std::size_t width
Final width size (pixels) of the SVG header. If 0, size is automatic.
Definition: io.hh:713
std::function< std::string(const double &)> color_fill
Color fill for default options.
Definition: io.hh:743
std::size_t block_size
size (pixels) of the deepst block/value of the matrix
Definition: io.hh:709
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:802
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:766
std::size_t interspace
size (pixels) of the interspace between blocks
Definition: io.hh:711
bool write_block_title
(Helper) Whether to write a title on the rectangle value
Definition: io.hh:758
std::size_t height
Final height size (pixels) of the SVG header. If 0, size is automatic.
Definition: io.hh:715
bool write_header
Whether to write the SVG header.
Definition: io.hh:717
std::string blockStyleClass(const RowPrefix &, const ColPrefix &) const
Helper function that returns an style class for a given prefix.
Definition: io.hh:751
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 & Uni Heidelberg  |  generated with Hugo v0.111.3 (Jan 8, 23:33, 2026)