Reference documentation for deal.II version 9.4.1
|
#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 |
Protected Types | |
enum | Inversion { gauss_jordan , householder , svd } |
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 |
void | initialize (const MatrixType &A, const AdditionalData parameters) |
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 () |
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 |
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 |
template<typename number2 > | |
void | inverse_vmult (size_type i, Vector< number2 > &dst, const Vector< number2 > &src) const |
template<typename number2 > | |
void | inverse_Tvmult (size_type i, Vector< number2 > &dst, const Vector< number2 > &src) const |
FullMatrix< number > & | inverse (size_type i) |
const FullMatrix< number > & | inverse (size_type i) const |
Householder< number > & | inverse_householder (size_type i) |
const Householder< number > & | inverse_householder (size_type i) const |
LAPACKFullMatrix< number > & | inverse_svd (size_type i) |
const LAPACKFullMatrix< number > & | inverse_svd (size_type i) const |
FullMatrix< number > & | diagonal (size_type i) |
const FullMatrix< number > & | diagonal (size_type i) const |
void | log_statistics () const |
Static Protected Member Functions | |
static ::ExceptionBase & | ExcWrongBlockSize (int arg1, int arg2) |
static ::ExceptionBase & | ExcInverseMatricesAlreadyExist () |
static ::ExceptionBase & | ExcDiagonalsNotStored () |
static ::ExceptionBase & | ExcInverseNotAvailable () |
Protected Attributes | |
size_type | blocksize |
SmartPointer< const MatrixType, PreconditionBlock< MatrixType, inverse_type > > | A |
double | relaxation |
std::vector< size_type > | permutation |
std::vector< size_type > | inverse_permutation |
Inversion | inversion |
Private Types | |
using | value_type = inverse_type |
Private Attributes | |
unsigned int | n_diagonal_blocks |
std::vector< FullMatrix< number > > | var_inverse_full |
std::vector< Householder< number > > | var_inverse_householder |
std::vector< LAPACKFullMatrix< number > > | var_inverse_svd |
std::vector< FullMatrix< number > > | var_diagonal |
bool | var_store_diagonals |
bool | var_same_diagonal |
bool | var_inverses_ready |
Subscriptor functionality | |
Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class. | |
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) |
using | map_value_type = decltype(counter_map)::value_type |
using | map_iterator = decltype(counter_map)::iterator |
std::atomic< unsigned int > | counter |
std::map< std::string, unsigned int > | counter_map |
std::vector< std::atomic< bool > * > | validity_pointers |
const std::type_info * | object_info |
static std::mutex | mutex |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
void | check_no_subscribers () const noexcept |
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 655 of file precondition_block.h.
|
inherited |
Choose a method for inverting the blocks, and thus the data type for the inverses.
Enumerator | |
---|---|
gauss_jordan | Use the standard Gauss-Jacobi method implemented in FullMatrix::inverse(). |
householder | Use QR decomposition of the Householder class. |
svd | Use the singular value decomposition of LAPACKFullMatrix. |
Definition at line 70 of file precondition_block_base.h.
|
inlineinherited |
Resize to this number of diagonal blocks with the given block size. If compress
is true, then only one block will be stored.
Definition at line 333 of file precondition_block_base.h.
|
inlineinherited |
Tell the class that inverses are computed.
Definition at line 598 of file precondition_block_base.h.
|
inlineinherited |
Does the matrix use only one diagonal block?
Definition at line 582 of file precondition_block_base.h.
|
inlineinherited |
Check, whether diagonal blocks (not their inverses) should be stored.
Definition at line 590 of file precondition_block_base.h.
|
inlineinherited |
Return true, if inverses are ready for use.
Definition at line 606 of file precondition_block_base.h.
|
inlineinherited |
The number of blocks.
Definition at line 411 of file precondition_block_base.h.
|
inlineinherited |
Multiply with the inverse block at position i
.
Definition at line 419 of file precondition_block_base.h.
|
inlineinherited |
Multiply with the transposed inverse block at position i
.
Definition at line 448 of file precondition_block_base.h.
|
inlineinherited |
Access to the inverse diagonal blocks if Inversion is gauss_jordan.
Definition at line 526 of file precondition_block_base.h.
|
inlineinherited |
Access to the inverse diagonal blocks.
Definition at line 476 of file precondition_block_base.h.
|
inlineinherited |
Access to the inverse diagonal blocks if Inversion is householder.
Definition at line 540 of file precondition_block_base.h.
|
inlineinherited |
Access to the inverse diagonal blocks if Inversion is householder.
Definition at line 488 of file precondition_block_base.h.
|
inlineinherited |
Access to the inverse diagonal blocks if Inversion is householder.
Definition at line 554 of file precondition_block_base.h.
|
inlineinherited |
Access to the inverse diagonal blocks if Inversion is householder.
Definition at line 500 of file precondition_block_base.h.
|
inlineinherited |
Access to the diagonal blocks.
Definition at line 568 of file precondition_block_base.h.
|
inlineinherited |
Access to the diagonal blocks.
Definition at line 512 of file precondition_block_base.h.
|
inlineinherited |
Print some statistics about the inverses to deallog
. Output depends on Inversion. It is richest for svd, where we obtain statistics on extremal singular values and condition numbers.
Definition at line 614 of file precondition_block_base.h.
|
protectedinherited |
The method used for inverting blocks.
Definition at line 243 of file precondition_block_base.h.
|
privateinherited |
The number of (inverse) diagonal blocks, if only one is stored.
Definition at line 249 of file precondition_block_base.h.
|
privateinherited |
Storage of the inverse matrices of the diagonal blocks matrices as FullMatrix<number>
matrices, if Inversion gauss_jordan is used. Using number=float
saves memory in comparison with number=double
, but may introduce numerical instability.
Definition at line 257 of file precondition_block_base.h.
|
privateinherited |
Storage of the inverse matrices of the diagonal blocks matrices as Householder
matrices if Inversion householder is used. Using number=float
saves memory in comparison with number=double
, but may introduce numerical instability.
Definition at line 265 of file precondition_block_base.h.
|
privateinherited |
Storage of the inverse matrices of the diagonal blocks matrices as LAPACKFullMatrix
matrices if Inversion svd is used. Using number=float
saves memory in comparison with number=double
, but may introduce numerical instability.
Definition at line 273 of file precondition_block_base.h.
|
privateinherited |
Storage of the original diagonal blocks.
Used by the blocked SSOR method.
Definition at line 280 of file precondition_block_base.h.
|
privateinherited |
This is true, if the field var_diagonal is to be used.
Definition at line 286 of file precondition_block_base.h.
|
privateinherited |
This is true, if only one inverse is stored.
Definition at line 291 of file precondition_block_base.h.
|
privateinherited |
The inverse matrices are usable. Set by the parent class via inverses_computed().
Definition at line 297 of file precondition_block_base.h.