Dune::BCRSMatrix< B, A > Class Template Reference
[Sparse Matrix and Vector classes]
#include <bcrsmatrix.hh>
Inheritance diagram for Dune::BCRSMatrix< B, A >:

Detailed Description
template<class B, class A = ISTLAllocator>
class Dune::BCRSMatrix< B, A >
A sparse block matrix with compressed row storage.
Implements a block compressed row storage scheme. The block type B can be any type implementing the matrix interface.
Different ways to build up a compressed row storage matrix are supported:
1. Row-wise scheme 2. Random scheme
Error checking: no error checking is provided normally. Setting the compile time switch DUNE_ISTL_WITH_CHECKING enables error checking.
Details:
1. Row-wise scheme
Rows are built up in sequential order. Size of the row and the column indices are defined. A row can be used as soon as it is initialized. With respect to memory there are two variants of this scheme: (a) number of non-zeroes known in advance (application finite difference schemes), (b) number of non-zeroes not known in advance (application: Sparse LU, ILU(n)).
2. Random scheme
For general finite element implementations the number of rows n is known, the number of non-zeroes might also be known (e.g. #edges + #nodes for P1) but the size of a row and the indices of a row can not be defined in sequential order.
#include<dune/common/fmatrix.hh> #include<dune/istl/bcrsmatrix.hh> ... typedef FieldMatrix<double,2,2> M; BCRSMatrix<M> B(4,4,BCRSMatrix<M>::random); // initially set row size for each row B.setrowsize(0,1); B.setrowsize(3,4); B.setrowsize(2,1); B.setrowsize(1,1); // increase row size for row 2 B.incrementrowsize(2) // finalize row setup phase B.endrowsizes(); // add column entries to rows B.addindex(0,0); B.addindex(3,1); B.addindex(2,2); B.addindex(1,1); B.addindex(2,0); B.addindex(3,2); B.addindex(3,0); B.addindex(3,3); // finalize column setup phase B.endindices(); // set entries using the random access operator B[0][0] = 1; B[1][1] = 2; B[2][0] = 3; B[2][2] = 4; B[3][1] = 5; B[3][2] = 6; B[3][0] = 7; B[3][3] = 8;
Public Types | |
enum | { blocklevel = B::blocklevel+1 } |
increment block level counter More... | |
enum | BuildMode { row_wise, random, unknown } |
we support two modes More... | |
typedef B::field_type | field_type |
export the type representing the field | |
typedef B | block_type |
export the type representing the components | |
typedef A | allocator_type |
export the allocator type | |
typedef CompressedBlockVectorWindow< B, A > | row_type |
implement row_type with compressed vector | |
typedef A::size_type | size_type |
The type for the index access and the size. | |
typedef RealRowIterator< row_type > | iterator |
The iterator over the (mutable matrix rows. | |
typedef Iterator | RowIterator |
rename the iterators for easier access | |
typedef row_type::Iterator | ColIterator |
Iterator for the entries of each row. | |
typedef RealRowIterator< const row_type > | const_iterator |
The const iterator over the matrix rows. | |
typedef ConstIterator | ConstRowIterator |
rename the const row iterator for easier access | |
typedef row_type::ConstIterator | ConstColIterator |
Const iterator to the entries of a row. | |
Public Member Functions | |
row_type & | operator[] (size_type i) |
random access to the rows | |
const row_type & | operator[] (size_type i) const |
same for read only access | |
Iterator | begin () |
Get iterator to first row. | |
Iterator | end () |
Get iterator to one beyond last row. | |
Iterator | rbegin () |
Get iterator to last row. | |
Iterator | rend () |
Get iterator to one before first row. | |
ConstIterator | begin () const |
Get const iterator to first row. | |
ConstIterator | end () const |
Get const iterator to one beyond last row. | |
ConstIterator | rbegin () const |
Get const iterator to last row. | |
ConstIterator | rend () const |
Get const iterator to one before first row. | |
BCRSMatrix () | |
an empty matrix | |
BCRSMatrix (size_type _n, size_type _m, size_type _nnz, BuildMode bm) | |
matrix with known number of nonzeroes | |
BCRSMatrix (size_type _n, size_type _m, BuildMode bm) | |
matrix with unknown number of nonzeroes | |
BCRSMatrix (const BCRSMatrix &Mat) | |
copy constructor | |
~BCRSMatrix () | |
destructor | |
void | setBuildMode (BuildMode bm) |
Sets the build mode of the matrix. | |
void | setSize (size_type rows, size_type columns, size_type nnz=0) |
Set the size of the matrix. | |
BCRSMatrix & | operator= (const BCRSMatrix &Mat) |
assignment | |
BCRSMatrix & | operator= (const field_type &k) |
Assignment from a scalar. | |
CreateIterator | createbegin () |
get initial create iterator | |
CreateIterator | createend () |
get create iterator pointing to one after the last block | |
void | setrowsize (size_type i, size_type s) |
set number of indices in row i to s | |
size_type | getrowsize (size_type i) const |
get current number of indices in row i | |
void | incrementrowsize (size_type i, size_type s=1) |
increment size of row i by s (1 by default) | |
void | endrowsizes () |
indicate that size of all rows is defined | |
void | addindex (size_type row, size_type col) |
add index (row,col) to the matrix | |
void | endindices () |
indicate that all indices are defined, check consistency | |
BCRSMatrix & | operator *= (const field_type &k) |
vector space multiplication with scalar | |
BCRSMatrix & | operator/= (const field_type &k) |
vector space division by scalar | |
BCRSMatrix & | operator+= (const BCRSMatrix &b) |
Add the entries of another matrix to this one. | |
BCRSMatrix & | operator-= (const BCRSMatrix &b) |
Substract the entries of another matrix to this one. | |
template<class X, class Y> | |
void | mv (const X &x, Y &y) const |
y = A x | |
template<class X, class Y> | |
void | umv (const X &x, Y &y) const |
y += A x | |
template<class X, class Y> | |
void | mmv (const X &x, Y &y) const |
y -= A x | |
template<class X, class Y> | |
void | usmv (const field_type &alpha, const X &x, Y &y) const |
y += alpha A x | |
template<class X, class Y> | |
void | umtv (const X &x, Y &y) const |
y += A^T x | |
template<class X, class Y> | |
void | mmtv (const X &x, Y &y) const |
y -= A^T x | |
template<class X, class Y> | |
void | usmtv (const field_type &alpha, const X &x, Y &y) const |
y += alpha A^T x | |
template<class X, class Y> | |
void | umhv (const X &x, Y &y) const |
y += A^H x | |
template<class X, class Y> | |
void | mmhv (const X &x, Y &y) const |
y -= A^H x | |
template<class X, class Y> | |
void | usmhv (const field_type &alpha, const X &x, Y &y) const |
y += alpha A^H x | |
double | frobenius_norm2 () const |
square of frobenius norm, need for block recursion | |
double | frobenius_norm () const |
frobenius norm: sqrt(sum over squared values of entries) | |
double | infinity_norm () const |
infinity norm (row sum norm, how to generalize for blocks?) | |
double | infinity_norm_real () const |
simplified infinity norm (uses Manhattan norm for complex values) | |
size_type | N () const |
number of rows (counted in blocks) | |
size_type | M () const |
number of columns (counted in blocks) | |
size_type | nonzeroes () const |
number of blocks that are stored (the number of blocks that possibly are nonzero) | |
size_type | rowdim (size_type i) const |
row dimension of block r | |
size_type | coldim (size_type c) const |
Column dimension of block c. | |
size_type | rowdim () const |
dimension of the destination vector space | |
size_type | coldim () const |
dimension of the source vector space | |
bool | exists (size_type i, size_type j) const |
return true if (i,j) is in pattern | |
Classes | |
class | CreateIterator |
Iterator class for sequential creation of blocks More... | |
class | RealRowIterator |
Iterator access to matrix rows More... |
Member Enumeration Documentation
anonymous enum |
enum Dune::BCRSMatrix::BuildMode |
we support two modes
- Enumerator:
Member Function Documentation
void Dune::BCRSMatrix< B, A >::setBuildMode | ( | BuildMode | bm | ) | [inline] |
Sets the build mode of the matrix.
- Parameters:
-
bm The build mode to use.
void Dune::BCRSMatrix< B, A >::setSize | ( | size_type | rows, | |
size_type | columns, | |||
size_type | nnz = 0 | |||
) | [inline] |
Set the size of the matrix.
Sets the number of rows and columns of the matrix and allocates the memory needed for the storage of the matrix entries.
- Warning:
- After calling this methods on an already allocated (and probably setup matrix) results in all the structure and data being deleted. I.~e. one has to setup the matrix again.
- Parameters:
-
rows The number of rows the matrix should contain. columns the number of columns the matrix should contain. nnz The number of nonzero entries the matrix should hold (if omitted defaults to 0).
void Dune::BCRSMatrix< B, A >::addindex | ( | size_type | row, | |
size_type | col | |||
) | [inline] |
add index (row,col) to the matrix
This method can only be used when building the BCRSMatrix in random mode.
addindex adds a new column entry to the row. If this column entry already exists, nothing is done.
Don't call addindex after the setup phase is finished (after endindices is called).
BCRSMatrix& Dune::BCRSMatrix< B, A >::operator+= | ( | const BCRSMatrix< B, A > & | b | ) | [inline] |
Add the entries of another matrix to this one.
- Parameters:
-
b The matrix to add to this one. Its sparsity pattern has to be subset of the sparsity pattern of this matrix.
BCRSMatrix& Dune::BCRSMatrix< B, A >::operator-= | ( | const BCRSMatrix< B, A > & | b | ) | [inline] |
Substract the entries of another matrix to this one.
- Parameters:
-
b The matrix to add to this one. Its sparsity pattern has to be subset of the sparsity pattern of this matrix.
size_type Dune::BCRSMatrix< B, A >::rowdim | ( | size_type | i | ) | const [inline] |
size_type Dune::BCRSMatrix< B, A >::coldim | ( | size_type | c | ) | const [inline] |
The documentation for this class was generated from the following file: