Dune Core Modules (2.10.0)

Sequential DILU preconditioner. More...

#include <dune/istl/preconditioners.hh>

Public Types

using matrix_type = M
 The matrix type the preconditioner is for.
 
using block_type = typename matrix_type::block_type
 block type of matrix
 
using domain_type = X
 The domain type of the preconditioner.
 
using range_type = Y
 The range type of the preconditioner.
 
using field_type = typename X::field_type
 The field type of the preconditioner.
 
using scalar_field_type = Simd::Scalar< field_type >
 scalar type underlying the field_type
 
using real_field_type = typename FieldTraits< scalar_field_type >::real_type
 real scalar type underlying the field_type
 

Public Member Functions

 SeqDILU (const M &A, real_field_type w)
 Constructor. More...
 
 SeqDILU (const std::shared_ptr< const AssembledLinearOperator< M, X, Y > > &A, const ParameterTree &configuration)
 Constructor. More...
 
 SeqDILU (const M &A, const ParameterTree &config)
 Constructor. More...
 
virtual void pre (X &x, Y &b)
 Prepare the preconditioner. More...
 
virtual void apply (X &v, const Y &d)
 Apply the preconditioner. More...
 
virtual void post (X &x)
 Clean up. More...
 
virtual SolverCategory::Category category () const
 Category of the preconditioner (see SolverCategory::Category)
 

Protected Attributes

const M & _A_
 The matrix we operate on.
 
const real_field_type _w
 The relaxation factor to use.
 
const bool wNotIdentity_
 true if w != 1.0
 

Detailed Description

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

Sequential DILU preconditioner.

Wraps the naked ISTL generic DILU preconditioner into the solver framework.

The Diagonal Incomplete LU factorization (DILU) is a simplified variant of the ILU(0) factorization, where only the diagonal elements are factorised. This preconditioner can be written as

M = (D + L_A) D^{-1} (D + U_A),

where L_A and U_A are the strictly lower and upper parts of A and D is the diagonal matrix containing the generated pivot values. The matrix M has the property

diag(A) = diag(M) = diag((D + L_A) D^{-1} (D + U_A)) = diag(D + L_A D^{-1} U_A)

such that the diagonal matrix D can be constructed:

D_11 = A_11 D_22 = A22 - L_A_{21} D_{11}^{-1} U_A_{12} and etc.

Note that when applying the preconditioner, M is never explicitly created but instead a lower and upper triangualr solve is performed using the values of A and D^-1. When applying the preconditioner, we are working with the residual d = b - Ax and update v = x_{n+1} - x_n, such that:

v = M^{-1} d = (D + U_A)^{-1} D (D + L_A)^{-1} d

define y = (D + L_A)^{-1} d

lower triangular solve: (D + L_A) y = d upper triangular solve: (D + U_A) v = D y

This means that the preconditioner only requires an additional storage of the diagonal matrix D^{-1}

For more details, see: R. Barrett et al., "Templates for the Solution of Linear Systems: Building Blocks for Iterative Methods", 1994. Available from: https://www.netlib.org/templates/templates.pdf

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

◆ SeqDILU() [1/3]

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

Constructor.

Constructor invoking DILU gets all parameters to operate the prec.

Parameters
AThe matrix to operate on.
wThe relaxation factor.

◆ SeqDILU() [2/3]

template<class M , class X , class Y , int l = 1>
Dune::SeqDILU< M, X, Y, l >::SeqDILU ( const std::shared_ptr< const AssembledLinearOperator< M, X, Y > > &  A,
const ParameterTree configuration 
)
inline

Constructor.

Parameters
AThe assembled linear operator to use.
configurationParameterTree containing preconditioner parameters.
ParameterTree Key Meaning
relaxation The relaxation factor. default=1.0

See ISTL_Factory for the ParameterTree layout and examples.

◆ SeqDILU() [3/3]

template<class M , class X , class Y , int l = 1>
Dune::SeqDILU< M, X, Y, l >::SeqDILU ( const M &  A,
const ParameterTree config 
)
inline

Constructor.

Parameters
AThe matrix to operate on.
configParameterTree containing preconditioner parameters.
ParameterTree Key Meaning
relaxation The relaxation factor. default=1.0

See ISTL_Factory for the ParameterTree layout and examples.

Member Function Documentation

◆ apply()

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

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

Implements Dune::Preconditioner< X, Y >.

◆ post()

template<class M , class X , class Y , int l = 1>
virtual void Dune::SeqDILU< 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 >.

◆ pre()

template<class M , class X , class Y , int l = 1>
virtual void Dune::SeqDILU< 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.

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.
X x(0.0);
Y b = ...; // rhs
prec.pre(x,b); // prepare the preconditioner
prec.apply(x,b); // can be called multiple times now...
prec.post(x); // cleanup internal state
Base class for matrix free definition of preconditioners.
Definition: preconditioner.hh:33
Parameters
xThe left hand side of the equation.
bThe right hand side of the equation.

Implements Dune::Preconditioner< X, Y >.


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 (Dec 26, 23:30, 2024)