dune-istl  2.1.1
Classes | Modules | Files | Typedefs | Enumerations | Enumerator | Functions | Variables
Sparse Matrix and Vector classes
Iterative Solvers Template Library (ISTL)

Matrix and Vector classes that support a block recursive structure capable of representing the natural structure from Finite Element discretisations. More...

Collaboration diagram for Sparse Matrix and Vector classes:

Classes

class  Dune::BCRSMatrix< B, A >
 A sparse block matrix with compressed row storage. More...
class  Dune::BDMatrix< B, A >
 A block-diagonal matrix. More...
class  Dune::BTDMatrix< B, A >
 A block-tridiagonal matrix. More...
class  Dune::BlockVector< B, A >
 A vector of blocks with memory management. More...
class  Dune::Matrix< T, A >
 A generic dynamic dense matrix. More...
struct  Dune::MatMultMatResult< M1, M2 >
 Helper TMP to get the result type of a sparse matrix matrix multiplication ( $C=A*B$) More...
struct  Dune::MatMultMatResult< FieldMatrix< T, n, k >, FieldMatrix< T, k, m > >
struct  Dune::MatMultMatResult< BCRSMatrix< FieldMatrix< T, n, k >, A >, BCRSMatrix< FieldMatrix< T, k, m >, A1 > >
struct  Dune::CheckIfDiagonalPresent< Matrix, blocklevel, l >
 Check whether the a matrix has diagonal values on blocklevel recursion levels. More...
struct  Dune::CheckIfDiagonalPresent< Matrix, 0, l >
struct  Dune::CheckIfDiagonalPresent< MultiTypeBlockMatrix< T1, T2, T3, T4, T5, T6, T7, T8, T9 >, blocklevel, l >
class  Dune::VariableBlockVector< B, A >
 A Vector of blocks with different blocksizes. More...

Modules

 Block Recursive Iterative Kernels

Files

file  io.hh
 

Some generic functions for pretty printing vectors and matrices.


file  matrixmatrix.hh
 

provides functions for sparse matrix matrix multiplication.


file  matrixutils.hh
 

Some handy generic functions for ISTL matrices.


Typedefs

typedef BCRSMatrix
< FieldMatrix< T, n, m >, A >
::CreateIterator 
Dune::SparsityPatternInitializer< T, A, n, m >::CreateIterator
typedef BCRSMatrix
< FieldMatrix< T, n, m >, A >
::size_type 
Dune::SparsityPatternInitializer< T, A, n, m >::size_type
typedef Dune::BCRSMatrix
< FieldMatrix< T, n, m >, TA > 
Dune::MatrixInitializer< transpose, T, TA, n, m >::Matrix
typedef Matrix::CreateIterator Dune::MatrixInitializer< transpose, T, TA, n, m >::CreateIterator
typedef Matrix::size_type Dune::MatrixInitializer< transpose, T, TA, n, m >::size_type
typedef Dune::BCRSMatrix
< Dune::FieldMatrix< T, n, m >
, TA > 
Dune::MatrixInitializer< 1, T, TA, n, m >::Matrix
typedef Matrix::CreateIterator Dune::MatrixInitializer< 1, T, TA, n, m >::CreateIterator
typedef Matrix::size_type Dune::MatrixInitializer< 1, T, TA, n, m >::size_type
typedef BCRSMatrix
< FieldMatrix< T, n, m >, A
Dune::EntryAccumulatorFather< T, A, n, m >::Matrix
typedef Matrix::RowIterator Dune::EntryAccumulatorFather< T, A, n, m >::Row
typedef Matrix::ColIterator Dune::EntryAccumulatorFather< T, A, n, m >::Col
typedef BCRSMatrix
< FieldMatrix< T, n, m >, A
Dune::EntryAccumulator< T, A, n, m, transpose >::Matrix
typedef Matrix::size_type Dune::EntryAccumulator< T, A, n, m, transpose >::size_type
typedef BCRSMatrix
< FieldMatrix< T, n, m >, A
Dune::EntryAccumulator< T, A, n, m, 0 >::Matrix
typedef Matrix::size_type Dune::EntryAccumulator< T, A, n, m, 0 >::size_type
typedef BCRSMatrix
< FieldMatrix< T, n, m >, A
Dune::EntryAccumulator< T, A, n, m, 1 >::Matrix
typedef Matrix::size_type Dune::EntryAccumulator< T, A, n, m, 1 >::size_type
typedef BCRSMatrix
< FieldMatrix< T, n, m >, A
Dune::EntryAccumulator< T, A, n, m, 2 >::Matrix
typedef Matrix::size_type Dune::EntryAccumulator< T, A, n, m, 2 >::size_type
typedef FieldMatrix< T, n, m > Dune::MatMultMatResult< FieldMatrix< T, n, k >, FieldMatrix< T, k, m > >::type
typedef BCRSMatrix< typename
MatMultMatResult< FieldMatrix
< T, n, k >, FieldMatrix< T, k,
m > >::type, A
Dune::MatMultMatResult< BCRSMatrix< FieldMatrix< T, n, k >, A >, BCRSMatrix< FieldMatrix< T, k, m >, A1 > >::type

Enumerations

enum  { Dune::SparsityPatternInitializer< T, A, n, m >::do_break = true }
enum  { Dune::MatrixInitializer< transpose, T, TA, n, m >::do_break = true }
enum  { Dune::MatrixInitializer< 1, T, TA, n, m >::do_break = false }
enum  { Dune::EntryAccumulatorFather< T, A, n, m >::do_break = false }

Functions

template<class V >
void Dune::recursive_printvector (std::ostream &s, const V &v, std::string rowtext, int &counter, int columns, int width, int precision)
template<class K , int n>
void Dune::recursive_printvector (std::ostream &s, const FieldVector< K, n > &v, std::string rowtext, int &counter, int columns, int width, int precision)
template<class V >
void Dune::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
void Dune::fill_row (std::ostream &s, int m, int width, int precision)
 print a row of zeros for a non-existing block
template<class M >
void Dune::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
template<class K , int n, int m>
void Dune::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)
 print one row of a matrix, specialization for FieldMatrix
template<class K >
void Dune::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)
 print one row of a matrix, specialization for FieldMatrix<K,1,1>
template<class M >
void Dune::printmatrix (std::ostream &s, const M &A, std::string title, std::string rowtext, int width=10, int precision=2)
 Prints a generic block matrix.
template<class B , int n, int m, class A >
void Dune::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.
template<class FieldType , int rows, int cols>
void Dune::writeMatrixToMatlabHelper (const FieldMatrix< FieldType, rows, cols > &matrix, int rowOffset, int colOffset, std::ostream &s)
 Helper method for the writeMatrixToMatlab routine.
template<class MatrixType >
void Dune::writeMatrixToMatlabHelper (const MatrixType &matrix, int externalRowOffset, int externalColOffset, std::ostream &s)
 Helper method for the writeMatrixToMatlab routine.
template<class MatrixType >
void Dune::writeMatrixToMatlab (const MatrixType &matrix, const std::string &filename)
 Writes sparse matrix in a Matlab-readable format.
template<class T , class A , class A1 , class A2 , int n, int m, int k>
void Dune::matMultTransposeMat (BCRSMatrix< FieldMatrix< T, n, k >, A > &res, const BCRSMatrix< FieldMatrix< T, n, m >, A1 > &mat, const BCRSMatrix< FieldMatrix< T, k, m >, A2 > &matt, bool tryHard=false)
 Calculate product of a sparse matrix with a transposed sparse matrices ( $C=A*B^T$).
template<class T , class A , class A1 , class A2 , int n, int m, int k>
void Dune::matMultMat (BCRSMatrix< FieldMatrix< T, n, m >, A > &res, const BCRSMatrix< FieldMatrix< T, n, k >, A1 > &mat, const BCRSMatrix< FieldMatrix< T, k, m >, A2 > &matt, bool tryHard=false)
 Calculate product of two sparse matrices ( $C=A*B$).
template<class T , class A , class A1 , class A2 , int n, int m, int k>
void Dune::transposeMatMultMat (BCRSMatrix< FieldMatrix< T, n, m >, A > &res, const BCRSMatrix< FieldMatrix< T, k, n >, A1 > &mat, const BCRSMatrix< FieldMatrix< T, k, m >, A2 > &matt, bool tryHard=false)
 Calculate product of a transposed sparse matrix with another sparse matrices ( $C=A^T*B$).
template<class M >
int Dune::countNonZeros (const M &matrix)
 Get the number of nonzero fields in the matrix.

Variables

Matrix &   mat
Col   col

Detailed Description

Matrix and Vector classes that support a block recursive structure capable of representing the natural structure from Finite Element discretisations.

The interface of our matrices is designed according to what they represent from a mathematical point of view. The vector classes are representations of vector spaces:

The matrix classes represent linear maps $A: V \mapsto W$ from vector space $V$ to vector space $W$ the recursive block structure of the matrix rows and columns immediately follows from the recursive block structure of the vectors representing the domain and range of the mapping, respectively:


Typedef Documentation

template<class T , class A , int n, int m>
typedef Matrix::ColIterator Dune::EntryAccumulatorFather< T, A, n, m >::Col
template<class T , class A , int n, int m>
typedef BCRSMatrix<FieldMatrix<T,n,m>,A>::CreateIterator Dune::SparsityPatternInitializer< T, A, n, m >::CreateIterator
template<int transpose, class T , class TA , int n, int m>
typedef Matrix::CreateIterator Dune::MatrixInitializer< transpose, T, TA, n, m >::CreateIterator
template<class T , class TA , int n, int m>
typedef Matrix::CreateIterator Dune::MatrixInitializer< 1, T, TA, n, m >::CreateIterator
template<int transpose, class T , class TA , int n, int m>
typedef Dune::BCRSMatrix<FieldMatrix<T,n,m>,TA> Dune::MatrixInitializer< transpose, T, TA, n, m >::Matrix
template<class T , class TA , int n, int m>
typedef Dune::BCRSMatrix<Dune::FieldMatrix<T,n,m>,TA> Dune::MatrixInitializer< 1, T, TA, n, m >::Matrix
template<class T , class A , int n, int m>
typedef BCRSMatrix<FieldMatrix<T,n,m>,A> Dune::EntryAccumulatorFather< T, A, n, m >::Matrix
template<class T , class A , int n, int m, int transpose>
typedef BCRSMatrix<FieldMatrix<T,n,m>,A> Dune::EntryAccumulator< T, A, n, m, transpose >::Matrix
template<class T , class A , int n, int m>
typedef BCRSMatrix<FieldMatrix<T,n,m>,A> Dune::EntryAccumulator< T, A, n, m, 0 >::Matrix
template<class T , class A , int n, int m>
typedef BCRSMatrix<FieldMatrix<T,n,m>,A> Dune::EntryAccumulator< T, A, n, m, 1 >::Matrix
template<class T , class A , int n, int m>
typedef BCRSMatrix<FieldMatrix<T,n,m>,A> Dune::EntryAccumulator< T, A, n, m, 2 >::Matrix
template<class T , class A , int n, int m>
typedef Matrix::RowIterator Dune::EntryAccumulatorFather< T, A, n, m >::Row
template<class T , class A , int n, int m>
typedef BCRSMatrix<FieldMatrix<T,n,m>,A>::size_type Dune::SparsityPatternInitializer< T, A, n, m >::size_type
template<int transpose, class T , class TA , int n, int m>
typedef Matrix::size_type Dune::MatrixInitializer< transpose, T, TA, n, m >::size_type
template<class T , class TA , int n, int m>
typedef Matrix::size_type Dune::MatrixInitializer< 1, T, TA, n, m >::size_type
template<class T , class A , int n, int m, int transpose>
typedef Matrix::size_type Dune::EntryAccumulator< T, A, n, m, transpose >::size_type
template<class T , class A , int n, int m>
typedef Matrix::size_type Dune::EntryAccumulator< T, A, n, m, 0 >::size_type
template<class T , class A , int n, int m>
typedef Matrix::size_type Dune::EntryAccumulator< T, A, n, m, 1 >::size_type
template<class T , class A , int n, int m>
typedef Matrix::size_type Dune::EntryAccumulator< T, A, n, m, 2 >::size_type
template<typename T , int n, int k, int m>
typedef FieldMatrix<T,n,m> Dune::MatMultMatResult< FieldMatrix< T, n, k >, FieldMatrix< T, k, m > >::type
template<typename T , typename A , typename A1 , int n, int k, int m>
typedef BCRSMatrix<typename MatMultMatResult<FieldMatrix<T,n,k>, FieldMatrix<T,k,m> >::type,A> Dune::MatMultMatResult< BCRSMatrix< FieldMatrix< T, n, k >, A >, BCRSMatrix< FieldMatrix< T, k, m >, A1 > >::type

Enumeration Type Documentation

template<class T , class A , int n, int m>
anonymous enum
Enumerator:
do_break 
template<int transpose, class T , class TA , int n, int m>
anonymous enum
Enumerator:
do_break 
template<class T , class TA , int n, int m>
anonymous enum
Enumerator:
do_break 
template<class T , class A , int n, int m>
anonymous enum
Enumerator:
do_break 

Function Documentation

template<class M >
int Dune::countNonZeros ( const M &  matrix) [inline]

Get the number of nonzero fields in the matrix.

This is not the number of nonzero blocks, but the number of non zero scalar entries (on blocklevel 1) if the matrix is viewed as a flat matrix.

For FieldMatrix this is simply the number of columns times the number of rows, for a BCRSMatrix<FieldMatrix<K,n,m>> this is the number of nonzero blocks time n*m.

References count.

Referenced by Dune::Amg::MatrixHierarchy< M, PI, A >::build().

void Dune::fill_row ( std::ostream &  s,
int  m,
int  width,
int  precision 
) [inline]

print a row of zeros for a non-existing block

#include <dune/istl/io.hh>

Referenced by Dune::print_row().

template<class T , class A , class A1 , class A2 , int n, int m, int k>
void Dune::matMultMat ( BCRSMatrix< FieldMatrix< T, n, m >, A > &  res,
const BCRSMatrix< FieldMatrix< T, n, k >, A1 > &  mat,
const BCRSMatrix< FieldMatrix< T, k, m >, A2 > &  matt,
bool  tryHard = false 
)

Calculate product of two sparse matrices ( $C=A*B$).

Parameters:
resMatrix for the result of the computation.
matMatrix A.
mattMatrix B.
tryHardignored

References mat.

template<class T , class A , class A1 , class A2 , int n, int m, int k>
void Dune::matMultTransposeMat ( BCRSMatrix< FieldMatrix< T, n, k >, A > &  res,
const BCRSMatrix< FieldMatrix< T, n, m >, A1 > &  mat,
const BCRSMatrix< FieldMatrix< T, k, m >, A2 > &  matt,
bool  tryHard = false 
)

Calculate product of a sparse matrix with a transposed sparse matrices ( $C=A*B^T$).

Parameters:
resMatrix for the result of the computation.
matMatrix A.
mattMatrix B, which will be transposed before the multiplication.
tryHardignored

References mat.

template<class M >
void Dune::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

#include <dune/istl/io.hh>

References Dune::fill_row().

Referenced by Dune::printmatrix().

template<class K , int n, int m>
void Dune::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 
)

print one row of a matrix, specialization for FieldMatrix

#include <dune/istl/io.hh>
template<class K >
void Dune::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 
)

print one row of a matrix, specialization for FieldMatrix<K,1,1>

#include <dune/istl/io.hh>

References A.

template<class M >
void Dune::printmatrix ( std::ostream &  s,
const M &  A,
std::string  title,
std::string  rowtext,
int  width = 10,
int  precision = 2 
)

Prints a generic block matrix.

#include <dune/istl/io.hh>
Bug:
Empty rows and columns are omitted by this method. (FlySpray #7)

References Dune::print_row().

Referenced by Dune::redistributeMatrixEntries(), and test_IO().

template<class B , int n, int m, class A >
void Dune::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.

#include <dune/istl/io.hh>

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.

Parameters:
sThe ostream to print to.
matThe matrix to print.
titleThe title for the matrix.
rowtextThe text to prepend to each print out of a matrix row.
widthThe number of nonzero blocks to print in one line.
precisionThe precision to use when printing the numbers.

References A, Dune::Matrix< T, A >::begin(), col, count, Dune::Matrix< T, A >::end(), Dune::Matrix< T, A >::M(), mat, Dune::Matrix< T, A >::N(), and row.

template<class V >
void Dune::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

#include <dune/istl/io.hh>

References Dune::recursive_printvector().

Referenced by test_IO().

template<class V >
void Dune::recursive_printvector ( std::ostream &  s,
const V &  v,
std::string  rowtext,
int &  counter,
int  columns,
int  width,
int  precision 
)
#include <dune/istl/io.hh>

Referenced by Dune::printvector().

template<class K , int n>
void Dune::recursive_printvector ( std::ostream &  s,
const FieldVector< K, n > &  v,
std::string  rowtext,
int &  counter,
int  columns,
int  width,
int  precision 
)
#include <dune/istl/io.hh>
template<class T , class A , class A1 , class A2 , int n, int m, int k>
void Dune::transposeMatMultMat ( BCRSMatrix< FieldMatrix< T, n, m >, A > &  res,
const BCRSMatrix< FieldMatrix< T, k, n >, A1 > &  mat,
const BCRSMatrix< FieldMatrix< T, k, m >, A2 > &  matt,
bool  tryHard = false 
)

Calculate product of a transposed sparse matrix with another sparse matrices ( $C=A^T*B$).

Parameters:
resMatrix for the result of the computation.
matMatrix A, which will be transposed before the multiplication.
mattMatrix B.
tryHardignored

References mat.

template<class MatrixType >
void Dune::writeMatrixToMatlab ( const MatrixType &  matrix,
const std::string &  filename 
)

Writes sparse matrix in a Matlab-readable format.

#include <dune/istl/io.hh>

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

new_mat = spconvert(load('filename'));

References Dune::writeMatrixToMatlabHelper().

template<class FieldType , int rows, int cols>
void Dune::writeMatrixToMatlabHelper ( const FieldMatrix< FieldType, rows, cols > &  matrix,
int  rowOffset,
int  colOffset,
std::ostream &  s 
)

Helper method for the writeMatrixToMatlab routine.

#include <dune/istl/io.hh>

This specialization for FieldMatrices ends the recursion

Referenced by Dune::writeMatrixToMatlab(), and Dune::writeMatrixToMatlabHelper().

template<class MatrixType >
void Dune::writeMatrixToMatlabHelper ( const MatrixType &  matrix,
int  externalRowOffset,
int  externalColOffset,
std::ostream &  s 
)

Helper method for the writeMatrixToMatlab routine.

#include <dune/istl/io.hh>

References row, and Dune::writeMatrixToMatlabHelper().


Variable Documentation

Matrix& A
Col col
template<class T , class A , int n, int m>
Col Dune::EntryAccumulatorFather< T, A, n, m >::col [protected]
std::size_t count
Matrix& mat
template<class T , class A , int n, int m>
Matrix& Dune::EntryAccumulatorFather< T, A, n, m >::mat [protected]
Row row
CreateIterator rowiter