42 #ifndef __Belos_SolverFactory_hpp 
   43 #define __Belos_SolverFactory_hpp 
   52 #include <Teuchos_Describable.hpp> 
   53 #include <Teuchos_StandardCatchMacros.hpp> 
   72                   const std::vector<std::string>& array);
 
   94 template<
class Scalar, 
class MV, 
class OP>
 
  162   getSolver (
const std::string& solverName,
 
  176   create (
const std::string& solverName,
 
  194   virtual bool isSupported (
const std::string& solverName) 
const;
 
  219       instance == Teuchos::null,
 
  220       std::invalid_argument, 
"Belos::SolverFactoryParent::registerSolver " 
  221       "was given a null solver to register.");
 
  223     get_solverManagers()[solverName] = instance;
 
  230     return (get_solverManagers().find(solverName) != get_solverManagers().end());
 
  252   static std::vector<Teuchos::RCP<custom_solver_factory_type> > factories_;
 
  254   static std::map<
const std::string, 
Teuchos::RCP<
typename 
  256     get_solverManagers() {
 
  257       static std::map<
const std::string, 
Teuchos::RCP<
typename 
  259       return solverManagers;
 
  263 template<
class Scalar, 
class MV, 
class OP>
 
  264 std::vector<Teuchos::RCP<typename SolverFactoryParent<Scalar, MV, OP>::custom_solver_factory_type> >
 
  265 SolverFactoryParent<Scalar, MV, OP>::factories_;
 
  267 template<
class SolverClass, 
class Scalar, 
class MV, 
class OP>
 
  283 template<
class SC, 
class MV, 
class OP>
 
  302 template<
class SC, 
class MV, 
class OP>
 
  303 using SolverFactory = typename ::Belos::Impl::SolverFactorySelector<SC, MV, OP>::type;
 
  307 template<
class Scalar, 
class MV, 
class OP>
 
  314   RCP<solver_base_type> solver = this->getSolver (solverName, solverParams);
 
  316     (solver.is_null (), std::invalid_argument,
 
  317      "Invalid or unsupported Belos solver name \"" << solverName << 
"\".");
 
  321 template<
class Scalar, 
class MV, 
class OP>
 
  330   for (std::size_t k = 0; k < factories_.size (); ++k) {
 
  331     RCP<CustomSolverFactory<Scalar, MV, OP> > factory = factories_[k];
 
  332     if (! factory.is_null ()) {
 
  333       RCP<SolverManager<Scalar, MV, OP> > solver =
 
  334         factory->getSolver (solverName, solverParams);
 
  335       if (! solver.is_null ()) {
 
  345   std::pair<std::string, bool> aliasResult =
 
  347   const std::string candidateCanonicalName = aliasResult.first;
 
  348   const bool isAnAlias = aliasResult.second;
 
  351   std::string standardized_name = isAnAlias ?
 
  352                                   candidateCanonicalName :
 
  360     solverParams.
is_null() ? Teuchos::parameterList() : solverParams;
 
  369     it = get_solverManagers().find (standardized_name);
 
  372     it == get_solverManagers().end(),
 
  373     std::invalid_argument, 
"Belos solver manager " << solverNameUC <<
 
  374     " with standardized name " << standardized_name << 
" has not been" 
  378     it->second == Teuchos::null,
 
  379     std::logic_error, 
"Belos::SolverFactoryParent: The registered " 
  380     "clone source for " << solverNameUC << 
" with standardized name " 
  381     << standardized_name << 
" is null which should never happen." 
  382     ".  Please report this bug to the Belos developers.");
 
  385   RCP<solver_base_type> solver = (it->second)->clone ();
 
  388     solver == Teuchos::null,
 
  389     std::logic_error, 
"Belos::SolverFactoryParent: Failed " 
  390     "to clone SolverManager with name " << solverNameUC << 
" with standardized" 
  391     " name" << standardized_name << 
"." 
  392     ".  Please report this bug to the Belos developers.");
 
  399     pl = Teuchos::parameterList (solver->getValidParameters ()->name ());
 
  403     pl.
is_null(), std::logic_error,
 
  404     "Belos::SolverFactory: ParameterList to pass to solver is null.  This " 
  405     "should never happen.  Please report this bug to the Belos developers.");
 
  406   solver->setParameters (pl);
 
  411 template<
class Scalar, 
class MV, 
class OP>
 
  416   factories_.push_back (factory);
 
  420 template<
class Scalar, 
class MV, 
class OP>
 
  427   std::ostringstream out;
 
  428   out << 
"\"Belos::SolverFactory\": {";
 
  429   if (this->getObjectLabel () != 
"") {
 
  430     out << 
"Label: " << this->getObjectLabel () << 
", ";
 
  432   out << 
"Scalar: \"" << TypeNameTraits<Scalar>::name ()
 
  433       << 
"\", MV: \"" << TypeNameTraits<MV>::name ()
 
  434       << 
"\", OP: \"" << TypeNameTraits<OP>::name ()
 
  440 template<
class Scalar, 
class MV, 
class OP>
 
  462   out << 
"\"Belos::SolverFactory\":" << endl;
 
  463   if (this->getObjectLabel () != 
"") {
 
  464     out << 
"Label: " << this->getObjectLabel () << endl;
 
  467     out << 
"Template parameters:" << endl;
 
  469     out << 
"Scalar: \"" << TypeNameTraits<Scalar>::name () << 
"\"" << endl
 
  470         << 
"MV: \"" << TypeNameTraits<MV>::name () << 
"\"" << endl
 
  471         << 
"OP: \"" << TypeNameTraits<OP>::name () << 
"\"" << endl;
 
  479     out << 
"Canonical solver names: ";
 
  483     out << 
"Aliases to canonical names: ";
 
  489 template<
class Scalar, 
class MV, 
class OP>
 
  494   int numSupported = 0;
 
  497   for (std::size_t k = 0; k < factories_.size (); ++k) {
 
  499     RCP<custom_solver_factory_type> factory = factories_[k];
 
  500     if (! factory.is_null ()) {
 
  501       numSupported += factory->numSupportedSolvers ();
 
  509 template<
class Scalar, 
class MV, 
class OP>
 
  514   typedef std::vector<std::string>::const_iterator iter_type;
 
  518   const std::size_t numFactories = factories_.
size ();
 
  519   for (std::size_t factInd = 0; factInd < numFactories; ++factInd) {
 
  522       std::vector<std::string> supportedSolvers =
 
  523         factory->supportedSolverNames ();
 
  524       const std::size_t numSolvers = supportedSolvers.size ();
 
  525       for (std::size_t solvInd = 0; solvInd < numSolvers; ++solvInd) {
 
  526         names.
push_back (supportedSolvers[solvInd]);
 
  533     for (iter_type iter = aliases.begin (); iter != aliases.end (); ++iter) {
 
  539     for (iter_type iter = canonicalNames.begin ();
 
  540          iter != canonicalNames.end (); ++iter) {
 
  547 template<
class Scalar, 
class MV, 
class OP>
 
  553   const std::size_t numFactories = factories_.size ();
 
  554   for (std::size_t factInd = 0; factInd < numFactories; ++factInd) {
 
  556     RCP<custom_solver_factory_type> factory = factories_[factInd];
 
  557     if (! factory.is_null ()) {
 
  558       if (factory->isSupported (solverName)) {
 
  569   std::pair<std::string, bool> aliasResult =
 
  571   const std::string candidateCanonicalName = aliasResult.first;
 
  572   const bool validCanonicalName =
 
  573     (get_solverManagers().find(candidateCanonicalName) != get_solverManagers().end());
 
  574   return validCanonicalName;
 
  587 #endif // __Belos_SolverFactory_hpp 
Interface for custom Belos solver factories. 
 
Class which manages the output and verbosity of the Belos solvers. 
 
SolverFactoryParent< SC, MV, OP > type
 
static void registerSolver(const std::string &solverName, Teuchos::RCP< SolverFactoryParent< Scalar, MV, OP >::solver_base_type > instance)
register a solver for Inverted Injection (DII). 
 
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
 
std::string upperCase(const std::string &s)
Return the upper-case version of s. 
 
virtual int numSupportedSolvers() const 
Number of supported solvers. 
 
void addFactory(const Teuchos::RCP< custom_solver_factory_type > &factory)
Add a custom solver factory. 
 
Declaration of alias functions for solver names. 
 
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const 
Describe this object. 
 
Pure virtual base class which describes the basic interface for a solver manager. ...
 
std::vector< std::string > solverNameAliases()
List of supported aliases (to canonical solver names). 
 
void printStringArray(std::ostream &out, const Teuchos::ArrayView< const std::string > &array)
Print the given array of strings, in YAML format, to out. 
 
int numSupportedSolvers()
Number of Belos solvers supported for any linear algebra implementation ("generically"). 
 
The Belos::SolverManager is a templated virtual base class that defines the basic interface that any ...
 
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 req...
 
static const EVerbosityLevel verbLevel_default
 
virtual Teuchos::Array< std::string > supportedSolverNames() const 
List of supported solver names. 
 
void push_back(const value_type &x)
 
virtual Teuchos::RCP< solver_base_type > create(const std::string &solverName, const Teuchos::RCP< Teuchos::ParameterList > &solverParams)
Create, configure, and return the specified solver. 
 
Specializations of Belos::SolverFactory may inherit from this class to get basic SolverFactory functi...
 
std::pair< std::string, bool > getCanonicalNameFromAlias(const std::string &candidateAlias)
Get the candidate canonical name for a given candidate alias. 
 
::Belos::SolverManager< Scalar, MV, OP > solver_base_type
The type of the solver returned by create(). 
 
void reviseParameterListForAlias(const std::string &aliasName, Teuchos::ParameterList &solverParams)
Modify the input ParameterList appropriately for the given solver alias. 
 
void registerSolverSubclassForTypes(const std::string &solverName)
 
virtual bool isSupported(const std::string &solverName) const 
Whether the given solver name names a supported solver. 
 
Belos header file which uses auto-configuration information to include necessary C++ headers...
 
typename::Belos::Impl::SolverFactorySelector< SC, MV, OP >::type SolverFactory
 
virtual std::string description() const 
A string description of this object. 
 
static bool isSolverRegistered(const std::string &solverName)
is solver registered for Inverted Injection (DII). 
 
CustomSolverFactory< Scalar, MV, OP > custom_solver_factory_type
The type of a solver factory that users may give to addFactory() (which see below) ...
 
std::vector< std::string > canonicalSolverNames()
List of canonical solver names.