dune-istl
2.4.1
|
A block-tridiagonal matrix. More...
#include <dune/istl/btdmatrix.hh>
Public Types | |
enum | { blocklevel = B::blocklevel+1 } |
increment block level counter More... | |
typedef B::field_type | field_type |
export the type representing the field More... | |
typedef B | block_type |
export the type representing the components More... | |
typedef A | allocator_type |
export the allocator type More... | |
typedef A::size_type | size_type |
implement row_type with compressed vector More... | |
enum | BuildStage { notbuilt =0, notAllocated =0, building =1, rowSizesBuilt =2, built =3 } |
enum | { blocklevel = B::blocklevel+1 } |
increment block level counter More... | |
enum | BuildMode { row_wise, random, implicit, unknown } |
we support two modes More... | |
typedef CompressedBlockVectorWindow< B, A > | row_type |
implement row_type with compressed vector More... | |
typedef ::Dune::CompressionStatistics < size_type > | CompressionStatistics |
The type for the statistics object returned by compress() More... | |
typedef RealRowIterator< row_type > | iterator |
The iterator over the (mutable matrix rows. More... | |
typedef RealRowIterator< row_type > | Iterator |
typedef Iterator | RowIterator |
rename the iterators for easier access More... | |
typedef row_type::Iterator | ColIterator |
Iterator for the entries of each row. More... | |
typedef RealRowIterator< const row_type > | const_iterator |
The const iterator over the matrix rows. More... | |
typedef RealRowIterator< const row_type > | ConstIterator |
typedef ConstIterator | ConstRowIterator |
rename the const row iterator for easier access More... | |
typedef row_type::ConstIterator | ConstColIterator |
Const iterator to the entries of a row. More... | |
Public Member Functions | |
BTDMatrix () | |
Default constructor. More... | |
BTDMatrix (int size) | |
BTDMatrix & | operator= (const BTDMatrix &other) |
assignment More... | |
BTDMatrix & | operator= (const field_type &k) |
assignment from scalar More... | |
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_type & | operator[] (size_type i) |
random access to the rows More... | |
const row_type & | operator[] (size_type i) const |
same for read only access More... | |
Iterator | begin () |
Get iterator to first row. More... | |
ConstIterator | begin () const |
Get const iterator to first row. More... | |
Iterator | end () |
Get iterator to one beyond last row. More... | |
ConstIterator | end () const |
Get const iterator to one beyond last row. More... | |
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 _overflow) |
Set parameters needed for creation in implicit build mode. More... | |
CreateIterator | createbegin () |
get initial create iterator More... | |
CreateIterator | createend () |
get create iterator pointing to one after the last block More... | |
size_type | getrowsize (size_type i) const |
get current number of indices in row i More... | |
void | incrementrowsize (size_type i, size_type s=1) |
increment size of row i by s (1 by default) More... | |
template<typename It > | |
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... | |
BCRSMatrix & | operator*= (const field_type &k) |
vector space multiplication with scalar More... | |
BCRSMatrix & | operator/= (const field_type &k) |
vector space division by scalar More... | |
BCRSMatrix & | operator+= (const BCRSMatrix &b) |
Add the entries of another matrix to this one. More... | |
BCRSMatrix & | operator-= (const BCRSMatrix &b) |
Substract the entries of another matrix to this one. More... | |
BCRSMatrix & | axpy (field_type alpha, const BCRSMatrix &b) |
Add the scaled entries of another matrix to this one. More... | |
template<class X , class Y > | |
void | mv (const X &x, Y &y) const |
y = A x More... | |
template<class X , class Y > | |
void | umv (const X &x, Y &y) const |
y += A x More... | |
template<class X , class Y > | |
void | mmv (const X &x, Y &y) const |
y -= A x More... | |
template<typename F , class X , class Y > | |
void | usmv (F &&alpha, const X &x, Y &y) const |
y += alpha A x More... | |
template<class X , class Y > | |
void | mtv (const X &x, Y &y) const |
y = A^T x More... | |
template<class X , class Y > | |
void | umtv (const X &x, Y &y) const |
y += A^T x More... | |
template<class X , class Y > | |
void | mmtv (const X &x, Y &y) const |
y -= A^T x More... | |
template<class X , class Y > | |
void | usmtv (const field_type &alpha, const X &x, Y &y) const |
y += alpha A^T x More... | |
template<class X , class Y > | |
void | umhv (const X &x, Y &y) const |
y += A^H x More... | |
template<class X , class Y > | |
void | mmhv (const X &x, Y &y) const |
y -= A^H x More... | |
template<class X , class Y > | |
void | usmhv (const field_type &alpha, const X &x, Y &y) const |
y += alpha A^H x More... | |
FieldTraits< field_type > ::real_type | frobenius_norm2 () const |
square of frobenius norm, need for block recursion More... | |
FieldTraits< field_type > ::real_type | frobenius_norm () const |
frobenius norm: sqrt(sum over squared values of entries) More... | |
FieldTraits< field_type > ::real_type | infinity_norm () const |
infinity norm (row sum norm, how to generalize for blocks?) More... | |
FieldTraits< field_type > ::real_type | infinity_norm_real () const |
simplified infinity norm (uses Manhattan norm for complex values) More... | |
size_type | N () const |
number of rows (counted in blocks) More... | |
size_type | M () const |
number of columns (counted in blocks) More... | |
size_type | nonzeroes () const |
number of blocks that are stored (the number of blocks that possibly are nonzero) More... | |
BuildStage | buildStage () const |
The current build stage of the matrix. More... | |
BuildMode | buildMode () const |
The currently selected build mode of the matrix. More... | |
bool | exists (size_type i, size_type j) const |
return true if (i,j) is in pattern More... | |
Protected Types | |
typedef std::map< std::pair < size_type, size_type >, B > | OverflowType |
Protected Member Functions | |
void | setWindowPointers (ConstRowIterator row) |
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. More... | |
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 | allocateData () |
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. More... | |
Protected Attributes | |
BuildMode | build_mode |
BuildStage | ready |
A::template rebind< B >::other | allocator_ |
A::template rebind< row_type > ::other | rowAllocator_ |
A::template rebind< size_type > ::other | sizeAllocator_ |
size_type | n |
size_type | m |
size_type | nnz_ |
size_type | allocationSize |
row_type * | r |
B * | a |
std::shared_ptr< size_type > | j_ |
size_type | avg |
double | overflowsize |
OverflowType | overflow |
A block-tridiagonal matrix.
typedef A Dune::BTDMatrix< B, A >::allocator_type |
export the allocator type
typedef B Dune::BTDMatrix< B, A >::block_type |
export the type representing the components
|
inherited |
Iterator for the entries of each row.
|
inherited |
The type for the statistics object returned by compress()
|
inherited |
The const iterator over the matrix rows.
|
inherited |
Const iterator to the entries of a row.
|
inherited |
|
inherited |
rename the const row iterator for easier access
typedef B::field_type Dune::BTDMatrix< B, A >::field_type |
export the type representing the field
|
inherited |
The iterator over the (mutable matrix rows.
|
inherited |
|
protectedinherited |
|
inherited |
implement row_type with compressed vector
|
inherited |
rename the iterators for easier access
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
|
inline |
Default constructor.
|
inlineexplicit |
|
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.
row | The number of rows the matrix should contain. |
columns | the number of columns the matrix should contain. |
allocationSize_ | The number of nonzero entries the matrix should hold (if omitted defaults to 0). |
allocateRow | Whether we have to allocate the row pointers, too. (Defaults to true) |
References Dune::BCRSMatrix< B, A >::allocateData(), and Dune::BCRSMatrix< B, A >::building.
Referenced by Dune::BCRSMatrix< B, A >::BCRSMatrix(), Dune::BCRSMatrix< B, A >::endrowsizes(), Dune::BCRSMatrix< B, A >::implicit_allocate(), Dune::BCRSMatrix< B, A >::operator=(), and Dune::BCRSMatrix< B, A >::setSize().
|
inlineprotectedinherited |
|
inlineinherited |
Add the scaled entries of another matrix to this one.
Matrix axpy operation: *this += alpha * b
alpha | Scaling factor. |
b | The 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 >::built, Dune::BCRSMatrix< B, A >::end(), Dune::BCRSMatrix< B, A >::M(), Dune::BCRSMatrix< B, A >::N(), and Dune::BCRSMatrix< B, A >::ready.
|
inlineinherited |
References Dune::BCRSMatrix< B, A >::r.
|
inlineinherited |
References Dune::BCRSMatrix< B, A >::r.
|
inlineinherited |
References Dune::BCRSMatrix< B, A >::n, and Dune::BCRSMatrix< B, A >::r.
|
inlineinherited |
References Dune::BCRSMatrix< B, A >::n, and Dune::BCRSMatrix< B, A >::r.
|
inlineinherited |
Get iterator to first row.
References Dune::BCRSMatrix< B, A >::r.
Referenced by Dune::SeqOverlappingSchwarzAssemblerHelper< S< BCRSMatrix< FieldMatrix< T, m, n >, A > >, true >::assembleLocalProblems(), Dune::BCRSMatrix< B, A >::axpy(), Dune::MatrixDimension< BCRSMatrix< B, TA > >::coldim(), Dune::BCRSMatrix< B, A >::compress(), Dune::BCRSMatrix< B, A >::copyWindowStructure(), Dune::BCRSMatrix< B, A >::endindices(), Dune::BCRSMatrix< B, A >::endrowsizes(), Dune::BCRSMatrix< B, A >::entry(), Dune::BCRSMatrix< B, A >::frobenius_norm2(), Dune::BCRSMatrix< B, A >::infinity_norm(), Dune::BCRSMatrix< B, A >::infinity_norm_real(), Dune::BCRSMatrix< B, A >::mmhv(), Dune::BCRSMatrix< B, A >::mmtv(), Dune::BCRSMatrix< B, A >::mmv(), Dune::BCRSMatrix< B, A >::mv(), Dune::BCRSMatrix< B, A >::operator*=(), Dune::BCRSMatrix< B, A >::operator+=(), Dune::BCRSMatrix< B, A >::operator-=(), Dune::BCRSMatrix< B, A >::operator/=(), test_IO(), test_Iter(), Dune::BCRSMatrix< B, A >::umhv(), Dune::BCRSMatrix< B, A >::umtv(), Dune::BCRSMatrix< B, A >::umv(), Dune::BCRSMatrix< B, A >::usmhv(), Dune::BCRSMatrix< B, A >::usmtv(), and Dune::BCRSMatrix< B, A >::usmv().
|
inlineinherited |
Get const iterator to first row.
References Dune::BCRSMatrix< B, A >::r.
|
inlineinherited |
The currently selected build mode of the matrix.
References Dune::BCRSMatrix< B, A >::build_mode.
|
inlineinherited |
The current build stage of the matrix.
References Dune::BCRSMatrix< B, A >::ready.
|
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.
References Dune::BCRSMatrix< B, A >::a, Dune::BCRSMatrix< B, A >::allocationSize, Dune::CompressionStatistics< size_type >::avg, Dune::BCRSMatrix< B, A >::begin(), Dune::BCRSMatrix< B, A >::build_mode, Dune::BCRSMatrix< B, A >::building, Dune::BCRSMatrix< B, A >::built, Dune::CompressedBlockVectorWindow< B, A >::getindexptr(), Dune::CompressedBlockVectorWindow< B, A >::getsize(), Dune::BCRSMatrix< B, A >::implicit, Dune::BCRSMatrix< B, A >::j_, Dune::CompressionStatistics< size_type >::maximum, Dune::CompressionStatistics< size_type >::mem_ratio, Dune::BCRSMatrix< B, A >::n, Dune::BCRSMatrix< B, A >::nnz_, Dune::BCRSMatrix< B, A >::notAllocated, Dune::BCRSMatrix< B, A >::overflow, Dune::CompressionStatistics< size_type >::overflow_total, Dune::BCRSMatrix< B, A >::r, Dune::BCRSMatrix< B, A >::ready, Dune::CompressedBlockVectorWindow< B, A >::setindexptr(), Dune::CompressedBlockVectorWindow< B, A >::setptr(), and Dune::CompressedBlockVectorWindow< B, A >::setsize().
|
inlineprotectedinherited |
Copy the window structure from another matrix.
References Dune::BCRSMatrix< B, A >::begin(), Dune::BCRSMatrix< B, A >::built, Dune::BCRSMatrix< B, A >::n, Dune::BCRSMatrix< B, A >::r, Dune::BCRSMatrix< B, A >::row_wise, and Dune::BCRSMatrix< B, A >::setWindowPointers().
Referenced by Dune::BCRSMatrix< B, A >::BCRSMatrix(), and Dune::BCRSMatrix< B, A >::operator=().
|
inlineinherited |
get initial create iterator
References Dune::BCRSMatrix< B, A >::CreateIterator.
Referenced by test_IO(), and test_Iter().
|
inlineinherited |
get create iterator pointing to one after the last block
References Dune::BCRSMatrix< B, A >::CreateIterator, and Dune::BCRSMatrix< B, A >::n.
Referenced by test_IO(), and test_Iter().
|
inlineprotectedinherited |
deallocate memory of the matrix.
deallocateRows | Whether we have to deallocate the row pointers, too. If false they will not be touched. (Defaults to true). |
References col, Dune::BCRSMatrix< B, A >::n, Dune::BCRSMatrix< B, A >::notAllocated, and Dune::CompressedBlockVectorWindow< B, A >::set().
Referenced by Dune::BCRSMatrix< B, A >::operator=(), Dune::BCRSMatrix< B, A >::setSize(), and Dune::BCRSMatrix< B, A >::~BCRSMatrix().
|
inlineinherited |
Get iterator to one beyond last row.
References Dune::BCRSMatrix< B, A >::n, and Dune::BCRSMatrix< B, A >::r.
Referenced by Dune::BCRSMatrix< B, A >::addindex(), Dune::BCRSMatrix< B, A >::axpy(), Dune::MatrixDimension< BCRSMatrix< B, TA > >::coldim(), Dune::BCRSMatrix< B, A >::endindices(), Dune::BCRSMatrix< B, A >::entry(), Dune::BCRSMatrix< B, A >::exists(), Dune::BCRSMatrix< B, A >::frobenius_norm2(), Dune::BCRSMatrix< B, A >::infinity_norm(), Dune::BCRSMatrix< B, A >::infinity_norm_real(), Dune::BCRSMatrix< B, A >::mmhv(), Dune::BCRSMatrix< B, A >::mmtv(), Dune::BCRSMatrix< B, A >::mmv(), Dune::BCRSMatrix< B, A >::mv(), Dune::BCRSMatrix< B, A >::operator*=(), Dune::BCRSMatrix< B, A >::operator+=(), Dune::BCRSMatrix< B, A >::operator-=(), Dune::BCRSMatrix< B, A >::operator/=(), Dune::DynamicMatrixSubdomainSolver< BCRSMatrix< FieldMatrix< K, n, n >, Al >, X, Y >::setSubMatrix(), test_IO(), test_Iter(), Dune::BCRSMatrix< B, A >::umhv(), Dune::BCRSMatrix< B, A >::umtv(), Dune::BCRSMatrix< B, A >::umv(), Dune::BCRSMatrix< B, A >::usmhv(), Dune::BCRSMatrix< B, A >::usmtv(), and Dune::BCRSMatrix< B, A >::usmv().
|
inlineinherited |
Get const iterator to one beyond last row.
References Dune::BCRSMatrix< B, A >::n, and Dune::BCRSMatrix< B, A >::r.
|
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.
References Dune::BCRSMatrix< B, A >::avg, Dune::BCRSMatrix< B, A >::begin(), Dune::BCRSMatrix< B, A >::build_mode, Dune::BCRSMatrix< B, A >::building, Dune::BCRSMatrix< B, A >::built, col, Dune::BCRSMatrix< B, A >::end(), Dune::CompressedBlockVectorWindow< B, A >::getindexptr(), Dune::CompressedBlockVectorWindow< B, A >::getptr(), Dune::CompressedBlockVectorWindow< B, A >::getsize(), Dune::BCRSMatrix< B, A >::implicit, Dune::BCRSMatrix< B, A >::m, Dune::BCRSMatrix< B, A >::n, Dune::BCRSMatrix< B, A >::notAllocated, Dune::BCRSMatrix< B, A >::overflow, Dune::BCRSMatrix< B, A >::r, Dune::BCRSMatrix< B, A >::ready, row, and Dune::CompressedBlockVectorWindow< B, A >::setsize().
|
inlineinherited |
return true if (i,j) is in pattern
References Dune::BCRSMatrix< B, A >::end(), Dune::BCRSMatrix< B, A >::m, Dune::BCRSMatrix< B, A >::n, and Dune::BCRSMatrix< B, A >::r.
|
inlineinherited |
frobenius norm: sqrt(sum over squared values of entries)
References Dune::BCRSMatrix< B, A >::frobenius_norm2().
|
inlineinherited |
square of frobenius norm, need for block recursion
References Dune::BCRSMatrix< B, A >::begin(), Dune::BCRSMatrix< B, A >::built, Dune::BCRSMatrix< B, A >::end(), and Dune::BCRSMatrix< B, A >::ready.
Referenced by Dune::BCRSMatrix< B, A >::frobenius_norm().
|
inlineinherited |
get current number of indices in row i
References Dune::CompressedBlockVectorWindow< B, A >::getsize(), Dune::BCRSMatrix< B, A >::n, and Dune::BCRSMatrix< B, A >::r.
|
inlineprotectedinherited |
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.
References Dune::BCRSMatrix< B, A >::allocate(), Dune::BCRSMatrix< B, A >::avg, Dune::BCRSMatrix< B, A >::building, Dune::BCRSMatrix< B, A >::implicit, Dune::BCRSMatrix< B, A >::n, Dune::BCRSMatrix< B, A >::notAllocated, and Dune::CompressedBlockVectorWindow< B, A >::set().
Referenced by Dune::BCRSMatrix< B, A >::BCRSMatrix(), and Dune::BCRSMatrix< B, A >::setSize().
|
inlineinherited |
increment size of row i by s (1 by default)
References Dune::BCRSMatrix< B, A >::build_mode, Dune::BCRSMatrix< B, A >::building, Dune::BCRSMatrix< B, A >::r, Dune::BCRSMatrix< B, A >::random, Dune::BCRSMatrix< B, A >::ready, and Dune::CompressedBlockVectorWindow< B, A >::setsize().
|
inlineinherited |
infinity norm (row sum norm, how to generalize for blocks?)
References Dune::BCRSMatrix< B, A >::begin(), Dune::BCRSMatrix< B, A >::built, Dune::BCRSMatrix< B, A >::end(), and Dune::BCRSMatrix< B, A >::ready.
|
inlineinherited |
simplified infinity norm (uses Manhattan norm for complex values)
References Dune::BCRSMatrix< B, A >::begin(), Dune::BCRSMatrix< B, A >::built, Dune::BCRSMatrix< B, A >::end(), and Dune::BCRSMatrix< B, A >::ready.
|
inlineinherited |
number of columns (counted in blocks)
References Dune::BCRSMatrix< B, A >::m.
Referenced by Dune::BCRSMatrix< B, A >::axpy(), Dune::MatrixDimension< BCRSMatrix< B, TA > >::coldim(), Dune::MatrixDimension< BCRSMatrix< FieldMatrix< B, n, m >, TA > >::coldim(), Dune::BCRSMatrix< B, A >::mmhv(), Dune::BCRSMatrix< B, A >::mmtv(), Dune::BCRSMatrix< B, A >::mmv(), Dune::BCRSMatrix< B, A >::mtv(), Dune::BCRSMatrix< B, A >::mv(), Dune::BCRSMatrix< B, A >::operator+=(), Dune::BCRSMatrix< B, A >::operator-=(), Dune::SuperLUMatrix< BCRSMatrix< FieldMatrix< B, n, m >, TA > >::setMatrix(), Dune::BCRSMatrix< B, A >::umhv(), Dune::BCRSMatrix< B, A >::umtv(), Dune::BCRSMatrix< B, A >::umv(), Dune::BCRSMatrix< B, A >::usmhv(), Dune::BCRSMatrix< B, A >::usmtv(), and Dune::BCRSMatrix< B, A >::usmv().
|
inlineinherited |
|
inlineinherited |
y -= A^T x
References Dune::BCRSMatrix< B, A >::begin(), Dune::BCRSMatrix< B, A >::end(), Dune::BCRSMatrix< B, A >::M(), and Dune::BCRSMatrix< B, A >::N().
|
inlineinherited |
y -= A x
References Dune::BCRSMatrix< B, A >::begin(), Dune::BCRSMatrix< B, A >::built, Dune::BCRSMatrix< B, A >::end(), Dune::BCRSMatrix< B, A >::M(), Dune::BCRSMatrix< B, A >::N(), and Dune::BCRSMatrix< B, A >::ready.
Referenced by test_Iter().
|
inlineinherited |
|
inlineinherited |
|
inlineinherited |
number of rows (counted in blocks)
References Dune::BCRSMatrix< B, A >::n.
Referenced by Dune::BCRSMatrix< B, A >::axpy(), Dune::MatrixDimension< BCRSMatrix< B, TA > >::coldim(), Dune::BDMatrix< B, A >::invert(), Dune::BCRSMatrix< B, A >::mmhv(), Dune::BCRSMatrix< B, A >::mmtv(), Dune::BCRSMatrix< B, A >::mmv(), Dune::BCRSMatrix< B, A >::mtv(), Dune::BCRSMatrix< B, A >::mv(), Dune::BCRSMatrix< B, A >::operator+=(), Dune::BCRSMatrix< B, A >::operator-=(), Dune::MatrixDimension< BCRSMatrix< B, TA > >::rowdim(), Dune::MatrixDimension< BCRSMatrix< FieldMatrix< B, n, m >, TA > >::rowdim(), Dune::SuperLUMatrix< BCRSMatrix< FieldMatrix< B, n, m >, TA > >::setMatrix(), Dune::BTDMatrix< B, A >::solve(), Dune::BCRSMatrix< B, A >::umhv(), Dune::BCRSMatrix< B, A >::umtv(), Dune::BCRSMatrix< B, A >::umv(), Dune::BCRSMatrix< B, A >::usmhv(), Dune::BCRSMatrix< B, A >::usmtv(), and Dune::BCRSMatrix< B, A >::usmv().
|
inlineinherited |
number of blocks that are stored (the number of blocks that possibly are nonzero)
References Dune::BCRSMatrix< B, A >::nnz_.
|
inlineinherited |
vector space multiplication with scalar
References Dune::BCRSMatrix< B, A >::a, Dune::BCRSMatrix< B, A >::begin(), Dune::BCRSMatrix< B, A >::built, Dune::BCRSMatrix< B, A >::end(), Dune::BCRSMatrix< B, A >::nnz_, and Dune::BCRSMatrix< B, A >::ready.
|
inlineinherited |
Add the entries of another matrix to this one.
b | The 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 >::built, Dune::BCRSMatrix< B, A >::end(), Dune::BCRSMatrix< B, A >::M(), Dune::BCRSMatrix< B, A >::N(), and Dune::BCRSMatrix< B, A >::ready.
|
inlineinherited |
Substract the entries of another matrix to this one.
b | The 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 >::built, Dune::BCRSMatrix< B, A >::end(), Dune::BCRSMatrix< B, A >::M(), Dune::BCRSMatrix< B, A >::N(), and Dune::BCRSMatrix< B, A >::ready.
|
inlineinherited |
vector space division by scalar
References Dune::BCRSMatrix< B, A >::a, Dune::BCRSMatrix< B, A >::begin(), Dune::BCRSMatrix< B, A >::built, Dune::BCRSMatrix< B, A >::end(), Dune::BCRSMatrix< B, A >::nnz_, and Dune::BCRSMatrix< B, A >::ready.
|
inline |
assignment
References Dune::BCRSMatrix< B, A >::operator=().
|
inline |
assignment from scalar
References Dune::BCRSMatrix< B, A >::operator=().
|
inlineinherited |
random access to the rows
References Dune::BCRSMatrix< B, A >::build_mode, Dune::BCRSMatrix< B, A >::built, Dune::BCRSMatrix< B, A >::implicit, Dune::BCRSMatrix< B, A >::n, Dune::BCRSMatrix< B, A >::r, and Dune::BCRSMatrix< B, A >::ready.
|
inlineinherited |
same for read only access
References Dune::BCRSMatrix< B, A >::build_mode, Dune::BCRSMatrix< B, A >::built, Dune::BCRSMatrix< B, A >::implicit, Dune::BCRSMatrix< B, A >::n, Dune::BCRSMatrix< B, A >::r, and Dune::BCRSMatrix< B, A >::ready.
|
inlineinherited |
Sets the build mode of the matrix.
bm | The build mode to use. |
References Dune::BCRSMatrix< B, A >::build_mode, Dune::BCRSMatrix< B, A >::building, Dune::BCRSMatrix< B, A >::notAllocated, Dune::BCRSMatrix< B, A >::random, Dune::BCRSMatrix< B, A >::ready, Dune::BCRSMatrix< B, A >::row_wise, and Dune::BCRSMatrix< B, A >::unknown.
|
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).
References Dune::BCRSMatrix< B, A >::n, row, Dune::CompressedBlockVectorWindow< B, A >::set(), Dune::CompressedBlockVectorWindow< B, A >::setindexptr(), and Dune::CompressedBlockVectorWindow< B, A >::setsize().
Referenced by Dune::BCRSMatrix< B, A >::endrowsizes().
|
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().
References Dune::BCRSMatrix< B, A >::a, Dune::CompressedBlockVectorWindow< B, A >::getsize(), Dune::BCRSMatrix< B, A >::n, Dune::CompressedBlockVectorWindow< B, A >::set(), and Dune::CompressedBlockVectorWindow< B, A >::setptr().
Referenced by Dune::BCRSMatrix< B, A >::endindices(), and Dune::BCRSMatrix< B, A >::CreateIterator::operator++().
|
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
_avg | expected average number of entries per row |
_overflowsize | fraction of _n*_avg which is expected to be needed for elements that exceed _avg entries per row. |
References Dune::BCRSMatrix< B, A >::avg, Dune::BCRSMatrix< B, A >::notAllocated, Dune::BCRSMatrix< B, A >::overflowsize, and Dune::BCRSMatrix< B, A >::ready.
|
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.
Calling this method overwrites any previously set column indices!
References Dune::CompressedBlockVectorWindow< B, A >::getindexptr(), Dune::BCRSMatrix< B, A >::r, row, and Dune::compressed_base_array_unmanaged< B, A >::size().
|
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.
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). Must be omitted in implicit mode. |
References Dune::BCRSMatrix< B, A >::allocate(), Dune::BCRSMatrix< B, A >::build_mode, Dune::BCRSMatrix< B, A >::deallocate(), Dune::BCRSMatrix< B, A >::implicit, and Dune::BCRSMatrix< B, A >::implicit_allocate().
|
inlineprotectedinherited |
References Dune::BCRSMatrix< B, A >::n, row, and Dune::CompressedBlockVectorWindow< B, A >::set().
Referenced by Dune::BCRSMatrix< B, A >::copyWindowStructure().
|
inline |
Use the Thomas algorithm to solve the system Ax=b in O(n) time.
ISTLError | if the matrix is singular |
References Dune::BCRSMatrix< B, A >::N().
|
inlineinherited |
|
inlineinherited |
|
inlineinherited |
y += A x
References Dune::BCRSMatrix< B, A >::begin(), Dune::BCRSMatrix< B, A >::built, Dune::BCRSMatrix< B, A >::end(), Dune::BCRSMatrix< B, A >::M(), Dune::BCRSMatrix< B, A >::N(), and Dune::BCRSMatrix< B, A >::ready.
Referenced by test_Iter().
|
inlineinherited |
|
inlineinherited |
|
inlineinherited |
|
protectedinherited |
|
protectedinherited |
Referenced by Dune::BCRSMatrix< B, A >::compress().
|
protectedinherited |
Referenced by Dune::BCRSMatrix< B, A >::CreateIterator::operator++().
|
protectedinherited |
|
protectedinherited |
Referenced by Dune::BCRSMatrix< B, A >::addindex(), Dune::BCRSMatrix< B, A >::buildMode(), Dune::BCRSMatrix< B, A >::compress(), Dune::BCRSMatrix< B, A >::CreateIterator::CreateIterator(), Dune::BCRSMatrix< B, A >::endindices(), Dune::BCRSMatrix< B, A >::endrowsizes(), Dune::BCRSMatrix< B, A >::entry(), Dune::BCRSMatrix< B, A >::incrementrowsize(), Dune::BCRSMatrix< B, A >::operator[](), Dune::BCRSMatrix< B, A >::setBuildMode(), Dune::BCRSMatrix< B, A >::setrowsize(), and Dune::BCRSMatrix< B, A >::setSize().
|
protectedinherited |
|
protectedinherited |
Referenced by Dune::BCRSMatrix< B, A >::addindex(), Dune::BCRSMatrix< B, A >::BCRSMatrix(), Dune::BCRSMatrix< B, A >::endindices(), Dune::BCRSMatrix< B, A >::endrowsizes(), Dune::BCRSMatrix< B, A >::entry(), Dune::BCRSMatrix< B, A >::exists(), Dune::BCRSMatrix< B, A >::M(), and Dune::BCRSMatrix< B, A >::operator=().
|
protectedinherited |
Referenced by Dune::BCRSMatrix< B, A >::BCRSMatrix(), Dune::BCRSMatrix< B, A >::beforeEnd(), Dune::BCRSMatrix< B, A >::compress(), Dune::BCRSMatrix< B, A >::copyWindowStructure(), Dune::BCRSMatrix< B, A >::createend(), Dune::BCRSMatrix< B, A >::deallocate(), Dune::BCRSMatrix< B, A >::end(), Dune::BCRSMatrix< B, A >::endrowsizes(), Dune::BCRSMatrix< B, A >::entry(), Dune::BCRSMatrix< B, A >::exists(), Dune::BCRSMatrix< B, A >::getrowsize(), Dune::BCRSMatrix< B, A >::implicit_allocate(), Dune::BCRSMatrix< B, A >::N(), Dune::BCRSMatrix< B, A >::CreateIterator::operator++(), Dune::BCRSMatrix< B, A >::operator=(), Dune::BCRSMatrix< B, A >::operator[](), Dune::BCRSMatrix< B, A >::setColumnPointers(), Dune::BCRSMatrix< B, A >::setDataPointers(), and Dune::BCRSMatrix< B, A >::setWindowPointers().
|
protectedinherited |
Referenced by Dune::BCRSMatrix< B, A >::BCRSMatrix(), Dune::MatrixDimension< BCRSMatrix< B, TA > >::coldim(), Dune::BCRSMatrix< B, A >::compress(), Dune::BCRSMatrix< B, A >::endrowsizes(), Dune::BCRSMatrix< B, A >::nonzeroes(), Dune::BCRSMatrix< B, A >::operator*=(), Dune::BCRSMatrix< B, A >::CreateIterator::operator++(), Dune::BCRSMatrix< B, A >::operator/=(), and Dune::BCRSMatrix< B, A >::operator=().
|
protectedinherited |
Referenced by Dune::BCRSMatrix< B, A >::compress(), and Dune::BCRSMatrix< B, A >::entry().
|
protectedinherited |
Referenced by Dune::BCRSMatrix< B, A >::setImplicitBuildModeParameters().
|
protectedinherited |
Referenced by Dune::BCRSMatrix< B, A >::addindex(), Dune::BCRSMatrix< B, A >::BCRSMatrix(), Dune::BCRSMatrix< B, A >::beforeBegin(), Dune::BCRSMatrix< B, A >::beforeEnd(), Dune::BCRSMatrix< B, A >::begin(), Dune::MatrixDimension< BCRSMatrix< B, TA > >::coldim(), Dune::BCRSMatrix< B, A >::compress(), Dune::BCRSMatrix< B, A >::copyWindowStructure(), Dune::BCRSMatrix< B, A >::end(), Dune::BCRSMatrix< B, A >::endindices(), Dune::BCRSMatrix< B, A >::endrowsizes(), Dune::BCRSMatrix< B, A >::entry(), Dune::BCRSMatrix< B, A >::exists(), Dune::BCRSMatrix< B, A >::getrowsize(), Dune::BCRSMatrix< B, A >::incrementrowsize(), Dune::BCRSMatrix< B, A >::CreateIterator::operator++(), Dune::BCRSMatrix< B, A >::operator=(), Dune::BCRSMatrix< B, A >::operator[](), Dune::MatrixDimension< BCRSMatrix< B, TA > >::rowdim(), Dune::BCRSMatrix< B, A >::setIndices(), and Dune::BCRSMatrix< B, A >::setrowsize().
|
protectedinherited |
Referenced by Dune::BCRSMatrix< B, A >::addindex(), Dune::BCRSMatrix< B, A >::axpy(), Dune::BCRSMatrix< B, A >::BCRSMatrix(), Dune::BCRSMatrix< B, A >::buildStage(), Dune::BCRSMatrix< B, A >::compress(), Dune::BCRSMatrix< B, A >::CreateIterator::CreateIterator(), Dune::BCRSMatrix< B, A >::endindices(), Dune::BCRSMatrix< B, A >::endrowsizes(), Dune::BCRSMatrix< B, A >::entry(), Dune::BCRSMatrix< B, A >::frobenius_norm2(), Dune::BCRSMatrix< B, A >::incrementrowsize(), Dune::BCRSMatrix< B, A >::infinity_norm(), Dune::BCRSMatrix< B, A >::infinity_norm_real(), Dune::BCRSMatrix< B, A >::mmhv(), Dune::BCRSMatrix< B, A >::mmv(), Dune::BCRSMatrix< B, A >::mtv(), Dune::BCRSMatrix< B, A >::mv(), Dune::BCRSMatrix< B, A >::operator*=(), Dune::BCRSMatrix< B, A >::CreateIterator::operator++(), Dune::BCRSMatrix< B, A >::operator+=(), Dune::BCRSMatrix< B, A >::operator-=(), Dune::BCRSMatrix< B, A >::operator/=(), Dune::BCRSMatrix< B, A >::operator=(), Dune::BCRSMatrix< B, A >::operator[](), Dune::BCRSMatrix< B, A >::setBuildMode(), Dune::BCRSMatrix< B, A >::setImplicitBuildModeParameters(), Dune::BCRSMatrix< B, A >::setrowsize(), Dune::BCRSMatrix< B, A >::umhv(), Dune::BCRSMatrix< B, A >::umtv(), Dune::BCRSMatrix< B, A >::umv(), Dune::BCRSMatrix< B, A >::usmhv(), Dune::BCRSMatrix< B, A >::usmtv(), and Dune::BCRSMatrix< B, A >::usmv().
|
protectedinherited |
Referenced by Dune::BCRSMatrix< B, A >::operator=().
|
protectedinherited |
Referenced by Dune::BCRSMatrix< B, A >::CreateIterator::operator++().