Dune::SeqOverlappingSchwarz< M, X, TM, TA > Class Template Reference
[ISTL]
#include <overlappingschwarz.hh>

Detailed Description
template<class M, class X, class TM = AdditiveSchwarzMode, class TA = std::allocator<X>>
class Dune::SeqOverlappingSchwarz< M, X, TM, TA >
Sequential overlapping Schwarz preconditioner. Public Types | |
enum | { category = SolverCategory::sequential } |
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. | |
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 < std::less< 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 SuperLU< matrix_type > | slu |
The type for the SuperLU solver in use. | |
typedef std::vector< slu, typename TA::template rebind < slu >::other > | slu_vector |
The vector type containing SuperLU solvers. | |
Public Member Functions | |
SeqOverlappingSchwarz (const matrix_type &mat, const subdomain_vector &subDomains, field_type relaxationFactor=1) | |
Construct the overlapping Schwarz method. | |
SeqOverlappingSchwarz (const matrix_type &mat, const rowtodomain_vector &rowToDomain, field_type relaxationFactor=1) | |
virtual void | pre (X &x, X &b) |
Prepare the preconditioner. | |
virtual void | apply (X &v, const X &d) |
Apply the precondtioner. | |
template<bool forward> | |
void | apply (X &v, const X &d) |
Apply one step of the preconditioner to the system A(v)=d. | |
virtual void | post (X &x) |
Clean up. |
Member Typedef Documentation
typedef TM Dune::SeqOverlappingSchwarz< M, X, TM, TA >::Mode |
The mode (additive or multiplicative) of the Schwarz method.
Either AdditiveSchwarzMode or MultiplicativeSchwarzMode
Member Enumeration Documentation
anonymous enum |
Constructor & Destructor Documentation
Dune::SeqOverlappingSchwarz< M, X, TM, TA >::SeqOverlappingSchwarz | ( | const matrix_type & | mat, | |
const rowtodomain_vector & | rowToDomain, | |||
field_type | relaxationFactor = 1 | |||
) |
Construct the overlapping Schwarz method
- Parameters:
-
mat The matrix to precondition. rowToDomain The mapping of the rows onto the domains.
Member Function Documentation
virtual void Dune::SeqOverlappingSchwarz< M, X, TM, TA >::pre | ( | X & | x, | |
X & | b | |||
) | [inline, virtual] |
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:
-
x The left hand side of the equation. b The right hand side of the equation.
Implements Dune::Preconditioner< X, X >.
void Dune::SeqOverlappingSchwarz< M, X, TM, TA >::apply | ( | X & | v, | |
const X & | d | |||
) | [inline, 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 where
is the approximate inverse of the operator
characterizing the preconditioner.
- Parameters:
-
[out] v The update to be computed d The current defect.
Implements Dune::Preconditioner< X, X >.
virtual void Dune::SeqOverlappingSchwarz< M, X, TM, TA >::post | ( | X & | x | ) | [inline, virtual] |
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:
-
x The right hand side of the equation.
Implements Dune::Preconditioner< X, X >.
The documentation for this class was generated from the following file: