10#ifndef __Belos_OrthoManagerFactory_hpp
11#define __Belos_OrthoManagerFactory_hpp
23#include <Teuchos_StandardCatchMacros.hpp>
49 template<
class Scalar,
class MV,
class OP>
53 std::vector<std::string> theList_;
73 return (
name ==
"TSQR");
83 theList_[index++] =
"ICGS";
84 theList_[index++] =
"IMGS";
85 theList_[index++] =
"DGKS";
87 theList_[index++] =
"TSQR";
89 theList_[index++] =
"Simple";
100 const std::vector<std::string>&
107 return (std::find (theList_.begin(), theList_.end(),
name) != theList_.end());
116 "Invalid number " <<
numValid <<
" of valid MatOrtho"
117 "Manager names. Please report this bug to the Belos "
121 out <<
"\"" << theList_[
k] <<
"\", ";
135 std::ostringstream
os;
157 Teuchos::RCP<const Teuchos::ParameterList>
160 if (
name ==
"DGKS") {
163#ifdef HAVE_BELOS_TSQR
164 else if (
name ==
"TSQR") {
166 return orthoMan.getValidParameters ();
169 else if (
name ==
"ICGS") {
172 else if (
name ==
"IMGS") {
175 else if (
name ==
"Simple") {
177 return orthoMan.getValidParameters ();
181 "Invalid orthogonalization manager name \"" <<
name
183 <<
". For many of the test executables, the "
184 "orthogonalization manager name often corresponds "
185 "to the \"ortho\" command-line argument.");
205 Teuchos::RCP<const Teuchos::ParameterList>
208 if (
name ==
"DGKS") {
211#ifdef HAVE_BELOS_TSQR
212 else if (
name ==
"TSQR") {
214 return orthoMan.getFastParameters ();
217 else if (
name ==
"ICGS") {
220 else if (
name ==
"IMGS") {
223 else if (
name ==
"Simple") {
225 return orthoMan.getFastParameters ();
229 "Invalid orthogonalization manager name \"" <<
name
231 <<
". For many of the test executables, the "
232 "orthogonalization manager name often corresponds "
233 "to the \"ortho\" command-line argument.");
259 Teuchos::RCP<Belos::MatOrthoManager<Scalar, MV, OP> >
261 const Teuchos::RCP<const OP>&
M,
263 const std::string&
label,
264 const Teuchos::RCP<Teuchos::ParameterList>&
params)
266#ifdef HAVE_BELOS_TSQR
275 if (
ortho ==
"DGKS") {
279#ifdef HAVE_BELOS_TSQR
280 else if (
ortho ==
"TSQR") {
285 else if (
ortho ==
"ICGS") {
289 else if (
ortho ==
"IMGS") {
293 else if (
ortho ==
"Simple") {
295 "SimpleOrthoManager does not yet support "
296 "the MatOrthoManager interface");
299 "Invalid orthogonalization manager name: Valid names"
301 "the test executables, the orthogonalization manager"
302 " name often corresponds to the \"ortho\" command-"
323 Teuchos::RCP<Belos::OrthoManager<Scalar, MV> >
325 const Teuchos::RCP<const OP>&
M,
327 const std::string&
label,
328 const Teuchos::RCP<Teuchos::ParameterList>&
params)
330#ifdef HAVE_BELOS_TSQR
335 if (
ortho ==
"Simple") {
337 "SimpleOrthoManager is not yet supported "
338 "when the operator M is nontrivial (i.e., "
342#ifdef HAVE_BELOS_TSQR
351 else if (
ortho ==
"TSQR" &&
M.is_null()) {
Belos header file which uses auto-configuration information to include necessary C++ headers.
Classical Gram-Schmidt (with DGKS correction) implementation of the Belos::OrthoManager class.
Iterated Classical Gram-Schmidt (ICGS) implementation of the Belos::OrthoManager class.
Iterated Modified Gram-Schmidt (IMGS) implementation of the Belos::OrthoManager class.
Class which manages the output and verbosity of the Belos solvers.
Simple OrthoManager implementation for benchmarks.
Orthogonalization manager based on Tall Skinny QR (TSQR)
An implementation of the Belos::MatOrthoManager that performs orthogonalization using (potentially) m...
An implementation of the Belos::MatOrthoManager that performs orthogonalization using multiple steps ...
An implementation of the Belos::MatOrthoManager that performs orthogonalization using multiple steps ...
Alternative run-time polymorphic interface for operators.
Enumeration of all valid Belos (Mat)OrthoManager classes.
Teuchos::RCP< Belos::OrthoManager< Scalar, MV > > makeOrthoManager(const std::string &ortho, const Teuchos::RCP< const OP > &M, const Teuchos::RCP< OutputManager< Scalar > > &outMan, const std::string &label, const Teuchos::RCP< Teuchos::ParameterList > ¶ms)
Return an instance of the specified OrthoManager subclass.
const std::vector< std::string > & validNames() const
List of MatOrthoManager subclasses this factory recognizes.
std::ostream & printValidNames(std::ostream &out) const
Print all recognized MatOrthoManager names to the given ostream.
OrthoManagerFactory()
Constructor.
const std::string & defaultName() const
Name of the "default" MatOrthoManager subclass.
bool isValidName(const std::string &name) const
Whether this factory recognizes the MatOrthoManager with the given name.
static int numOrthoManagers()
Number of MatOrthoManager subclasses this factory recognizes.
std::string validNamesString() const
List (as a string) of recognized MatOrthoManager names.
Teuchos::RCP< Belos::MatOrthoManager< Scalar, MV, OP > > makeMatOrthoManager(const std::string &ortho, const Teuchos::RCP< const OP > &M, const Teuchos::RCP< OutputManager< Scalar > > &, const std::string &label, const Teuchos::RCP< Teuchos::ParameterList > ¶ms)
Return an instance of the specified MatOrthoManager subclass.
Teuchos::RCP< const Teuchos::ParameterList > getDefaultParameters(const std::string &name) const
Default parameters for the given MatOrthoManager subclass.
Teuchos::RCP< const Teuchos::ParameterList > getFastParameters(const std::string &name) const
"Fast" parameters for the given MatOrthoManager subclass.
static bool isRankRevealing(const std::string &name)
Is the given MatOrthoManager subclass rank-reealing?
Simple OrthoManager implementation for benchmarks.
MatOrthoManager subclass using TSQR or DGKS.
TSQR-based OrthoManager subclass.