![]() |
Reference documentation for deal.II version 9.0.0
|
#include <deal.II/lac/precondition_block.h>
Classes | |
class | const_iterator |
Public Types | |
typedef types::global_dof_index | size_type |
Public Member Functions | |
template<typename number2 > | |
void | vmult (Vector< number2 > &, const Vector< number2 > &) const |
template<typename number2 > | |
void | Tvmult (Vector< number2 > &, const Vector< number2 > &) const |
template<typename number2 > | |
void | vmult_add (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 |
const_iterator | begin () const |
const_iterator | end () const |
const_iterator | begin (const size_type r) const |
const_iterator | end (const size_type r) const |
![]() | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) noexcept | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) noexcept |
void | subscribe (const char *identifier=nullptr) const |
void | unsubscribe (const char *identifier=nullptr) const |
unsigned int | n_subscriptions () const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Private Types | |
typedef MatrixType::value_type | number |
![]() | |
typedef types::global_dof_index | size_type |
Private Member Functions | |
template<typename number2 > | |
void | do_vmult (Vector< number2 > &, const Vector< number2 > &, bool adding) const |
![]() | |
PreconditionBlock (bool store_diagonals=false) | |
~PreconditionBlock ()=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 |
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 () |
Additional Inherited Members | |
![]() | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
![]() | |
static ::ExceptionBase & | ExcWrongBlockSize (int arg1, int arg2) |
static ::ExceptionBase & | ExcInverseMatricesAlreadyExist () |
![]() | |
size_type | blocksize |
SmartPointer< const MatrixType, PreconditionBlock< MatrixType, inverse_type > > | A |
double | relaxation |
std::vector< size_type > | permutation |
std::vector< size_type > | inverse_permutation |
Block Jacobi preconditioning. See PreconditionBlock for requirements on the matrix. This class satisfies the relaxation concept.
<float> and <double>
; others can be generated in application programs (see the section on Template instantiations in the manual).Definition at line 367 of file precondition_block.h.
|
private |
Define number type of matrix.
Definition at line 374 of file precondition_block.h.
typedef types::global_dof_index PreconditionBlockJacobi< MatrixType, inverse_type >::size_type |
Declare type for container size.
Definition at line 380 of file precondition_block.h.
template void PreconditionBlockJacobi< MatrixType, inverse_type >::vmult< double > | ( | Vector< number2 > & | , |
const Vector< number2 > & | |||
) | const |
Execute block Jacobi preconditioning.
This function will automatically use the inverse matrices if they exist, if not then BlockJacobi will need much time inverting the diagonal block matrices in each preconditioning step.
template void PreconditionBlockJacobi< MatrixType, inverse_type >::Tvmult< double > | ( | Vector< number2 > & | , |
const Vector< number2 > & | |||
) | const |
Same as vmult
, since Jacobi is symmetric.
template void PreconditionBlockJacobi< MatrixType, inverse_type >::vmult_add< double > | ( | Vector< number2 > & | , |
const Vector< number2 > & | |||
) | const |
Execute block Jacobi preconditioning, adding to dst
.
This function will automatically use the inverse matrices if they exist, if not then BlockJacobi will need much time inverting the diagonal block matrices in each preconditioning step.
template void PreconditionBlockJacobi< MatrixType, inverse_type >::Tvmult_add< double > | ( | Vector< number2 > & | , |
const Vector< number2 > & | |||
) | const |
Same as vmult_add
, since Jacobi is symmetric.
void PreconditionBlockJacobi< MatrixType, inverse_type >::step | ( | Vector< number2 > & | dst, |
const Vector< number2 > & | rhs | ||
) | const |
Perform one step of the Jacobi iteration.
void PreconditionBlockJacobi< MatrixType, inverse_type >::Tstep | ( | Vector< number2 > & | dst, |
const Vector< number2 > & | rhs | ||
) | const |
Perform one step of the Jacobi iteration.
const_iterator PreconditionBlockJacobi< MatrixType, inverse_type >::begin | ( | ) | const |
Iterator starting at the first entry.
const_iterator PreconditionBlockJacobi< MatrixType, inverse_type >::end | ( | ) | const |
Final iterator.
const_iterator PreconditionBlockJacobi< MatrixType, inverse_type >::begin | ( | const size_type | r | ) | const |
Iterator starting at the first entry of row r
.
const_iterator PreconditionBlockJacobi< MatrixType, inverse_type >::end | ( | const size_type | r | ) | const |
Final iterator of row r
.
|
private |
Actual implementation of the preconditioner.
Depending on adding
, the result of preconditioning is added to the destination vector.