Belos  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
BelosOrthoManagerFactory.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Belos: Block Linear Solvers Package
4 //
5 // Copyright 2004-2016 NTESS and the Belos contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef __Belos_OrthoManagerFactory_hpp
11 #define __Belos_OrthoManagerFactory_hpp
12 
13 #include <BelosConfigDefs.hpp>
14 #ifdef HAVE_BELOS_TSQR
15 # include <BelosTsqrOrthoManager.hpp>
16 #endif // HAVE_BELOS_TSQR
21 #include <BelosOutputManager.hpp>
22 
23 #include <Teuchos_StandardCatchMacros.hpp>
24 
25 #include <algorithm>
26 #include <sstream>
27 #include <stdexcept> // #include <string>
28 #include <vector>
29 
31 
32 namespace Belos {
33 
49  template<class Scalar, class MV, class OP>
51  private:
53  std::vector<std::string> theList_;
54 
55  public:
57  static int numOrthoManagers () {
58 #ifdef HAVE_BELOS_TSQR
59  return 5;
60 #else
61  return 4;
62 #endif // HAVE_BELOS_TSQR
63  }
64 
70  static bool isRankRevealing (const std::string& name) {
71 #ifdef HAVE_BELOS_TSQR
72  // Currently only TSQR has a full rank-revealing capability.
73  return (name == "TSQR");
74 #else
75  return false;
76 #endif // HAVE_BELOS_TSQR
77  }
78 
81  {
82  int index = 0;
83  theList_[index++] = "ICGS";
84  theList_[index++] = "IMGS";
85  theList_[index++] = "DGKS";
86 #ifdef HAVE_BELOS_TSQR
87  theList_[index++] = "TSQR";
88 #endif // HAVE_BELOS_TSQR
89  theList_[index++] = "Simple";
90  }
91 
100  const std::vector<std::string>&
101  validNames () const { return theList_; }
102 
104  bool
105  isValidName (const std::string& name) const
106  {
107  return (std::find (theList_.begin(), theList_.end(), name) != theList_.end());
108  }
109 
111  std::ostream&
112  printValidNames (std::ostream& out) const
113  {
114  const int numValid = numOrthoManagers();
115  TEUCHOS_TEST_FOR_EXCEPTION(numValid <= 0, std::logic_error,
116  "Invalid number " << numValid << " of valid MatOrtho"
117  "Manager names. Please report this bug to the Belos "
118  "developers." );
119  if (numValid > 1) {
120  for (int k = 0; k < numValid - 1; ++k)
121  out << "\"" << theList_[k] << "\", ";
122  out << "or ";
123  }
124  out << "\"" << theList_[numValid-1] << "\"";
125  return out;
126  }
127 
132  std::string
134  {
135  std::ostringstream os;
136  (void) printValidNames (os);
137  return os.str();
138  }
139 
146  const std::string& defaultName () const { return theList_[0]; }
147 
158  getDefaultParameters (const std::string& name) const
159  {
160  if (name == "DGKS") {
161  return Belos::getDGKSDefaultParameters<Scalar, MV, OP> ();
162  }
163 #ifdef HAVE_BELOS_TSQR
164  else if (name == "TSQR") {
166  return orthoMan.getValidParameters ();
167  }
168 #endif // HAVE_BELOS_TSQR
169  else if (name == "ICGS") {
170  return Belos::getICGSDefaultParameters<Scalar, MV, OP> ();
171  }
172  else if (name == "IMGS") {
173  return Belos::getIMGSDefaultParameters<Scalar, MV, OP> ();
174  }
175  else if (name == "Simple") {
177  return orthoMan.getValidParameters ();
178  }
179  else {
180  TEUCHOS_TEST_FOR_EXCEPTION(true, std::invalid_argument,
181  "Invalid orthogonalization manager name \"" << name
182  << "\": Valid names are " << validNamesString()
183  << ". For many of the test executables, the "
184  "orthogonalization manager name often corresponds "
185  "to the \"ortho\" command-line argument.");
186  // Placate the compiler if necessary; we should never reach
187  // this point.
188  TEUCHOS_UNREACHABLE_RETURN(Teuchos::null);
189  }
190  }
191 
206  getFastParameters (const std::string& name) const
207  {
208  if (name == "DGKS") {
209  return Belos::getDGKSFastParameters<Scalar, MV, OP> ();
210  }
211 #ifdef HAVE_BELOS_TSQR
212  else if (name == "TSQR") {
214  return orthoMan.getFastParameters ();
215  }
216 #endif // HAVE_BELOS_TSQR
217  else if (name == "ICGS") {
218  return Belos::getICGSFastParameters<Scalar, MV, OP> ();
219  }
220  else if (name == "IMGS") {
221  return Belos::getIMGSFastParameters<Scalar, MV, OP> ();
222  }
223  else if (name == "Simple") {
225  return orthoMan.getFastParameters ();
226  }
227  else {
228  TEUCHOS_TEST_FOR_EXCEPTION(true, std::invalid_argument,
229  "Invalid orthogonalization manager name \"" << name
230  << "\": Valid names are " << validNamesString()
231  << ". For many of the test executables, the "
232  "orthogonalization manager name often corresponds "
233  "to the \"ortho\" command-line argument.");
234  // Placate the compiler if necessary; we should never reach
235  // this point.
236  TEUCHOS_UNREACHABLE_RETURN(Teuchos::null);
237  }
238  }
239 
260  makeMatOrthoManager (const std::string& ortho,
261  const Teuchos::RCP<const OP>& M,
262  const Teuchos::RCP<OutputManager<Scalar> >& /* outMan */,
263  const std::string& label,
265  {
266 #ifdef HAVE_BELOS_TSQR
268 #endif // HAVE_BELOS_TSQR
273  using Teuchos::rcp;
274 
275  if (ortho == "DGKS") {
276  typedef DGKSOrthoManager<Scalar, MV, OP> ortho_type;
277  return rcp (new ortho_type (params, label, M));
278  }
279 #ifdef HAVE_BELOS_TSQR
280  else if (ortho == "TSQR") {
281  typedef TsqrMatOrthoManager<Scalar, MV, OP> ortho_type;
282  return rcp (new ortho_type (params, label, M));
283  }
284 #endif // HAVE_BELOS_TSQR
285  else if (ortho == "ICGS") {
286  typedef ICGSOrthoManager<Scalar, MV, OP> ortho_type;
287  return rcp (new ortho_type (params, label, M));
288  }
289  else if (ortho == "IMGS") {
290  typedef IMGSOrthoManager<Scalar, MV, OP> ortho_type;
291  return rcp (new ortho_type (params, label, M));
292  }
293  else if (ortho == "Simple") {
294  TEUCHOS_TEST_FOR_EXCEPTION(ortho == "Simple", std::logic_error,
295  "SimpleOrthoManager does not yet support "
296  "the MatOrthoManager interface");
297  }
298  TEUCHOS_TEST_FOR_EXCEPTION(true, std::invalid_argument,
299  "Invalid orthogonalization manager name: Valid names"
300  " are " << validNamesString() << ". For many of "
301  "the test executables, the orthogonalization manager"
302  " name often corresponds to the \"ortho\" command-"
303  "line argument.");
304  TEUCHOS_UNREACHABLE_RETURN(Teuchos::null); // Guard to avoid compiler warnings.
305  }
306 
324  makeOrthoManager (const std::string& ortho,
325  const Teuchos::RCP<const OP>& M,
326  const Teuchos::RCP<OutputManager<Scalar> >& outMan,
327  const std::string& label,
329  {
330 #ifdef HAVE_BELOS_TSQR
332 #endif // HAVE_BELOS_TSQR
333  using Teuchos::rcp;
334 
335  if (ortho == "Simple") {
336  TEUCHOS_TEST_FOR_EXCEPTION(! M.is_null(), std::logic_error,
337  "SimpleOrthoManager is not yet supported "
338  "when the operator M is nontrivial (i.e., "
339  "M != null).");
340  return rcp (new SimpleOrthoManager<Scalar, MV> (outMan, label, params));
341  }
342 #ifdef HAVE_BELOS_TSQR
343  // TsqrMatOrthoManager has to store more things and do more work
344  // than TsqrOrthoManager, in order for the former to be correct
345  // for the case of a nondefault (non-Euclidean) inner product.
346  // Thus, it's better to create a TsqrOrthoManager, when we know
347  // the operator is the default operator (M is null). Of course,
348  // a MatOrthoManager is-an OrthoManager, so returning a
349  // TsqrMatOrthoManager would still be correct; this is just an
350  // optimization.
351  else if (ortho == "TSQR" && M.is_null()) {
352  return rcp (new TsqrOrthoManager<Scalar, MV> (params, label));
353  }
354 #endif // HAVE_BELOS_TSQR
355  else {
356  // A MatOrthoManager is-an OrthoManager.
357  return makeMatOrthoManager (ortho, M, outMan, label, params);
358  }
359  }
360  };
361 
362 } // namespace Belos
363 
364 #endif // __Belos_OrthoManagerFactory_hpp
365 
Belos&#39;s basic output manager for sending information of select verbosity levels to the appropriate ou...
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 > &params)
Return an instance of the specified OrthoManager subclass.
Class which manages the output and verbosity of the Belos solvers.
static bool isRankRevealing(const std::string &name)
Is the given MatOrthoManager subclass rank-reealing?
Teuchos::RCP< const Teuchos::ParameterList > getFastParameters()
Get &quot;fast&quot; parameters for TsqrMatOrthoManager.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Teuchos::RCP< const Teuchos::ParameterList > getDefaultParameters(const std::string &name) const
Default parameters for the given MatOrthoManager subclass.
static int numOrthoManagers()
Number of MatOrthoManager subclasses this factory recognizes.
std::string validNamesString() const
List (as a string) of recognized MatOrthoManager names.
Iterated Modified Gram-Schmidt (IMGS) implementation of the Belos::OrthoManager class.
Simple OrthoManager implementation for benchmarks.
Classical Gram-Schmidt (with DGKS correction) implementation of the Belos::OrthoManager class...
Teuchos::RCP< const Teuchos::ParameterList > getFastParameters()
Get a &quot;fast&quot; list of parameters.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Simple OrthoManager implementation for benchmarks.
bool isValidName(const std::string &name) const
Whether this factory recognizes the MatOrthoManager with the given name.
std::ostream & printValidNames(std::ostream &out) const
Print all recognized MatOrthoManager names to the given ostream.
const std::string & defaultName() const
Name of the &quot;default&quot; MatOrthoManager subclass.
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Get a default list of parameters.
Enumeration of all valid Belos (Mat)OrthoManager classes.
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Get default parameters for TsqrMatOrthoManager.
Iterated Classical Gram-Schmidt (ICGS) implementation of the Belos::OrthoManager class.
An implementation of the Belos::MatOrthoManager that performs orthogonalization using multiple steps ...
An implementation of the Belos::MatOrthoManager that performs orthogonalization using multiple steps ...
#define TEUCHOS_UNREACHABLE_RETURN(dummyReturnVal)
An implementation of the Belos::MatOrthoManager that performs orthogonalization using (potentially) m...
Orthogonalization manager based on Tall Skinny QR (TSQR)
Teuchos::RCP< const Teuchos::ParameterList > getFastParameters(const std::string &name) const
&quot;Fast&quot; parameters for the given MatOrthoManager subclass.
const std::vector< std::string > & validNames() const
List of MatOrthoManager subclasses this factory recognizes.
MatOrthoManager subclass using TSQR or DGKS.
Belos header file which uses auto-configuration information to include necessary C++ headers...
TSQR-based OrthoManager subclass.
bool is_null() const
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 > &params)
Return an instance of the specified MatOrthoManager subclass.

Generated on Fri Nov 22 2024 09:23:06 for Belos by doxygen 1.8.5