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

Factory for all solvers which Belos supports. More...

#include <BelosSolverFactory.hpp>

Inheritance diagram for Belos::SolverFactory< Scalar, MV, OP >:
Inheritance graph
[legend]

Public Types

typedef
parent_type::solver_base_type 
solver_base_type
 The type of the solver returned by create(). More...
 
typedef
parent_type::custom_solver_factory_type 
custom_solver_factory_type
 The type of a solver factory that users may give to addFactory() (which see). More...
 
- Public Types inherited from Belos::Impl::SolverFactoryParent< Scalar, MV, OP >
typedef ::Belos::SolverManager
< Scalar, MV, OP > 
solver_base_type
 The type of the solver returned by create(). More...
 
typedef CustomSolverFactory
< Scalar, MV, OP > 
custom_solver_factory_type
 The type of a solver factory that users may give to addFactory() (which see below) More...
 

Additional Inherited Members

- Public Member Functions inherited from Belos::Impl::SolverFactoryParent< Scalar, MV, OP >
virtual Teuchos::RCP
< solver_base_type
create (const std::string &solverName, const Teuchos::RCP< Teuchos::ParameterList > &solverParams)
 Create, configure, and return the specified solver. More...
 
virtual int numSupportedSolvers () const
 Number of supported solvers. More...
 
virtual Teuchos::Array
< std::string > 
supportedSolverNames () const
 List of supported solver names. More...
 
virtual bool isSupported (const std::string &solverName) const
 Whether the given solver name names a supported solver. More...
 
void addFactory (const Teuchos::RCP< custom_solver_factory_type > &factory)
 Add a custom solver factory. More...
 
virtual std::string description () const
 A string description of this object. More...
 
virtual void describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
 Describe this object. 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)
 
- Static Public Attributes inherited from Teuchos::Describable
static const EVerbosityLevel verbLevel_default
 
- Protected Member Functions inherited from Belos::Impl::SolverFactoryParent< Scalar, MV, OP >
virtual Teuchos::RCP
< solver_base_type
getSolver (const std::string &solverName, const Teuchos::RCP< Teuchos::ParameterList > &solverParams)
 Return an instance of the specified solver, or Teuchos::null if this factory does not provide the requested solver. More...
 

Detailed Description

template<class Scalar, class MV, class OP>
class Belos::SolverFactory< Scalar, MV, OP >

Factory for all solvers which Belos supports.

Author
Mark Hoemmen

New Belos users should start by creating an instance of this class, and using it to create the solver they want.

Belos implements several different iterative solvers. The usual way in which users interact with these solvers is through appropriately named subclasses of SolverManager. This factory class tells users which solvers are supported. It can initialize and return any supported subclass of SolverManager, given a short name of the subclass (such as "GMRES" or "CG").

Users ask for the solver they want by a string name, and supply an optional (but recommended) list of parameters (Teuchos::ParameterList) for the solver. The solver may fill in the parameter list with all the valid parameters and their default values, which users may later inspect and modify. Valid solver names include both "canonical names" (each maps one-to-one to a specific SolverManager subclass) and "aliases." Some aliases are short nicknames for canonical names, like "GMRES" for "Pseudoblock GMRES". Other aliases refer to a canonical solver name, but also modify the user's parameter list. For example, "Flexible GMRES" is an alias for "Block GMRES", and also sets the "Flexible Gmres" parameter to true in the input parameter list.

Mapping of solver names and aliases to Belos classes
Solver name Aliases SolverManager subclass
Pseudoblock GMRES GMRES, Pseudo Block GMRES, PseudoBlockGMRES, PseudoBlockGmres PseudoBlockGmresSolMgr
Block GMRES Flexible GMRES BlockGmresSolMgr
Block CG Block CG BlockCGSolMgr
Pseudoblock CG PseudoBlockCG, Pseudo Block CG PseudoBlockCGSolMgr
Pseudoblock Stochastic CG Stochastic CG PseudoBlockStochasticCGSolMgr
GCRODR Recycling GMRES GCRODRSolMgr
RCG Recycling CG RCGSolMgr
MINRES MINRES MinresSolMgr
LSQR LSQR LSQRSolMgr
TFQMR TFQMR, Transpose-Free QMR TFQMRSolMgr
Pseudoblock TFQMR Pseudoblock TFQMR, Pseudo Block Transpose-Free QMR PseudoBlockTFQMRSolMgr
Hybrid Block GMRES GmresPoly, Seed GMRES GmresPolySolMgr
PCPG CGPoly, Seed CG PCPGSolMgr

This class' template parameters are the same as those of Belos::SolverManager. Scalar is the scalar type (of entries in the multivector), MV is the multivector type, and OP is the operator type. For example: Scalar=double, MV=Epetra_MultiVector, and OP=Epetra_Operator will access the Epetra specialization of the Belos solvers.

Here is a simple example of how to use SolverFactory to create a GMRES solver for your linear system. Your code needs to include BelosSolverFactory.hpp and whatever linear algebra library header files you would normally use. Suppose that Scalar, MV, and OP have been previously typedef'd to the scalar resp. multivector resp. operator type in your application.

using Teuchos::parameterList;
using Teuchos::rcp; // Save some typing
// The ellipses represent the code you would normally use to create
// the sparse matrix, preconditioner, right-hand side, and initial
// guess for the linear system AX=B you want to solve.
RCP<OP> A = ...; // The sparse matrix / operator A
RCP<OP> M = ...; // The (right) preconditioner M
RCP<MV> B = ...; // Right-hand side of AX=B
RCP<MV> X = ...; // Initial guess for the solution
// Make an empty new parameter list.
RCP<ParameterList> solverParams = parameterList();
// Set some GMRES parameters.
//
// "Num Blocks" = Maximum number of Krylov vectors to store. This
// is also the restart length. "Block" here refers to the ability
// of this particular solver (and many other Belos solvers) to solve
// multiple linear systems at a time, even though we are only solving
// one linear system in this example.
solverParams->set ("Num Blocks", 40);
solverParams->set ("Maximum Iterations", 400);
solverParams->set ("Convergence Tolerance", 1.0e-8);
// Create the GMRES solver.
RCP<Belos::SolverManager<Scalar, MV, OP> > solver =
factory.create ("GMRES", solverParams);
// Create a LinearProblem struct with the problem to solve.
// A, X, B, and M are passed by (smart) pointer, not copied.
RCP<Belos::LinearProblem<Scalar, MV, OP> > problem =
problem->setRightPrec (M);
// Tell the solver what problem you want to solve.
solver->setProblem (problem);
// Attempt to solve the linear system. result == Belos::Converged
// means that it was solved to the desired tolerance. This call
// overwrites X with the computed approximate solution.
Belos::ReturnType result = solver->solve();
// Ask the solver how many iterations the last solve() took.
const int numIters = solver->getNumIters();

Belos developers who have implemented a new solver (i.e., a new subclass of SolverManager) and who want to make the solver available through the factory should do the following:

  1. Add a new symbol corresponding to their solver to the Details::EBelosSolverType enum.
  2. If necessary, specialize Details::makeSolverManagerTmpl for their SolverManager subclass. In most cases, the default implementation suffices.
  3. Add a case for their enum symbol that instantiates their solver to the long switch-case statement in Details::makeSolverManagerFromEnum.
  4. In the SolverFactory constructor, define a canonical string name for their solver and its mapping to the corresponding enum value, following the examples and comments there. (This takes one line of code.)

Definition at line 373 of file BelosSolverFactory.hpp.

Member Typedef Documentation

template<class Scalar, class MV, class OP>
typedef parent_type::solver_base_type Belos::SolverFactory< Scalar, MV, OP >::solver_base_type

The type of the solver returned by create().

This is a specialization of SolverManager for the same scalar, multivector, and operator types as the template parameters of this factory.

Definition at line 384 of file BelosSolverFactory.hpp.

template<class Scalar, class MV, class OP>
typedef parent_type::custom_solver_factory_type Belos::SolverFactory< Scalar, MV, OP >::custom_solver_factory_type

The type of a solver factory that users may give to addFactory() (which see).

Definition at line 389 of file BelosSolverFactory.hpp.


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

Generated on Wed Jun 20 2018 14:44:43 for Belos by doxygen 1.8.6