dune-istl 2.8.0
|
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. | |
SeqILDL (const matrix_type &A, const ParameterTree &config) | |
Constructor. | |
SeqILDL (const matrix_type &A, real_field_type relax=real_field_type(1)) | |
constructor | |
void | pre (X &x, Y &b) override |
Prepare the preconditioner. | |
void | apply (X &v, const Y &d) override |
Apply one step of the preconditioner to the system A(v)=d. | |
void | post (X &x) override |
Clean up. | |
SolverCategory::Category | category () const override |
Category of the preconditioner (see SolverCategory::Category) | |
sequential ILDL preconditioner
Wraps the naked ISTL generic ILDL preconditioner into the solver framework.
M | type of matrix to operate on |
X | type of update |
Y | type of defect |
typedef X Dune::SeqILDL< M, X, Y >::domain_type |
domain type of the preconditioner
typedef X::field_type Dune::SeqILDL< M, X, Y >::field_type |
field type of the preconditioner
typedef std::remove_const_t< M > Dune::SeqILDL< M, X, Y >::matrix_type |
type of matrix the preconditioner is for
typedef Y Dune::SeqILDL< M, X, Y >::range_type |
range type of the preconditioner
typedef FieldTraits<scalar_field_type>::real_type Dune::SeqILDL< M, X, Y >::real_field_type |
real scalar type underlying the field_type
typedef Simd::Scalar<field_type> Dune::SeqILDL< M, X, Y >::scalar_field_type |
scalar type underlying the field_type
|
inline |
Constructor.
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.
|
inline |
Constructor.
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.
|
inlineexplicit |
constructor
The constructor copies the matrix A and computes its ILDL decomposition.
[in] | A | matrix to operate on |
[in] | relax | relaxation factor |
|
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
[out] | v | The update to be computed |
d | The current defect. |
Implements Dune::Preconditioner< X, Y >.
|
inlineoverridevirtual |
Category of the preconditioner (see SolverCategory::Category)
Implements Dune::Preconditioner< X, Y >.
|
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.
x | The right hand side of the equation. |
Implements Dune::Preconditioner< X, Y >.
|
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.
x | The left hand side of the equation. |
b | The right hand side of the equation. |
Implements Dune::Preconditioner< X, Y >.