dune-istl  2.2.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"
16 #include <dune/common/fvector.hh>
17 #include <dune/common/fmatrix.hh>
18 
19 #include <dune/istl/matrix.hh>
20 #include "bcrsmatrix.hh"
21 
22 
23 namespace Dune {
24 
36 
37  //
38  // pretty printing of vectors
39  //
40 
41  // recursively print all the blocks
47  template<class V>
48  void recursive_printvector (std::ostream& s, const V& v, std::string rowtext,
49  int& counter, int columns, int width,
50  int precision)
51  {
52  for (typename V::ConstIterator i=v.begin(); i!=v.end(); ++i)
53  recursive_printvector(s,*i,rowtext,counter,columns,width,precision);
54  }
55 
56  // recursively print all the blocks -- specialization for FieldVector
62  template<class K, int n>
63  void recursive_printvector (std::ostream& s, const FieldVector<K,n>& v,
64  std::string rowtext, int& counter, int columns,
65  int width, int precision)
66  {
67  // we now can print n numbers
68  for (int i=0; i<n; i++)
69  {
70  if (counter%columns==0)
71  {
72  s << rowtext; // start a new row
73  s << " "; // space in front of each entry
74  s.width(4); // set width for counter
75  s << counter; // number of first entry in a line
76  }
77  s << " "; // space in front of each entry
78  s.width(width); // set width for each entry anew
79  s << v[i]; // yeah, the number !
80  counter++; // increment the counter
81  if (counter%columns==0)
82  s << std::endl; // start a new line
83  }
84  }
85 
86 
88 
93  template<class V>
94  void printvector (std::ostream& s, const V& v, std::string title,
95  std::string rowtext, int columns=1, int width=10,
96  int precision=2)
97  {
98  // count the numbers printed to make columns
99  int counter=0;
100 
101  // remember old flags
102  std::ios_base::fmtflags oldflags = s.flags();
103 
104  // set the output format
105  s.setf(std::ios_base::scientific, std::ios_base::floatfield);
106  int oldprec = s.precision();
107  s.precision(precision);
108 
109  // print title
110  s << title << " [blocks=" << v.N() << ",dimension=" << v.dim() << "]"
111  << std::endl;
112 
113  // print data from all blocks
114  recursive_printvector(s,v,rowtext,counter,columns,width,precision);
115 
116  // check if new line is required
117  if (counter%columns!=0)
118  s << std::endl;
119 
120  // reset the output format
121  s.flags(oldflags);
122  s.precision(oldprec);
123  }
124 
125 
127  //
128  // pretty printing of matrices
129  //
130 
132 
137  inline void fill_row (std::ostream& s, int m, int width, int precision)
138  {
139  for (int j=0; j<m; j++)
140  {
141  s << " "; // space in front of each entry
142  s.width(width); // set width for each entry anew
143  s << "."; // yeah, the number !
144  }
145  }
146 
148 
153  template<class M>
154  void print_row (std::ostream& s, const M& A, typename M::size_type I,
155  typename M::size_type J, typename M::size_type therow,
156  int width, int precision)
157  {
158  typename M::size_type i0=I;
159  for (typename M::size_type i=0; i<A.N(); i++)
160  {
161  if (therow>=i0 && therow<i0+MatrixDimension<M>::rowdim(A,i))
162  {
163  // the row is in this block row !
164  typename M::size_type j0=J;
165  for (typename M::size_type j=0; j<A.M(); j++)
166  {
167  // find this block
168  typename M::ConstColIterator it = A[i].find(j);
169 
170  // print row or filler
171  if (it!=A[i].end())
172  print_row(s,*it,i0,j0,therow,width,precision);
173  else
174  fill_row(s,MatrixDimension<M>::coldim(A,j),width,precision);
175 
176  // advance columns
177  j0 += MatrixDimension<M>::coldim(A,j);
178  }
179  }
180  // advance rows
181  i0 += MatrixDimension<M>::rowdim(A,i);
182  }
183  }
184 
186 
191  template<class K, int n, int m>
192  void print_row (std::ostream& s, const FieldMatrix<K,n,m>& A,
193  typename FieldMatrix<K,n,m>::size_type I,
194  typename FieldMatrix<K,n,m>::size_type J,
195  typename FieldMatrix<K,n,m>::size_type therow, int width,
196  int precision)
197  {
198  typedef typename FieldMatrix<K,n,m>::size_type size_type;
199 
200  for (size_type i=0; i<n; i++)
201  if (I+i==therow)
202  for (int j=0; j<m; j++)
203  {
204  s << " "; // space in front of each entry
205  s.width(width); // set width for each entry anew
206  s << A[i][j]; // yeah, the number !
207  }
208  }
209 
211 
216  template<class K>
217  void print_row (std::ostream& s, const FieldMatrix<K,1,1>& A,
218  typename FieldMatrix<K,1,1>::size_type I,
219  typename FieldMatrix<K,1,1>::size_type J,
220  typename FieldMatrix<K,1,1>::size_type therow,
221  int width, int precision)
222  {
223  if (I==therow)
224  {
225  s << " "; // space in front of each entry
226  s.width(width); // set width for each entry anew
227  s << static_cast<K>(A); // yeah, the number !
228  }
229  }
230 
232 
238  template<class M>
239  void printmatrix (std::ostream& s, const M& A, std::string title,
240  std::string rowtext, int width=10, int precision=2)
241  {
242 
243  // remember old flags
244  std::ios_base::fmtflags oldflags = s.flags();
245 
246  // set the output format
247  s.setf(std::ios_base::scientific, std::ios_base::floatfield);
248  int oldprec = s.precision();
249  s.precision(precision);
250 
251  // print title
252  s << title
253  << " [n=" << A.N()
254  << ",m=" << A.M()
255  << ",rowdim=" << MatrixDimension<M>::rowdim(A)
256  << ",coldim=" << MatrixDimension<M>::coldim(A)
257  << "]" << std::endl;
258 
259  // print all rows
260  for (typename M::size_type i=0; i<MatrixDimension<M>::rowdim(A); i++)
261  {
262  s << rowtext; // start a new row
263  s << " "; // space in front of each entry
264  s.width(4); // set width for counter
265  s << i; // number of first entry in a line
266  print_row(s,A,0,0,i,width,precision); // generic print
267  s << std::endl;// start a new line
268  }
269 
270  // reset the output format
271  s.flags(oldflags);
272  s.precision(oldprec);
273  }
274 
276 
295  template<class B, int n, int m, class A>
296  void printSparseMatrix(std::ostream& s,
297  const BCRSMatrix<FieldMatrix<B,n,m>,A>& mat,
298  std::string title, std::string rowtext,
299  int width=3, int precision=2)
300  {
302  // remember old flags
303  std::ios_base::fmtflags oldflags = s.flags();
304  // set the output format
305  s.setf(std::ios_base::scientific, std::ios_base::floatfield);
306  int oldprec = s.precision();
307  s.precision(precision);
308  // print title
309  s << title
310  << " [n=" << mat.N()
311  << ",m=" << mat.M()
312  << ",rowdim=" << MatrixDimension<Matrix>::rowdim(mat)
313  << ",coldim=" << MatrixDimension<Matrix>::coldim(mat)
314  << "]" << std::endl;
315 
316  typedef typename Matrix::ConstRowIterator Row;
317 
318  for(Row row=mat.begin(); row != mat.end();++row) {
319  int skipcols=0;
320  bool reachedEnd=false;
321 
322  while(!reachedEnd) {
323  for(int innerrow=0; innerrow<n; ++innerrow) {
324  int count=0;
325  typedef typename Matrix::ConstColIterator Col;
326  Col col=row->begin();
327  for(; col != row->end(); ++col,++count) {
328  if(count<skipcols)
329  continue;
330  if(count>=skipcols+width)
331  break;
332  if(innerrow==0) {
333  if(count==skipcols) {
334  s << rowtext; // start a new row
335  s << " "; // space in front of each entry
336  s.width(4); // set width for counter
337  s << row.index()<<": "; // number of first entry in a line
338  }
339  s.width(4);
340  s<<col.index()<<": |";
341  } else {
342  if(count==skipcols){
343  for(typename std::string::size_type i=0; i < rowtext.length(); i++)
344  s<<" ";
345  s<<" ";
346  }
347  s<<" |";
348  }
349  for(int innercol=0; innercol < m; ++innercol) {
350  s.width(9);
351  s<<(*col)[innerrow][innercol]<<" ";
352  }
353 
354  s<<"|";
355  }
356  if(innerrow==n-1 && col==row->end())
357  reachedEnd = true;
358  else
359  s << std::endl;
360  }
361  skipcols += width;
362  s << std::endl;
363  }
364  s << std::endl;
365  }
366 
367  // reset the output format
368  s.flags(oldflags);
369  s.precision(oldprec);
370  }
371 
372  namespace
373  {
374  template<typename T>
375  struct MatlabPODWriter
376  {
377  static std::ostream& write(const T& t, std::ostream& s)
378  {
379  s << t;
380  return s;
381  }
382  };
383  template<typename T>
384  struct MatlabPODWriter<std::complex<T> >
385  {
386  static std::ostream& write(const std::complex<T>& t, std::ostream& s)
387  {
388  s << t.real() << " " << t.imag();
389  return s;
390  }
391  };
392  } // anonymous namespace
393 
395 
402  template <class FieldType, int rows, int cols>
404  ( const FieldMatrix<FieldType,rows,cols>& matrix, int rowOffset,
405  int colOffset, std::ostream& s)
406  {
407  for (int i=0; i<rows; i++)
408  for (int j=0; j<cols; j++){
409  //+1 for Matlab numbering
410  s << rowOffset + i + 1 << " " << colOffset + j + 1 << " ";
411  MatlabPODWriter<FieldType>::write(matrix[i][j], s)<< std::endl;
412  }
413  }
414 
416 
421  template <class MatrixType>
422  void writeMatrixToMatlabHelper(const MatrixType& matrix,
423  int externalRowOffset, int externalColOffset,
424  std::ostream& s)
425  {
426  // Precompute the accumulated sizes of the columns
427  std::vector<typename MatrixType::size_type> colOffset(matrix.M());
428  if (colOffset.size() > 0)
429  colOffset[0] = 0;
430 
431  for (typename MatrixType::size_type i=0; i<matrix.M()-1; i++)
432  colOffset[i+1] = colOffset[i] +
434 
435  typename MatrixType::size_type rowOffset = 0;
436 
437  // Loop over all matrix rows
438  for (typename MatrixType::size_type rowIdx=0; rowIdx<matrix.N(); rowIdx++)
439  {
440 
441  const typename MatrixType::row_type& row = matrix[rowIdx];
442 
443  typename MatrixType::row_type::ConstIterator cIt = row.begin();
444  typename MatrixType::row_type::ConstIterator cEndIt = row.end();
445 
446  // Loop over all columns in this row
447  for (; cIt!=cEndIt; ++cIt)
449  externalRowOffset+rowOffset,
450  externalColOffset + colOffset[cIt.index()],
451  s);
452 
453  rowOffset += MatrixDimension<MatrixType>::rowdim(matrix, rowIdx);
454  }
455 
456  }
457 
459 
476  template <class MatrixType>
477  void writeMatrixToMatlab(const MatrixType& matrix,
478  const std::string& filename, int outputPrecision = 18)
479  {
480  std::ofstream outStream(filename.c_str());
481  int oldPrecision = outStream.precision();
482  outStream.precision(outputPrecision);
483 
484  writeMatrixToMatlabHelper(matrix, 0, 0, outStream);
485  outStream.precision(oldPrecision);
486  }
487 
490 } // namespace Dune
491 
492 #endif