Reference documentation for deal.II version 9.2.0
|
#include <deal.II/lac/precondition_block.h>
Public Types | |
using | size_type = types::global_dof_index |
using | number = typename MatrixType::value_type |
Public Member Functions | |
PreconditionBlockSOR () | |
template<typename number2 > | |
void | vmult (Vector< number2 > &, const Vector< number2 > &) const |
template<typename number2 > | |
void | vmult_add (Vector< number2 > &, const Vector< number2 > &) const |
template<typename number2 > | |
void | Tvmult (Vector< number2 > &, const Vector< number2 > &) const |
template<typename number2 > | |
void | Tvmult_add (Vector< number2 > &, const Vector< number2 > &) const |
template<typename number2 > | |
void | step (Vector< number2 > &dst, const Vector< number2 > &rhs) const |
template<typename number2 > | |
void | Tstep (Vector< number2 > &dst, const Vector< number2 > &rhs) const |
Public Member Functions inherited from Subscriptor | |
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) |
Protected Member Functions | |
PreconditionBlockSOR (bool store) | |
template<typename number2 > | |
void | forward (Vector< number2 > &, const Vector< number2 > &, const bool transpose_diagonal, const bool adding) const |
template<typename number2 > | |
void | backward (Vector< number2 > &, const Vector< number2 > &, const bool transpose_diagonal, const bool adding) const |
Protected Member Functions inherited from PreconditionBlock< MatrixType, typename MatrixType::value_type > | |
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 () |
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 () |
void | forward_step (Vector< number2 > &dst, const Vector< number2 > &prev, const Vector< number2 > &src, const bool transpose_diagonal) const |
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 |
Protected Member Functions inherited from PreconditionBlockBase< typename MatrixType::value_type > | |
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 |
Block SOR preconditioning. This class satisfies the relaxation concept.
The functions vmult
and Tvmult
execute a (transposed) block-SOR step, based on the blocks in PreconditionBlock. The elements outside the diagonal blocks may be distributed arbitrarily.
See PreconditionBlock for requirements on the matrix. The blocks used in this class must be contiguous and non-overlapping. An overlapping Schwarz relaxation method can be found in RelaxationBlockSOR; that class does not offer preconditioning, though.
Optionally, the entries of the source vector can be treated in the order of the indices in the permutation vector set by set_permutation (or the opposite order for Tvmult()). The inverse permutation is used for storing elements back into this vector. This functionality is automatically enabled after a call to set_permutation() with vectors of nonzero size.
<float> and <double>
; others can be generated in application programs (see the section on Template instantiations in the manual).Definition at line 659 of file precondition_block.h.