dune-istl  2.2.1
Public Types | Public Member Functions | List of all members
Dune::BTDMatrix< B, A > Class Template Reference

A block-tridiagonal matrix. More...

#include <dune/istl/btdmatrix.hh>

Inheritance diagram for Dune::BTDMatrix< B, A >:
Dune::BCRSMatrix< B, A >

Public Types

enum  { blocklevel = B::blocklevel+1 }
 increment block level counter 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 A::size_type size_type
 implement row_type with compressed vector
enum  { blocklevel = B::blocklevel+1 }
 increment block level counter More...
enum  BuildMode { row_wise, random, unknown }
 we support two modes More...
typedef
CompressedBlockVectorWindow< B,
A
row_type
 implement row_type with compressed vector
typedef RealRowIterator< row_typeiterator
 The iterator over the (mutable matrix rows.
typedef RealRowIterator< row_typeIterator
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 RealRowIterator< const
row_type
ConstIterator
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

 BTDMatrix ()
 Default constructor.
 BTDMatrix (int size)
BTDMatrixoperator= (const BTDMatrix &other)
 assignment
BTDMatrixoperator= (const field_type &k)
 assignment from scalar
template<class V >
void solve (V &x, const V &rhs) const
 Use the Thomas algorithm to solve the system Ax=b in O(n) time.
row_typeoperator[] (size_type i)
 random access to the rows
const row_typeoperator[] (size_type i) const
 same for read only access
Iterator begin ()
 Get iterator to first row.
ConstIterator begin () const
 Get const iterator to first row.
Iterator end ()
 Get iterator to one beyond last row.
ConstIterator end () const
 Get const iterator to one beyond last row.
Iterator beforeEnd ()
ConstIterator beforeEnd () const
Iterator beforeBegin ()
ConstIterator beforeBegin () const
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.
CreateIterator createbegin ()
 get initial create iterator
CreateIterator createend ()
 get create iterator pointing to one after the last block
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)
BCRSMatrixoperator*= (const field_type &k)
 vector space multiplication with scalar
BCRSMatrixoperator/= (const field_type &k)
 vector space division by scalar
BCRSMatrixoperator+= (const BCRSMatrix &b)
 Add the entries of another matrix to this one.
BCRSMatrixoperator-= (const BCRSMatrix &b)
 Substract the entries of another matrix to this one.
BCRSMatrixaxpy (field_type alpha, const BCRSMatrix &b)
 Add the scaled 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 mtv (const X &x, Y &y) const
 y = A^T 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)
bool exists (size_type i, size_type j) const
 return true if (i,j) is in pattern

Detailed Description

template<class B, class A = std::allocator<B>>
class Dune::BTDMatrix< B, A >

A block-tridiagonal matrix.

Todo:
It would be safer and more efficient to have a real implementation of a block-tridiagonal matrix and not just subclassing from BCRSMatrix. But that's quite a lot of work for that little advantage.

Member Typedef Documentation

template<class B , class A = std::allocator<B>>
typedef A Dune::BTDMatrix< B, A >::allocator_type

export the allocator type

template<class B , class A = std::allocator<B>>
typedef B Dune::BTDMatrix< B, A >::block_type

export the type representing the components

template<class B, class A = std::allocator<B>>
typedef row_type::Iterator Dune::BCRSMatrix< B, A >::ColIterator
inherited

Iterator for the entries of each row.

template<class B, class A = std::allocator<B>>
typedef RealRowIterator<const row_type> Dune::BCRSMatrix< B, A >::const_iterator
inherited

The const iterator over the matrix rows.

template<class B, class A = std::allocator<B>>
typedef row_type::ConstIterator Dune::BCRSMatrix< B, A >::ConstColIterator
inherited

Const iterator to the entries of a row.

template<class B, class A = std::allocator<B>>
typedef RealRowIterator<const row_type> Dune::BCRSMatrix< B, A >::ConstIterator
inherited
template<class B, class A = std::allocator<B>>
typedef ConstIterator Dune::BCRSMatrix< B, A >::ConstRowIterator
inherited

rename the const row iterator for easier access

template<class B , class A = std::allocator<B>>
typedef B::field_type Dune::BTDMatrix< B, A >::field_type

export the type representing the field

template<class B, class A = std::allocator<B>>
typedef RealRowIterator<row_type> Dune::BCRSMatrix< B, A >::iterator
inherited

The iterator over the (mutable matrix rows.

template<class B, class A = std::allocator<B>>
typedef RealRowIterator<row_type> Dune::BCRSMatrix< B, A >::Iterator
inherited
template<class B, class A = std::allocator<B>>
typedef CompressedBlockVectorWindow<B,A> Dune::BCRSMatrix< B, A >::row_type
inherited

implement row_type with compressed vector

template<class B, class A = std::allocator<B>>
typedef Iterator Dune::BCRSMatrix< B, A >::RowIterator
inherited

rename the iterators for easier access

template<class B , class A = std::allocator<B>>
typedef A::size_type Dune::BTDMatrix< B, A >::size_type

implement row_type with compressed vector

The type for the index access and the size

Member Enumeration Documentation

template<class B, class A = std::allocator<B>>
anonymous enum
inherited

increment block level counter

Enumerator:
blocklevel 

The number of blocklevels the matrix contains.

template<class B , class A = std::allocator<B>>
anonymous enum

increment block level counter

Enumerator:
blocklevel 
template<class B, class A = std::allocator<B>>
enum Dune::BCRSMatrix::BuildMode
inherited

we support two modes

Enumerator:
row_wise 

Build in a row-wise manner.

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)).

random 

Build entries randomly.

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.

unknown 

Build mode not set!

Constructor & Destructor Documentation

template<class B , class A = std::allocator<B>>
Dune::BTDMatrix< B, A >::BTDMatrix ( )
inline

Default constructor.

template<class B , class A = std::allocator<B>>
Dune::BTDMatrix< B, A >::BTDMatrix ( int  size)
inlineexplicit

Member Function Documentation

template<class B, class A = std::allocator<B>>
BCRSMatrix& Dune::BCRSMatrix< B, A >::axpy ( field_type  alpha,
const BCRSMatrix< B, A > &  b 
)
inlineinherited

Add the scaled entries of another matrix to this one.

Matrix axpy operation: *this += alpha * b

Parameters
alphaScaling factor.
bThe matrix to add to this one. Its sparsity pattern has to be subset of the sparsity pattern of this matrix.

References Dune::BCRSMatrix< B, A >::begin(), Dune::BCRSMatrix< B, A >::end(), Dune::BCRSMatrix< B, A >::M(), and Dune::BCRSMatrix< B, A >::N().

template<class B, class A = std::allocator<B>>
Iterator Dune::BCRSMatrix< B, A >::beforeBegin ( )
inlineinherited
Returns
an iterator that is positioned before the first row of the matrix.
template<class B, class A = std::allocator<B>>
ConstIterator Dune::BCRSMatrix< B, A >::beforeBegin ( ) const
inlineinherited
Returns
an iterator that is positioned before the first row of the matrix.
template<class B, class A = std::allocator<B>>
Iterator Dune::BCRSMatrix< B, A >::beforeEnd ( )
inlineinherited
Returns
an iterator that is positioned before the end iterator of the rows, i.e. at the last row.
template<class B, class A = std::allocator<B>>
ConstIterator Dune::BCRSMatrix< B, A >::beforeEnd ( ) const
inlineinherited
Returns
an iterator that is positioned before the end iterator of the rows. i.e. at the last row.
template<class B, class A = std::allocator<B>>
Iterator Dune::BCRSMatrix< B, A >::begin ( )
inlineinherited
template<class B, class A = std::allocator<B>>
ConstIterator Dune::BCRSMatrix< B, A >::begin ( ) const
inlineinherited

Get const iterator to first row.

template<class B, class A = std::allocator<B>>
CreateIterator Dune::BCRSMatrix< B, A >::createbegin ( )
inlineinherited

get initial create iterator

References Dune::BCRSMatrix< B, A >::CreateIterator.

Referenced by test_IO(), and test_Iter().

template<class B, class A = std::allocator<B>>
CreateIterator Dune::BCRSMatrix< B, A >::createend ( )
inlineinherited

get create iterator pointing to one after the last block

References Dune::BCRSMatrix< B, A >::CreateIterator.

Referenced by test_IO(), and test_Iter().

template<class B, class A = std::allocator<B>>
Iterator Dune::BCRSMatrix< B, A >::end ( )
inlineinherited
template<class B, class A = std::allocator<B>>
ConstIterator Dune::BCRSMatrix< B, A >::end ( ) const
inlineinherited

Get const iterator to one beyond last row.

template<class B, class A = std::allocator<B>>
bool Dune::BCRSMatrix< B, A >::exists ( size_type  i,
size_type  j 
) const
inlineinherited

return true if (i,j) is in pattern

References Dune::BCRSMatrix< B, A >::end().

template<class B, class A = std::allocator<B>>
double Dune::BCRSMatrix< B, A >::frobenius_norm ( ) const
inlineinherited

frobenius norm: sqrt(sum over squared values of entries)

References Dune::BCRSMatrix< B, A >::frobenius_norm2().

template<class B, class A = std::allocator<B>>
double Dune::BCRSMatrix< B, A >::frobenius_norm2 ( ) const
inlineinherited

square of frobenius norm, need for block recursion

References Dune::BCRSMatrix< B, A >::begin(), and Dune::BCRSMatrix< B, A >::end().

Referenced by Dune::BCRSMatrix< B, A >::frobenius_norm().

template<class B, class A = std::allocator<B>>
size_type Dune::BCRSMatrix< B, A >::getrowsize ( size_type  i) const
inlineinherited

get current number of indices in row i

References Dune::CompressedBlockVectorWindow< B, A >::getsize().

template<class B, class A = std::allocator<B>>
void Dune::BCRSMatrix< B, A >::incrementrowsize ( size_type  i,
size_type  s = 1 
)
inlineinherited

increment size of row i by s (1 by default)

References Dune::BCRSMatrix< B, A >::random, and Dune::CompressedBlockVectorWindow< B, A >::setsize().

template<class B, class A = std::allocator<B>>
double Dune::BCRSMatrix< B, A >::infinity_norm ( ) const
inlineinherited

infinity norm (row sum norm, how to generalize for blocks?)

References Dune::BCRSMatrix< B, A >::begin(), and Dune::BCRSMatrix< B, A >::end().

template<class B, class A = std::allocator<B>>
double Dune::BCRSMatrix< B, A >::infinity_norm_real ( ) const
inlineinherited

simplified infinity norm (uses Manhattan norm for complex values)

References Dune::BCRSMatrix< B, A >::begin(), and Dune::BCRSMatrix< B, A >::end().

template<class B, class A = std::allocator<B>>
size_type Dune::BCRSMatrix< B, A >::M ( ) const
inlineinherited
template<class B, class A = std::allocator<B>>
template<class X , class Y >
void Dune::BCRSMatrix< B, A >::mmhv ( const X &  x,
Y &  y 
) const
inlineinherited
template<class B, class A = std::allocator<B>>
template<class X , class Y >
void Dune::BCRSMatrix< B, A >::mmtv ( const X &  x,
Y &  y 
) const
inlineinherited
template<class B, class A = std::allocator<B>>
template<class X , class Y >
void Dune::BCRSMatrix< B, A >::mmv ( const X &  x,
Y &  y 
) const
inlineinherited
template<class B, class A = std::allocator<B>>
template<class X , class Y >
void Dune::BCRSMatrix< B, A >::mtv ( const X &  x,
Y &  y 
) const
inlineinherited
template<class B, class A = std::allocator<B>>
template<class X , class Y >
void Dune::BCRSMatrix< B, A >::mv ( const X &  x,
Y &  y 
) const
inlineinherited
template<class B, class A = std::allocator<B>>
size_type Dune::BCRSMatrix< B, A >::N ( ) const
inlineinherited
template<class B, class A = std::allocator<B>>
size_type Dune::BCRSMatrix< B, A >::nonzeroes ( ) const
inlineinherited

number of blocks that are stored (the number of blocks that possibly are nonzero)

template<class B, class A = std::allocator<B>>
BCRSMatrix& Dune::BCRSMatrix< B, A >::operator*= ( const field_type k)
inlineinherited

vector space multiplication with scalar

References Dune::BCRSMatrix< B, A >::begin(), and Dune::BCRSMatrix< B, A >::end().

template<class B, class A = std::allocator<B>>
BCRSMatrix& Dune::BCRSMatrix< B, A >::operator+= ( const BCRSMatrix< B, A > &  b)
inlineinherited

Add the entries of another matrix to this one.

Parameters
bThe matrix to add to this one. Its sparsity pattern has to be subset of the sparsity pattern of this matrix.

References Dune::BCRSMatrix< B, A >::begin(), Dune::BCRSMatrix< B, A >::end(), Dune::BCRSMatrix< B, A >::M(), and Dune::BCRSMatrix< B, A >::N().

template<class B, class A = std::allocator<B>>
BCRSMatrix& Dune::BCRSMatrix< B, A >::operator-= ( const BCRSMatrix< B, A > &  b)
inlineinherited

Substract the entries of another matrix to this one.

Parameters
bThe matrix to add to this one. Its sparsity pattern has to be subset of the sparsity pattern of this matrix.

References Dune::BCRSMatrix< B, A >::begin(), Dune::BCRSMatrix< B, A >::end(), Dune::BCRSMatrix< B, A >::M(), and Dune::BCRSMatrix< B, A >::N().

template<class B, class A = std::allocator<B>>
BCRSMatrix& Dune::BCRSMatrix< B, A >::operator/= ( const field_type k)
inlineinherited

vector space division by scalar

References Dune::BCRSMatrix< B, A >::begin(), and Dune::BCRSMatrix< B, A >::end().

template<class B , class A = std::allocator<B>>
BTDMatrix& Dune::BTDMatrix< B, A >::operator= ( const BTDMatrix< B, A > &  other)
inline

assignment

Referenced by Dune::BTDMatrix< B, A >::operator=().

template<class B , class A = std::allocator<B>>
BTDMatrix& Dune::BTDMatrix< B, A >::operator= ( const field_type k)
inline

assignment from scalar

Reimplemented from Dune::BCRSMatrix< B, A >.

References Dune::BTDMatrix< B, A >::operator=().

template<class B, class A = std::allocator<B>>
row_type& Dune::BCRSMatrix< B, A >::operator[] ( size_type  i)
inlineinherited

random access to the rows

template<class B, class A = std::allocator<B>>
const row_type& Dune::BCRSMatrix< B, A >::operator[] ( size_type  i) const
inlineinherited

same for read only access

template<class B, class A = std::allocator<B>>
void Dune::BCRSMatrix< B, A >::setBuildMode ( BuildMode  bm)
inlineinherited

Sets the build mode of the matrix.

Parameters
bmThe build mode to use.
template<class B, class A = std::allocator<B>>
void Dune::BCRSMatrix< B, A >::setSize ( size_type  rows,
size_type  columns,
size_type  nnz = 0 
)
inlineinherited

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
rowsThe number of rows the matrix should contain.
columnsthe number of columns the matrix should contain.
nnzThe number of nonzero entries the matrix should hold (if omitted defaults to 0).
template<class B , class A = std::allocator<B>>
template<class V >
void Dune::BTDMatrix< B, A >::solve ( V &  x,
const V &  rhs 
) const
inline

Use the Thomas algorithm to solve the system Ax=b in O(n) time.

Exceptions
ISTLErrorif the matrix is singular

References Dune::BCRSMatrix< B, A >::N().

template<class B, class A = std::allocator<B>>
template<class X , class Y >
void Dune::BCRSMatrix< B, A >::umhv ( const X &  x,
Y &  y 
) const
inlineinherited
template<class B, class A = std::allocator<B>>
template<class X , class Y >
void Dune::BCRSMatrix< B, A >::umtv ( const X &  x,
Y &  y 
) const
inlineinherited
template<class B, class A = std::allocator<B>>
template<class X , class Y >
void Dune::BCRSMatrix< B, A >::umv ( const X &  x,
Y &  y 
) const
inlineinherited
template<class B, class A = std::allocator<B>>
template<class X , class Y >
void Dune::BCRSMatrix< B, A >::usmhv ( const field_type alpha,
const X &  x,
Y &  y 
) const
inlineinherited
template<class B, class A = std::allocator<B>>
template<class X , class Y >
void Dune::BCRSMatrix< B, A >::usmtv ( const field_type alpha,
const X &  x,
Y &  y 
) const
inlineinherited
template<class B, class A = std::allocator<B>>
template<class X , class Y >
void Dune::BCRSMatrix< B, A >::usmv ( const field_type alpha,
const X &  x,
Y &  y 
) const
inlineinherited

The documentation for this class was generated from the following file: