|
| PreconditionBlock (bool store_diagonals=false) |
|
| ~PreconditionBlock () override=default |
|
void | initialize (const MatrixType &A, const AdditionalData parameters) |
|
void | clear () |
|
bool | empty () const |
|
value_type | el (size_type i, size_type j) const |
|
void | invert_diagblocks () |
|
template<typename number2 > |
void | forward_step (Vector< number2 > &dst, const Vector< number2 > &prev, const Vector< number2 > &src, const bool transpose_diagonal) const |
|
template<typename number2 > |
void | backward_step (Vector< number2 > &dst, const Vector< number2 > &prev, const Vector< number2 > &src, const bool transpose_diagonal) const |
|
size_type | block_size () const |
|
std::size_t | memory_consumption () const |
|
| Subscriptor () |
|
| Subscriptor (const Subscriptor &) |
|
| Subscriptor (Subscriptor &&) noexcept |
|
virtual | ~Subscriptor () |
|
Subscriptor & | operator= (const Subscriptor &) |
|
Subscriptor & | operator= (Subscriptor &&) noexcept |
|
void | subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
|
void | unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
|
unsigned int | n_subscriptions () const |
|
template<typename StreamType > |
void | list_subscribers (StreamType &stream) const |
|
void | list_subscribers () const |
|
template<class Archive > |
void | serialize (Archive &ar, const unsigned int version) |
|
|
void | initialize (const MatrixType &A, const std::vector< size_type > &permutation, const std::vector< size_type > &inverse_permutation, const AdditionalData parameters) |
|
void | set_permutation (const std::vector< size_type > &permutation, const std::vector< size_type > &inverse_permutation) |
|
void | invert_permuted_diagblocks () |
|
| PreconditionBlockBase (bool store_diagonals=false, Inversion method=gauss_jordan) |
|
| ~PreconditionBlockBase ()=default |
|
void | clear () |
|
void | reinit (unsigned int nblocks, size_type blocksize, bool compress, Inversion method=gauss_jordan) |
|
void | inverses_computed (bool are_they) |
|
bool | same_diagonal () const |
|
bool | store_diagonals () const |
|
bool | inverses_ready () const |
|
unsigned int | size () const |
|
void | inverse_vmult (size_type i, Vector< number2 > &dst, const Vector< number2 > &src) const |
|
void | inverse_Tvmult (size_type i, Vector< number2 > &dst, const Vector< number2 > &src) const |
|
FullMatrix< typename MatrixType::value_type > & | inverse (size_type i) |
|
const FullMatrix< typename MatrixType::value_type > & | inverse (size_type i) const |
|
Householder< typename MatrixType::value_type > & | inverse_householder (size_type i) |
|
const Householder< typename MatrixType::value_type > & | inverse_householder (size_type i) const |
|
LAPACKFullMatrix< typename MatrixType::value_type > & | inverse_svd (size_type i) |
|
const LAPACKFullMatrix< typename MatrixType::value_type > & | inverse_svd (size_type i) const |
|
FullMatrix< typename MatrixType::value_type > & | diagonal (size_type i) |
|
const FullMatrix< typename MatrixType::value_type > & | diagonal (size_type i) const |
|
void | log_statistics () const |
|
std::size_t | memory_consumption () const |
|
template<typename MatrixType, typename inverse_type = typename MatrixType::value_type>
class PreconditionBlock< MatrixType, inverse_type >
Base class for actual block preconditioners. This class assumes the MatrixType
consisting of invertible blocks of blocksize
on the diagonal and provides the inversion of the diagonal blocks of the matrix. It is not necessary for this class that the matrix be block diagonal; rather, it applies to matrices of arbitrary structure with the minimal property of having invertible blocks on the diagonal. Still the matrix must have access to single matrix entries. Therefore, BlockMatrixArray and similar classes are not a possible matrix class template arguments.
The block matrix structure used by this class is given, e.g., for the DG method for the transport equation. For a downstream numbering the matrices even have got a block lower left matrix structure, i.e. the matrices are empty above the diagonal blocks.
- Note
- This class is intended to be used for matrices whose structure is given by local contributions from disjoint cells, such as for DG methods. It is not intended for problems where the block structure results from different physical variables such as in the Stokes equations considered in step-22.
For all matrices that are empty above and below the diagonal blocks (i.e. for all block diagonal matrices) the BlockJacobi
preconditioner is a direct solver. For all matrices that are empty only above the diagonal blocks (e.g. the matrices one gets by the DG method with downstream numbering) BlockSOR
is a direct solver.
This first implementation of the PreconditionBlock
assumes the matrix has blocks each of the same block size. Varying block sizes within the matrix must still be implemented if needed.
The first template parameter denotes the type of number representation in the sparse matrix, the second denotes the type of number representation in which the inverted diagonal block matrices are stored within this class by invert_diagblocks()
. If you don't want to use the block inversion as an exact solver, but rather as a preconditioner, you may probably want to store the inverted blocks with less accuracy than the original matrix; for example, number==double, inverse_type=float
might be a viable choice.
- See also
- Block (linear algebra)
- Author
- Ralf Hartmann, Guido Kanschat
- Date
- 1999, 2000, 2010
Definition at line 84 of file precondition_block.h.