Dune Core Modules (2.9.0)

Sequential SOR preconditioner. More...

#include <dune/istl/preconditioners.hh>

Public Types

typedef M 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 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

 SeqSOR (const M &A, int n, real_field_type w)
 Constructor. More...
 
 SeqSOR (const std::shared_ptr< const AssembledLinearOperator< M, X, Y >> &A, const ParameterTree &configuration)
 Constructor. More...
 
 SeqSOR (const M &A, const ParameterTree &configuration)
 Constructor. More...
 
virtual void pre ([[maybe_unused]] X &x, [[maybe_unused]] Y &b)
 Prepare the preconditioner. More...
 
virtual void apply (X &v, const Y &d)
 Apply the preconditioner. More...
 
template<bool forward>
void apply (X &v, const Y &d)
 Apply the preconditioner in a special direction. More...
 
virtual void post ([[maybe_unused]] X &x)
 Clean up. More...
 
virtual SolverCategory::Category category () const
 Category of the preconditioner (see SolverCategory::Category)
 
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, int l = 1>
class Dune::SeqSOR< M, X, Y, l >

Sequential SOR preconditioner.

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

Template Parameters
MThe matrix type to operate on
XType of the update
YType of the defect
lThe block level to invert. Default is 1

Constructor & Destructor Documentation

◆ SeqSOR() [1/3]

template<class M , class X , class Y , int l = 1>
Dune::SeqSOR< M, X, Y, l >::SeqSOR ( const M &  A,
int  n,
real_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.

◆ SeqSOR() [2/3]

template<class M , class X , class Y , int l = 1>
Dune::SeqSOR< M, X, Y, l >::SeqSOR ( 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
iterations The number of iterations to perform. default=1
relaxation The relaxation factor. default=1.0

See ISTL_Factory for the ParameterTree layout and examples.

◆ SeqSOR() [3/3]

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

Constructor.

Parameters
AThe matrix to operate on.
configurationParameterTree containing preconditioner parameters.
ParameterTree Key Meaning
iterations The number of iterations to perform. default=1
relaxation The relaxation factor. default=1.0

See ISTL_Factory for the ParameterTree layout and examples.

Member Function Documentation

◆ apply() [1/2]

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

◆ apply() [2/2]

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

Apply the preconditioner in a special direction.

The template parameter forward indications the direction the smoother is applied. If true The application is started at the lowest index in the vector v, if false at the highest index of vector v.

Implements Dune::Preconditioner< X, Y >.

References Dune::bsorb(), and Dune::bsorf().

◆ post() [1/2]

template<class M , class X , class Y , int l = 1>
virtual void Dune::SeqSOR< M, X, Y, l >::post ( [[maybe_unused] ] 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.

◆ 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 , int l = 1>
virtual void Dune::SeqSOR< M, X, Y, l >::pre ( [[maybe_unused] ] X &  x,
[[maybe_unused] ] 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 mangement.
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 mangement.
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 >, 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 >.


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 (Mar 27, 23:31, 2024)