Dune Core Modules (2.9.1)

io.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (C) 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 void print_row (std::ostream& s, const K& value,
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,
155 int width,
156 [[maybe_unused]] int precision,
157 typename std::enable_if_t<Dune::IsNumber<K>::value>* sfinae = nullptr)
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 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,
175 typename std::enable_if_t<!Dune::IsNumber<M>::value>* sfinae = nullptr)
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
270 template<class B, int n, int m, class A>
271 void printSparseMatrix(std::ostream& s,
272 const BCRSMatrix<FieldMatrix<B,n,m>,A>& mat,
273 std::string title, std::string rowtext,
274 int width=3, int precision=2)
275 {
277 // remember old flags
278 std::ios_base::fmtflags oldflags = s.flags();
279 // set the output format
280 s.setf(std::ios_base::scientific, std::ios_base::floatfield);
281 int oldprec = s.precision();
282 s.precision(precision);
283 // print title
284 s << title
285 << " [n=" << mat.N()
286 << ",m=" << mat.M()
287 << ",rowdim=" << MatrixDimension<Matrix>::rowdim(mat)
288 << ",coldim=" << MatrixDimension<Matrix>::coldim(mat)
289 << "]" << std::endl;
290
291 typedef typename Matrix::ConstRowIterator Row;
292
293 for(Row row=mat.begin(); row != mat.end(); ++row) {
294 int skipcols=0;
295 bool reachedEnd=false;
296
297 while(!reachedEnd) {
298 for(int innerrow=0; innerrow<n; ++innerrow) {
299 int count=0;
300 typedef typename Matrix::ConstColIterator Col;
301 Col col=row->begin();
302 for(; col != row->end(); ++col,++count) {
303 if(count<skipcols)
304 continue;
305 if(count>=skipcols+width)
306 break;
307 if(innerrow==0) {
308 if(count==skipcols) {
309 s << rowtext; // start a new row
310 s << " "; // space in front of each entry
311 s.width(4); // set width for counter
312 s << row.index()<<": "; // number of first entry in a line
313 }
314 s.width(4);
315 s<<col.index()<<": |";
316 } else {
317 if(count==skipcols) {
318 for(typename std::string::size_type i=0; i < rowtext.length(); i++)
319 s<<" ";
320 s<<" ";
321 }
322 s<<" |";
323 }
324 for(int innercol=0; innercol < m; ++innercol) {
325 s.width(9);
326 s<<(*col)[innerrow][innercol]<<" ";
327 }
328
329 s<<"|";
330 }
331 if(innerrow==n-1 && col==row->end())
332 reachedEnd = true;
333 else
334 s << std::endl;
335 }
336 skipcols += width;
337 s << std::endl;
338 }
339 s << std::endl;
340 }
341
342 // reset the output format
343 s.flags(oldflags);
344 s.precision(oldprec);
345 }
346
347 namespace
348 {
349 template<typename T>
350 struct MatlabPODWriter
351 {
352 static std::ostream& write(const T& t, std::ostream& s)
353 {
354 s << t;
355 return s;
356 }
357 };
358 template<typename T>
359 struct MatlabPODWriter<std::complex<T> >
360 {
361 static std::ostream& write(const std::complex<T>& t, std::ostream& s)
362 {
363 s << t.real() << " " << t.imag();
364 return s;
365 }
366 };
367 } // anonymous namespace
368
378 template <class FieldType>
379 void writeMatrixToMatlabHelper(const FieldType& value,
380 int rowOffset, int colOffset,
381 std::ostream& s,
382 typename std::enable_if_t<Dune::IsNumber<FieldType>::value>* sfinae = nullptr)
383 {
384 //+1 for Matlab numbering
385 s << rowOffset + 1 << " " << colOffset + 1 << " ";
386 MatlabPODWriter<FieldType>::write(value, s)<< std::endl;
387 }
388
396 template <class MatrixType>
397 void writeMatrixToMatlabHelper(const MatrixType& matrix,
398 int externalRowOffset, int externalColOffset,
399 std::ostream& s,
400 typename std::enable_if_t<!Dune::IsNumber<MatrixType>::value>* sfinae = nullptr)
401 {
402 // Precompute the accumulated sizes of the columns
403 std::vector<typename MatrixType::size_type> colOffset(matrix.M());
404 if (colOffset.size() > 0)
405 colOffset[0] = 0;
406
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);
410
411 typename MatrixType::size_type rowOffset = 0;
412
413 // Loop over all matrix rows
414 for (typename MatrixType::size_type rowIdx=0; rowIdx<matrix.N(); rowIdx++)
415 {
416 auto cIt = matrix[rowIdx].begin();
417 auto cEndIt = matrix[rowIdx].end();
418
419 // Loop over all columns in this row
420 for (; cIt!=cEndIt; ++cIt)
422 externalRowOffset+rowOffset,
423 externalColOffset + colOffset[cIt.index()],
424 s);
425
426 rowOffset += MatrixDimension<MatrixType>::rowdim(matrix, rowIdx);
427 }
428
429 }
430
450 template <class MatrixType>
451 void writeMatrixToMatlab(const MatrixType& matrix,
452 const std::string& filename, int outputPrecision = 18)
453 {
454 std::ofstream outStream(filename.c_str());
455 int oldPrecision = outStream.precision();
456 outStream.precision(outputPrecision);
457
458 writeMatrixToMatlabHelper(matrix, 0, 0, outStream);
459 outStream.precision(oldPrecision);
460 }
461
462 // Recursively write vector entries to a stream
463 template<class V>
464 void writeVectorToMatlabHelper (const V& v, std::ostream& stream)
465 {
466 if constexpr (IsNumber<V>()) {
467 stream << v << std::endl;
468 } else {
469 for (const auto& entry : v)
470 writeVectorToMatlabHelper(entry, stream);
471 }
472 }
473
491 template <class VectorType>
492 void writeVectorToMatlab(const VectorType& vector,
493 const std::string& filename, int outputPrecision = 18)
494 {
495 std::ofstream outStream(filename.c_str());
496 int oldPrecision = outStream.precision();
497 outStream.precision(outputPrecision);
498
499 writeVectorToMatlabHelper(vector, outStream);
500 outStream.precision(oldPrecision);
501 }
502
503 namespace Impl {
504
506 struct NullStream {
507 template <class Any>
508 friend NullStream &operator<<(NullStream &dev0, Any &&) {
509 return dev0;
510 }
511 };
512
514 // svg shall be closed with a group and an svg. i.e. "</g></svg>"
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;
521 // if empty, we try to figure out a sensible value of width and height
522 if (opts.width == 0 and opts.height == 0)
523 width = height = 500;
524 if (opts.width == 0)
525 width = opts.height * (double(col_offset) / row_offset);
526 if (opts.height == 0)
527 height = opts.width * (double(row_offset) / col_offset);
528
529 // scale group w.r.t final offsets
530 double scale_width = width / col_offset;
531 double scale_height = height / row_offset;
532
533 // write the header text
534 out << "<svg xmlns='http://www.w3.org/2000/svg' width='" << std::ceil(width)
535 << "' height='" << std::ceil(height) << "' version='1.1'>\n"
536 << "<style>\n"
537 << opts.style << "</style>\n"
538 << "<g transform='scale(" << scale_width << " " << scale_height
539 << ")'>\n";
540 }
541
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) {
548 // get values to fill the offests
549 const auto& block_size = opts.block_size;
550 const auto& interspace = opts.interspace;
551
552 const std::size_t rows = mat.N();
553 const std::size_t cols = mat.M();
554
555 // disable header write for recursive calls
556 const bool write_header = opts.write_header;
557 opts.write_header = false;
558
559 // counter of offsets for every block
560 std::size_t row_offset = interspace;
561 std::size_t col_offset = interspace;
562
563 // lambda helper: for-each value
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);
568 }
569 }
570 };
571
572 // accumulate content in another stream so that we write in correct order
573 std::stringstream ss;
574
575 // we need to append current row and col values to the prefixes
576 row_prefix.push_back(0);
577 col_prefix.push_back(0);
578
579 // do we need to write nested matrix blocks?
580 if constexpr (Dune::blockLevel<typename Mat::block_type>() == 0) {
581 // simple case: write svg block content to stream for each value
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});
589 });
590 col_offset += cols * (block_size + interspace);
591 row_offset += rows * (block_size + interspace);
592 } else {
593 // before we write anything, we need to calculate the
594 // offset for every {row,col} index
595 const auto null_offset = std::numeric_limits<std::size_t>::max();
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) {
599 NullStream dev0;
600 // get size of sub-block
601 auto sub_size =
602 writeSVGMatrix(val, dev0, opts, row_prefix, col_prefix);
603
604 // if we didn't see col size before
605 if (col_offsets[col + 1] == null_offset) // write it in the offset vector
606 col_offsets[col + 1] = sub_size.first;
607
608 // repeat proces for row sizes
609 if (row_offsets[row + 1] == null_offset)
610 row_offsets[row + 1] = sub_size.second;
611 });
612
613 // if some rows/cols were not visited, make an educated guess with the minimum offset
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);
618
619 // we have sizes for every block: to get offsets we make a partial sum
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;
626
627 for_each_entry([&](const auto &row, const auto &col, const auto &val) {
628 // calculate svg view from offsets
629 std::size_t width =
630 col_offsets[col + 1] - col_offsets[col] - interspace;
631 std::size_t height =
632 row_offsets[row + 1] - row_offsets[row] - interspace;
633 row_prefix.back() = row;
634 col_prefix.back() = col;
635 // content of the sub-block has origin at {0,0}: shift it to the correct place
636 ss << "<svg x='" << col_offsets[col] << "' y='" << row_offsets[row]
637 << "' width='" << width << "' height='" << height << "'>\n";
638 // write a nested svg with the contents of the sub-block
639 writeSVGMatrix(val, ss, opts, row_prefix, col_prefix);
640 ss << "</svg>\n";
641 });
642 col_offset = col_offsets.back();
643 row_offset = row_offsets.back();
644 }
645
646 // write content in order!
647 // (i) if required, first header
648 if (write_header)
649 writeSVGMatrixHeader(out, opts, {col_offset, row_offset});
650
651 col_prefix.pop_back();
652 row_prefix.pop_back();
653 // (ii) an svg block for this level
654 opts.writeSVGBlock(out, row_prefix, col_prefix, mat,
655 {0, 0, col_offset, row_offset});
656 // (iii) the content of the matrix
657 out << ss.str();
658 // (iv) if required, close the header
659 if (write_header)
660 out << "</g>\n</svg>\n";
661
662 // return the total required for this block
663 return {col_offset, row_offset};
664 }
665 } // namespace Impl
666
667
676 std::size_t block_size = 10;
678 std::size_t interspace = 5;
680 std::size_t width = 500;
682 std::size_t height = 0;
684 bool write_header = true;
686 std::string style = " .matrix-block {\n"
687 " fill: cornflowerblue;\n"
688 " fill-opacity: 0.4;\n"
689 " stroke-width: 2;\n"
690 " stroke: black;\n"
691 " stroke-opacity: 0.5;\n"
692 " }\n"
693 " .matrix-block:hover {\n"
694 " fill: lightcoral;\n"
695 " fill-opacity: 0.4;\n"
696 " stroke-opacity: 1;\n"
697 " }\n";
698
710 std::function<std::string(const double&)> color_fill;
711
717 template <class RowPrefix, class ColPrefix>
718 std::string blockStyleClass(const RowPrefix &row_prefix,
719 const ColPrefix &col_prefix) const {
720 // here, you can potentially give a different style to each block
721 return "matrix-block";
722 }
723
725 bool write_block_title = true;
726
732 template <class Stream, class RowPrefix, class ColPrefix, class Block>
733 void writeBlockTitle(Stream& out, const RowPrefix &row_prefix,
734 const ColPrefix &col_prefix,
735 const Block &block) const {
736 if (this->write_block_title) {
737 out << "<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;
743 out << "</title>\n";
744 }
745 }
746
768 template <class Stream, class RowPrefix, class ColPrefix, class Block>
769 void writeSVGBlock(Stream &out,
770 const RowPrefix &row_prefix,
771 const ColPrefix &col_prefix, const Block block,
772 const std::array<std::size_t, 4> &svg_box) const {
773 // get bounding box values
774 auto &[x_off, y_off, width, height] = svg_box;
775 // get style class
776 std::string block_class = this->blockStyleClass(row_prefix, col_prefix);
777 // write a rectangle on the bounding box
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>{})
781 if (color_fill)
782 out << " style='fill-opacity: 1;fill:" << color_fill(double(block)) << "'";
783
784 out << ">\n";
785 // give the rectangle a title (in html this shows info about the block)
786 this->writeBlockTitle(out,row_prefix, col_prefix, block);
787 // close rectangle
788 out << "</rect>\n";
789 }
790 };
791
806 template <class Mat, class SVGOptions = DefaultSVGMatrixOptions>
807 void writeSVGMatrix(const Mat &mat, std::ostream &out,
808 SVGOptions opts = {}) {
809 // We need a vector that can fit all the multi-indices for rows and colums
811 // Call overload for Mat type
812 Impl::writeSVGMatrix(mat, out, opts, IndexPrefix{}, IndexPrefix{});
813 }
814
817} // namespace Dune
818
819#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
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
STL namespace.
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
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 21, 23:30, 2024)