|
Belos
Version of the Day
|
A solver manager for the Block GCRO-DR (Block Recycling GMRES) linear solver. More...
#include <BelosBlockGCRODRSolMgr.hpp>

Public Member Functions | |
Constructors/Destructor | |
| BlockGCRODRSolMgr () | |
| Default constructor. | |
| BlockGCRODRSolMgr (const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem, const Teuchos::RCP< Teuchos::ParameterList > &pl) | |
| Basic constructor for GCRODRSolMgr. | |
| virtual | ~BlockGCRODRSolMgr () |
| Destructor. | |
Implementation of the Teuchos::Describable interface | |
| std::string | description () const |
| A description of the Block GCRODR solver manager. | |
Accessor methods | |
| const LinearProblem < ScalarType, MV, OP > & | getProblem () const |
| Get current linear problem being solved for in this object. | |
| Teuchos::RCP< const Teuchos::ParameterList > | getValidParameters () const |
| Get a parameter list containing the valid parameters for this object. | |
| Teuchos::RCP< const Teuchos::ParameterList > | getCurrentParameters () const |
| Get a parameter list containing the current parameters for this object. | |
| MagnitudeType | achievedTol () const |
| Get the residual for the most recent call to solve(). | |
| int | getNumIters () const |
Get the iteration count for the most recent call to solve(). | |
| bool | isLOADetected () const |
| Whether a loss of accuracy was detected during the most recent solve. | |
Set methods | |
| void | setProblem (const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem) |
Set the linear problem to solve on the next call to solve(). | |
| void | setParameters (const Teuchos::RCP< Teuchos::ParameterList > ¶ms) |
| Set the parameters the solver should use to solve the linear problem. | |
Reset methods | |
| void | reset (const ResetType type) |
Performs a reset of the solver manager specified by the ResetType. | |
Solver application methods | |
| ReturnType | solve () |
| Solve the current linear problem. | |
A solver manager for the Block GCRO-DR (Block Recycling GMRES) linear solver.
GCRO-DR (also called Recycling GMRES) is a variant of GMRES that can more efficiently solve sequences linear systems. It does so by "recycling" Krylov basis information from previous solves.
The original GCRO-DR algorithm can only solve one right-hand side at a time. Block GCRO-DR extends GCRO-DR so that it can solve multiple right-hand sides at a time; thus, it can solve sequences of block systems.
Definition at line 127 of file BelosBlockGCRODRSolMgr.hpp.
| Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::BlockGCRODRSolMgr | ( | ) |
Default constructor.
This constructor sets up the solver with default parameters. The linear problem must be passed in using setProblem() before solve() is called on this object. The solver values can be changed using setParameters().
Definition at line 439 of file BelosBlockGCRODRSolMgr.hpp.
| Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::BlockGCRODRSolMgr | ( | const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > & | problem, |
| const Teuchos::RCP< Teuchos::ParameterList > & | pl | ||
| ) |
Basic constructor for GCRODRSolMgr.
This constructor accepts the LinearProblem to be solved in addition to a parameter list of options for the solver manager. Some of the more important options include the following:
int specifying the number of blocks allocated for the Krylov basis. Default: 50.int specifying the number of right hand sides being solved at a time.int specifying the number of blocks allocated for the Krylov basis. Default: 5.int specifying the maximum number of iterations the underlying solver is allowed to perform. Default: 5000.int specifying the maximum number of restarts the underlying solver is allowed to perform. Default: 100.std::string specifying the desired orthogonalization. Currently supported values: "DGKS", "ICGS", "IMGS", and "TSQR" (if Belos was built with TSQR support). Default: "DGKS".MagnitudeType specifying the level that residual norms must reach to decide convergence. Default: 1e-8.Other supported options:
MagnitudeType corresponding to the "depTol" parameter of DGKS orthogonalization. Ignored unless DGKS orthogonalization is used. DGKS decides the default value. Definition at line 446 of file BelosBlockGCRODRSolMgr.hpp.
| virtual Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::~BlockGCRODRSolMgr | ( | ) | [inline, virtual] |
Destructor.
Definition at line 181 of file BelosBlockGCRODRSolMgr.hpp.
| std::string Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::description | ( | ) | const [virtual] |
A description of the Block GCRODR solver manager.
Reimplemented from Teuchos::Describable.
Definition at line 542 of file BelosBlockGCRODRSolMgr.hpp.
| const LinearProblem<ScalarType,MV,OP>& Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::getProblem | ( | ) | const [inline, virtual] |
Get current linear problem being solved for in this object.
Implements Belos::SolverManager< ScalarType, MV, OP >.
Definition at line 197 of file BelosBlockGCRODRSolMgr.hpp.
| Teuchos::RCP< const Teuchos::ParameterList > Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::getValidParameters | ( | ) | const [virtual] |
Get a parameter list containing the valid parameters for this object.
Implements Belos::SolverManager< ScalarType, MV, OP >.
Definition at line 557 of file BelosBlockGCRODRSolMgr.hpp.
| Teuchos::RCP<const Teuchos::ParameterList> Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::getCurrentParameters | ( | ) | const [inline, virtual] |
Get a parameter list containing the current parameters for this object.
Implements Belos::SolverManager< ScalarType, MV, OP >.
Definition at line 205 of file BelosBlockGCRODRSolMgr.hpp.
| MagnitudeType Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::achievedTol | ( | ) | const [inline, virtual] |
Get the residual for the most recent call to solve().
Reimplemented from Belos::SolverManager< ScalarType, MV, OP >.
Definition at line 208 of file BelosBlockGCRODRSolMgr.hpp.
| int Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::getNumIters | ( | ) | const [inline, virtual] |
Get the iteration count for the most recent call to solve().
Implements Belos::SolverManager< ScalarType, MV, OP >.
Definition at line 213 of file BelosBlockGCRODRSolMgr.hpp.
| bool Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::isLOADetected | ( | ) | const [inline, virtual] |
Whether a loss of accuracy was detected during the most recent solve.
Implements Belos::SolverManager< ScalarType, MV, OP >.
Definition at line 216 of file BelosBlockGCRODRSolMgr.hpp.
| void Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::setProblem | ( | const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > & | problem | ) | [inline] |
Set the linear problem to solve on the next call to solve().
Definition at line 224 of file BelosBlockGCRODRSolMgr.hpp.
| void Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::setParameters | ( | const Teuchos::RCP< Teuchos::ParameterList > & | params | ) |
Set the parameters the solver should use to solve the linear problem.
Definition at line 645 of file BelosBlockGCRODRSolMgr.hpp.
| void Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::reset | ( | const ResetType | type | ) | [inline, virtual] |
Performs a reset of the solver manager specified by the ResetType.
This informs the solver manager that the solver should prepare for the next call to solve by resetting certain elements of the iterative solver strategy.
Implements Belos::SolverManager< ScalarType, MV, OP >.
Definition at line 253 of file BelosBlockGCRODRSolMgr.hpp.
| ReturnType Belos::BlockGCRODRSolMgr< ScalarType, MV, OP >::solve | ( | ) | [virtual] |
Solve the current linear problem.
This method performs possibly repeated calls to the underlying linear solver's iterate() routine until the problem has been solved (as decided by the solver manager) or the solver manager decides to quit.
This method calls BlockGCRODRIter::iterate(), which will return either because a specially constructed status test evaluates to Passed or an exception is thrown.
A return from BlockGCRODRIter::iterate() signifies one of the following scenarios:
Implements Belos::SolverManager< ScalarType, MV, OP >.
Definition at line 1797 of file BelosBlockGCRODRSolMgr.hpp.
1.7.6.1