DUNE-FEM (unstable)
Sequential ILU preconditioner. More...
#include <dune/istl/preconditioners.hh>
Public Types | |
typedef std::remove_const< M >::type | matrix_type |
The matrix type the preconditioner is for. | |
typedef matrix_type::block_type | block_type |
block type of matrix | |
typedef X | domain_type |
The domain type of the preconditioner. | |
typedef Y | range_type |
The range type of the preconditioner. | |
typedef X::field_type | field_type |
The 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 | |
typedef ILU::CRS< block_type, typename M::allocator_type > | CRS |
type of ILU storage | |
Public Member Functions | |
SeqILU (const M &A, real_field_type w, const bool resort=false) | |
Constructor. More... | |
SeqILU (const std::shared_ptr< const AssembledLinearOperator< M, X, Y > > &A, const ParameterTree &configuration) | |
Constructor. More... | |
SeqILU (const M &A, const ParameterTree &config) | |
Constructor. More... | |
SeqILU (const M &A, int n, real_field_type w, const bool resort=false) | |
Constructor. More... | |
void | pre (X &x, Y &b) override |
Prepare the preconditioner. More... | |
void | apply (X &v, const Y &d) override |
Apply the preconditioner. More... | |
void | post (X &x) override |
Clean up. More... | |
SolverCategory::Category | category () const override |
Category of the preconditioner (see SolverCategory::Category) | |
Protected Attributes | |
std::unique_ptr< matrix_type > | ILU_ |
The ILU(n) decomposition of the matrix. As storage a BCRSMatrix is used. | |
CRS | lower_ |
The ILU(n) decomposition of the matrix. As storage a CRS structure is used. | |
const real_field_type | w_ |
The relaxation factor to use. | |
const bool | wNotIdentity_ |
true if w != 1.0 | |
Detailed Description
class Dune::SeqILU< M, X, Y, l >
Sequential ILU preconditioner.
Wraps the naked ISTL generic ILU preconditioner into the solver framework.
- Template Parameters
-
M The matrix type to operate on X Type of the update Y Type of the defect l Ignored. Just there to have the same number of template arguments as other preconditioners.
Constructor & Destructor Documentation
◆ SeqILU() [1/4]
|
inline |
Constructor.
Constructor invoking ILU(0) gets all parameters to operate the prec.
- Parameters
-
A The matrix to operate on. w The relaxation factor. resort true if a resort of the computed ILU for improved performance should be done.
◆ SeqILU() [2/4]
|
inline |
Constructor.
- Parameters
-
A The assembled linear operator to use. configuration ParameterTree containing preconditioner parameters.
ParameterTree Key | Meaning |
---|---|
n | The order of the ILU decomposition. default=0 |
relaxation | The relaxation factor. default=1.0 |
resort | True if a resort of the computed ILU for improved performance should be done. default=false |
See ISTL_Factory for the ParameterTree layout and examples.
◆ SeqILU() [3/4]
|
inline |
Constructor.
- Parameters
-
A The matrix to operate on. config ParameterTree containing preconditioner parameters.
ParameterTree Key | Meaning |
---|---|
n | The order of the ILU decomposition. default=0 |
relaxation | The relaxation factor. default=1.0 |
resort | True if a resort of the computed ILU for improved performance should be done. default=false |
See ISTL_Factory for the ParameterTree layout and examples.
◆ SeqILU() [4/4]
|
inline |
Constructor.
Constructor invoking ILU(n).
- Parameters
-
A The matrix to operate on. n The order of the ILU decomposition. w The relaxation factor. resort true if a resort of the computed ILU for improved performance should be done.
Member Function Documentation
◆ apply()
|
inlineoverridevirtual |
Apply the preconditioner.
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 >.
References Dune::SeqILU< M, X, Y, l >::ILU_, Dune::SeqILU< M, X, Y, l >::lower_, Dune::SeqILU< M, X, Y, l >::w_, and Dune::SeqILU< M, X, Y, l >::wNotIdentity_.
◆ post()
|
inlineoverridevirtual |
Clean up.
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.
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 management.
- 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