Reference documentation for deal.II version 9.1.1
|
#include <deal.II/lac/block_matrix_array.h>
Public Types | |
using | size_type = types::global_dof_index |
Public Member Functions | |
BlockTrianglePrecondition () | |
BlockTrianglePrecondition (const unsigned int n_blocks) | |
void | reinit (const unsigned int n_block_rows) |
template<typename MatrixType > | |
void | enter (const MatrixType &matrix, const size_type row, const size_type col, const number prefix=1., const bool transpose=false) |
void | vmult (BlockVectorType &dst, const BlockVectorType &src) const |
void | vmult_add (BlockVectorType &dst, const BlockVectorType &src) const |
void | Tvmult (BlockVectorType &dst, const BlockVectorType &src) const |
void | Tvmult_add (BlockVectorType &dst, const BlockVectorType &src) const |
void | subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
void | unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
Static Public Member Functions | |
static ::ExceptionBase & | ExcNoDiagonal (size_type arg1) |
static ::ExceptionBase & | ExcMultipleDiagonal (size_type arg1) |
Private Member Functions | |
void | do_row (BlockVectorType &dst, size_type row_num) const |
Private Member Functions inherited from BlockMatrixArray< number, BlockVectorType > | |
BlockMatrixArray () | |
BlockMatrixArray (const unsigned int n_block_rows, const unsigned int n_block_cols) | |
void | initialize (const unsigned int n_block_rows, const unsigned int n_block_cols) |
void | reinit (const unsigned int n_block_rows, const unsigned int n_block_cols) |
template<typename MatrixType > | |
void | enter (const MatrixType &matrix, const unsigned int row, const unsigned int col, const number prefix=1., const bool transpose=false) |
void | clear () |
unsigned int | n_block_rows () const |
unsigned int | n_block_cols () const |
void | vmult (BlockVectorType &dst, const BlockVectorType &src) const |
void | vmult_add (BlockVectorType &dst, const BlockVectorType &src) const |
void | Tvmult (BlockVectorType &dst, const BlockVectorType &src) const |
void | Tvmult_add (BlockVectorType &dst, const BlockVectorType &src) const |
number | matrix_scalar_product (const BlockVectorType &u, const BlockVectorType &v) const |
number | matrix_norm_square (const BlockVectorType &u) const |
template<class StreamType > | |
void | print_latex (StreamType &out) const |
Private 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) |
Private Attributes | |
bool | backward |
Private Attributes inherited from BlockMatrixArray< number, BlockVectorType > | |
std::vector< Entry > | entries |
Additional Inherited Members | |
Private Types inherited from BlockMatrixArray< number, BlockVectorType > | |
using | size_type = types::global_dof_index |
Static Private Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
Inversion of a block-triangular matrix.
In this block matrix, the inverses of the diagonal blocks are stored together with the off-diagonal blocks of a block matrix. Then, forward or backward insertion is performed block-wise. The diagonal blocks are NOT inverted for this purpose!
Like for all preconditioners, the preconditioning operation is performed by the vmult() member function.
The implementation may be a little clumsy, but it should be sufficient as long as the block sizes are much larger than the number of blocks.
Here, we document the second part of examples/doxygen/block_matrix_array.cc
. For the beginning of this file, see BlockMatrixArray.
In order to set up the preconditioner, we have to compute the inverses of the diagonal blocks ourselves. Since we used FullMatrix objects, this is fairly easy.
After creating a 2x2 BlockTrianglePrecondition object, we only fill its diagonals. The scaling factor 1/2 used for A
is the reciprocal of the scaling factor used for the matrix
itself. Remember, this preconditioner actually multiplies with the diagonal blocks.
Now, we have a block Jacobi preconditioner, which is still symmetric, since the blocks are symmetric. Therefore, we can still use the preconditioned conjugate gradient method.
Now, we enter the subdiagonal block. This is the same as in matrix
.
Since the preconditioner is not symmetric anymore, we use the GMRES method for solving.
Definition at line 421 of file block_matrix_array.h.
using BlockTrianglePrecondition< number, BlockVectorType >::size_type = types::global_dof_index |
Declare type for container size.
Definition at line 428 of file block_matrix_array.h.
BlockTrianglePrecondition< number, BlockVectorType >::BlockTrianglePrecondition | ( | ) |
Default constructor creating a useless object. initialize() must be called before using it.
Definition at line 258 of file block_matrix_array.cc.
BlockTrianglePrecondition< number, BlockVectorType >::BlockTrianglePrecondition | ( | const unsigned int | n_blocks | ) |
Constructor. This matrix must be block-quadratic, and n_blocks
is the number of blocks in each direction.
Definition at line 265 of file block_matrix_array.cc.
void BlockTrianglePrecondition< number, BlockVectorType >::reinit | ( | const unsigned int | n_block_rows | ) |
Resize preconditioner to a new size and clear all blocks.
Definition at line 274 of file block_matrix_array.cc.
void BlockTrianglePrecondition< number, BlockVectorType >::enter | ( | const MatrixType & | matrix, |
const size_type | row, | ||
const size_type | col, | ||
const number | prefix = 1. , |
||
const bool | transpose = false |
||
) |
Enter a block. This calls BlockMatrixArray::enter(). Remember that the diagonal blocks should actually be inverse matrices or preconditioners.
void BlockTrianglePrecondition< number, BlockVectorType >::vmult | ( | BlockVectorType & | dst, |
const BlockVectorType & | src | ||
) | const |
Preconditioning.
Definition at line 387 of file block_matrix_array.cc.
void BlockTrianglePrecondition< number, BlockVectorType >::vmult_add | ( | BlockVectorType & | dst, |
const BlockVectorType & | src | ||
) | const |
Preconditioning adding to dst
.
Definition at line 368 of file block_matrix_array.cc.
void BlockTrianglePrecondition< number, BlockVectorType >::Tvmult | ( | BlockVectorType & | dst, |
const BlockVectorType & | src | ||
) | const |
Transposed preconditioning
Definition at line 412 of file block_matrix_array.cc.
void BlockTrianglePrecondition< number, BlockVectorType >::Tvmult_add | ( | BlockVectorType & | dst, |
const BlockVectorType & | src | ||
) | const |
Transposed preconditioning adding to dst
.
Definition at line 422 of file block_matrix_array.cc.
|
private |
Add all off-diagonal contributions and return the entry of the diagonal element for one row.
Definition at line 282 of file block_matrix_array.cc.
void Subscriptor::subscribe |
Subscribes a user of the object by storing the pointer validity
. The subscriber may be identified by text supplied as identifier
.
Definition at line 130 of file subscriptor.cc.
void Subscriptor::unsubscribe |
Unsubscribes a user from the object.
identifier
and the validity
pointer must be the same as the one supplied to subscribe(). Definition at line 150 of file subscriptor.cc.
|
private |
Flag for backward insertion.
Definition at line 536 of file block_matrix_array.h.