/usr/include/dune/istl/io.hh is in libdune-istl-dev 2.2.1-2.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 | // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
// vi: set et ts=4 sw=2 sts=2:
#ifndef DUNE_ISTLIO_HH
#define DUNE_ISTLIO_HH
#include <cmath>
#include <complex>
#include <limits>
#include <ios>
#include <iomanip>
#include <fstream>
#include <string>
#include "matrixutils.hh"
#include "istlexception.hh"
#include <dune/common/fvector.hh>
#include <dune/common/fmatrix.hh>
#include <dune/istl/matrix.hh>
#include "bcrsmatrix.hh"
namespace Dune {
/**
@addtogroup ISTL_IO
@{
*/
/** \file
\brief Some generic functions for pretty printing vectors and matrices
*/
////////////////////////////////////////////////////////////////////////
//
// pretty printing of vectors
//
// recursively print all the blocks
/**
* \code
#include <dune/istl/io.hh>
* \endcode
*/
template<class V>
void recursive_printvector (std::ostream& s, const V& v, std::string rowtext,
int& counter, int columns, int width,
int precision)
{
for (typename V::ConstIterator i=v.begin(); i!=v.end(); ++i)
recursive_printvector(s,*i,rowtext,counter,columns,width,precision);
}
// recursively print all the blocks -- specialization for FieldVector
/**
* \code
#include <dune/istl/io.hh>
* \endcode
*/
template<class K, int n>
void recursive_printvector (std::ostream& s, const FieldVector<K,n>& v,
std::string rowtext, int& counter, int columns,
int width, int precision)
{
// we now can print n numbers
for (int i=0; i<n; i++)
{
if (counter%columns==0)
{
s << rowtext; // start a new row
s << " "; // space in front of each entry
s.width(4); // set width for counter
s << counter; // number of first entry in a line
}
s << " "; // space in front of each entry
s.width(width); // set width for each entry anew
s << v[i]; // yeah, the number !
counter++; // increment the counter
if (counter%columns==0)
s << std::endl; // start a new line
}
}
//! print an ISTL vector
/**
* \code
#include <dune/istl/io.hh>
* \endcode
*/
template<class V>
void printvector (std::ostream& s, const V& v, std::string title,
std::string rowtext, int columns=1, int width=10,
int precision=2)
{
// count the numbers printed to make columns
int counter=0;
// remember old flags
std::ios_base::fmtflags oldflags = s.flags();
// set the output format
s.setf(std::ios_base::scientific, std::ios_base::floatfield);
int oldprec = s.precision();
s.precision(precision);
// print title
s << title << " [blocks=" << v.N() << ",dimension=" << v.dim() << "]"
<< std::endl;
// print data from all blocks
recursive_printvector(s,v,rowtext,counter,columns,width,precision);
// check if new line is required
if (counter%columns!=0)
s << std::endl;
// reset the output format
s.flags(oldflags);
s.precision(oldprec);
}
////////////////////////////////////////////////////////////////////////
//
// pretty printing of matrices
//
//! print a row of zeros for a non-existing block
/**
* \code
#include <dune/istl/io.hh>
* \endcode
*/
inline void fill_row (std::ostream& s, int m, int width, int precision)
{
for (int j=0; j<m; j++)
{
s << " "; // space in front of each entry
s.width(width); // set width for each entry anew
s << "."; // yeah, the number !
}
}
//! print one row of a matrix
/**
* \code
#include <dune/istl/io.hh>
* \endcode
*/
template<class M>
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)
{
typename M::size_type i0=I;
for (typename M::size_type i=0; i<A.N(); i++)
{
if (therow>=i0 && therow<i0+MatrixDimension<M>::rowdim(A,i))
{
// the row is in this block row !
typename M::size_type j0=J;
for (typename M::size_type j=0; j<A.M(); j++)
{
// find this block
typename M::ConstColIterator it = A[i].find(j);
// print row or filler
if (it!=A[i].end())
print_row(s,*it,i0,j0,therow,width,precision);
else
fill_row(s,MatrixDimension<M>::coldim(A,j),width,precision);
// advance columns
j0 += MatrixDimension<M>::coldim(A,j);
}
}
// advance rows
i0 += MatrixDimension<M>::rowdim(A,i);
}
}
//! print one row of a matrix, specialization for FieldMatrix
/**
* \code
#include <dune/istl/io.hh>
* \endcode
*/
template<class K, int n, int m>
void print_row (std::ostream& s, const FieldMatrix<K,n,m>& A,
typename FieldMatrix<K,n,m>::size_type I,
typename FieldMatrix<K,n,m>::size_type J,
typename FieldMatrix<K,n,m>::size_type therow, int width,
int precision)
{
typedef typename FieldMatrix<K,n,m>::size_type size_type;
for (size_type i=0; i<n; i++)
if (I+i==therow)
for (int j=0; j<m; j++)
{
s << " "; // space in front of each entry
s.width(width); // set width for each entry anew
s << A[i][j]; // yeah, the number !
}
}
//! print one row of a matrix, specialization for FieldMatrix<K,1,1>
/**
* \code
#include <dune/istl/io.hh>
* \endcode
*/
template<class K>
void print_row (std::ostream& s, const FieldMatrix<K,1,1>& A,
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)
{
if (I==therow)
{
s << " "; // space in front of each entry
s.width(width); // set width for each entry anew
s << static_cast<K>(A); // yeah, the number !
}
}
//! Prints a generic block matrix
/**
* \code
#include <dune/istl/io.hh>
* \endcode
* \bug Empty rows and columns are omitted by this method. (FlySpray #7)
*/
template<class M>
void printmatrix (std::ostream& s, const M& A, std::string title,
std::string rowtext, int width=10, int precision=2)
{
// remember old flags
std::ios_base::fmtflags oldflags = s.flags();
// set the output format
s.setf(std::ios_base::scientific, std::ios_base::floatfield);
int oldprec = s.precision();
s.precision(precision);
// print title
s << title
<< " [n=" << A.N()
<< ",m=" << A.M()
<< ",rowdim=" << MatrixDimension<M>::rowdim(A)
<< ",coldim=" << MatrixDimension<M>::coldim(A)
<< "]" << std::endl;
// print all rows
for (typename M::size_type i=0; i<MatrixDimension<M>::rowdim(A); i++)
{
s << rowtext; // start a new row
s << " "; // space in front of each entry
s.width(4); // set width for counter
s << i; // number of first entry in a line
print_row(s,A,0,0,i,width,precision); // generic print
s << std::endl;// start a new line
}
// reset the output format
s.flags(oldflags);
s.precision(oldprec);
}
//! Prints a BCRSMatrix with fixed sized blocks.
/**
* \code
#include <dune/istl/io.hh>
* \endcode
*
* Only the nonzero entries will be printed as matrix blocks
* together with their
* corresponding column index and all others will be omitted.
*
* This might be preferable over printmatrix in the case of big
* sparse matrices with nonscalar blocks.
*
* @param s The ostream to print to.
* @param mat The matrix to print.
* @param title The title for the matrix.
* @param rowtext The text to prepend to each print out of a matrix row.
* @param width The number of nonzero blocks to print in one line.
* @param precision The precision to use when printing the numbers.
*/
template<class B, int n, int m, class A>
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)
{
typedef BCRSMatrix<FieldMatrix<B,n,m>,A> Matrix;
// remember old flags
std::ios_base::fmtflags oldflags = s.flags();
// set the output format
s.setf(std::ios_base::scientific, std::ios_base::floatfield);
int oldprec = s.precision();
s.precision(precision);
// print title
s << title
<< " [n=" << mat.N()
<< ",m=" << mat.M()
<< ",rowdim=" << MatrixDimension<Matrix>::rowdim(mat)
<< ",coldim=" << MatrixDimension<Matrix>::coldim(mat)
<< "]" << std::endl;
typedef typename Matrix::ConstRowIterator Row;
for(Row row=mat.begin(); row != mat.end();++row) {
int skipcols=0;
bool reachedEnd=false;
while(!reachedEnd) {
for(int innerrow=0; innerrow<n; ++innerrow) {
int count=0;
typedef typename Matrix::ConstColIterator Col;
Col col=row->begin();
for(; col != row->end(); ++col,++count) {
if(count<skipcols)
continue;
if(count>=skipcols+width)
break;
if(innerrow==0) {
if(count==skipcols) {
s << rowtext; // start a new row
s << " "; // space in front of each entry
s.width(4); // set width for counter
s << row.index()<<": "; // number of first entry in a line
}
s.width(4);
s<<col.index()<<": |";
} else {
if(count==skipcols){
for(typename std::string::size_type i=0; i < rowtext.length(); i++)
s<<" ";
s<<" ";
}
s<<" |";
}
for(int innercol=0; innercol < m; ++innercol) {
s.width(9);
s<<(*col)[innerrow][innercol]<<" ";
}
s<<"|";
}
if(innerrow==n-1 && col==row->end())
reachedEnd = true;
else
s << std::endl;
}
skipcols += width;
s << std::endl;
}
s << std::endl;
}
// reset the output format
s.flags(oldflags);
s.precision(oldprec);
}
namespace
{
template<typename T>
struct MatlabPODWriter
{
static std::ostream& write(const T& t, std::ostream& s)
{
s << t;
return s;
}
};
template<typename T>
struct MatlabPODWriter<std::complex<T> >
{
static std::ostream& write(const std::complex<T>& t, std::ostream& s)
{
s << t.real() << " " << t.imag();
return s;
}
};
} // anonymous namespace
//! Helper method for the writeMatrixToMatlab routine.
/**
* \code
#include <dune/istl/io.hh>
* \endcode
*
* This specialization for FieldMatrices ends the recursion
*/
template <class FieldType, int rows, int cols>
void writeMatrixToMatlabHelper
( const FieldMatrix<FieldType,rows,cols>& matrix, int rowOffset,
int colOffset, std::ostream& s)
{
for (int i=0; i<rows; i++)
for (int j=0; j<cols; j++){
//+1 for Matlab numbering
s << rowOffset + i + 1 << " " << colOffset + j + 1 << " ";
MatlabPODWriter<FieldType>::write(matrix[i][j], s)<< std::endl;
}
}
//! Helper method for the writeMatrixToMatlab routine.
/**
* \code
#include <dune/istl/io.hh>
* \endcode
*/
template <class MatrixType>
void writeMatrixToMatlabHelper(const MatrixType& matrix,
int externalRowOffset, int externalColOffset,
std::ostream& s)
{
// Precompute the accumulated sizes of the columns
std::vector<typename MatrixType::size_type> colOffset(matrix.M());
if (colOffset.size() > 0)
colOffset[0] = 0;
for (typename MatrixType::size_type i=0; i<matrix.M()-1; i++)
colOffset[i+1] = colOffset[i] +
MatrixDimension<MatrixType>::coldim(matrix,i);
typename MatrixType::size_type rowOffset = 0;
// Loop over all matrix rows
for (typename MatrixType::size_type rowIdx=0; rowIdx<matrix.N(); rowIdx++)
{
const typename MatrixType::row_type& row = matrix[rowIdx];
typename MatrixType::row_type::ConstIterator cIt = row.begin();
typename MatrixType::row_type::ConstIterator cEndIt = row.end();
// Loop over all columns in this row
for (; cIt!=cEndIt; ++cIt)
writeMatrixToMatlabHelper(*cIt,
externalRowOffset+rowOffset,
externalColOffset + colOffset[cIt.index()],
s);
rowOffset += MatrixDimension<MatrixType>::rowdim(matrix, rowIdx);
}
}
//! Writes sparse matrix in a Matlab-readable format
/**
* \code
#include <dune/istl/io.hh>
* \endcode
* This routine writes the argument BCRSMatrix to a file with the name given
* by the filename argument. The file format is ASCII, with no header, and
* three data columns. Each row describes a scalar matrix entry and
* consists of the matrix row and column numbers (both counted starting from
* 1), and the matrix entry. Such a file can be read from Matlab using the
* command
* \code
new_mat = spconvert(load('filename'));
* \endcode
* @param matrix reference to matrix
* @param filename
* @param outputPrecision (number of digits) which is used to write the output file
*/
template <class MatrixType>
void writeMatrixToMatlab(const MatrixType& matrix,
const std::string& filename, int outputPrecision = 18)
{
std::ofstream outStream(filename.c_str());
int oldPrecision = outStream.precision();
outStream.precision(outputPrecision);
writeMatrixToMatlabHelper(matrix, 0, 0, outStream);
outStream.precision(oldPrecision);
}
/** @} end documentation */
} // namespace Dune
#endif
|