Dune Core Modules (2.6.0)

Sequential overlapping Schwarz preconditioner. More...

#include <dune/istl/overlappingschwarz.hh>

Public Types

typedef M matrix_type
 The type of the matrix to precondition.
 
typedef X domain_type
 The domain type of the preconditioner.
 
typedef X range_type
 The range type of the preconditioner.
 
typedef TM Mode
 The mode (additive or multiplicative) of the Schwarz method. More...
 
typedef X::field_type field_type
 The field type of the preconditioner.
 
typedef matrix_type::size_type size_type
 The return type of the size method.
 
typedef TA allocator
 The allocator to use.
 
typedef std::set< size_type, std::less< size_type >, typename TA::template rebind< size_type >::other > subdomain_type
 The type for the subdomain to row index mapping.
 
typedef std::vector< subdomain_type, typename TA::template rebind< subdomain_type >::other > subdomain_vector
 The vector type containing the subdomain to row index mapping.
 
typedef SLList< size_type, typename TA::template rebind< size_type >::other > subdomain_list
 The type for the row to subdomain mapping.
 
typedef std::vector< subdomain_list, typename TA::template rebind< subdomain_list >::other > rowtodomain_vector
 The vector type containing the row index to subdomain mapping.
 
typedef TD slu
 The type for the subdomain solver in use.
 
typedef std::vector< slu, typename TA::template rebind< slu >::other > slu_vector
 The vector type containing subdomain solvers.
 

Public Member Functions

 SeqOverlappingSchwarz (const matrix_type &mat, const subdomain_vector &subDomains, field_type relaxationFactor=1, bool onTheFly_=true)
 Construct the overlapping Schwarz method. More...
 
 SeqOverlappingSchwarz (const matrix_type &mat, const rowtodomain_vector &rowToDomain, field_type relaxationFactor=1, bool onTheFly_=true)
 
virtual void pre (X &x, X &b)
 Prepare the preconditioner. More...
 
virtual void apply (X &v, const X &d)
 Apply the precondtioner. More...
 
virtual void post (X &x)
 Postprocess the preconditioner. More...
 
template<bool forward>
void apply (X &v, const X &d)
 Apply one step of the preconditioner to the system A(v)=d. More...
 
virtual SolverCategory::Category category () const
 Category of the preconditioner (see SolverCategory::Category)
 

Detailed Description

template<class M, class X, class TM = AdditiveSchwarzMode, class TD = ILU0SubdomainSolver<M,X,X>, class TA = std::allocator<X>>
class Dune::SeqOverlappingSchwarz< M, X, TM, TD, TA >

Sequential overlapping Schwarz preconditioner.

Template Parameters
MThe matrix type.
XThe range and domain type.
TMThe Schwarz mode. Currently supported modes are AdditiveSchwarzMode, MultiplicativeSchwarzMode, and SymmetricMultiplicativeSchwarzMode. (Default values is AdditiveSchwarzMode)
TDThe type of the local subdomain solver to be used.
TAThe type of the allocator to use.

Member Typedef Documentation

◆ Mode

template<class M , class X , class TM = AdditiveSchwarzMode, class TD = ILU0SubdomainSolver<M,X,X>, class TA = std::allocator<X>>
typedef TM Dune::SeqOverlappingSchwarz< M, X, TM, TD, TA >::Mode

The mode (additive or multiplicative) of the Schwarz method.

Either AdditiveSchwarzMode or MultiplicativeSchwarzMode

Member Function Documentation

◆ apply()

template<class M , class X , class TM = AdditiveSchwarzMode, class TD = ILU0SubdomainSolver<M,X,X>, class TA = std::allocator<X>>
template<bool forward>
void Dune::SeqOverlappingSchwarz< M, X, TM, TD, TA >::apply ( X &  v,
const X &  d 
)
virtual

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, X >.

◆ post()

template<class M , class X , class TM = AdditiveSchwarzMode, class TD = ILU0SubdomainSolver<M,X,X>, class TA = std::allocator<X>>
virtual void Dune::SeqOverlappingSchwarz< M, X, TM, TD, TA >::post ( X &  x)
inlinevirtual

Postprocess the preconditioner.

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, X >.

References DUNE_UNUSED_PARAMETER.

◆ pre()

template<class M , class X , class TM = AdditiveSchwarzMode, class TD = ILU0SubdomainSolver<M,X,X>, class TA = std::allocator<X>>
virtual void Dune::SeqOverlappingSchwarz< M, X, TM, TD, TA >::pre ( X &  x,
X &  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, X >.

References DUNE_UNUSED_PARAMETER.


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