Dune Core Modules (2.3.1)

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_ISTLIO_HH
4#define DUNE_ISTLIO_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"
19#include <dune/common/unused.hh>
20
21#include <dune/istl/matrix.hh>
23#include "bcrsmatrix.hh"
24
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 int precision)
55 {
56 for (typename V::ConstIterator i=v.begin(); i!=v.end(); ++i)
57 recursive_printvector(s,*i,rowtext,counter,columns,width,precision);
58 }
59
67 template<class K, int n>
68 void recursive_printvector (std::ostream& s, const FieldVector<K,n>& v,
69 std::string rowtext, int& counter, int columns,
70 int width, int precision)
71 {
72 DUNE_UNUSED_PARAMETER(precision);
73 // we now can print n numbers
74 for (int i=0; i<n; i++)
75 {
76 if (counter%columns==0)
77 {
78 s << rowtext; // start a new row
79 s << " "; // space in front of each entry
80 s.width(4); // set width for counter
81 s << counter; // number of first entry in a line
82 }
83 s << " "; // space in front of each entry
84 s.width(width); // set width for each entry anew
85 s << v[i]; // yeah, the number !
86 counter++; // increment the counter
87 if (counter%columns==0)
88 s << std::endl; // start a new line
89 }
90 }
91
92
100 template<class V>
101 void printvector (std::ostream& s, const V& v, std::string title,
102 std::string rowtext, int columns=1, int width=10,
103 int precision=2)
104 {
105 // count the numbers printed to make columns
106 int counter=0;
107
108 // remember old flags
109 std::ios_base::fmtflags oldflags = s.flags();
110
111 // set the output format
112 s.setf(std::ios_base::scientific, std::ios_base::floatfield);
113 int oldprec = s.precision();
114 s.precision(precision);
115
116 // print title
117 s << title << " [blocks=" << v.N() << ",dimension=" << v.dim() << "]"
118 << std::endl;
119
120 // print data from all blocks
121 recursive_printvector(s,v,rowtext,counter,columns,width,precision);
122
123 // check if new line is required
124 if (counter%columns!=0)
125 s << std::endl;
126
127 // reset the output format
128 s.flags(oldflags);
129 s.precision(oldprec);
130 }
131
132
134 //
135 // pretty printing of matrices
136 //
137
145 inline void fill_row (std::ostream& s, int m, int width, int precision)
146 {
147 DUNE_UNUSED_PARAMETER(precision);
148 for (int j=0; j<m; j++)
149 {
150 s << " "; // space in front of each entry
151 s.width(width); // set width for each entry anew
152 s << "."; // yeah, the number !
153 }
154 }
155
163 template<class M>
164 void print_row (std::ostream& s, const M& A, typename M::size_type I,
165 typename M::size_type J, typename M::size_type therow,
166 int width, int precision)
167 {
168 typename M::size_type i0=I;
169 for (typename M::size_type i=0; i<A.N(); i++)
170 {
171 if (therow>=i0 && therow<i0+MatrixDimension<M>::rowdim(A,i))
172 {
173 // the row is in this block row !
174 typename M::size_type j0=J;
175 for (typename M::size_type j=0; j<A.M(); j++)
176 {
177 // find this block
178 typename M::ConstColIterator it = A[i].find(j);
179
180 // print row or filler
181 if (it!=A[i].end())
182 print_row(s,*it,i0,j0,therow,width,precision);
183 else
184 fill_row(s,MatrixDimension<M>::coldim(A,j),width,precision);
185
186 // advance columns
187 j0 += MatrixDimension<M>::coldim(A,j);
188 }
189 }
190 // advance rows
191 i0 += MatrixDimension<M>::rowdim(A,i);
192 }
193 }
194
202 template<class K, int n, int m>
203 void print_row (std::ostream& s, const FieldMatrix<K,n,m>& A,
204 typename FieldMatrix<K,n,m>::size_type I,
205 typename FieldMatrix<K,n,m>::size_type J,
206 typename FieldMatrix<K,n,m>::size_type therow, int width,
207 int precision)
208 {
210 DUNE_UNUSED_PARAMETER(precision);
211
212 typedef typename FieldMatrix<K,n,m>::size_type size_type;
213
214 for (size_type i=0; i<n; i++)
215 if (I+i==therow)
216 for (int j=0; j<m; j++)
217 {
218 s << " "; // space in front of each entry
219 s.width(width); // set width for each entry anew
220 s << A[i][j]; // yeah, the number !
221 }
222 }
223
231 template<class K>
232 void print_row (std::ostream& s, const FieldMatrix<K,1,1>& A,
233 typename FieldMatrix<K,1,1>::size_type I,
234 typename FieldMatrix<K,1,1>::size_type J,
235 typename FieldMatrix<K,1,1>::size_type therow,
236 int width, int precision)
237 {
239 DUNE_UNUSED_PARAMETER(precision);
240
241 if (I==therow)
242 {
243 s << " "; // space in front of each entry
244 s.width(width); // set width for each entry anew
245 s << static_cast<K>(A); // yeah, the number !
246 }
247 }
248
257 template<class M>
258 void printmatrix (std::ostream& s, const M& A, std::string title,
259 std::string rowtext, int width=10, int precision=2)
260 {
261
262 // remember old flags
263 std::ios_base::fmtflags oldflags = s.flags();
264
265 // set the output format
266 s.setf(std::ios_base::scientific, std::ios_base::floatfield);
267 int oldprec = s.precision();
268 s.precision(precision);
269
270 // print title
271 s << title
272 << " [n=" << A.N()
273 << ",m=" << A.M()
274 << ",rowdim=" << MatrixDimension<M>::rowdim(A)
275 << ",coldim=" << MatrixDimension<M>::coldim(A)
276 << "]" << std::endl;
277
278 // print all rows
279 for (typename M::size_type i=0; i<MatrixDimension<M>::rowdim(A); i++)
280 {
281 s << rowtext; // start a new row
282 s << " "; // space in front of each entry
283 s.width(4); // set width for counter
284 s << i; // number of first entry in a line
285 print_row(s,A,0,0,i,width,precision); // generic print
286 s << std::endl; // start a new line
287 }
288
289 // reset the output format
290 s.flags(oldflags);
291 s.precision(oldprec);
292 }
293
315 template<class B, int n, int m, class A>
316 void printSparseMatrix(std::ostream& s,
317 const BCRSMatrix<FieldMatrix<B,n,m>,A>& mat,
318 std::string title, std::string rowtext,
319 int width=3, int precision=2)
320 {
322 // remember old flags
323 std::ios_base::fmtflags oldflags = s.flags();
324 // set the output format
325 s.setf(std::ios_base::scientific, std::ios_base::floatfield);
326 int oldprec = s.precision();
327 s.precision(precision);
328 // print title
329 s << title
330 << " [n=" << mat.N()
331 << ",m=" << mat.M()
332 << ",rowdim=" << MatrixDimension<Matrix>::rowdim(mat)
333 << ",coldim=" << MatrixDimension<Matrix>::coldim(mat)
334 << "]" << std::endl;
335
336 typedef typename Matrix::ConstRowIterator Row;
337
338 for(Row row=mat.begin(); row != mat.end(); ++row) {
339 int skipcols=0;
340 bool reachedEnd=false;
341
342 while(!reachedEnd) {
343 for(int innerrow=0; innerrow<n; ++innerrow) {
344 int count=0;
345 typedef typename Matrix::ConstColIterator Col;
346 Col col=row->begin();
347 for(; col != row->end(); ++col,++count) {
348 if(count<skipcols)
349 continue;
350 if(count>=skipcols+width)
351 break;
352 if(innerrow==0) {
353 if(count==skipcols) {
354 s << rowtext; // start a new row
355 s << " "; // space in front of each entry
356 s.width(4); // set width for counter
357 s << row.index()<<": "; // number of first entry in a line
358 }
359 s.width(4);
360 s<<col.index()<<": |";
361 } else {
362 if(count==skipcols) {
363 for(typename std::string::size_type i=0; i < rowtext.length(); i++)
364 s<<" ";
365 s<<" ";
366 }
367 s<<" |";
368 }
369 for(int innercol=0; innercol < m; ++innercol) {
370 s.width(9);
371 s<<(*col)[innerrow][innercol]<<" ";
372 }
373
374 s<<"|";
375 }
376 if(innerrow==n-1 && col==row->end())
377 reachedEnd = true;
378 else
379 s << std::endl;
380 }
381 skipcols += width;
382 s << std::endl;
383 }
384 s << std::endl;
385 }
386
387 // reset the output format
388 s.flags(oldflags);
389 s.precision(oldprec);
390 }
391
392 namespace
393 {
394 template<typename T>
395 struct MatlabPODWriter
396 {
397 static std::ostream& write(const T& t, std::ostream& s)
398 {
399 s << t;
400 return s;
401 }
402 };
403 template<typename T>
404 struct MatlabPODWriter<std::complex<T> >
405 {
406 static std::ostream& write(const std::complex<T>& t, std::ostream& s)
407 {
408 s << t.real() << " " << t.imag();
409 return s;
410 }
411 };
412 } // anonymous namespace
413
423 template <class FieldType, int dim>
424 void writeMatrixToMatlabHelper(const ScaledIdentityMatrix<FieldType,dim>& matrix, int rowOffset, int colOffset, std::ostream& s)
425 {
426 for (int i=0; i<dim; i++)
427 {
428 //+1 for Matlab numbering
429 s << rowOffset + i + 1 << " " << colOffset + i + 1 << " ";
430 MatlabPODWriter<FieldType>::write(matrix.scalar(), s)<< std::endl;
431 }
432 }
433
443 template <class FieldType, int dim>
444 void writeMatrixToMatlabHelper(const DiagonalMatrix<FieldType,dim>& matrix, int rowOffset, int colOffset, std::ostream& s)
445 {
446 for (int i=0; i<dim; i++)
447 {
448 //+1 for Matlab numbering
449 s << rowOffset + i + 1 << " " << colOffset + i + 1 << " ";
450 MatlabPODWriter<FieldType>::write(matrix.diagonal(i), s)<< std::endl;
451 }
452 }
453
463 template <class FieldType, int rows, int cols>
465 ( const FieldMatrix<FieldType,rows,cols>& matrix, int rowOffset,
466 int colOffset, std::ostream& s)
467 {
468 for (int i=0; i<rows; i++)
469 for (int j=0; j<cols; j++) {
470 //+1 for Matlab numbering
471 s << rowOffset + i + 1 << " " << colOffset + j + 1 << " ";
472 MatlabPODWriter<FieldType>::write(matrix[i][j], s)<< std::endl;
473 }
474 }
475
483 template <class MatrixType>
484 void writeMatrixToMatlabHelper(const MatrixType& matrix,
485 int externalRowOffset, int externalColOffset,
486 std::ostream& s)
487 {
488 // Precompute the accumulated sizes of the columns
489 std::vector<typename MatrixType::size_type> colOffset(matrix.M());
490 if (colOffset.size() > 0)
491 colOffset[0] = 0;
492
493 for (typename MatrixType::size_type i=0; i<matrix.M()-1; i++)
494 colOffset[i+1] = colOffset[i] +
495 MatrixDimension<MatrixType>::coldim(matrix,i);
496
497 typename MatrixType::size_type rowOffset = 0;
498
499 // Loop over all matrix rows
500 for (typename MatrixType::size_type rowIdx=0; rowIdx<matrix.N(); rowIdx++)
501 {
502
503 const typename MatrixType::row_type& row = matrix[rowIdx];
504
505 typename MatrixType::row_type::ConstIterator cIt = row.begin();
506 typename MatrixType::row_type::ConstIterator cEndIt = row.end();
507
508 // Loop over all columns in this row
509 for (; cIt!=cEndIt; ++cIt)
511 externalRowOffset+rowOffset,
512 externalColOffset + colOffset[cIt.index()],
513 s);
514
515 rowOffset += MatrixDimension<MatrixType>::rowdim(matrix, rowIdx);
516 }
517
518 }
519
539 template <class MatrixType>
540 void writeMatrixToMatlab(const MatrixType& matrix,
541 const std::string& filename, int outputPrecision = 18)
542 {
543 std::ofstream outStream(filename.c_str());
544 int oldPrecision = outStream.precision();
545 outStream.precision(outputPrecision);
546
547 writeMatrixToMatlabHelper(matrix, 0, 0, outStream);
548 outStream.precision(oldPrecision);
549 }
550
553} // namespace Dune
554
555#endif
Implementation of the BCRSMatrix class.
A sparse block matrix with compressed row storage.
Definition: bcrsmatrix.hh:414
A diagonal matrix of static size.
Definition: diagonalmatrix.hh:49
A dense n x m matrix.
Definition: fmatrix.hh:67
A generic dynamic dense matrix.
Definition: matrix.hh:25
RowIterator end()
Get iterator to one beyond last row.
Definition: matrix.hh:85
RowIterator begin()
Get iterator to first row.
Definition: matrix.hh:79
row_type::const_iterator ConstColIterator
Const iterator for the entries of each row.
Definition: matrix.hh:53
size_type M() const
Return the number of columns.
Definition: matrix.hh:165
size_type N() const
Return the number of rows.
Definition: matrix.hh:160
A multiple of the identity matrix of static size.
Definition: scaledidmatrix.hh:28
const K & scalar() const
Get const reference to the scalar diagonal value.
Definition: scaledidmatrix.hh:459
ConstIterator class for sequential access.
Definition: vbvector.hh:647
This file implements a quadratic diagonal matrix of fixed size.
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.
const K & diagonal(size_type i) const
Get const reference to diagonal entry.
Definition: diagonalmatrix.hh:495
void writeMatrixToMatlab(const MatrixType &matrix, const std::string &filename, int outputPrecision=18)
Writes sparse matrix in a Matlab-readable format.
Definition: io.hh:540
void recursive_printvector(std::ostream &s, const V &v, std::string rowtext, int &counter, int columns, int width, int precision)
Recursively print all the blocks.
Definition: io.hh:52
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:258
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:101
void print_row(std::ostream &s, const M &A, typename M::size_type I, typename M::size_type J, typename M::size_type therow, int width, int precision)
Print one row of a matrix.
Definition: io.hh:164
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:316
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:145
void writeMatrixToMatlabHelper(const ScaledIdentityMatrix< FieldType, dim > &matrix, int rowOffset, int colOffset, std::ostream &s)
Helper method for the writeMatrixToMatlab routine.
Definition: io.hh:424
A dynamic dense block matrix class.
Some handy generic functions for ISTL matrices.
Dune namespace.
Definition: alignment.hh:14
STL namespace.
This file implements a quadratic matrix of fixed size which is a multiple of the identity.
Definition of the DUNE_UNUSED macro for the case that config.h is not available.
#define DUNE_UNUSED_PARAMETER(parm)
A macro to mark intentional unused function parameters with.
Definition: unused.hh:18
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 12, 23:30, 2024)