Sequential overlapping Schwarz preconditioner. More...
#include <dune/istl/overlappingschwarz.hh>
Classes | |
struct | Assigner |
struct | Assigner< BlockVector< FieldVector< T, n >, A > > |
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 < 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. |
Sequential overlapping Schwarz preconditioner.
M | The matrix type. | |
X | The range and domain type. | |
TM | The Schwarz mode. Currently supported modes are AdditiveSchwarzMode, MultiplicativeSchwarzMode, and SymmetricMultiplicativeSchwarzMode. (Default values is AdditiveSchwarzMode) | |
onTheFly | If true the decomposition of the exact local solvers is computed on the fly for each subdomain and iteration step. If false all decompositions are computed in pre and only forward and backward substitution takes place in the iteration steps. | |
TA | The type of the allocator to use. |
typedef TA Dune::SeqOverlappingSchwarz< M, X, TM, onTheFly, TA >::allocator |
The allocator to use.
typedef X Dune::SeqOverlappingSchwarz< M, X, TM, onTheFly, TA >::domain_type |
The domain type of the preconditioner.
Reimplemented from Dune::Preconditioner< X, X >.
typedef X::field_type Dune::SeqOverlappingSchwarz< M, X, TM, onTheFly, TA >::field_type |
The field type of the preconditioner.
Reimplemented from Dune::Preconditioner< X, X >.
typedef M Dune::SeqOverlappingSchwarz< M, X, TM, onTheFly, TA >::matrix_type |
The type of the matrix to precondition.
typedef TM Dune::SeqOverlappingSchwarz< M, X, TM, onTheFly, TA >::Mode |
The mode (additive or multiplicative) of the Schwarz method.
Either AdditiveSchwarzMode or MultiplicativeSchwarzMode
typedef X Dune::SeqOverlappingSchwarz< M, X, TM, onTheFly, TA >::range_type |
The range type of the preconditioner.
Reimplemented from Dune::Preconditioner< X, X >.
typedef std::vector<subdomain_list, typename TA::template rebind<subdomain_list>::other > Dune::SeqOverlappingSchwarz< M, X, TM, onTheFly, TA >::rowtodomain_vector |
The vector type containing the row index to subdomain mapping.
typedef matrix_type::size_type Dune::SeqOverlappingSchwarz< M, X, TM, onTheFly, TA >::size_type |
The return type of the size method.
typedef SuperLU<matrix_type> Dune::SeqOverlappingSchwarz< M, X, TM, onTheFly, TA >::slu |
The type for the SuperLU solver in use.
typedef std::vector<slu, typename TA::template rebind<slu>::other> Dune::SeqOverlappingSchwarz< M, X, TM, onTheFly, TA >::slu_vector |
The vector type containing SuperLU solvers.
typedef SLList<size_type, typename TA::template rebind<size_type>::other> Dune::SeqOverlappingSchwarz< M, X, TM, onTheFly, TA >::subdomain_list |
The type for the row to subdomain mapping.
typedef std::set<size_type, std::less<size_type>, typename TA::template rebind<size_type>::other> Dune::SeqOverlappingSchwarz< M, X, TM, onTheFly, TA >::subdomain_type |
The type for the subdomain to row index mapping.
typedef std::vector<subdomain_type, typename TA::template rebind<subdomain_type>::other> Dune::SeqOverlappingSchwarz< M, X, TM, onTheFly, TA >::subdomain_vector |
The vector type containing the subdomain to row index mapping.
anonymous enum |
void Dune::SeqOverlappingSchwarz< M, X, TM, onTheFly, 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.
[out] | v | The update to be computed |
d | The current defect. |
Implements Dune::Preconditioner< X, X >.
virtual void Dune::SeqOverlappingSchwarz< M, X, TM, onTheFly, 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.
x | The right hand side of the equation. |
Implements Dune::Preconditioner< X, X >.
virtual void Dune::SeqOverlappingSchwarz< M, X, TM, onTheFly, 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.
x | The left hand side of the equation. | |
b | The right hand side of the equation. |
Implements Dune::Preconditioner< X, X >.