Reference documentation for deal.II version 9.2.0
|
Modules | |
Exceptions and assertions | |
This module contains classes that are used in the exception mechanism of deal.II. | |
Enumerations | |
enum | SparseBlockVanka< number >::BlockingStrategy { SparseBlockVanka< number >::index_intervals, SparseBlockVanka< number >::adaptive } |
Friends | |
template<typename T > | |
class | SparseVanka< number >::SparseBlockVanka |
Preconditioners are used to accelerate the iterative solution of linear systems. Typical preconditioners are Jacobi, Gauss-Seidel, or SSOR, but the library also supports more complex ones such as Vanka or incomplete LU decompositions (ILU). In addition, sparse direct solvers can be used as preconditioners when available.
Broadly speaking, preconditioners are operators, which are multiplied with a matrix to improve conditioning. The idea is, that the preconditioned system P-1Ax = P-1b is much easier to solve than the original system Ax = b. What this means exactly depends on the structure of the matrix and cannot be discussed here in generality. For symmetric, positive definite matrices A and P, it means that the spectral condition number (the quotient of greatest and smallest eigenvalue) of P-1A is much smaller than the one of A.
At hand of the simplest example, Richardson iteration, implemented in SolverRichardson, the preconditioned iteration looks like
\[ x^{k+1} = x^k - P^{-1} \bigl(A x^k - b\bigr). \]
Accordingly, preconditioning amounts to applying a linear operator to the residual, and consequently, the action of the preconditioner P-1 is implemented as vmult()
. Templates in deal.II that require a preconditioner indicate the requirement with the PreconditionerType concept. In practice, one can usually treat any matrix-like object which defines vmult()
and Tvmult()
as a preconditioner. All preconditioner classes in this module implement this interface.
When used in Krylov space methods, it is up to the method, whether it simply replaces multiplications with A by those with P-1A (for instance SolverBicgstab), or does more sophisticated things. SolverCG for instance uses P-1 to define an inner product, which is the reason why it requires a symmetric, positive definite operator P.
Many preconditioners rely on an additive splitting A = P - N into two matrices. In this case, the iteration step of the Richardson method above can be simplified to
\[ x^{k+1} = P^{-1} \bigl(N x^k + b\bigr), \]
thus avoiding multiplication with A completely. We call operators mapping the previous iterate xk to the next iterate in this way relaxation operators. Their generic interface is given by the RelaxationType concept. The classes with names starting with Relaxation
in this module implement this interface, as well as the preconditioners PreconditionJacobi, PreconditionSOR, PreconditionBlockJacobi, PreconditionBlockSOR, and PreconditionBlockSSOR.
In this section, we discuss the interface preconditioners usually have to provide to work inside the deal.II library.
In order to be able to be stored in containers, all preconditioners have a constructor with no arguments. Since this will typically produce a useless object, all preconditioners have a function
This function receives the matrix to be preconditioned as well as additional required parameters and sets up the internal structures of the preconditioner.
Some preconditioners, like SOR and Jacobi, were used as iterative solvers long before they were used as preconditioners. Thus, they satisfy both MatrixType and RelaxationType concepts.
Declare type for container size.
Definition at line 86 of file precondition.h.
Declare type for container size.
Definition at line 205 of file precondition.h.
using PreconditionUseMatrix< MatrixType, VectorType >::function_ptr = void (MatrixType::*)(VectorType &, const VectorType &) const |
Type of the preconditioning function of the matrix.
Definition at line 372 of file precondition.h.
using PreconditionRelaxation< MatrixType >::size_type = typename MatrixType::size_type |
Declare type for container size.
Definition at line 417 of file precondition.h.
using PreconditionJacobi< MatrixType >::AdditionalData = typename PreconditionRelaxation<MatrixType>::AdditionalData |
An alias to the base class AdditionalData.
Definition at line 515 of file precondition.h.
using PreconditionSOR< MatrixType >::AdditionalData = typename PreconditionRelaxation<MatrixType>::AdditionalData |
An alias to the base class AdditionalData.
Definition at line 603 of file precondition.h.
using PreconditionSSOR< MatrixType >::AdditionalData = typename PreconditionRelaxation<MatrixType>::AdditionalData |
An alias to the base class AdditionalData.
Definition at line 672 of file precondition.h.
using PreconditionSSOR< MatrixType >::size_type = typename MatrixType::size_type |
Declare type for container size.
Definition at line 677 of file precondition.h.
using PreconditionSSOR< MatrixType >::BaseClass = PreconditionRelaxation<MatrixType> |
An alias to the base class.
Definition at line 682 of file precondition.h.
using PreconditionPSOR< MatrixType >::size_type = typename MatrixType::size_type |
Declare type for container size.
Definition at line 773 of file precondition.h.
using PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::size_type = types::global_dof_index |
Declare type for container size.
Definition at line 1018 of file precondition.h.
|
private |
Define number type of matrix.
Definition at line 91 of file precondition_block.h.
|
private |
Value type for inverse matrices.
Definition at line 96 of file precondition_block.h.
using PreconditionBlock< MatrixType, inverse_type >::size_type = types::global_dof_index |
Declare type for container size.
Definition at line 102 of file precondition_block.h.
using PreconditionSelector< MatrixType, VectorType >::size_type = typename MatrixType::size_type |
Declare type for container size.
Definition at line 109 of file precondition_selector.h.
using SparseLUDecomposition< number >::size_type = typename SparseMatrix<number>::size_type |
Declare type for container size.
Definition at line 126 of file sparse_decomposition.h.
using SparseILU< number >::size_type = typename SparseLUDecomposition<number>::size_type |
Declare type for container size.
Definition at line 67 of file sparse_ilu.h.
using SparseILU< number >::AdditionalData = typename SparseLUDecomposition<number>::AdditionalData |
Make SparseLUDecomposition::AdditionalData accessible to this class as well.
Definition at line 81 of file sparse_ilu.h.
using SparseMIC< number >::size_type = types::global_dof_index |
Declare type for container size.
Definition at line 55 of file sparse_mic.h.
using SparseMIC< number >::AdditionalData = typename SparseLUDecomposition<number>::AdditionalData |
Make the AdditionalData
type in the base class accessible to this class as well.
Definition at line 79 of file sparse_mic.h.
using SparseVanka< number >::size_type = types::global_dof_index |
Declare type for container size.
Definition at line 146 of file sparse_vanka.h.
using SparseBlockVanka< number >::size_type = types::global_dof_index |
Declare type for container size.
Definition at line 527 of file sparse_vanka.h.
enum SparseBlockVanka::BlockingStrategy |
Enumeration of the different methods by which the DoFs are distributed to the blocks on which we are to work.
Enumerator | |
---|---|
index_intervals | Block by index intervals. |
adaptive | Block with an adaptive strategy. |
Definition at line 533 of file sparse_vanka.h.
|
default |
Constructor.
PreconditionIdentity::PreconditionIdentity | ( | ) |
Constructor, sets the domain and range sizes to their defaults.
void PreconditionIdentity::initialize | ( | const MatrixType & | matrix, |
const AdditionalData & | additional_data = AdditionalData() |
||
) |
The matrix argument is ignored and here just for compatibility with more complex preconditioners.
void PreconditionIdentity::vmult | ( | VectorType & | , |
const VectorType & | |||
) | const |
Apply preconditioner.
void PreconditionIdentity::Tvmult | ( | VectorType & | , |
const VectorType & | |||
) | const |
Apply transpose preconditioner. Since this is the identity, this function is the same as vmult().
void PreconditionIdentity::vmult_add | ( | VectorType & | , |
const VectorType & | |||
) | const |
Apply preconditioner, adding to the previous value.
void PreconditionIdentity::Tvmult_add | ( | VectorType & | , |
const VectorType & | |||
) | const |
Apply transpose preconditioner, adding. Since this is the identity, this function is the same as vmult_add().
|
inline |
This function is only present to provide the interface of a preconditioner to be handed to a smoother. This does nothing.
Definition at line 149 of file precondition.h.
size_type PreconditionIdentity::m | ( | ) | const |
Return the dimension of the codomain (or range) space. Note that the matrix is of dimension \(m \times n\).
size_type PreconditionIdentity::n | ( | ) | const |
Return the dimension of the domain space. Note that the matrix is of dimension \(m \times n\).
PreconditionRichardson::AdditionalData::AdditionalData | ( | const double | relaxation = 1. | ) |
Constructor. Block size must be given since there is no reasonable default parameter.
PreconditionRichardson::PreconditionRichardson | ( | ) |
Constructor, sets the relaxation parameter, domain and range sizes to their default.
void PreconditionRichardson::initialize | ( | const AdditionalData & | parameters | ) |
Change the relaxation parameter.
void PreconditionRichardson::initialize | ( | const MatrixType & | matrix, |
const AdditionalData & | parameters | ||
) |
Change the relaxation parameter in a way consistent with other preconditioners. The matrix argument is ignored and here just for compatibility with more complex preconditioners.
void PreconditionRichardson::vmult | ( | VectorType & | , |
const VectorType & | |||
) | const |
Apply preconditioner.
void PreconditionRichardson::Tvmult | ( | VectorType & | , |
const VectorType & | |||
) | const |
Apply transpose preconditioner. Since this is the identity, this function is the same as vmult().
void PreconditionRichardson::vmult_add | ( | VectorType & | , |
const VectorType & | |||
) | const |
Apply preconditioner, adding to the previous value.
void PreconditionRichardson::Tvmult_add | ( | VectorType & | , |
const VectorType & | |||
) | const |
Apply transpose preconditioner, adding. Since this is the identity, this function is the same as vmult_add().
|
inline |
This function is only present to provide the interface of a preconditioner to be handed to a smoother. This does nothing.
Definition at line 280 of file precondition.h.
size_type PreconditionRichardson::m | ( | ) | const |
Return the dimension of the codomain (or range) space. Note that the matrix is of dimension \(m \times n\).
size_type PreconditionRichardson::n | ( | ) | const |
Return the dimension of the domain space. Note that the matrix is of dimension \(m \times n\).
PreconditionUseMatrix< MatrixType, VectorType >::PreconditionUseMatrix | ( | const MatrixType & | M, |
const function_ptr | method | ||
) |
Constructor. This constructor stores a reference to the matrix object for later use and selects a preconditioning method, which must be a member function of that matrix.
void PreconditionUseMatrix< MatrixType, VectorType >::vmult | ( | VectorType & | dst, |
const VectorType & | src | ||
) | const |
Execute preconditioning. Calls the function passed to the constructor of this object with the two arguments given here.
PreconditionRelaxation< MatrixType >::AdditionalData::AdditionalData | ( | const double | relaxation = 1. | ) |
Constructor.
void PreconditionRelaxation< MatrixType >::initialize | ( | const MatrixType & | A, |
const AdditionalData & | parameters = AdditionalData() |
||
) |
Initialize matrix and relaxation parameter. The matrix is just stored in the preconditioner object. The relaxation parameter should be larger than zero and smaller than 2 for numerical reasons. It defaults to 1.
void PreconditionRelaxation< MatrixType >::clear | ( | ) |
Release the matrix and reset its pointer.
size_type PreconditionRelaxation< MatrixType >::m | ( | ) | const |
Return the dimension of the codomain (or range) space. Note that the matrix is of dimension \(m \times n\).
size_type PreconditionRelaxation< MatrixType >::n | ( | ) | const |
Return the dimension of the domain space. Note that the matrix is of dimension \(m \times n\).
void PreconditionJacobi< MatrixType >::vmult | ( | VectorType & | , |
const VectorType & | |||
) | const |
Apply preconditioner.
void PreconditionJacobi< MatrixType >::Tvmult | ( | VectorType & | , |
const VectorType & | |||
) | const |
Apply transpose preconditioner. Since this is a symmetric preconditioner, this function is the same as vmult().
void PreconditionJacobi< MatrixType >::step | ( | VectorType & | x, |
const VectorType & | rhs | ||
) | const |
Perform one step of the preconditioned Richardson iteration.
void PreconditionJacobi< MatrixType >::Tstep | ( | VectorType & | x, |
const VectorType & | rhs | ||
) | const |
Perform one transposed step of the preconditioned Richardson iteration.
void PreconditionSOR< MatrixType >::vmult | ( | VectorType & | , |
const VectorType & | |||
) | const |
Apply preconditioner.
void PreconditionSOR< MatrixType >::Tvmult | ( | VectorType & | , |
const VectorType & | |||
) | const |
Apply transpose preconditioner.
void PreconditionSOR< MatrixType >::step | ( | VectorType & | x, |
const VectorType & | rhs | ||
) | const |
Perform one step of the preconditioned Richardson iteration.
void PreconditionSOR< MatrixType >::Tstep | ( | VectorType & | x, |
const VectorType & | rhs | ||
) | const |
Perform one transposed step of the preconditioned Richardson iteration.
void PreconditionSSOR< MatrixType >::initialize | ( | const MatrixType & | A, |
const typename BaseClass::AdditionalData & | parameters = typename BaseClass::AdditionalData() |
||
) |
Initialize matrix and relaxation parameter. The matrix is just stored in the preconditioner object. The relaxation parameter should be larger than zero and smaller than 2 for numerical reasons. It defaults to 1.
void PreconditionSSOR< MatrixType >::vmult | ( | VectorType & | , |
const VectorType & | |||
) | const |
Apply preconditioner.
void PreconditionSSOR< MatrixType >::Tvmult | ( | VectorType & | , |
const VectorType & | |||
) | const |
Apply transpose preconditioner. Since this is a symmetric preconditioner, this function is the same as vmult().
void PreconditionSSOR< MatrixType >::step | ( | VectorType & | x, |
const VectorType & | rhs | ||
) | const |
Perform one step of the preconditioned Richardson iteration
void PreconditionSSOR< MatrixType >::Tstep | ( | VectorType & | x, |
const VectorType & | rhs | ||
) | const |
Perform one transposed step of the preconditioned Richardson iteration.
PreconditionPSOR< MatrixType >::AdditionalData::AdditionalData | ( | const std::vector< size_type > & | permutation, |
const std::vector< size_type > & | inverse_permutation, | ||
const typename PreconditionRelaxation< MatrixType >::AdditionalData & | parameters = typename PreconditionRelaxation< MatrixType >::AdditionalData() |
||
) |
Constructor. For the parameters' description, see below.
The permutation vectors are stored as a reference. Therefore, it has to be assured that the lifetime of the vector exceeds the lifetime of the preconditioner.
The relaxation parameter should be larger than zero and smaller than 2 for numerical reasons. It defaults to 1.
void PreconditionPSOR< MatrixType >::initialize | ( | const MatrixType & | A, |
const std::vector< size_type > & | permutation, | ||
const std::vector< size_type > & | inverse_permutation, | ||
const typename PreconditionRelaxation< MatrixType >::AdditionalData & | parameters = typename PreconditionRelaxation< MatrixType >::AdditionalData() |
||
) |
Initialize matrix and relaxation parameter. The matrix is just stored in the preconditioner object.
The permutation vector is stored as a pointer. Therefore, it has to be assured that the lifetime of the vector exceeds the lifetime of the preconditioner.
The relaxation parameter should be larger than zero and smaller than 2 for numerical reasons. It defaults to 1.
void PreconditionPSOR< MatrixType >::initialize | ( | const MatrixType & | A, |
const AdditionalData & | additional_data | ||
) |
Initialize matrix and relaxation parameter. The matrix is just stored in the preconditioner object.
For more detail about possible parameters, see the class documentation and the documentation of the PreconditionPSOR::AdditionalData class.
After this function is called the preconditioner is ready to be used (using the vmult
function of derived classes).
void PreconditionPSOR< MatrixType >::vmult | ( | VectorType & | , |
const VectorType & | |||
) | const |
Apply preconditioner.
void PreconditionPSOR< MatrixType >::Tvmult | ( | VectorType & | , |
const VectorType & | |||
) | const |
Apply transpose preconditioner.
PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::AdditionalData | ( | const unsigned int | degree = 1 , |
const double | smoothing_range = 0. , |
||
const unsigned int | eig_cg_n_iterations = 8 , |
||
const double | eig_cg_residual = 1e-2 , |
||
const double | max_eigenvalue = 1 |
||
) |
Constructor.
AdditionalData& PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::operator= | ( | const AdditionalData & | other_data | ) |
Copy assignment operator.
PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::PreconditionChebyshev | ( | ) |
Constructor.
void PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::initialize | ( | const MatrixType & | matrix, |
const AdditionalData & | additional_data = AdditionalData() |
||
) |
Initialize function. Takes the matrix which is used to form the preconditioner, and additional flags if there are any. This function works only if the input matrix has an operator el(i,i)
for accessing all the elements in the diagonal. Alternatively, the diagonal can be supplied with the help of the AdditionalData field.
This function calculates an estimate of the eigenvalue range of the matrix weighted by its diagonal using a modified CG iteration in case the given number of iterations is positive.
void PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::vmult | ( | VectorType & | dst, |
const VectorType & | src | ||
) | const |
Compute the action of the preconditioner on src
, storing the result in dst
.
void PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::Tvmult | ( | VectorType & | dst, |
const VectorType & | src | ||
) | const |
Compute the action of the transposed preconditioner on src
, storing the result in dst
.
void PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::step | ( | VectorType & | dst, |
const VectorType & | src | ||
) | const |
Perform one step of the preconditioned Richardson iteration.
void PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::Tstep | ( | VectorType & | dst, |
const VectorType & | src | ||
) | const |
Perform one transposed step of the preconditioned Richardson iteration.
void PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::clear | ( | ) |
Resets the preconditioner.
size_type PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::m | ( | ) | const |
Return the dimension of the codomain (or range) space. Note that the matrix is of dimension \(m \times n\).
size_type PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::n | ( | ) | const |
Return the dimension of the domain space. Note that the matrix is of dimension \(m \times n\).
|
inline |
Constructor initializing with invalid values.
Definition at line 1193 of file precondition.h.
EigenvalueInformation PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::estimate_eigenvalues | ( | const VectorType & | src | ) | const |
Compute eigenvalue estimates required for the preconditioner.
This function is called automatically on first use of the preconditioner if it is not called by the user. The layout of the vector src
is used to create internal temporary vectors and its content does not matter.
Initializes the factors theta and delta based on an eigenvalue computation. If the user set provided values for the largest eigenvalue in AdditionalData, no computation is performed and the information given by the user is used.
PreconditionBlock< MatrixType, inverse_type >::PreconditionBlock | ( | bool | store_diagonals = false | ) |
Constructor.
|
overridedefault |
Destructor.
void PreconditionBlock< MatrixType, inverse_type >::initialize | ( | const MatrixType & | A, |
const AdditionalData | parameters | ||
) |
Initialize matrix and block size. We store the matrix and the block size in the preconditioner object. In a second step, the inverses of the diagonal blocks may be computed.
Additionally, a relaxation parameter for derived classes may be provided.
|
protected |
Initialize matrix and block size for permuted preconditioning. Additionally to the parameters of the other initialize() function, we hand over two index vectors with the permutation and its inverse. For the meaning of these vectors see PreconditionBlockSOR.
In a second step, the inverses of the diagonal blocks may be computed. Make sure you use invert_permuted_diagblocks() to yield consistent data.
Additionally, a relaxation parameter for derived classes may be provided.
|
protected |
Set either the permutation of rows or the permutation of blocks, depending on the size of the vector.
If the size of the permutation vectors is equal to the dimension of the linear system, it is assumed that rows are permuted individually. In this case, set_permutation() must be called before initialize(), since the diagonal blocks are built from the permuted entries of the matrix.
If the size of the permutation vector is not equal to the dimension of the system, the diagonal blocks are computed from the unpermuted entries. Instead, the relaxation methods step() and Tstep() are executed applying the blocks in the order given by the permutation vector. They will throw an exception if length of this vector is not equal to the number of blocks.
|
protected |
Replacement of invert_diagblocks() for permuted preconditioning.
void PreconditionBlock< MatrixType, inverse_type >::clear | ( | ) |
Deletes the inverse diagonal block matrices if existent, sets the blocksize to 0, hence leaves the class in the state that it had directly after calling the constructor.
bool PreconditionBlock< MatrixType, inverse_type >::empty | ( | ) | const |
Check whether the object is empty.
value_type PreconditionBlock< MatrixType, inverse_type >::el | ( | size_type | i, |
size_type | j | ||
) | const |
Read-only access to entries. This function is only possible if the inverse diagonal blocks are stored.
void PreconditionBlock< MatrixType, inverse_type >::invert_diagblocks | ( | ) |
Stores the inverse of the diagonal blocks in inverse
. This costs some additional memory - for DG methods about 1/3 (for double inverses) or 1/6 (for float inverses) of that used for the matrix - but it makes the preconditioning much faster.
It is not allowed to call this function twice (will produce an error) before a call of clear(...)
because at the second time there already exist the inverse matrices.
After this function is called, the lock on the matrix given through the use_matrix
function is released, i.e. you may overwrite of delete it. You may want to do this in case you use this matrix to precondition another matrix.
void PreconditionBlock< MatrixType, inverse_type >::forward_step | ( | Vector< number2 > & | dst, |
const Vector< number2 > & | prev, | ||
const Vector< number2 > & | src, | ||
const bool | transpose_diagonal | ||
) | const |
Perform one block relaxation step in forward numbering.
Depending on the arguments dst
and pref
, this performs an SOR step (both reference the same vector) of a Jacobi step (both different vectors). For the Jacobi step, the calling function must copy dst
to pref
after this.
void PreconditionBlock< MatrixType, inverse_type >::backward_step | ( | Vector< number2 > & | dst, |
const Vector< number2 > & | prev, | ||
const Vector< number2 > & | src, | ||
const bool | transpose_diagonal | ||
) | const |
Perform one block relaxation step in backward numbering.
Depending on the arguments dst
and pref
, this performs an SOR step (both reference the same vector) of a Jacobi step (both different vectors). For the Jacobi step, the calling function must copy dst
to pref
after this.
size_type PreconditionBlock< MatrixType, inverse_type >::block_size | ( | ) | const |
Return the size of the blocks.
std::size_t PreconditionBlock< MatrixType, inverse_type >::memory_consumption | ( | ) | const |
Determine an estimate for the memory consumption (in bytes) of this object.
PreconditionSelector< MatrixType, VectorType >::PreconditionSelector | ( | const std::string & | preconditioning, |
const typename VectorType::value_type & | omega = 1. |
||
) |
Constructor. omega
denotes the damping parameter of the preconditioning.
Definition at line 208 of file precondition_selector.h.
|
overridevirtual |
Destructor.
Definition at line 217 of file precondition_selector.h.
void PreconditionSelector< MatrixType, VectorType >::use_matrix | ( | const MatrixType & | M | ) |
Takes the matrix that is needed for preconditionings that involves a matrix. e.g. for precondition_jacobi
, ~_sor
, ~_ssor
.
Definition at line 226 of file precondition_selector.h.
|
inline |
Return the dimension of the codomain (or range) space. Note that the matrix is of dimension \(m \times n\).
Definition at line 234 of file precondition_selector.h.
|
inline |
Return the dimension of the domain space. Note that the matrix is of dimension \(m \times n\).
Definition at line 243 of file precondition_selector.h.
|
virtual |
Precondition procedure. Calls the preconditioning that was specified in the constructor.
Definition at line 253 of file precondition_selector.h.
|
virtual |
Transpose precondition procedure. Calls the preconditioning that was specified in the constructor.
Definition at line 284 of file precondition_selector.h.
|
static |
Get the names of all implemented preconditionings. The list of possible options includes:
Definition at line 316 of file precondition_selector.h.
|
protected |
Constructor. Does nothing.
Call the initialize() function before using this object as preconditioner (vmult()).
|
overridepure virtual |
Destruction. Mark the destructor pure to ensure that this class isn't used directly, but only its derived classes.
|
overridevirtual |
Deletes all member variables. Leaves the class in the state that it had directly after calling the constructor
Reimplemented from SparseMatrix< number >.
Reimplemented in SparseMIC< number >.
template void SparseLUDecomposition< number >::initialize< float > | ( | const SparseMatrix< somenumber > & | matrix, |
const AdditionalData | parameters | ||
) |
This function needs to be called before an object of this class is used as preconditioner.
For more detail about possible parameters, see the class documentation and the documentation of the SparseLUDecomposition::AdditionalData class.
According to the parameters
, this function creates a new SparsityPattern or keeps the previous sparsity or takes the sparsity given by the user to data
. Then, this function performs the LU decomposition.
After this function is called the preconditioner is ready to be used (using the vmult
function of derived classes).
bool SparseLUDecomposition< number >::empty | ( | ) | const |
Return whether the object is empty. It calls the inherited SparseMatrix::empty() function.
size_type SparseLUDecomposition< number >::m | ( | ) | const |
Return the dimension of the codomain (or range) space. It calls the inherited SparseMatrix::m() function. Note that the matrix is of dimension \(m \times n\).
size_type SparseLUDecomposition< number >::n | ( | ) | const |
Return the dimension of the domain space. It calls the inherited SparseMatrix::n() function. Note that the matrix is of dimension \(m \times n\).
void SparseLUDecomposition< number >::vmult_add | ( | OutVector & | dst, |
const InVector & | src | ||
) | const |
Adding Matrix-vector multiplication. Add M*src on dst with M being this matrix.
Source and destination must not be the same vector.
void SparseLUDecomposition< number >::Tvmult_add | ( | OutVector & | dst, |
const InVector & | src | ||
) | const |
Adding Matrix-vector multiplication. Add MT*src to dst with M being this matrix. This function does the same as vmult_add() but takes the transposed matrix.
Source and destination must not be the same vector.
|
virtual |
Determine an estimate for the memory consumption (in bytes) of this object.
Reimplemented in SparseILU< number >, and SparseMIC< number >.
Constructor. Does nothing.
Call the initialize
function before using this object as preconditioner.
template void SparseILU< number >::initialize< float > | ( | const SparseMatrix< somenumber > & | matrix, |
const AdditionalData & | parameters = AdditionalData() |
||
) |
Perform the incomplete LU factorization of the given matrix.
This function needs to be called before an object of this class is used as preconditioner.
For more details about possible parameters, see the class documentation of SparseLUDecomposition and the documentation of the SparseLUDecomposition::AdditionalData
class.
According to the parameters
, this function creates a new SparsityPattern or keeps the previous sparsity or takes the sparsity given by the user to data
. Then, this function performs the LU decomposition.
After this function is called the preconditioner is ready to be used.
template void SparseILU< number >::vmult< float > | ( | Vector< somenumber > & | dst, |
const Vector< somenumber > & | src | ||
) | const |
Apply the incomplete decomposition, i.e. do one forward-backward step \(dst=(LU)^{-1}src\).
The initialize() function needs to be called before.
template void SparseILU< number >::Tvmult< float > | ( | Vector< somenumber > & | dst, |
const Vector< somenumber > & | src | ||
) | const |
Apply the transpose of the incomplete decomposition, i.e. do one forward- backward step \(dst=(LU)^{-T}src\).
The initialize() function needs to be called before.
|
overridevirtual |
Determine an estimate for the memory consumption (in bytes) of this object.
Reimplemented from SparseLUDecomposition< number >.
Constructor. Does nothing, so you have to call decompose
sometimes afterwards.
|
overridevirtual |
Deletes all member variables. Leaves the class in the state that it had directly after calling the constructor
Reimplemented from SparseLUDecomposition< number >.
template void SparseMIC< number >::initialize< float > | ( | const SparseMatrix< somenumber > & | matrix, |
const AdditionalData & | parameters = AdditionalData() |
||
) |
Perform the incomplete LU factorization of the given matrix.
This function needs to be called before an object of this class is used as preconditioner.
For more details about possible parameters, see the class documentation of SparseLUDecomposition and the documentation of the SparseLUDecomposition::AdditionalData
class.
According to the parameters
, this function creates a new SparsityPattern or keeps the previous sparsity or takes the sparsity given by the user to data
. Then, this function performs the MIC decomposition.
After this function is called the preconditioner is ready to be used.
template void SparseMIC< number >::vmult< float > | ( | Vector< somenumber > & | dst, |
const Vector< somenumber > & | src | ||
) | const |
Apply the incomplete decomposition, i.e. do one forward-backward step \(dst=(LU)^{-1}src\).
Call initialize
before calling this function.
template void SparseMIC< number >::Tvmult< float > | ( | Vector< somenumber > & | dst, |
const Vector< somenumber > & | src | ||
) | const |
Apply the transpose of the incomplete decomposition, i.e. do one forward- backward step \(dst=(LU)^{-1}src\).
Call initialize
before calling this function.
|
overridevirtual |
Determine an estimate for the memory consumption (in bytes) of this object.
Reimplemented from SparseLUDecomposition< number >.
SparseVanka< number >::SparseVanka | ( | ) |
Constructor. Does nothing.
Call the initialize() function before using this object as preconditioner (vmult()).
SparseVanka< number >::SparseVanka | ( | const SparseMatrix< number > & | M, |
const std::vector< bool > & | selected, | ||
const bool | conserve_memory = false , |
||
const unsigned int | n_threads = MultithreadInfo::n_threads() |
||
) |
Constructor. Gets the matrix for preconditioning and a bit vector with entries true
for all rows to be updated. A reference to this vector will be stored, so it must persist longer than the Vanka object. The same is true for the matrix.
The matrix M
which is passed here may or may not be the same matrix for which this object shall act as preconditioner. In particular, it is conceivable that the preconditioner is build up for one matrix once, but is used for subsequent steps in a nonlinear process as well, where the matrix changes in each step slightly.
If conserve_mem
is false
, then the inverses of the local systems are computed now; if the flag is true
, then they are computed every time the preconditioner is applied. This saves some memory, but makes preconditioning very slow. Note also, that if the flag is false
, then the contents of the matrix M
at the time of calling this constructor are used, while if the flag is true
, then the values in M
at the time of preconditioning are used. This may lead to different results, obviously, of M
changes.
The parameter n_threads
determines how many threads shall be used in parallel when building the inverses of the diagonal blocks. This parameter is ignored if not in multithreaded mode.
SparseVanka< number >::~SparseVanka | ( | ) |
Destructor. Delete all allocated matrices.
SparseVanka< number >::AdditionalData::AdditionalData | ( | const std::vector< bool > & | selected, |
const bool | conserve_memory = false , |
||
const unsigned int | n_threads = MultithreadInfo::n_threads() |
||
) |
Constructor. For the parameters' description, see below.
void SparseVanka< number >::initialize | ( | const SparseMatrix< number > & | M, |
const AdditionalData & | additional_data | ||
) |
If the default constructor is used then this function needs to be called before an object of this class is used as preconditioner.
For more detail about possible parameters, see the class documentation and the documentation of the SparseVanka::AdditionalData class.
After this function is called the preconditioner is ready to be used (using the vmult
function of derived classes).
template void SparseVanka< number >::vmult< double > | ( | Vector< number2 > & | dst, |
const Vector< number2 > & | src | ||
) | const |
Do the preconditioning. This function takes the residual in src
and returns the resulting update vector in dst
.
void SparseVanka< number >::Tvmult | ( | Vector< number2 > & | dst, |
const Vector< number2 > & | src | ||
) | const |
Apply transpose preconditioner. This function takes the residual in src
and returns the resulting update vector in dst
.
size_type SparseVanka< number >::m | ( | ) | const |
Return the dimension of the codomain (or range) space. Note that the matrix is of dimension \(m \times n\).
size_type SparseVanka< number >::n | ( | ) | const |
Return the dimension of the domain space. Note that the matrix is of dimension \(m \times n\).
|
protected |
Apply the inverses corresponding to those degrees of freedom that have a true
value in dof_mask
, to the src
vector and move the result into dst
. Actually, only values for allowed indices are written to dst
, so the application of this function only does what is announced in the general documentation if the given mask sets all values to zero
The reason for providing the mask anyway is that in derived classes we may want to apply the preconditioner to parts of the matrix only, in order to parallelize the application. Then, it is important to only write to some slices of dst
, in order to eliminate the dependencies of threads of each other.
If a null pointer is passed instead of a pointer to the dof_mask
(as is the default value), then it is assumed that we shall work on all degrees of freedom. This is then equivalent to calling the function with a vector<bool>(n_dofs,true)
.
The vmult
of this class of course calls this function with a null pointer
|
protected |
Determine an estimate for the memory consumption (in bytes) of this object.
|
private |
Compute the inverses of all selected diagonal elements.
|
private |
Compute the inverses at positions in the range [begin,end)
. In non-multithreaded mode, compute_inverses()
calls this function with the whole range, but in multithreaded mode, several copies of this function are spawned.
|
private |
Compute the inverse of the block located at position row
. Since the vector is used quite often, it is generated only once in the caller of this function and passed to this function which first clears it. Reusing the vector makes the process significantly faster than in the case where this function re-creates it each time.
SparseBlockVanka< number >::SparseBlockVanka | ( | const SparseMatrix< number > & | M, |
const std::vector< bool > & | selected, | ||
const unsigned int | n_blocks, | ||
const BlockingStrategy | blocking_strategy, | ||
const bool | conserve_memory = false , |
||
const unsigned int | n_threads = MultithreadInfo::n_threads() |
||
) |
Constructor. Pass all arguments except for n_blocks
to the base class.
template void SparseBlockVanka< number >::vmult< double > | ( | Vector< number2 > & | dst, |
const Vector< number2 > & | src | ||
) | const |
Apply the preconditioner.
std::size_t SparseBlockVanka< number >::memory_consumption | ( | ) | const |
Determine an estimate for the memory consumption (in bytes) of this object.
|
private |
Compute the contents of the field dof_masks
. This function is called from the constructor.
|
private |
The dimension of the range space.
Definition at line 176 of file precondition.h.
|
private |
The dimension of the domain space.
Definition at line 181 of file precondition.h.
double PreconditionRichardson::AdditionalData::relaxation |
Relaxation parameter.
Definition at line 222 of file precondition.h.
|
private |
The relaxation parameter multiplied with the vectors.
Definition at line 307 of file precondition.h.
|
private |
The dimension of the range space.
Definition at line 312 of file precondition.h.
|
private |
The dimension of the domain space.
Definition at line 317 of file precondition.h.
|
private |
Pointer to the matrix in use.
Definition at line 392 of file precondition.h.
|
private |
Pointer to the preconditioning function.
Definition at line 397 of file precondition.h.
double PreconditionRelaxation< MatrixType >::AdditionalData::relaxation |
Relaxation parameter.
Definition at line 433 of file precondition.h.
|
protected |
Pointer to the matrix object.
Definition at line 469 of file precondition.h.
|
protected |
Relaxation parameter.
Definition at line 474 of file precondition.h.
|
private |
An array that stores for each matrix row where the first position after the diagonal is located.
Definition at line 730 of file precondition.h.
const std::vector<size_type>& PreconditionPSOR< MatrixType >::AdditionalData::permutation |
Storage for the permutation vector.
Definition at line 801 of file precondition.h.
const std::vector<size_type>& PreconditionPSOR< MatrixType >::AdditionalData::inverse_permutation |
Storage for the inverse permutation vector.
Definition at line 805 of file precondition.h.
PreconditionRelaxation<MatrixType>::AdditionalData PreconditionPSOR< MatrixType >::AdditionalData::parameters |
Relaxation parameters
Definition at line 809 of file precondition.h.
|
private |
Storage for the permutation vector.
Definition at line 862 of file precondition.h.
|
private |
Storage for the inverse permutation vector.
Definition at line 866 of file precondition.h.
unsigned int PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::degree |
This determines the degree of the Chebyshev polynomial. The degree of the polynomial gives the number of matrix-vector products to be performed for one application of the vmult() operation. Degree one corresponds to a damped Jacobi method.
If the degree is set to numbers::invalid_unsigned_int, the algorithm will automatically determine the number of necessary iterations based on the usual Chebyshev error formula as mentioned in the discussion of the main class.
Definition at line 1052 of file precondition.h.
double PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::smoothing_range |
This sets the range between the largest eigenvalue in the matrix and the smallest eigenvalue to be treated. If the parameter is set to a number less than 1, an estimate for the largest and for the smallest eigenvalue will be calculated internally. For a smoothing range larger than one, the Chebyshev polynomial will act in the interval \([\lambda_\mathrm{max}/ \tt{smoothing\_range}, \lambda_\mathrm{max}]\), where \(\lambda_\mathrm{max}\) is an estimate of the maximum eigenvalue of the matrix. A choice of smoothing_range
between 5 and 20 is useful in case the preconditioner is used as a smoother in multigrid.
Definition at line 1065 of file precondition.h.
unsigned int PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::eig_cg_n_iterations |
Maximum number of CG iterations performed for finding the maximum eigenvalue. If set to zero, no computations are performed. Instead, the user must supply a largest eigenvalue via the variable PreconditionChebyshev::AdditionalData::max_eigenvalue.
Definition at line 1073 of file precondition.h.
double PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::eig_cg_residual |
Tolerance for CG iterations performed for finding the maximum eigenvalue.
Definition at line 1079 of file precondition.h.
double PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::max_eigenvalue |
Maximum eigenvalue to work with. Only in effect if eig_cg_n_iterations
is set to zero, otherwise this parameter is ignored.
Definition at line 1086 of file precondition.h.
AffineConstraints<double> PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::constraints |
Constraints to be used for the operator given. This variable is used to zero out the correct entries when creating an initial guess.
Definition at line 1092 of file precondition.h.
std::shared_ptr<PreconditionerType> PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::preconditioner |
Stores the preconditioner object that the Chebyshev is wrapped around.
Definition at line 1097 of file precondition.h.
double PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::EigenvalueInformation::min_eigenvalue_estimate |
Estimate for the smallest eigenvalue.
Definition at line 1176 of file precondition.h.
double PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::EigenvalueInformation::max_eigenvalue_estimate |
Estimate for the largest eigenvalue.
Definition at line 1180 of file precondition.h.
unsigned int PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::EigenvalueInformation::cg_iterations |
Number of CG iterations performed or 0.
Definition at line 1184 of file precondition.h.
unsigned int PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::EigenvalueInformation::degree |
The degree of the Chebyshev polynomial (either as set using AdditionalData::degree or estimated as described there).
Definition at line 1189 of file precondition.h.
|
private |
A pointer to the underlying matrix.
Definition at line 1223 of file precondition.h.
|
mutableprivate |
Internal vector used for the vmult
operation.
Definition at line 1228 of file precondition.h.
|
mutableprivate |
Internal vector used for the vmult
operation.
Definition at line 1233 of file precondition.h.
|
mutableprivate |
Internal vector used for the vmult
operation.
Definition at line 1238 of file precondition.h.
|
private |
Stores the additional data passed to the initialize function, obtained through a copy operation.
Definition at line 1244 of file precondition.h.
|
private |
Average of the largest and smallest eigenvalue under consideration.
Definition at line 1249 of file precondition.h.
|
private |
Half the interval length between the largest and smallest eigenvalue under consideration.
Definition at line 1255 of file precondition.h.
|
private |
Stores whether the preconditioner has been set up and eigenvalues have been computed.
Definition at line 1261 of file precondition.h.
|
mutableprivate |
A mutex to avoid that multiple vmult() invocations by different threads overwrite the temporary vectors.
Definition at line 1267 of file precondition.h.
const std::vector<bool>& SparseVanka< number >::AdditionalData::selected |
Indices of those degrees of freedom that we shall work on.
Definition at line 207 of file sparse_vanka.h.
const bool SparseVanka< number >::AdditionalData::conserve_mem |
Conserve memory flag.
Definition at line 212 of file sparse_vanka.h.
const unsigned int SparseVanka< number >::AdditionalData::n_threads |
Number of threads to be used when building the inverses. Only relevant in multithreaded mode.
Definition at line 218 of file sparse_vanka.h.
|
private |
Pointer to the matrix.
Definition at line 311 of file sparse_vanka.h.
|
private |
Conserve memory flag.
Definition at line 316 of file sparse_vanka.h.
|
private |
Indices of those degrees of freedom that we shall work on.
Definition at line 321 of file sparse_vanka.h.
|
private |
Number of threads to be used when building the inverses. Only relevant in multithreaded mode.
Definition at line 327 of file sparse_vanka.h.
|
mutableprivate |
Array of inverse matrices, one for each degree of freedom. Only those elements will be used that are tagged in selected
.
Definition at line 334 of file sparse_vanka.h.
|
private |
The dimension of the range space.
Definition at line 339 of file sparse_vanka.h.
|
private |
The dimension of the domain space.
Definition at line 344 of file sparse_vanka.h.
|
private |
Store the number of blocks.
Definition at line 573 of file sparse_vanka.h.
|
private |
In this field, we precompute for each block which degrees of freedom belong to it. Thus, if dof_masks[i][j]==true
, then DoF j
belongs to block i
. Of course, no other dof_masks[l][j]
may be true
for l!=i
. This computation is done in the constructor, to avoid recomputing each time the preconditioner is called.
Definition at line 582 of file sparse_vanka.h.
Definition at line 382 of file sparse_vanka.h.