Dune Core Modules (2.6.0)

Sequential ILU(n) 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 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 SimdScalar< field_typescalar_field_type
 scalar type underlying the field_type
 

Public Member Functions

 SeqILUn (const M &A, int n, scalar_field_type w)
 Constructor. More...
 
virtual void pre (X &x, Y &b)
 Prepare the preconditioner. More...
 
virtual void apply (X &v, const Y &d)
 Apply the precondioner. More...
 
virtual void post (X &x)
 Clean up. More...
 
virtual SolverCategory::Category category () const
 Category of the preconditioner (see SolverCategory::Category)
 

Detailed Description

template<class M, class X, class Y, int l = 1>
class Dune::SeqILUn< M, X, Y, l >

Sequential ILU(n) preconditioner.

Wraps the naked ISTL generic ILU(n) preconditioner into the solver framework.

Template Parameters
MThe matrix type to operate on
XType of the update
YType of the defect
lIgnored. Just there to have the same number of template arguments as other preconditioners.

Constructor & Destructor Documentation

◆ SeqILUn()

template<class M , class X , class Y , int l = 1>
Dune::SeqILUn< M, X, Y, l >::SeqILUn ( const M &  A,
int  n,
scalar_field_type  w 
)
inline

Constructor.

Constructor gets all parameters to operate the prec.

Parameters
AThe matrix to operate on.
nThe number of iterations to perform.
wThe relaxation factor.

Member Function Documentation

◆ apply()

template<class M , class X , class Y , int l = 1>
virtual void Dune::SeqILUn< M, X, Y, l >::apply ( X &  v,
const Y &  d 
)
inlinevirtual

Apply the precondioner.

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]vThe update to be computed
dThe current defect.

Implements Dune::Preconditioner< X, Y >.

References Dune::bilu_backsolve().

◆ post()

template<class M , class X , class Y , int l = 1>
virtual void Dune::SeqILUn< M, X, Y, l >::post ( X &  x)
inlinevirtual

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
xThe right hand side of the equation.

Implements Dune::Preconditioner< X, Y >.

References DUNE_UNUSED_PARAMETER.

◆ pre()

template<class M , class X , class Y , int l = 1>
virtual void Dune::SeqILUn< M, X, Y, l >::pre ( X &  x,
Y &  b 
)
inlinevirtual

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.

Parameters
xThe left hand side of the equation.
bThe right hand side of the equation.

Implements Dune::Preconditioner< X, Y >.

References DUNE_UNUSED_PARAMETER.


The documentation for this class was generated from the following file:
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 24, 23:30, 2024)