Belos Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Private Types | Private Attributes | Static Private Attributes | List of all members
Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP > Class Template Reference

Interface to standard and "pseudoblock" GMRES. More...

#include <BelosPseudoBlockGmresSolMgr.hpp>

Inheritance diagram for Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >:
Inheritance graph
[legend]

Private Types

typedef MultiVecTraits
< ScalarType, MV > 
MVT
 
typedef OperatorTraits
< ScalarType, MV, OP > 
OPT
 
typedef Teuchos::ScalarTraits
< ScalarType > 
SCT
 
typedef Teuchos::ScalarTraits
< ScalarType >::magnitudeType 
MagnitudeType
 
typedef Teuchos::ScalarTraits
< MagnitudeType
MT
 

Private Attributes

Teuchos::RCP< LinearProblem
< ScalarType, MV, OP > > 
problem_
 The current linear problem to solve. More...
 
Teuchos::RCP< OutputManager
< ScalarType > > 
printer_
 
Teuchos::RCP< std::ostream > outputStream_
 
Teuchos::RCP< StatusTest
< ScalarType, MV, OP > > 
userConvStatusTest_
 
Teuchos::RCP< StatusTest
< ScalarType, MV, OP > > 
debugStatusTest_
 
Teuchos::RCP< StatusTest
< ScalarType, MV, OP > > 
sTest_
 
Teuchos::RCP
< StatusTestMaxIters
< ScalarType, MV, OP > > 
maxIterTest_
 
Teuchos::RCP< StatusTest
< ScalarType, MV, OP > > 
convTest_
 
Teuchos::RCP
< StatusTestResNorm
< ScalarType, MV, OP > > 
impConvTest_
 
Teuchos::RCP
< StatusTestResNorm
< ScalarType, MV, OP > > 
expConvTest_
 
Teuchos::RCP< StatusTestOutput
< ScalarType, MV, OP > > 
outputTest_
 
StatusTestCombo< ScalarType,
MV, OP >::ComboType 
comboType_
 
Teuchos::RCP< std::map
< std::string, Teuchos::RCP
< StatusTest< ScalarType, MV,
OP > > > > 
taggedTests_
 
Teuchos::RCP< MatOrthoManager
< ScalarType, MV, OP > > 
ortho_
 
Teuchos::RCP
< Teuchos::ParameterList
params_
 
MagnitudeType convtol_
 
MagnitudeType orthoKappa_
 
MagnitudeType achievedTol_
 
int maxRestarts_
 
int maxIters_
 
int numIters_
 
int blockSize_
 
int numBlocks_
 
int verbosity_
 
int outputStyle_
 
int outputFreq_
 
int defQuorum_
 
bool showMaxResNormOnly_
 
std::string orthoType_
 
std::string impResScale_
 
std::string expResScale_
 
MagnitudeType resScaleFactor_
 
std::string label_
 
Teuchos::RCP< Teuchos::TimetimerSolve_
 
bool isSet_
 
bool isSTSet_
 
bool expResTest_
 
bool loaDetected_
 

Static Private Attributes

static constexpr int maxRestarts_default_ = 20
 
static constexpr int maxIters_default_ = 1000
 
static constexpr bool showMaxResNormOnly_default_ = false
 
static constexpr int blockSize_default_ = 1
 
static constexpr int numBlocks_default_ = 300
 
static constexpr int verbosity_default_ = Belos::Errors
 
static constexpr int outputStyle_default_ = Belos::General
 
static constexpr int outputFreq_default_ = -1
 
static constexpr int defQuorum_default_ = 1
 
static constexpr const char * impResScale_default_ = "Norm of Preconditioned Initial Residual"
 
static constexpr const char * expResScale_default_ = "Norm of Initial Residual"
 
static constexpr const char * label_default_ = "Belos"
 
static constexpr const char * orthoType_default_ = "ICGS"
 
static constexpr std::ostream * outputStream_default_ = &std::cout
 

PseudoBlockGmresSolMgr Exceptions

bool checkStatusTest ()
 Check current status tests against current linear problem. More...
 

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 > &params) 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...
 

Additional Inherited Members

- 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
DescribableStreamManipulatorState describe (const Describable &describable, const EVerbosityLevel verbLevel=Describable::verbLevel_default)
 
std::ostream & operator<< (std::ostream &os, const DescribableStreamManipulatorState &d)
 
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
 
- Static Public Attributes inherited from Teuchos::Describable
static const EVerbosityLevel verbLevel_default
 

Detailed Description

template<class ScalarType, class MV, class OP>
class Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >

Interface to standard and "pseudoblock" GMRES.

Author
Heidi Thornquist, Chris Baker, and Teri Barth

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.

Member Typedef Documentation

template<class ScalarType , class MV , class OP >
typedef MultiVecTraits<ScalarType,MV> Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::MVT
private

Definition at line 124 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
typedef OperatorTraits<ScalarType,MV,OP> Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::OPT
private

Definition at line 125 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
typedef Teuchos::ScalarTraits<ScalarType> Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::SCT
private

Definition at line 126 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
typedef Teuchos::ScalarTraits<ScalarType>::magnitudeType Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::MagnitudeType
private

Definition at line 127 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
typedef Teuchos::ScalarTraits<MagnitudeType> Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::MT
private

Definition at line 128 of file BelosPseudoBlockGmresSolMgr.hpp.

Constructor & Destructor Documentation

template<class ScalarType , class MV , class OP >
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.

template<class ScalarType , class MV , class OP >
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.

Parameters
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:

  • "Num Blocks" (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.
  • "Maximum Iterations" (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.
  • "Maximum Restarts" (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:

  • "Convergence Tolerance" (MagnitudeType): The level that residual norms must reach in order for the solver to stop iterating.
  • "Implicit Residual Scaling" (std::string): How to scale the implicit residual norm. The default is the norm of the preconditioned initial residual.
  • "Explicit Residual Scaling" (std::string): How to scale the explicit residual norm. The default is the norm of the (unpreconditioned) initial residual.
  • "Deflation Quorum" (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.
  • "Orthogonalization" (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:

  • "Norm of Initial Residual"
  • "Norm of Preconditioned Initial Residual"
  • "Norm of RHS" (RHS stands for "right-hand side")
  • "None" (no scaling factor)

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.

  • "Verbosity": a sum of MsgType enum values specifying the verbosity. Default: Belos::Errors.
  • "Output Frequency" (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.
  • "Output Style" (OutputType): The style of output. Accepted values are General and Brief. Default: General.
  • "Output Stream" (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.
  • "Show Maximum Residual Norm Only": When solving for multiple right-hand sides, this controls whether output shows residual norms for all the right-hand sides, or just the current maximum residual norm over all right-hand sides.
    Parameters
    pl[in] ParameterList with construction information

Definition at line 558 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
virtual Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::~PseudoBlockGmresSolMgr ( )
inlinevirtual

Destructor.

Definition at line 258 of file BelosPseudoBlockGmresSolMgr.hpp.

Member Function Documentation

template<class ScalarType , class MV , class OP >
Teuchos::RCP<SolverManager<ScalarType, MV, OP> > Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::clone ( ) const
inlineoverridevirtual

clone for Inverted Injection (DII)

Implements Belos::SolverManager< ScalarType, MV, OP >.

Definition at line 261 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
const LinearProblem<ScalarType,MV,OP>& Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::getProblem ( ) const
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.

template<class ScalarType , class MV , class OP >
Teuchos::RCP< const Teuchos::ParameterList > Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::getValidParameters ( ) const
overridevirtual

A list of valid default parameters for this solver.

Implements Belos::SolverManager< ScalarType, MV, OP >.

Definition at line 1015 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<const Teuchos::ParameterList> Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::getCurrentParameters ( ) const
inlineoverridevirtual

The current parameters for this solver.

Implements Belos::SolverManager< ScalarType, MV, OP >.

Definition at line 277 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::Array<Teuchos::RCP<Teuchos::Time> > Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::getTimers ( ) const
inline

Return the timers for this object.

The timers are ordered as follows:

Definition at line 284 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
MagnitudeType Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::achievedTol ( ) const
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.

Warning
This result may not be meaningful if there was a loss of accuracy during the solve. You should first call 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.

template<class ScalarType , class MV , class OP >
int Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::getNumIters ( ) const
inlineoverridevirtual

Iteration count for the most recent call to solve().

Implements Belos::SolverManager< ScalarType, MV, OP >.

Definition at line 303 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
bool Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::isLOADetected ( ) const
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 $X_k$, applying the global operator $A$ to $X_k$, 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 $M_L^{-1} A X = M_L^{-1} B$ instead of the original problem $A X = B$. 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.

Note
Calling 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.

template<class ScalarType , class MV , class OP >
const StatusTestResNorm<ScalarType,MV,OP>* Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::getResidualStatusTest ( ) const
inline

Return the residual status test.

Definition at line 366 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
void Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::setProblem ( const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &  problem)
inlineoverridevirtual

Set the linear problem to solve.

Implements Belos::SolverManager< ScalarType, MV, OP >.

Definition at line 373 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
void Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::setParameters ( const Teuchos::RCP< Teuchos::ParameterList > &  params)
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.

template<class ScalarType , class MV , class OP >
void Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::setUserConvStatusTest ( const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &  userConvStatusTest,
const typename StatusTestCombo< ScalarType, MV, OP >::ComboType &  comboType = StatusTestCombo<ScalarType,MV,OP>::SEQ 
)
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.

template<class ScalarType , class MV , class OP >
void Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::setDebugStatusTest ( const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &  debugStatusTest)
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.

template<class ScalarType , class MV , class OP >
void Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::reset ( const ResetType  type)
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.

template<class ScalarType , class MV , class OP >
ReturnType Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::solve ( )
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:

  • the maximum number of restarts has been exceeded. In this scenario, the current solutions to the linear system will be placed in the linear problem and return Unconverged.
  • global convergence has been met. In this case, the current solutions to the linear system will be placed in the linear problem and the solver manager will return Converged
Returns
ReturnType specifying:
  • Converged: the linear problem was solved to the specification required by the solver manager.
  • Unconverged: the linear problem was not solved to the specification desired by the solver manager.

Implements Belos::SolverManager< ScalarType, MV, OP >.

Definition at line 1191 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
void Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::describe ( Teuchos::FancyOStream out,
const Teuchos::EVerbosityLevel  verbLevel = Teuchos::Describable::verbLevel_default 
) const
overridevirtual

Print the object with the given verbosity level to a FancyOStream.

Parameters
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.

template<class ScalarType , class MV , class OP >
std::string Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::description ( ) const
overridevirtual

Return a one-line description of this object.

Reimplemented from Teuchos::Describable.

Definition at line 1596 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
bool Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::checkStatusTest ( )
private

Check current status tests against current linear problem.

(Re)create all the status tests, based on the current solve parameters and the current linear problem to solve. This is necessary whenever the linear problem is set or changed via setProblem(), because the residual norm test to use depends on whether or not the (new) linear problem defines a left preconditioner. Furthermore, include the user's custom convergence test if they set one via setUserConvStatusTest().

Returns
False if we were able to (re)create all the status tests correctly, else true. The solve() routine may call this method. If it does, it checks the return value.

Definition at line 1079 of file BelosPseudoBlockGmresSolMgr.hpp.

Member Data Documentation

template<class ScalarType , class MV , class OP >
Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::problem_
private

The current linear problem to solve.

Definition at line 468 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<OutputManager<ScalarType> > Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::printer_
private

Definition at line 471 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<std::ostream> Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::outputStream_
private

Definition at line 472 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<StatusTest<ScalarType,MV,OP> > Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::userConvStatusTest_
private

Definition at line 475 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<StatusTest<ScalarType,MV,OP> > Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::debugStatusTest_
private

Definition at line 476 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<StatusTest<ScalarType,MV,OP> > Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::sTest_
private

Definition at line 477 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<StatusTestMaxIters<ScalarType,MV,OP> > Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::maxIterTest_
private

Definition at line 478 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<StatusTest<ScalarType,MV,OP> > Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::convTest_
private

Definition at line 479 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<StatusTestResNorm<ScalarType,MV,OP> > Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::impConvTest_
private

Definition at line 480 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<StatusTestResNorm<ScalarType,MV,OP> > Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::expConvTest_
private

Definition at line 480 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<StatusTestOutput<ScalarType,MV,OP> > Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::outputTest_
private

Definition at line 481 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
StatusTestCombo<ScalarType,MV,OP>::ComboType Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::comboType_
private

Definition at line 482 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<std::map<std::string, Teuchos::RCP<StatusTest<ScalarType, MV, OP> > > > Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::taggedTests_
private

Definition at line 483 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<MatOrthoManager<ScalarType,MV,OP> > Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::ortho_
private

Definition at line 486 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<Teuchos::ParameterList> Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::params_
private

Definition at line 489 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
constexpr int Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::maxRestarts_default_ = 20
staticprivate

Definition at line 492 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
constexpr int Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::maxIters_default_ = 1000
staticprivate

Definition at line 493 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
constexpr bool Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::showMaxResNormOnly_default_ = false
staticprivate

Definition at line 494 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
constexpr int Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::blockSize_default_ = 1
staticprivate

Definition at line 495 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
constexpr int Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::numBlocks_default_ = 300
staticprivate

Definition at line 496 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
constexpr int Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::verbosity_default_ = Belos::Errors
staticprivate

Definition at line 497 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
constexpr int Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::outputStyle_default_ = Belos::General
staticprivate

Definition at line 498 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
constexpr int Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::outputFreq_default_ = -1
staticprivate

Definition at line 499 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
constexpr int Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::defQuorum_default_ = 1
staticprivate

Definition at line 500 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
constexpr const char* Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::impResScale_default_ = "Norm of Preconditioned Initial Residual"
staticprivate

Definition at line 501 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
constexpr const char* Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::expResScale_default_ = "Norm of Initial Residual"
staticprivate

Definition at line 502 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
constexpr const char* Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::label_default_ = "Belos"
staticprivate

Definition at line 503 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
constexpr const char* Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::orthoType_default_ = "ICGS"
staticprivate

Definition at line 504 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
constexpr std::ostream* Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::outputStream_default_ = &std::cout
staticprivate

Definition at line 505 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
MagnitudeType Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::convtol_
private

Definition at line 508 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
MagnitudeType Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::orthoKappa_
private

Definition at line 508 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
MagnitudeType Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::achievedTol_
private

Definition at line 508 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
int Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::maxRestarts_
private

Definition at line 509 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
int Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::maxIters_
private

Definition at line 509 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
int Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::numIters_
private

Definition at line 509 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
int Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::blockSize_
private

Definition at line 510 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
int Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::numBlocks_
private

Definition at line 510 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
int Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::verbosity_
private

Definition at line 510 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
int Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::outputStyle_
private

Definition at line 510 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
int Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::outputFreq_
private

Definition at line 510 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
int Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::defQuorum_
private

Definition at line 510 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
bool Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::showMaxResNormOnly_
private

Definition at line 511 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
std::string Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::orthoType_
private

Definition at line 512 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
std::string Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::impResScale_
private

Definition at line 513 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
std::string Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::expResScale_
private

Definition at line 513 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
MagnitudeType Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::resScaleFactor_
private

Definition at line 514 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
std::string Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::label_
private

Definition at line 517 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<Teuchos::Time> Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::timerSolve_
private

Definition at line 518 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
bool Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::isSet_
private

Definition at line 521 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
bool Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::isSTSet_
private

Definition at line 521 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
bool Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::expResTest_
private

Definition at line 521 of file BelosPseudoBlockGmresSolMgr.hpp.

template<class ScalarType , class MV , class OP >
bool Belos::PseudoBlockGmresSolMgr< ScalarType, MV, OP >::loaDetected_
private

Definition at line 522 of file BelosPseudoBlockGmresSolMgr.hpp.


The documentation for this class was generated from the following file: