BEM++
2.0
|
Default iterative solver for boundary integral equations. More...
#include </home/wojtek/Projects/BEM/bempp-sven/bempp/lib/linalg/default_iterative_solver.hpp>
Public Types | |
typedef Solver < BasisFunctionType, ResultType > | Base |
Public Member Functions | |
DefaultIterativeSolver (const BoundaryOperator< BasisFunctionType, ResultType > &boundaryOp, ConvergenceTestMode::Mode mode=ConvergenceTestMode::TEST_CONVERGENCE_IN_DUAL_TO_RANGE) | |
Constructor of the DefaultIterativeSolver class. More... | |
DefaultIterativeSolver (const BlockedBoundaryOperator< BasisFunctionType, ResultType > &boundaryOp, ConvergenceTestMode::Mode mode=ConvergenceTestMode::TEST_CONVERGENCE_IN_DUAL_TO_RANGE) | |
Constructor of the DefaultIterativeSolver class. More... | |
BEMPP_DEPRECATED void | setPreconditioner (const Preconditioner< ResultType > &preconditioner) |
Define a preconditioner. More... | |
void | initializeSolver (const Teuchos::RCP< Teuchos::ParameterList > ¶mList) |
Initialize a Belos iterative solver. More... | |
void | initializeSolver (const Teuchos::RCP< Teuchos::ParameterList > ¶mList, const Preconditioner< ResultType > &preconditioner) |
Initialize a preconditioned Belos iterative solver. More... | |
![]() | |
Solution< BasisFunctionType, ResultType > | solve (const GridFunction< BasisFunctionType, ResultType > &rhs) const |
Solve a standard (non-blocked) boundary integral equation. More... | |
BlockedSolution < BasisFunctionType, ResultType > | solve (const std::vector< GridFunction< BasisFunctionType, ResultType > > &rhs) const |
Solve a block-operator system of boundary integral equations. More... | |
Private Member Functions | |
virtual Solution < BasisFunctionType, ResultType > | solveImplNonblocked (const GridFunction< BasisFunctionType, ResultType > &rhs) const |
virtual BlockedSolution < BasisFunctionType, ResultType > | solveImplBlocked (const std::vector< GridFunction< BasisFunctionType, ResultType > > &rhs) const |
Private Attributes | |
boost::scoped_ptr< Impl > | m_impl |
Additional Inherited Members | |
![]() | |
static void | checkConsistency (const BoundaryOperator< BasisFunctionType, ResultType > &boundaryOp, const GridFunction< BasisFunctionType, ResultType > &rhs, ConvergenceTestMode::Mode mode) |
static void | checkConsistency (const BlockedBoundaryOperator< BasisFunctionType, ResultType > &boundaryOp, const std::vector< GridFunction< BasisFunctionType, ResultType > > &rhs, ConvergenceTestMode::Mode mode) |
static std::vector < GridFunction < BasisFunctionType, ResultType > > | canonicalizeBlockedRhs (const BlockedBoundaryOperator< BasisFunctionType, ResultType > &boundaryOp, const std::vector< GridFunction< BasisFunctionType, ResultType > > &rhs, ConvergenceTestMode::Mode mode) |
static void | constructBlockedGridFunction (const arma::Col< ResultType > &solution, const BlockedBoundaryOperator< BasisFunctionType, ResultType > &boundaryOp, std::vector< GridFunction< BasisFunctionType, ResultType > > &solutionFunctions) |
Default iterative solver for boundary integral equations.
This class provides an interface to various iterative solvers available via the Stratimikos interface to the Belos solver family from Trilinos (see Stratimikos documentation). Convergence can be tested either in range space or in the dual space to the range space. A standard Galerkin discretisation of the form , maps into the dual space of the range of the operator. By choosing to test in the range space the equation
is solved, where
is the mass matrix, mapping from the range space into its dual and
is its pseudoinverse.
Bempp::DefaultIterativeSolver< BasisFunctionType, ResultType >::DefaultIterativeSolver | ( | const BoundaryOperator< BasisFunctionType, ResultType > & | boundaryOp, |
ConvergenceTestMode::Mode | mode = ConvergenceTestMode::TEST_CONVERGENCE_IN_DUAL_TO_RANGE |
||
) |
Constructor of the DefaultIterativeSolver
class.
[in] | boundaryOp | Non-blocked boundary operator. |
[in] | mode | Convergence test mode. Default: TEST_CONVERGENCE_IN_DUAL_TO_RANGE . |
Bempp::DefaultIterativeSolver< BasisFunctionType, ResultType >::DefaultIterativeSolver | ( | const BlockedBoundaryOperator< BasisFunctionType, ResultType > & | boundaryOp, |
ConvergenceTestMode::Mode | mode = ConvergenceTestMode::TEST_CONVERGENCE_IN_DUAL_TO_RANGE |
||
) |
Constructor of the DefaultIterativeSolver
class.
[in] | boundaryOp | Blocked boundary operator. |
[in] | mode | Convergence test mode. Default: TEST_CONVERGENCE_IN_DUAL_TO_RANGE . |
void Bempp::DefaultIterativeSolver< BasisFunctionType, ResultType >::initializeSolver | ( | const Teuchos::RCP< Teuchos::ParameterList > & | paramList | ) |
Initialize a Belos iterative solver.
[in] | paramList | Parameter lists can be read in from XML files or defined in code. Use defaultGmresParameterList() and defaultCgParameterList() to construct default parameter lists for the GMRES and CG solvers. |
void Bempp::DefaultIterativeSolver< BasisFunctionType, ResultType >::initializeSolver | ( | const Teuchos::RCP< Teuchos::ParameterList > & | paramList, |
const Preconditioner< ResultType > & | preconditioner | ||
) |
Initialize a preconditioned Belos iterative solver.
[in] | paramList | Parameter lists can be read in from XML files or defined in code. Use defaultGmresParameterList() and defaultCgParameterList() to construct default parameter lists for the GMRES and CG solvers. |
[in] | preconditioner | Preconditioner to be used by the solver. |
void Bempp::DefaultIterativeSolver< BasisFunctionType, ResultType >::setPreconditioner | ( | const Preconditioner< ResultType > & | preconditioner | ) |
Define a preconditioner.
The preconditioner is passed on to the Belos solver.
[in] | preconditioner |