Dune Core Modules (2.7.0)

io.hh
Go to the documentation of this file.
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3#ifndef DUNE_ISTL_IO_HH
4#define DUNE_ISTL_IO_HH
5
6#include <cmath>
7#include <complex>
8#include <limits>
9#include <ios>
10#include <iomanip>
11#include <fstream>
12#include <string>
13
14#include "matrixutils.hh"
15#include "istlexception.hh"
18#include <dune/common/hybridutilities.hh>
19#include <dune/common/unused.hh>
20
22
23
24namespace Dune {
25
38 //
39 // pretty printing of vectors
40 //
41
49 template<class V>
50 void recursive_printvector (std::ostream& s, const V& v, std::string rowtext,
51 int& counter, int columns, int width)
52 {
54 [&](auto id) {
55 // Print one number
56 if (counter%columns==0)
57 {
58 s << rowtext; // start a new row
59 s << " "; // space in front of each entry
60 s.width(4); // set width for counter
61 s << counter; // number of first entry in a line
62 }
63 s << " "; // space in front of each entry
64 s.width(width); // set width for each entry anew
65 s << v; // yeah, the number !
66 counter++; // increment the counter
67 if (counter%columns==0)
68 s << std::endl; // start a new line
69 },
70 [&](auto id) {
71 // Recursively print a vector
72 for (const auto& entry : id(v))
73 recursive_printvector(s,id(entry),rowtext,counter,columns,width);
74 });
75 }
76
77
85 template<class V>
86 void printvector (std::ostream& s, const V& v, std::string title,
87 std::string rowtext, int columns=1, int width=10,
88 int precision=2)
89 {
90 // count the numbers printed to make columns
91 int counter=0;
92
93 // remember old flags
94 std::ios_base::fmtflags oldflags = s.flags();
95
96 // set the output format
97 s.setf(std::ios_base::scientific, std::ios_base::floatfield);
98 int oldprec = s.precision();
99 s.precision(precision);
100
101 // print title
102 s << title << " [blocks=" << v.N() << ",dimension=" << v.dim() << "]"
103 << std::endl;
104
105 // print data from all blocks
106 recursive_printvector(s,v,rowtext,counter,columns,width);
107
108 // check if new line is required
109 if (counter%columns!=0)
110 s << std::endl;
111
112 // reset the output format
113 s.flags(oldflags);
114 s.precision(oldprec);
115 }
116
117
119 //
120 // pretty printing of matrices
121 //
122
130 inline void fill_row (std::ostream& s, int m, int width, int precision)
131 {
132 DUNE_UNUSED_PARAMETER(precision);
133 for (int j=0; j<m; j++)
134 {
135 s << " "; // space in front of each entry
136 s.width(width); // set width for each entry anew
137 s << "."; // yeah, the number !
138 }
139 }
140
148 template<class K>
149 void print_row (std::ostream& s, const K& value,
150 typename FieldMatrix<K,1,1>::size_type I,
151 typename FieldMatrix<K,1,1>::size_type J,
152 typename FieldMatrix<K,1,1>::size_type therow,
153 int width, int precision,
154 typename std::enable_if_t<Dune::IsNumber<K>::value>* sfinae = nullptr)
155 {
158 DUNE_UNUSED_PARAMETER(therow);
159 DUNE_UNUSED_PARAMETER(precision);
160
161 s << " "; // space in front of each entry
162 s.width(width); // set width for each entry anew
163 s << value;
164 }
165
173 template<class M>
174 void print_row (std::ostream& s, const M& A, typename M::size_type I,
175 typename M::size_type J, typename M::size_type therow,
176 int width, int precision,
177 typename std::enable_if_t<!Dune::IsNumber<M>::value>* sfinae = nullptr)
178 {
179 typename M::size_type i0=I;
180 for (typename M::size_type i=0; i<A.N(); i++)
181 {
182 if (therow>=i0 && therow<i0+MatrixDimension<M>::rowdim(A,i))
183 {
184 // the row is in this block row !
185 typename M::size_type j0=J;
186 for (typename M::size_type j=0; j<A.M(); j++)
187 {
188 // find this block
189 typename M::ConstColIterator it = A[i].find(j);
190
191 // print row or filler
192 if (it!=A[i].end())
193 print_row(s,*it,i0,j0,therow,width,precision);
194 else
195 fill_row(s,MatrixDimension<M>::coldim(A,j),width,precision);
196
197 // advance columns
198 j0 += MatrixDimension<M>::coldim(A,j);
199 }
200 }
201 // advance rows
202 i0 += MatrixDimension<M>::rowdim(A,i);
203 }
204 }
205
214 template<class M>
215 void printmatrix (std::ostream& s, const M& A, std::string title,
216 std::string rowtext, int width=10, int precision=2)
217 {
218
219 // remember old flags
220 std::ios_base::fmtflags oldflags = s.flags();
221
222 // set the output format
223 s.setf(std::ios_base::scientific, std::ios_base::floatfield);
224 int oldprec = s.precision();
225 s.precision(precision);
226
227 // print title
228 s << title
229 << " [n=" << A.N()
230 << ",m=" << A.M()
231 << ",rowdim=" << MatrixDimension<M>::rowdim(A)
232 << ",coldim=" << MatrixDimension<M>::coldim(A)
233 << "]" << std::endl;
234
235 // print all rows
236 for (typename M::size_type i=0; i<MatrixDimension<M>::rowdim(A); i++)
237 {
238 s << rowtext; // start a new row
239 s << " "; // space in front of each entry
240 s.width(4); // set width for counter
241 s << i; // number of first entry in a line
242 print_row(s,A,0,0,i,width,precision); // generic print
243 s << std::endl; // start a new line
244 }
245
246 // reset the output format
247 s.flags(oldflags);
248 s.precision(oldprec);
249 }
250
272 template<class B, int n, int m, class A>
273 void printSparseMatrix(std::ostream& s,
274 const BCRSMatrix<FieldMatrix<B,n,m>,A>& mat,
275 std::string title, std::string rowtext,
276 int width=3, int precision=2)
277 {
279 // remember old flags
280 std::ios_base::fmtflags oldflags = s.flags();
281 // set the output format
282 s.setf(std::ios_base::scientific, std::ios_base::floatfield);
283 int oldprec = s.precision();
284 s.precision(precision);
285 // print title
286 s << title
287 << " [n=" << mat.N()
288 << ",m=" << mat.M()
289 << ",rowdim=" << MatrixDimension<Matrix>::rowdim(mat)
290 << ",coldim=" << MatrixDimension<Matrix>::coldim(mat)
291 << "]" << std::endl;
292
293 typedef typename Matrix::ConstRowIterator Row;
294
295 for(Row row=mat.begin(); row != mat.end(); ++row) {
296 int skipcols=0;
297 bool reachedEnd=false;
298
299 while(!reachedEnd) {
300 for(int innerrow=0; innerrow<n; ++innerrow) {
301 int count=0;
302 typedef typename Matrix::ConstColIterator Col;
303 Col col=row->begin();
304 for(; col != row->end(); ++col,++count) {
305 if(count<skipcols)
306 continue;
307 if(count>=skipcols+width)
308 break;
309 if(innerrow==0) {
310 if(count==skipcols) {
311 s << rowtext; // start a new row
312 s << " "; // space in front of each entry
313 s.width(4); // set width for counter
314 s << row.index()<<": "; // number of first entry in a line
315 }
316 s.width(4);
317 s<<col.index()<<": |";
318 } else {
319 if(count==skipcols) {
320 for(typename std::string::size_type i=0; i < rowtext.length(); i++)
321 s<<" ";
322 s<<" ";
323 }
324 s<<" |";
325 }
326 for(int innercol=0; innercol < m; ++innercol) {
327 s.width(9);
328 s<<(*col)[innerrow][innercol]<<" ";
329 }
330
331 s<<"|";
332 }
333 if(innerrow==n-1 && col==row->end())
334 reachedEnd = true;
335 else
336 s << std::endl;
337 }
338 skipcols += width;
339 s << std::endl;
340 }
341 s << std::endl;
342 }
343
344 // reset the output format
345 s.flags(oldflags);
346 s.precision(oldprec);
347 }
348
349 namespace
350 {
351 template<typename T>
352 struct MatlabPODWriter
353 {
354 static std::ostream& write(const T& t, std::ostream& s)
355 {
356 s << t;
357 return s;
358 }
359 };
360 template<typename T>
361 struct MatlabPODWriter<std::complex<T> >
362 {
363 static std::ostream& write(const std::complex<T>& t, std::ostream& s)
364 {
365 s << t.real() << " " << t.imag();
366 return s;
367 }
368 };
369 } // anonymous namespace
370
380 template <class FieldType>
381 void writeMatrixToMatlabHelper(const FieldType& value,
382 int rowOffset, int colOffset,
383 std::ostream& s,
384 typename std::enable_if_t<Dune::IsNumber<FieldType>::value>* sfinae = nullptr)
385 {
386 //+1 for Matlab numbering
387 s << rowOffset + 1 << " " << colOffset + 1 << " ";
388 MatlabPODWriter<FieldType>::write(value, s)<< std::endl;
389 }
390
398 template <class MatrixType>
399 void writeMatrixToMatlabHelper(const MatrixType& matrix,
400 int externalRowOffset, int externalColOffset,
401 std::ostream& s,
402 typename std::enable_if_t<!Dune::IsNumber<MatrixType>::value>* sfinae = nullptr)
403 {
404 // Precompute the accumulated sizes of the columns
405 std::vector<typename MatrixType::size_type> colOffset(matrix.M());
406 if (colOffset.size() > 0)
407 colOffset[0] = 0;
408
409 for (typename MatrixType::size_type i=0; i<matrix.M()-1; i++)
410 colOffset[i+1] = colOffset[i] +
411 MatrixDimension<MatrixType>::coldim(matrix,i);
412
413 typename MatrixType::size_type rowOffset = 0;
414
415 // Loop over all matrix rows
416 for (typename MatrixType::size_type rowIdx=0; rowIdx<matrix.N(); rowIdx++)
417 {
418 auto cIt = matrix[rowIdx].begin();
419 auto cEndIt = matrix[rowIdx].end();
420
421 // Loop over all columns in this row
422 for (; cIt!=cEndIt; ++cIt)
424 externalRowOffset+rowOffset,
425 externalColOffset + colOffset[cIt.index()],
426 s);
427
428 rowOffset += MatrixDimension<MatrixType>::rowdim(matrix, rowIdx);
429 }
430
431 }
432
452 template <class MatrixType>
453 void writeMatrixToMatlab(const MatrixType& matrix,
454 const std::string& filename, int outputPrecision = 18)
455 {
456 std::ofstream outStream(filename.c_str());
457 int oldPrecision = outStream.precision();
458 outStream.precision(outputPrecision);
459
460 writeMatrixToMatlabHelper(matrix, 0, 0, outStream);
461 outStream.precision(oldPrecision);
462 }
463
464 // Recursively write vector entries to a stream
465 template<class V>
466 void writeVectorToMatlabHelper (const V& v, std::ostream& stream)
467 {
468 Hybrid::ifElse(IsNumber<V>(),
469 [&](auto id) {
470 stream << id(v) << std::endl;
471 },
472 [&](auto id) {
473 for (const auto& entry : id(v))
474 writeVectorToMatlabHelper(entry, stream);
475 });
476 }
477
495 template <class VectorType>
496 void writeVectorToMatlab(const VectorType& vector,
497 const std::string& filename, int outputPrecision = 18)
498 {
499 std::ofstream outStream(filename.c_str());
500 int oldPrecision = outStream.precision();
501 outStream.precision(outputPrecision);
502
503 writeVectorToMatlabHelper(vector, outStream);
504 outStream.precision(oldPrecision);
505 }
506
509} // namespace Dune
510
511#endif
Implementation of the BCRSMatrix class.
A sparse block matrix with compressed row storage.
Definition: bcrsmatrix.hh:425
A dense n x m matrix.
Definition: fmatrix.hh:69
ConstIterator class for sequential access.
Definition: matrix.hh:401
A generic dynamic dense matrix.
Definition: matrix.hh:558
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:586
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
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.
#define DUNE_UNUSED_PARAMETER(parm)
A macro to mark intentionally unused function parameters with.
Definition: unused.hh:25
decltype(auto) ifElse(const Condition &condition, IfFunc &&ifFunc, ElseFunc &&elseFunc)
A conditional expression.
Definition: hybridutilities.hh:355
void writeMatrixToMatlab(const MatrixType &matrix, const std::string &filename, int outputPrecision=18)
Writes sparse matrix in a Matlab-readable format.
Definition: io.hh:453
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:149
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:215
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:86
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:381
void writeVectorToMatlab(const VectorType &vector, const std::string &filename, int outputPrecision=18)
Writes vectors in a Matlab-readable format.
Definition: io.hh:496
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:50
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:273
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:130
Some handy generic functions for ISTL matrices.
Dune namespace.
Definition: alignedallocator.hh:14
STL namespace.
Whether this type acts as a scalar in the context of (hierarchically blocked) containers.
Definition: typetraits.hh:163
Definition of the DUNE_UNUSED macro for the case that config.h is not available.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)