Belos  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
List of all members
Belos::GmresPolySolMgr< ScalarType, MV, OP > Class Template Reference

The GMRES polynomial can be created in conjunction with any standard preconditioner. More...

#include <BelosGmresPolySolMgr.hpp>

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

Public Member Functions

Constructors/Destructor
 GmresPolySolMgr ()
 Empty constructor for GmresPolySolMgr. This constructor takes no arguments and sets the default values for the solver. The linear problem must be passed in using setProblem() before solve() is called on this object. The solver values can be changed using setParameters(). More...
 
 GmresPolySolMgr (const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem, const Teuchos::RCP< Teuchos::ParameterList > &pl)
 Basic constructor for GmresPolySolMgr. More...
 
virtual ~GmresPolySolMgr ()
 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
 Get current linear problem being solved for in this object. More...
 
Teuchos::RCP< const
Teuchos::ParameterList
getValidParameters () const override
 Get a parameter list containing the valid parameters for this object. More...
 
Teuchos::RCP< const
Teuchos::ParameterList
getCurrentParameters () const override
 Get a parameter list containing the current parameters for this object. More...
 
MagnitudeType achievedTol () const override
 Tolerance achieved by the last solve() invocation. More...
 
Teuchos::Array< Teuchos::RCP
< Teuchos::Time > > 
getTimers () const
 Return the timers for this object. More...
 
int getNumIters () const override
 Get the iteration count for the most recent call to solve(). More...
 
bool isLOADetected () const override
 Return whether a loss of accuracy was detected by this solver during the most current solve. More...
 
Set methods
void setProblem (const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem) override
 Set the linear problem that needs to be solved. More...
 
void setParameters (const Teuchos::RCP< Teuchos::ParameterList > &params) override
 Set the parameters the solver manager should use to solve the linear problem. More...
 
Reset methods
void reset (const ResetType type) override
 Reset the solver. 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
std::string description () const override
 Method to return description of the hybrid block GMRES solver manager. More...
 
- Public Member Functions inherited from Belos::SolverManager< ScalarType, MV, OP >
 SolverManager ()
 Empty constructor. More...
 
virtual ~SolverManager ()
 Destructor. More...
 
virtual void setUserConvStatusTest (const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &, const typename StatusTestCombo< ScalarType, MV, OP >::ComboType &=StatusTestCombo< ScalarType, MV, OP >::SEQ)
 Set user-defined convergence status test. More...
 
virtual void setDebugStatusTest (const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &)
 Set user-defined debug status test. More...
 
- Public Member Functions inherited from Teuchos::Describable
virtual void describe (FancyOStream &out, const EVerbosityLevel verbLevel=verbLevel_default) const
 
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
 

Detailed Description

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

The GMRES polynomial can be created in conjunction with any standard preconditioner.

Hybrid block GMRES iterative linear solver.

Simply pass the preconditioner to the LinearProblem before calling the GmresPolySolMgr and your preconditioner will be combined with the polynomial automatically.

Here is a list of all the parameters that this solver accepts:

This solver manager provides three different implementations of the same polynomial preconditioner. The polynomial is the minimum residual polynomial from GMRES. The "Roots" version is default. It is the only implementation which provides the option of added roots for stability. These added roots can allow for high-degree polynomials. The "Arnoldi" version typically gives similar results to the "Roots" version but is slightly more expensive to apply. Both of these polynomials can be "damped", which is sometimes useful for indefinite or other ill-conditioned problems. The "Gmres" version is based on a power-basis implementation and is only stable for well-conditioned problems and low-degree polynomials. For more information on the implementation and formulas, see the following references: "Roots" version: https://arxiv.org/abs/1806.08020 (Includes explanation of root-adding and damping.) "Arnoldi" version: https://scholarship.rice.edu/handle/1911/17630 "Gmres" version: https://epubs.siam.org/doi/pdf/10.1137/140968276

Like all Belos solvers, parameters have relative or "delta" semantics. This means the following:

Author
Heidi Thornquist and Jennifer Loe
Examples:
BlockGmres/BlockGmresPolyEpetraExFile.cpp.

Definition at line 152 of file BelosGmresPolySolMgr.hpp.

Constructor & Destructor Documentation

template<class ScalarType , class MV , class OP >
Belos::GmresPolySolMgr< ScalarType, MV, OP >::GmresPolySolMgr ( )

Empty constructor for GmresPolySolMgr. This constructor takes no arguments and sets the default values for the solver. 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 376 of file BelosGmresPolySolMgr.hpp.

template<class ScalarType , class MV , class OP >
Belos::GmresPolySolMgr< ScalarType, MV, OP >::GmresPolySolMgr ( const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &  problem,
const Teuchos::RCP< Teuchos::ParameterList > &  pl 
)

Basic constructor for GmresPolySolMgr.

This constructor accepts the LinearProblem to be solved in addition to a parameter list of options for the solver manager. These options include the following:

  • "Polynomial Type" -a std::string specifying the type of polynomial: Roots, Arnoldi, or Gmres. Default: "Roots"
  • "Maximum Degree" - a int specifying the maximum degree of the polynomial. Default: 25
  • "Random RHS" - a bool indicates whether to generate polynomial using a random vector. Default: true
  • "Add Roots" - a bool to add roots to the polynomial as needed for stability. Default: true
  • "Damp Poly" - a bool to damp polynomial. Default: false
  • "Orthogonalization" - a std::string specifying the desired orthogonalization to create the polynomial: DGKS, ICGS, and IMGS. Default: "ICGS"
  • "Verbosity" - a sum of MsgType specifying the verbosity. Default: Belos::Errors
  • "Polynomial Tolerance" - a MagnitudeType specifying the polynomial tolerance (sometimes) used to generate polynomial. Default: 1e-8
  • "Outer Solver" -a std::string specifying name of outer solver in Belos solver factory. Default: ""
  • "Outer Solver Params" -a Teuchos::ParameterList giving parameters for the outer solver.
  • "Timer Label" -a std::string specifying the label on polynomial solve timers.

Definition at line 399 of file BelosGmresPolySolMgr.hpp.

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

Destructor.

Definition at line 195 of file BelosGmresPolySolMgr.hpp.

Member Function Documentation

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

clone for Inverted Injection (DII)

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

Definition at line 198 of file BelosGmresPolySolMgr.hpp.

template<class ScalarType , class MV , class OP >
const LinearProblem<ScalarType,MV,OP>& Belos::GmresPolySolMgr< ScalarType, MV, OP >::getProblem ( ) const
inlineoverridevirtual

Get current linear problem being solved for in this object.

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

Definition at line 208 of file BelosGmresPolySolMgr.hpp.

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

Get a parameter list containing the valid parameters for this object.

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

Definition at line 435 of file BelosGmresPolySolMgr.hpp.

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

Get a parameter list containing the current parameters for this object.

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

Definition at line 218 of file BelosGmresPolySolMgr.hpp.

template<class ScalarType , class MV , class OP >
MagnitudeType Belos::GmresPolySolMgr< 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 solver manager may be use as either a polynomial preconditioned iterative method or a polynomial preconditioner. In the later case, where a static polynomial is being applied through each call to solve(), there is not an outer solve that can provide the achieved tolerance.
This result may not be meaningful if there was a loss of accuracy during the outer 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 236 of file BelosGmresPolySolMgr.hpp.

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

Return the timers for this object.

The timers are ordered as follows:

Definition at line 245 of file BelosGmresPolySolMgr.hpp.

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

Get the iteration count for the most recent call to solve().

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

Definition at line 250 of file BelosGmresPolySolMgr.hpp.

template<class ScalarType , class MV , class OP >
bool Belos::GmresPolySolMgr< ScalarType, MV, OP >::isLOADetected ( ) const
inlineoverridevirtual

Return whether a loss of accuracy was detected by this solver during the most current solve.

Note
This flag will be reset the next time solve() is called.

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

Definition at line 257 of file BelosGmresPolySolMgr.hpp.

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

Set the linear problem that needs to be solved.

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

Definition at line 265 of file BelosGmresPolySolMgr.hpp.

template<class ScalarType , class MV , class OP >
void Belos::GmresPolySolMgr< 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 476 of file BelosGmresPolySolMgr.hpp.

template<class ScalarType , class MV , class OP >
void Belos::GmresPolySolMgr< ScalarType, MV, OP >::reset ( const ResetType  type)
inlineoverridevirtual

Reset the solver.

Parameters
type[in] How to reset the solver.

If type includes Belos::Problem, then reset the solver's state. This clears out the stored coefficients, so that the next call to solve() actually computes a full block GMRES solve, instead of just reusing the coefficients from the first solve.

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

Definition at line 282 of file BelosGmresPolySolMgr.hpp.

template<class ScalarType , class MV , class OP >
ReturnType Belos::GmresPolySolMgr< 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 BlockGmresIter::iterate(), which will return either because a specially constructed status test evaluates to Passed or an std::exception is thrown.

A return from BlockGmresIter::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 630 of file BelosGmresPolySolMgr.hpp.

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

Method to return description of the hybrid block GMRES solver manager.

Reimplemented from Teuchos::Describable.

Definition at line 740 of file BelosGmresPolySolMgr.hpp.


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

Generated on Fri Apr 26 2024 09:25:24 for Belos by doxygen 1.8.5