Dune Core Modules (2.9.1)
sequential ILDL preconditioner More...
#include <dune/istl/preconditioners.hh>
Public Types | |
typedef std::remove_const_t< M > | matrix_type |
type of matrix the preconditioner is for | |
typedef X | domain_type |
domain type of the preconditioner | |
typedef Y | range_type |
range type of the preconditioner | |
typedef X::field_type | field_type |
field type of the preconditioner | |
typedef Simd::Scalar< field_type > | scalar_field_type |
scalar type underlying the field_type | |
typedef FieldTraits< scalar_field_type >::real_type | real_field_type |
real scalar type underlying the field_type | |
Public Member Functions | |
SeqILDL (const std::shared_ptr< const AssembledLinearOperator< M, X, Y > > &A, const ParameterTree &configuration) | |
Constructor. More... | |
SeqILDL (const matrix_type &A, const ParameterTree &config) | |
Constructor. More... | |
SeqILDL (const matrix_type &A, real_field_type relax=real_field_type(1)) | |
constructor More... | |
void | pre (X &x, Y &b) override |
Prepare the preconditioner. More... | |
void | apply (X &v, const Y &d) override |
Apply one step of the preconditioner to the system A(v)=d. More... | |
void | post (X &x) override |
Clean up. More... | |
SolverCategory::Category | category () const override |
Category of the preconditioner (see SolverCategory::Category) More... | |
Detailed Description
class Dune::SeqILDL< M, X, Y >
sequential ILDL preconditioner
Wraps the naked ISTL generic ILDL preconditioner into the solver framework.
- Template Parameters
-
M type of matrix to operate on X type of update Y type of defect
Constructor & Destructor Documentation
◆ SeqILDL() [1/3]
|
inline |
Constructor.
- Parameters
-
A The linear operator to use. configuration ParameterTree containing preconditioner parameters.
ParameterTree Key | Meaning |
---|---|
relaxation | relaxation factor |
See ISTL_Factory for the ParameterTree layout and examples.
◆ SeqILDL() [2/3]
|
inline |
Constructor.
- Parameters
-
A The matrix to operate on. configuration ParameterTree containing preconditioner parameters.
ParameterTree Key | Meaning |
---|---|
relaxation | relaxation factor. default=1.0 |
See ISTL_Factory for the ParameterTree layout and examples.
◆ SeqILDL() [3/3]
|
inlineexplicit |
constructor
The constructor copies the matrix A and computes its ILDL decomposition.
- Parameters
-
[in] A matrix to operate on [in] relax relaxation factor
Member Function Documentation
◆ apply()
|
inlineoverridevirtual |
Apply one step of the preconditioner to the system A(v)=d.
On entry v=0 and d=b-A(x) (although this might not be computed in that way. On exit v contains the update, i.e one step computes \( v = M^{-1} d \) where \( M \) is the approximate inverse of the operator \( A \) characterizing the preconditioner.
- Parameters
-
[out] v The update to be computed d The current defect.
Implements Dune::Preconditioner< X, Y >.
◆ category()
|
inlineoverridevirtual |
Category of the preconditioner (see SolverCategory::Category)
Implements Dune::Preconditioner< X, Y >.
References Dune::SolverCategory::sequential.
◆ post()
|
inlineoverridevirtual |
Clean up.
This method is called after the last apply call for the linear system to be solved. Memory may be deallocated safely here. x is the solution of the linear equation.
- Parameters
-
x The right hand side of the equation.
Implements Dune::Preconditioner< X, Y >.
◆ pre()
|
inlineoverridevirtual |
Prepare the preconditioner.
A solver solves a linear operator equation A(x)=b by applying one or several steps of the preconditioner. The method pre() is called before the first apply operation. b and x are right hand side and solution vector of the linear system respectively. It may. e.g., scale the system, allocate memory or compute a (I)LU decomposition. Note: The ILU decomposition could also be computed in the constructor or with a separate method of the derived method if several linear systems with the same matrix are to be solved.
- Note
- if a preconditioner is copied (e.g. for a second thread) again the pre() method has to be called to ensure proper memory mangement.
- Parameters
-
x The left hand side of the equation. b The right hand side of the equation.
Implements Dune::Preconditioner< X, Y >.
The documentation for this class was generated from the following file:
- dune/istl/preconditioners.hh