Belos
Version of the Day
|
Interface to standard and "pseudoblock" GMRES. More...
#include <BelosPseudoBlockGmresSolMgr.hpp>
Public Member Functions | |
Constructors and destructor | |
PseudoBlockGmresSolMgr () | |
Empty constructor. More... | |
PseudoBlockGmresSolMgr (const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem, const Teuchos::RCP< Teuchos::ParameterList > &pl) | |
Constructor that takes the problem to solve, and a list of solver options. More... | |
virtual | ~PseudoBlockGmresSolMgr () |
Destructor. More... | |
Teuchos::RCP< SolverManager < ScalarType, MV, OP > > | clone () const override |
clone for Inverted Injection (DII) More... | |
Accessor methods | |
const LinearProblem < ScalarType, MV, OP > & | getProblem () const override |
Return a reference to the linear problem being solved by this solver manager. More... | |
Teuchos::RCP< const Teuchos::ParameterList > | getValidParameters () const override |
A list of valid default parameters for this solver. More... | |
Teuchos::RCP< const Teuchos::ParameterList > | getCurrentParameters () const override |
The current parameters for this solver. More... | |
Teuchos::Array< Teuchos::RCP < Teuchos::Time > > | getTimers () const |
Return the timers for this object. More... | |
MagnitudeType | achievedTol () const override |
Tolerance achieved by the last solve() invocation. More... | |
int | getNumIters () const override |
Iteration count for the most recent call to solve() . More... | |
bool | isLOADetected () const override |
Whether a "loss of accuracy" was detected during the last solve(). More... | |
const StatusTestResNorm < ScalarType, MV, OP > * | getResidualStatusTest () const |
Return the residual status test. More... | |
Set methods | |
void | setProblem (const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem) override |
Set the linear problem to solve. More... | |
void | setParameters (const Teuchos::RCP< Teuchos::ParameterList > ¶ms) override |
Set the parameters the solver manager should use to solve the linear problem. More... | |
virtual void | setUserConvStatusTest (const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &userConvStatusTest, const typename StatusTestCombo< ScalarType, MV, OP >::ComboType &comboType=StatusTestCombo< ScalarType, MV, OP >::SEQ) override |
Set a custom status test. More... | |
void | setDebugStatusTest (const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &debugStatusTest) override |
Set a debug status test that will be checked at the same time as the top-level status test. More... | |
Reset methods | |
void | reset (const ResetType type) override |
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. More... | |
Solver application methods | |
ReturnType | solve () override |
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. More... | |
Overridden from Teuchos::Describable | |
void | describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const override |
Print the object with the given verbosity level to a FancyOStream. More... | |
std::string | description () const override |
Return a one-line description of this object. More... | |
Public Member Functions inherited from Belos::SolverManager< ScalarType, MV, OP > | |
SolverManager () | |
Empty constructor. More... | |
virtual | ~SolverManager () |
Destructor. More... | |
Public Member Functions inherited from Teuchos::Describable | |
void | describe (std::ostream &out, const EVerbosityLevel verbLevel=verbLevel_default) const |
virtual | ~Describable () |
LabeledObject () | |
virtual | ~LabeledObject () |
virtual void | setObjectLabel (const std::string &objectLabel) |
virtual std::string | getObjectLabel () const |
DescribableStreamManipulatorState | describe (const Describable &describable, const EVerbosityLevel verbLevel=Describable::verbLevel_default) |
std::ostream & | operator<< (std::ostream &os, const DescribableStreamManipulatorState &d) |
Additional Inherited Members | |
Static Public Attributes inherited from Teuchos::Describable | |
static const EVerbosityLevel | verbLevel_default |
Interface to standard and "pseudoblock" GMRES.
This class provides an interface to the following iterative solvers:
If you are a new Belos user and just want standard GMRES, use this class. If you want Flexible GMRES, use BlockGmresSolMgr
with the appropriate option set.
"Pseudoblock" GMRES is a way to improve performance when solving systems with multiple right-hand sides, without changing the convergence characteristics. It is equivalent in terms of convergence to running a separate instance of (standard) GMRES for each right-hand side, but should often be faster. When solving for multiple right-hand sides, "Block GMRES" (as implemented by BlockGmresSolMgr
) is a different algorithm with different convergence characteristics than Pseudoblock GMRES.
Definition at line 121 of file BelosPseudoBlockGmresSolMgr.hpp.
Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::PseudoBlockGmresSolMgr | ( | ) |
Empty constructor.
This constructor takes no arguments. It sets default solver parameters, which you may change by calling setParameters(). Before you may call solve(), you must first give the solver a linear problem to solve, by calling setProblem().
Definition at line 528 of file BelosPseudoBlockGmresSolMgr.hpp.
Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::PseudoBlockGmresSolMgr | ( | const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > & | problem, |
const Teuchos::RCP< Teuchos::ParameterList > & | pl | ||
) |
Constructor that takes the problem to solve, and a list of solver options.
problem | [in/out] The linear problem to be solved. |
pl | [in/out] A list of solver options. |
Belos' solvers accept many different options. You may accept their default values, or set any of them yourself. We will explain the options by category.
The following options govern the number of iterations and restarts:
int
): The restart length. The number of vectors (or blocks, in the case of multiple right-hand sides) allocated for the Krylov basis. Its default value is 300.int
): The maximum number of iterations the solver is allowed to perform. This does not include computing the initial residual, but it does include iterations before and after any restarts. Its default value is 1000.int
): The maximum number of restarts. This does not include the first "Num Blocks" iterations (before the first restart). Its default value is 20.We do not currently perform any sanity checks for these options. This may affect you if you set some of them but let others keep their default values. For example, if you set "Num Blocks" to 2 and "Maximum Iterations" to 100, but don't set "Maximum Restarts", you will only get 40 = 20*2 total iterations, rather than 100. Thus, if you set one of these parameters, you should always set them all.
When solving with multiple right-hand sides, the "Block Size" (int
) parameter controls the number of right-hand sides for which the solver solves at once. This setting controls both performance and total memory use. Doubling it (approximately) doubles the total amount of memory used by the solver, but might make the solves faster by reducing synchronization overhead and improving memory bandwidth utilization. The gain from increasing this tends to level off quickly. Making this setting too large may actually hurt performance.
These options govern convergence and the numerical algorithm:
MagnitudeType
): The level that residual norms must reach in order for the solver to stop iterating.std::string
): How to scale the implicit residual norm. The default is the norm of the preconditioned initial residual.std::string
): How to scale the explicit residual norm. The default is the norm of the (unpreconditioned) initial residual.int
): When solving with multiple right-hand sides: the number of right-hand sides that must have converged to the given tolerance, before the solver will consider all the systems converged. If -1, then the solver will require that all the right-hand sides have converged before declaring all the systems converged. This must be no bigger than the "Block Size" parameter.std::string
): The desired orthogonalization method. Currently accepted values are "DGKS", "ICGS", "IMGS", and optionally "TSQR" (depending on build settings). Please refer to Belos' documentation for more details.For an explanation of "implicit" vs. "explicit" residuals, please see the documentation of isLOADetected(). The difference matters if using left preconditioning. Otherwise, it is not so important to most users.
The residual scaling parameters ("Implicit Residual Scaling" and "Explicit Residual Scaling") accept the following values:
GMRES always uses the 2 norm (square root of sum of squares of magnitudes of entries) to measure convergence.
Belos' solvers let users control intermediate "status" output. This output tells you the current iteration and the values of current convergence criteria. The following parameters control output. The default values are fine for users who only care about the final result and don't want to see status output.
MsgType
enum values specifying the verbosity. Default: Belos::Errors.int
): How often (in terms of number of iterations) to print intermediate status output. The default (-1) means not to print intermediate status output at all.OutputType
): The style of output. Accepted values are General and Brief. Default: General.Teuchos::RCP<std::ostream>
): A pointer to an output stream to which the solver will write status output. The default is a pointer to std::cout
. Currently, if Trilinos was built with MPI support, only the MPI process with rank 0 in MPI_COMM_WORLD will print to this output stream.pl | [in] ParameterList with construction information
|
Definition at line 558 of file BelosPseudoBlockGmresSolMgr.hpp.
|
inlinevirtual |
Destructor.
Definition at line 258 of file BelosPseudoBlockGmresSolMgr.hpp.
|
inlineoverridevirtual |
clone for Inverted Injection (DII)
Implements Belos::SolverManager< ScalarType, MV, OP >.
Definition at line 261 of file BelosPseudoBlockGmresSolMgr.hpp.
|
inlineoverridevirtual |
Return a reference to the linear problem being solved by this solver manager.
Implements Belos::SolverManager< ScalarType, MV, OP >.
Definition at line 269 of file BelosPseudoBlockGmresSolMgr.hpp.
|
overridevirtual |
A list of valid default parameters for this solver.
Implements Belos::SolverManager< ScalarType, MV, OP >.
Definition at line 1015 of file BelosPseudoBlockGmresSolMgr.hpp.
|
inlineoverridevirtual |
The current parameters for this solver.
Implements Belos::SolverManager< ScalarType, MV, OP >.
Definition at line 277 of file BelosPseudoBlockGmresSolMgr.hpp.
|
inline |
Return the timers for this object.
The timers are ordered as follows:
Definition at line 284 of file BelosPseudoBlockGmresSolMgr.hpp.
|
inlineoverridevirtual |
Tolerance achieved by the last solve()
invocation.
This is the maximum over all right-hand sides' achieved convergence tolerances, and is set whether or not the solve actually managed to achieve the desired convergence tolerance.
isLOADetected()
to check for a loss of accuracy during the last solve. Reimplemented from Belos::SolverManager< ScalarType, MV, OP >.
Definition at line 298 of file BelosPseudoBlockGmresSolMgr.hpp.
|
inlineoverridevirtual |
Iteration count for the most recent call to solve()
.
Implements Belos::SolverManager< ScalarType, MV, OP >.
Definition at line 303 of file BelosPseudoBlockGmresSolMgr.hpp.
|
inlineoverridevirtual |
Whether a "loss of accuracy" was detected during the last solve().
This solver uses two different residual norms to predict convergence: "implicit" (also called "native") and "explicit" (also called "exact," not to be confused with "exact arithmetic"). The "implicit" residuals are computed by the solver via a recurrence relation (the Arnoldi relation, in the case of GMRES). The "explicit" residuals are computed directly as $B - A X_k$. Implicit residuals are much cheaper to compute, since they are available almost "for free" from the recurrence relation. In contrast, computing exact residuals requires computing the current approximate solution , applying the global operator to , and then computing the norm of the resulting vector(s) via a global reduction. Thus, GMRES favors using the cheaper implicit residuals to predict convergence. Users typically want convergence with respect to explicit residuals, though.
Implicit and explicit residuals may differ due to rounding error. However, the difference between implicit and explicit residuals matters most when using a left (or split) preconditioner. In that case, the implicit residuals are those of the left-preconditioned problem instead of the original problem . The implicit residual norms may thus differ significantly from the explicit residual norms, even if one could compute without rounding error.
When using a left preconditioner, this solver tries to detect if the implicit residuals have converged but the explicit residuals have not. In that case, it will reduce the convergence tolerance and iterate a little while longer to attempt to reduce the explicit residual norm. However, if that doesn't work, it declares a "loss of accuracy" for the affected right-hand side(s), and stops iterating on them. (Not all right-hand sides may have experienced a loss of accuracy.) Thus, the affected right-hand sides may or may not have converged to the desired residual norm tolerance. Calling this method tells you whether a "loss of accuracy" (LOA) occurred during the last solve()
invocation.
When not using a left preconditioner, this solver will iterate until both the implicit and explicit residuals converge. (It does not start testing the explicit residuals until the implicit residuals have converged. This avoids whenever possible the cost of computing explicit residuals.) Implicit and explicit residuals may differ due to rounding error, even though they are identical when no rounding error occurs. In this case, the algorithm does not report a "loss of accuracy," since it continues iterating until the explicit residuals converge.
solve()
again resets the flag that reports whether a loss of accuracy was detected. Thus, you should call this method immediately after calling solve()
. Implements Belos::SolverManager< ScalarType, MV, OP >.
Definition at line 362 of file BelosPseudoBlockGmresSolMgr.hpp.
|
inline |
Return the residual status test.
Definition at line 366 of file BelosPseudoBlockGmresSolMgr.hpp.
|
inlineoverridevirtual |
Set the linear problem to solve.
Implements Belos::SolverManager< ScalarType, MV, OP >.
Definition at line 373 of file BelosPseudoBlockGmresSolMgr.hpp.
|
overridevirtual |
Set the parameters the solver manager should use to solve the linear problem.
Implements Belos::SolverManager< ScalarType, MV, OP >.
Definition at line 598 of file BelosPseudoBlockGmresSolMgr.hpp.
|
overridevirtual |
Set a custom status test.
A custom status test is not required. If you decide to set one, the current implementation will apply it sequentially (short-circuiting OR, like the || operator in C++) after Pseudoblock GMRES' standard convergence test.
Reimplemented from Belos::SolverManager< ScalarType, MV, OP >.
Definition at line 993 of file BelosPseudoBlockGmresSolMgr.hpp.
|
overridevirtual |
Set a debug status test that will be checked at the same time as the top-level status test.
Reimplemented from Belos::SolverManager< ScalarType, MV, OP >.
Definition at line 1004 of file BelosPseudoBlockGmresSolMgr.hpp.
|
inlineoverridevirtual |
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 403 of file BelosPseudoBlockGmresSolMgr.hpp.
|
overridevirtual |
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 PseudoBlockGmresIter::iterate(), which will return either because a specially constructed status test evaluates to Passed or an std::exception is thrown.
A return from PseudoBlockGmresIter::iterate() signifies one of the following scenarios:
Implements Belos::SolverManager< ScalarType, MV, OP >.
Definition at line 1191 of file BelosPseudoBlockGmresSolMgr.hpp.
|
overridevirtual |
Print the object with the given verbosity level to a FancyOStream.
out | [out] Output stream to which to print. |
verbLevel | [in] Verbosity level. The default verbosity (verbLevel=Teuchos::VERB_DEFAULT) is Teuchos::VERB_LOW. |
Reimplemented from Teuchos::Describable.
Definition at line 1616 of file BelosPseudoBlockGmresSolMgr.hpp.
|
overridevirtual |
Return a one-line description of this object.
Reimplemented from Teuchos::Describable.
Definition at line 1596 of file BelosPseudoBlockGmresSolMgr.hpp.