Dune Core Modules (unstable)

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_typescalar_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 ([[maybe_unused]] X &x, [[maybe_unused]] 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 ([[maybe_unused]] X &x) override
 Clean up. More...
 
SolverCategory::Category category () const override
 Category of the preconditioner (see SolverCategory::Category) More...
 
virtual void pre (X &x, Y &b)=0
 Prepare the preconditioner. More...
 
virtual void post (X &x)=0
 Clean up. More...
 

Detailed Description

template<class M, class X, class Y>
class Dune::SeqILDL< M, X, Y >

sequential ILDL preconditioner

Author
Martin Nolte

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

Template Parameters
Mtype of matrix to operate on
Xtype of update
Ytype of defect

Constructor & Destructor Documentation

◆ SeqILDL() [1/3]

template<class M , class X , class Y >
Dune::SeqILDL< M, X, Y >::SeqILDL ( const std::shared_ptr< const AssembledLinearOperator< M, X, Y >> &  A,
const ParameterTree configuration 
)
inline

Constructor.

Parameters
AThe linear operator to use.
configurationParameterTree containing preconditioner parameters.
ParameterTree Key Meaning
relaxation relaxation factor

See ISTL_Factory for the ParameterTree layout and examples.

◆ SeqILDL() [2/3]

template<class M , class X , class Y >
Dune::SeqILDL< M, X, Y >::SeqILDL ( const matrix_type A,
const ParameterTree config 
)
inline

Constructor.

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

See ISTL_Factory for the ParameterTree layout and examples.

◆ SeqILDL() [3/3]

template<class M , class X , class Y >
Dune::SeqILDL< M, X, Y >::SeqILDL ( const matrix_type A,
real_field_type  relax = real_field_type( 1 ) 
)
inlineexplicit

constructor

The constructor copies the matrix A and computes its ILDL decomposition.

Parameters
[in]Amatrix to operate on
[in]relaxrelaxation factor

Member Function Documentation

◆ apply()

template<class M , class X , class Y >
void Dune::SeqILDL< M, X, Y >::apply ( X &  v,
const Y &  d 
)
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]vThe update to be computed
dThe current defect.

Implements Dune::Preconditioner< X, Y >.

◆ category()

template<class M , class X , class Y >
SolverCategory::Category Dune::SeqILDL< M, X, Y >::category ( ) const
inlineoverridevirtual

Category of the preconditioner (see SolverCategory::Category)

Implements Dune::Preconditioner< X, Y >.

References Dune::SolverCategory::sequential.

◆ post() [1/2]

template<class M , class X , class Y >
void Dune::SeqILDL< M, X, Y >::post ( [[maybe_unused] ] X &  x)
inlineoverride

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.

◆ post() [2/2]

template<class X , class Y >
virtual void Dune::Preconditioner< X, Y >::post ( X &  x)
pure virtualinherited

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.

Implemented in Dune::BlockPreconditioner< X, Y, C, P >, Dune::NonoverlappingBlockPreconditioner< C, P >, Dune::InverseOperator2Preconditioner< O, c >, Dune::Amg::KAMG< M, X, S, PI, K, A >, Dune::Amg::FastAMG< M, X, PI, A >, Dune::Amg::AMG< M, X, S, PI, A >, Dune::Amg::AMG< M, X, S, SequentialInformation, std::allocator< X > >, and Dune::Amg::AMG< Operator, X, Smoother >.

◆ pre() [1/2]

template<class M , class X , class Y >
void Dune::SeqILDL< M, X, Y >::pre ( [[maybe_unused] ] X &  x,
[[maybe_unused] ] Y &  b 
)
inlineoverride

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
Preconditioner<X,Y> prec(...);
prec.pre(x,b); // prepare the preconditioner
prec.apply(x,b); // can be called multiple times now...
prec.post(x); // cleanup internal state
Parameters
xThe left hand side of the equation.
bThe right hand side of the equation.

◆ pre() [2/2]

template<class X , class Y >
virtual void Dune::Preconditioner< X, Y >::pre ( X &  x,
Y &  b 
)
pure virtualinherited

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
Preconditioner<X,Y> prec(...);
prec.pre(x,b); // prepare the preconditioner
prec.apply(x,b); // can be called multiple times now...
prec.post(x); // cleanup internal state
Parameters
xThe left hand side of the equation.
bThe right hand side of the equation.

Implemented in Dune::BlockPreconditioner< X, Y, C, P >, Dune::NonoverlappingBlockPreconditioner< C, P >, and Dune::InverseOperator2Preconditioner< O, c >.


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.80.0 (Apr 17, 22:29, 2024)