DUNE-FEM (unstable)

A block-tridiagonal matrix. More...

#include <dune/istl/btdmatrix.hh>

Public Types

using field_type = typename Imp::BlockTraits< B >::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 More...
 
enum  BuildMode
 we support two modes
 
typedef Imp::CompressedBlockVectorWindow< B, size_typerow_type
 implement row_type with compressed vector
 
typedef ::Dune::CompressionStatistics< size_typeCompressionStatistics
 The type for the statistics object returned by compress()
 
typedef RealRowIterator< row_typeiterator
 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_typeconst_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

 BTDMatrix ()
 Default constructor.
 
void setSize (size_type size)
 Resize the matrix. Invalidates the content!
 
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. More...
 
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. More...
 
void setSize (size_type rows, size_type columns, size_type nnz=0)
 Set the size of the matrix. More...
 
void setImplicitBuildModeParameters (size_type _avg, double compressionBufferSize)
 Set parameters needed for creation in implicit build mode. More...
 
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)
 
void setIndicesNoSort (size_type row, It begin, It end)
 Set all column indices for row from the given iterator range. More...
 
void setIndices (size_type row, It begin, It end)
 Set all column indices for row from the given iterator range. More...
 
B & entry (size_type row, size_type col)
 Returns reference to entry (row,col) of the matrix. More...
 
CompressionStatistics compress ()
 Finishes the buildstage in implicit mode. More...
 
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. More...
 
BCRSMatrixoperator-= (const BCRSMatrix &b)
 Subtract the entries of another matrix from this one. More...
 
BCRSMatrixaxpy (field_type alpha, const BCRSMatrix &b)
 Add the scaled entries of another matrix to this one. More...
 
void mv (const X &x, Y &y) const
 y = A x
 
void umv (const X &x, Y &y) const
 y += A x
 
void mmv (const X &x, Y &y) const
 y -= A x
 
void usmv (F &&alpha, const X &x, Y &y) const
 y += alpha A x
 
void mtv (const X &x, Y &y) const
 y = A^T x
 
void umtv (const X &x, Y &y) const
 y += A^T x
 
void mmtv (const X &x, Y &y) const
 y -= A^T x
 
void usmtv (const field_type &alpha, const X &x, Y &y) const
 y += alpha A^T x
 
void umhv (const X &x, Y &y) const
 y += A^H x
 
void mmhv (const X &x, Y &y) const
 y -= A^H x
 
void usmhv (const field_type &alpha, const X &x, Y &y) const
 y += alpha A^H x
 
FieldTraits< field_type >::real_type frobenius_norm2 () const
 square of frobenius norm, need for block recursion
 
FieldTraits< field_type >::real_type frobenius_norm () const
 frobenius norm: sqrt(sum over squared values of entries)
 
FieldTraits< ft >::real_type infinity_norm () const
 infinity norm (row sum norm, how to generalize for blocks?)
 
FieldTraits< ft >::real_type infinity_norm () const
 infinity norm (row sum norm, how to generalize for blocks?)
 
FieldTraits< ft >::real_type infinity_norm_real () const
 simplified infinity norm (uses Manhattan norm for complex values)
 
FieldTraits< ft >::real_type 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)
 
BuildStage buildStage () const
 The current build stage of the matrix.
 
BuildMode buildMode () const
 The currently selected build mode of the matrix.
 
bool exists (size_type i, size_type j) const
 return true if (i,j) is in pattern
 

Protected Member Functions

void setColumnPointers (ConstRowIterator row)
 Copy row sizes from iterator range starting at row and set column index pointers for all rows. More...
 
void setDataPointers ()
 Set data pointers for all rows. More...
 
void copyWindowStructure (const BCRSMatrix &Mat)
 Copy the window structure from another matrix.
 
void deallocate (bool deallocateRows=true)
 deallocate memory of the matrix. More...
 
void allocate (size_type rows, size_type columns, size_type allocationSize, bool allocateRows, bool allocate_data)
 Allocate memory for the matrix structure. More...
 
void implicit_allocate (size_type _n, size_type _m)
 organizes allocation implicit mode calculates correct array size to be allocated and sets the the window pointers to their correct positions for insertion. internally uses allocate() for the real allocation.
 

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

◆ size_type

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 Function Documentation

◆ allocate()

void Dune::BCRSMatrix< B, std::allocator< B > >::allocate ( size_type  rows,
size_type  columns,
size_type  allocationSize,
bool  allocateRows,
bool  allocate_data 
)
inlineprotectedinherited

Allocate memory for the matrix structure.

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 lost. Please call deallocate() before calling allocate in this case.
Parameters
rowsThe number of rows the matrix should contain.
columnsthe number of columns the matrix should contain.
allocationSizeThe number of nonzero entries the matrix should hold (if omitted defaults to 0).
allocateRowsWhether we have to allocate the row pointers, too. (Defaults to true)
allocate_dataWhether data array needs to be allocated

◆ axpy()

BCRSMatrix & Dune::BCRSMatrix< B, std::allocator< B > >::axpy ( field_type  alpha,
const BCRSMatrix< B, std::allocator< B > > &  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.

◆ beforeBegin() [1/2]

Iterator Dune::BCRSMatrix< B, std::allocator< B > >::beforeBegin ( )
inlineinherited
Returns
an iterator that is positioned before the first row of the matrix.

◆ beforeBegin() [2/2]

ConstIterator Dune::BCRSMatrix< B, std::allocator< B > >::beforeBegin ( ) const
inlineinherited
Returns
an iterator that is positioned before the first row of the matrix.

◆ beforeEnd() [1/2]

Iterator Dune::BCRSMatrix< B, std::allocator< B > >::beforeEnd ( )
inlineinherited
Returns
an iterator that is positioned before the end iterator of the rows, i.e. at the last row.

◆ beforeEnd() [2/2]

ConstIterator Dune::BCRSMatrix< B, std::allocator< B > >::beforeEnd ( ) const
inlineinherited
Returns
an iterator that is positioned before the end iterator of the rows. i.e. at the last row.

◆ compress()

CompressionStatistics Dune::BCRSMatrix< B, std::allocator< B > >::compress ( )
inlineinherited

Finishes the buildstage in implicit mode.

Performs compression of index and data arrays with linear complexity in the number of nonzeroes.

After calling this method, the matrix is in the built state and no more entries can be added.

Returns
An object with some statistics about the compression for future optimization.

◆ deallocate()

void Dune::BCRSMatrix< B, std::allocator< B > >::deallocate ( bool  deallocateRows = true)
inlineprotectedinherited

deallocate memory of the matrix.

Parameters
deallocateRowsWhether we have to deallocate the row pointers, too. If false they will not be touched. (Defaults to true).

◆ entry()

B & Dune::BCRSMatrix< B, std::allocator< B > >::entry ( size_type  row,
size_type  col 
)
inlineinherited

Returns reference to entry (row,col) of the matrix.

This method can only be used when the matrix is in implicit building mode.

A reference to entry (row, col) of the matrix is returned. If entry (row, col) is accessed for the first time, it is created on the fly.

This method can only be used while building the matrix, after compression operator[] gives a much better performance.

◆ operator+=()

BCRSMatrix & Dune::BCRSMatrix< B, std::allocator< B > >::operator+= ( const BCRSMatrix< B, std::allocator< B > > &  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.

◆ operator-=()

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

Subtract the entries of another matrix from this one.

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

◆ setBuildMode()

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

Sets the build mode of the matrix.

Parameters
bmThe build mode to use.

◆ setColumnPointers()

void Dune::BCRSMatrix< B, std::allocator< B > >::setColumnPointers ( ConstRowIterator  row)
inlineprotectedinherited

Copy row sizes from iterator range starting at row and set column index pointers for all rows.

This method does not modify the data pointers, as those are set only after building the pattern (to allow for a delayed allocation).

◆ setDataPointers()

void Dune::BCRSMatrix< B, std::allocator< B > >::setDataPointers ( )
inlineprotectedinherited

Set data pointers for all rows.

This method assumes that column pointers and row sizes have been correctly set up by a prior call to setColumnPointers().

◆ setImplicitBuildModeParameters()

void Dune::BCRSMatrix< B, std::allocator< B > >::setImplicitBuildModeParameters ( size_type  _avg,
double  compressionBufferSize 
)
inlineinherited

Set parameters needed for creation in implicit build mode.

Use this method before setSize() to define storage behaviour of a matrix in implicit build mode

Parameters
_avgexpected average number of entries per row
compressionBufferSizefraction of _n*_avg which is expected to be needed for elements that exceed _avg entries per row.

◆ setIndices()

void Dune::BCRSMatrix< B, std::allocator< B > >::setIndices ( size_type  row,
It  begin,
It  end 
)
inlineinherited

Set all column indices for row from the given iterator range.

The iterator range has to be of the same length as the previously set row size. The entries in the iterator range do not have to be in any particular order, but must not contain duplicate values. This method will insert the indices and sort them afterwards.

Calling this method overwrites any previously set column indices!

◆ setIndicesNoSort()

void Dune::BCRSMatrix< B, std::allocator< B > >::setIndicesNoSort ( size_type  row,
It  begin,
It  end 
)
inlineinherited

Set all column indices for row from the given iterator range.

The iterator range has to be of the same length as the previously set row size. The entries in the iterator range must be sorted and must not contain duplicate values. This method will insert the indices in the given order.

Calling this method overwrites any previously set column indices!

◆ setSize()

void Dune::BCRSMatrix< B, std::allocator< B > >::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). Must be omitted in implicit mode.

◆ solve()

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, std::allocator< B > >::N().


The documentation for this class was generated from the following file:
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 20, 23:30, 2024)