Belos Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
BelosGmresPolySolMgr.hpp
Go to the documentation of this file.
1 //@HEADER
2 // ************************************************************************
3 //
4 // Belos: Block Linear Solvers Package
5 // Copyright 2004 Sandia Corporation
6 //
7 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
8 // the U.S. Government retains certain rights in this software.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ************************************************************************
40 //@HEADER
41 //
42 
43 #ifndef BELOS_GMRES_POLY_SOLMGR_HPP
44 #define BELOS_GMRES_POLY_SOLMGR_HPP
45 
49 
50 #include "BelosConfigDefs.hpp"
51 #include "BelosTypes.hpp"
52 
53 #include "BelosLinearProblem.hpp"
54 #include "BelosSolverManager.hpp"
55 #include "BelosGmresPolyOp.hpp"
57 #include "Teuchos_as.hpp"
58 #ifdef BELOS_TEUCHOS_TIME_MONITOR
59 #include "Teuchos_TimeMonitor.hpp"
60 #endif
61 
62 
63 namespace Belos {
64 
66 
67 
75  GmresPolySolMgrLinearProblemFailure(const std::string& what_arg) : BelosError(what_arg)
76  {}};
77 
85  GmresPolySolMgrPolynomialFailure(const std::string& what_arg) : BelosError(what_arg)
86  {}};
87 
102 //
135 //
147 
148 template<class ScalarType, class MV, class OP>
149 class GmresPolySolMgr : public SolverManager<ScalarType,MV,OP> {
150 private:
151 
157 
158 public:
159 
161 
162 
168  GmresPolySolMgr();
169 
190 
192  virtual ~GmresPolySolMgr() {};
193 
197  }
199 
201 
202 
205  const LinearProblem<ScalarType,MV,OP>& getProblem() const override {
206  return *problem_;
207  }
208 
212 
216 
223  return Teuchos::tuple(timerPoly_);
224  }
225 
227  int getNumIters() const override {
228  return numIters_;
229  }
230 
234  bool isLOADetected() const override { return loaDetected_; }
235 
237 
239 
240 
242  void setProblem( const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem ) override { problem_ = problem; }
243 
245  void setParameters( const Teuchos::RCP<Teuchos::ParameterList> &params ) override;
246 
248 
250 
259  void reset( const ResetType type ) override {
260  if ((type & Belos::Problem) && ! problem_.is_null ()) {
261  problem_->setProblem ();
263  poly_dim_ = 0; // Rebuild the GMRES polynomial
264  }
265  }
266 
268 
270 
288  ReturnType solve() override;
289 
291 
294 
296  std::string description() const override;
297 
299 
300 private:
301 
302  // Linear problem.
304 
305  // Output manager.
307 
308  // Current parameter list.
311 
312  // Default solver values.
313  static constexpr int maxDegree_default_ = 25;
314  static constexpr int verbosity_default_ = Belos::Errors;
315  static constexpr const char * label_default_ = "Belos";
316  static constexpr const char * outerSolverType_default_ = "";
317  static constexpr const char * polyType_default_ = "Arnoldi";
318  static constexpr const char * orthoType_default_ = "DGKS";
319  static constexpr bool addRoots_default_ = true;
320  static constexpr bool dampPoly_default_ = false;
321  static constexpr bool randomRHS_default_ = true;
322  static constexpr std::ostream * outputStream_default_ = &std::cout;
323 
324  // Current solver values.
330  bool damp_;
331  bool addRoots_;
332  std::string polyType_;
333  std::string outerSolverType_;
334  std::string orthoType_;
335 
336  // Polynomial storage
339 
340  // Timers.
341  std::string label_;
343 
344  // Internal state variables.
345  bool isSet_;
347 
350 };
351 
352 
353 template<class ScalarType, class MV, class OP>
355  outputStream_ (Teuchos::rcp(outputStream_default_,false)),
356  polyTol_ (DefaultSolverParameters::polyTol),
357  maxDegree_ (maxDegree_default_),
358  numIters_ (0),
359  verbosity_ (verbosity_default_),
360  hasOuterSolver_ (false),
361  randomRHS_ (randomRHS_default_),
362  damp_ (dampPoly_default_),
363  addRoots_ (addRoots_default_),
364  polyType_ (polyType_default_),
365  outerSolverType_ (outerSolverType_default_),
366  orthoType_ (orthoType_default_),
367  poly_dim_ (0),
368  label_ (label_default_),
369  isSet_ (false),
370  loaDetected_ (false)
371 {}
372 
373 
374 template<class ScalarType, class MV, class OP>
378  problem_ (problem),
379  outputStream_ (Teuchos::rcp(outputStream_default_,false)),
380  polyTol_ (DefaultSolverParameters::polyTol),
381  maxDegree_ (maxDegree_default_),
382  numIters_ (0),
383  verbosity_ (verbosity_default_),
384  hasOuterSolver_ (false),
385  randomRHS_ (randomRHS_default_),
386  damp_ (dampPoly_default_),
387  addRoots_ (addRoots_default_),
388  polyType_ (polyType_default_),
389  outerSolverType_ (outerSolverType_default_),
390  orthoType_ (orthoType_default_),
391  poly_dim_ (0),
392  label_ (label_default_),
393  isSet_ (false),
394  loaDetected_ (false)
395 {
397  problem_.is_null (), std::invalid_argument,
398  "Belos::GmresPolySolMgr: The given linear problem is null. "
399  "Please call this constructor with a nonnull LinearProblem argument, "
400  "or call the constructor that does not take a LinearProblem.");
401 
402  // If the input parameter list is null, then the parameters take
403  // default values.
404  if (! pl.is_null ()) {
405  setParameters (pl);
406  }
407 }
408 
409 
410 template<class ScalarType, class MV, class OP>
413 {
414  if (validPL_.is_null ()) {
415  Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList ();
416 
417  // The static_cast is to resolve an issue with older clang versions which
418  // would cause the constexpr to link fail. With c++17 the problem is resolved.
419  pl->set("Polynomial Type", static_cast<const char *>(polyType_default_),
420  "The type of GMRES polynomial that is used as a preconditioner: Roots, Arnoldi, or Gmres.");
421  pl->set("Polynomial Tolerance", static_cast<MagnitudeType>(DefaultSolverParameters::polyTol),
422  "The relative residual tolerance that used to construct the GMRES polynomial.");
423  pl->set("Maximum Degree", static_cast<int>(maxDegree_default_),
424  "The maximum degree allowed for any GMRES polynomial.");
425  pl->set("Outer Solver", static_cast<const char *>(outerSolverType_default_),
426  "The outer solver that this polynomial is used to precondition.");
427  pl->set("Verbosity", static_cast<int>(verbosity_default_),
428  "What type(s) of solver information should be outputted\n"
429  "to the output stream.");
430  pl->set("Output Stream", Teuchos::rcp(outputStream_default_,false),
431  "A reference-counted pointer to the output stream where all\n"
432  "solver output is sent.");
433  pl->set("Timer Label", static_cast<const char *>(label_default_),
434  "The string to use as a prefix for the timer labels.");
435  pl->set("Orthogonalization", static_cast<const char *>(orthoType_default_),
436  "The type of orthogonalization to use to generate polynomial: DGKS, ICGS, or IMGS.");
437  pl->set("Random RHS", static_cast<bool>(randomRHS_default_),
438  "Add roots to polynomial for stability.");
439  pl->set("Add Roots", static_cast<bool>(addRoots_default_),
440  "Add roots to polynomial for stability.");
441  pl->set("Damp Poly", static_cast<bool>(dampPoly_default_),
442  "Damp polynomial for ill-conditioned problems.");
443  validPL_ = pl;
444  }
445  return validPL_;
446 }
447 
448 
449 template<class ScalarType, class MV, class OP>
452 {
453  // Create the internal parameter list if ones doesn't already exist.
454  if (params_.is_null ()) {
455  params_ = Teuchos::parameterList (*getValidParameters ());
456  }
457  else {
458  params->validateParameters (*getValidParameters ());
459  }
460 
461  // Check which Gmres polynomial to use
462  if (params->isParameter("Polynomial Type")) {
463  polyType_ = params->get("Polynomial Type", polyType_default_);
464  }
465 
466  // Update the outer solver in our list.
467  params_->set("Polynomial Type", polyType_);
468 
469  // Check if there is an outer solver for this Gmres Polynomial
470  if (params->isParameter("Outer Solver")) {
471  outerSolverType_ = params->get("Outer Solver", outerSolverType_default_);
472  }
473 
474  // Update the outer solver in our list.
475  params_->set("Outer Solver", outerSolverType_);
476 
477  // Check if there is a parameter list for the outer solver
478  if (params->isSublist("Outer Solver Params")) {
479  outerParams_ = Teuchos::parameterList( params->get<Teuchos::ParameterList>("Outer Solver Params") );
480  }
481 
482  // Check for maximum polynomial degree
483  if (params->isParameter("Maximum Degree")) {
484  maxDegree_ = params->get("Maximum Degree",maxDegree_default_);
485  }
486 
487  // Update parameter in our list.
488  params_->set("Maximum Degree", maxDegree_);
489 
490  // Check to see if the timer label changed.
491  if (params->isParameter("Timer Label")) {
492  std::string tempLabel = params->get("Timer Label", label_default_);
493 
494  // Update parameter in our list and solver timer
495  if (tempLabel != label_) {
496  label_ = tempLabel;
497 #ifdef BELOS_TEUCHOS_TIME_MONITOR
498  std::string polyLabel = label_ + ": GmresPolyOp creation time";
499  timerPoly_ = Teuchos::TimeMonitor::getNewCounter(polyLabel);
500 #endif
501  }
502  }
503 
504  // Update timer label
505  params_->set("Timer Label", label_);
506 
507  // Check if the orthogonalization changed.
508  if (params->isParameter("Orthogonalization")) {
509  std::string tempOrthoType = params->get("Orthogonalization",orthoType_default_);
510  TEUCHOS_TEST_FOR_EXCEPTION( tempOrthoType != "DGKS" && tempOrthoType != "ICGS" && tempOrthoType != "IMGS",
511  std::invalid_argument,
512  "Belos::GmresPolySolMgr: \"Orthogonalization\" must be either \"DGKS\", \"ICGS\", or \"IMGS\".");
513  if (tempOrthoType != orthoType_) {
514  orthoType_ = tempOrthoType;
515  }
516  }
517 
518  params_->set("Orthogonalization", orthoType_);
519 
520  // Check for a change in verbosity level
521  if (params->isParameter("Verbosity")) {
522  if (Teuchos::isParameterType<int>(*params,"Verbosity")) {
523  verbosity_ = params->get("Verbosity", verbosity_default_);
524  } else {
525  verbosity_ = (int)Teuchos::getParameter<Belos::MsgType>(*params,"Verbosity");
526  }
527  }
528 
529  // Update parameter in our list.
530  params_->set("Verbosity", verbosity_);
531 
532  // output stream
533  if (params->isParameter("Output Stream")) {
534  outputStream_ = Teuchos::getParameter<Teuchos::RCP<std::ostream> >(*params,"Output Stream");
535  }
536 
537  // Update parameter in our list.
538  params_->set("Output Stream", outputStream_);
539 
540  // Convergence
541  // Check for polynomial convergence tolerance
542  if (params->isParameter("Polynomial Tolerance")) {
543  if (params->isType<MagnitudeType> ("Polynomial Tolerance")) {
544  polyTol_ = params->get ("Polynomial Tolerance",
545  static_cast<MagnitudeType> (DefaultSolverParameters::polyTol));
546  }
547  else {
548  polyTol_ = params->get ("Polynomial Tolerance", DefaultSolverParameters::polyTol);
549  }
550  }
551 
552  // Update parameter in our list and residual tests.
553  params_->set("Polynomial Tolerance", polyTol_);
554 
555  // Check for maximum polynomial degree
556  if (params->isParameter("Random RHS")) {
557  randomRHS_ = params->get("Random RHS",randomRHS_default_);
558  }
559 
560  // Update parameter in our list.
561  params_->set("Random RHS", randomRHS_);
562 
563 
564  // Check for polynomial damping
565  if (params->isParameter("Damped Poly")) {
566  damp_ = params->get("Damped Poly",dampPoly_default_);
567  }
568  // Update parameter in our list.
569  params_->set("Damped Poly", damp_);
570 
571  // Check: Should we add roots for stability if needed?
572  if (params->isParameter("Add Roots")) {
573  addRoots_ = params->get("Add Roots",addRoots_default_);
574  }
575 
576  // Update parameter in our list.
577  params_->set("Add Roots", addRoots_);
578 
579  // Create the timers if we need to.
580 #ifdef BELOS_TEUCHOS_TIME_MONITOR
581  if (timerPoly_ == Teuchos::null) {
582  std::string polyLabel = label_ + ": GmresPolyOp creation time";
583  timerPoly_ = Teuchos::TimeMonitor::getNewCounter(polyLabel);
584  }
585 #endif
586 
587  // Check if we are going to perform an outer solve.
588  if (outerSolverType_ != "") {
589  hasOuterSolver_ = true;
590  }
591 
592  // Inform the solver manager that the current parameters were set.
593  isSet_ = true;
594 }
595 
596 
597 template<class ScalarType, class MV, class OP>
599 {
600  using Teuchos::RCP;
601  using Teuchos::rcp;
602  using Teuchos::rcp_const_cast;
603 
604  // Assume convergence is achieved if user does not require strict convergence.
606 
607  // Set the current parameters if they were not set before. NOTE:
608  // This may occur if the user generated the solver manager with the
609  // default constructor and then didn't set any parameters using
610  // setParameters().
611  if (! isSet_) {
612  setParameters (Teuchos::parameterList (*getValidParameters ()));
613  }
614 
616  problem_.is_null (), GmresPolySolMgrLinearProblemFailure,
617  "Belos::GmresPolySolMgr::solve: The linear problem has not been set yet, "
618  "or was set to null. Please call setProblem() with a nonnull input before "
619  "calling solve().");
620 
622  ! problem_->isProblemSet (), GmresPolySolMgrLinearProblemFailure,
623  "Belos::GmresPolySolMgr::solve: The linear problem is not ready. Please "
624  "call setProblem() on the LinearProblem object before calling solve().");
625 
626  // If the GMRES polynomial has not been constructed for this
627  // (nmatrix, preconditioner) pair, generate it.
628  if (!poly_dim_ && maxDegree_) {
629 #ifdef BELOS_TEUCHOS_TIME_MONITOR
630  Teuchos::TimeMonitor slvtimer(*timerPoly_);
631 #endif
632  poly_Op_ = Teuchos::rcp( new gmres_poly_t( problem_, params_ ) );
633  poly_dim_ = poly_Op_->polyDegree();
634 
636  "Belos::GmresPolyOp: Failed to generate polynomial that satisfied requirements.");
637  }
638 
639 
640  // Solve the linear system using the polynomial
641  if (hasOuterSolver_ && maxDegree_) {
642 
643  // Then the polynomial will be used as an operator for an outer solver.
644  // Use outer solver parameter list passed in a sublist.
646  RCP<SolverManager<ScalarType, MultiVec<ScalarType>, Operator<ScalarType> > > solver = factory.create( outerSolverType_, outerParams_ );
647  TEUCHOS_TEST_FOR_EXCEPTION( solver == Teuchos::null, std::invalid_argument,
648  "Belos::GmresPolySolMgr::solve(): Selected solver is not valid.");
649 
650  // Create a copy of the linear problem that uses the polynomial as a preconditioner.
651  // The original initial solution and right-hand side are thinly wrapped in the gmres_poly_mv_t
652  RCP<gmres_poly_mv_t> new_lhs = rcp( new gmres_poly_mv_t( problem_->getLHS() ) );
653  RCP<gmres_poly_mv_t> new_rhs = rcp( new gmres_poly_mv_t( rcp_const_cast<MV>( problem_->getRHS() ) ) );
654  RCP<gmres_poly_t> A = rcp( new gmres_poly_t( problem_ ) ); // This just performs problem_->applyOp
655  RCP<LinearProblem<ScalarType,MultiVec<ScalarType>,Operator<ScalarType> > > newProblem =
656  rcp( new LinearProblem<ScalarType,MultiVec<ScalarType>,Operator<ScalarType> >( A, new_lhs, new_rhs ) );
657  std::string solverLabel = label_ + ": Hybrid Gmres";
658  newProblem->setLabel(solverLabel);
659 
660  // If the preconditioner is left preconditioner, use Gmres poly as a left preconditioner.
661  if (problem_->getLeftPrec() != Teuchos::null)
662  newProblem->setLeftPrec( poly_Op_ );
663  else
664  newProblem->setRightPrec( poly_Op_ );
665  // Set the initial residual vector, if it has already been set in the original problem.
666  // Don't set the preconditioned residual vector, because it is not the GmresPoly preconditioned residual vector.
667  if (problem_->getInitResVec() != Teuchos::null)
668  newProblem->setInitResVec( rcp( new gmres_poly_mv_t( rcp_const_cast<MV>( problem_->getInitResVec() ) ) ) );
669  newProblem->setProblem();
670 
671  solver->setProblem( newProblem );
672 
673  ret = solver->solve();
674  numIters_ = solver->getNumIters();
675  loaDetected_ = solver->isLOADetected();
676 
677  } // if (hasOuterSolver_ && maxDegree_)
678  else if (hasOuterSolver_) {
679 
680  // There is no polynomial, just create the outer solver with the outerSolverType_ and outerParams_.
682  RCP<SolverManager<ScalarType, MV, OP> > solver = factory.create( outerSolverType_, outerParams_ );
683  TEUCHOS_TEST_FOR_EXCEPTION( solver == Teuchos::null, std::invalid_argument,
684  "Belos::GmresPolySolMgr::solve(): Selected solver is not valid.");
685 
686  solver->setProblem( problem_ );
687 
688  ret = solver->solve();
689  numIters_ = solver->getNumIters();
690  loaDetected_ = solver->isLOADetected();
691 
692  }
693  else if (maxDegree_) {
694 
695  // Apply the polynomial to the current linear system
696  poly_Op_->ApplyPoly( *problem_->getRHS(), *problem_->getLHS() );
697 
698  }
699 
700  return ret;
701 }
702 
703 
704 template<class ScalarType, class MV, class OP>
706 {
707  std::ostringstream out;
708 
709  out << "\"Belos::GmresPolySolMgr\": {"
710  << "ScalarType: " << Teuchos::TypeNameTraits<ScalarType>::name ()
711  << ", Poly Degree: " << poly_dim_
712  << ", Poly Max Degree: " << maxDegree_
713  << ", Poly Tol: " << polyTol_;
714  out << "}";
715  return out.str ();
716 }
717 
718 } // namespace Belos
719 
720 #endif // BELOS_GMRES_POLY_SOLMGR_HPP
Teuchos::RCP< const Teuchos::ParameterList > getCurrentParameters() const override
Get a parameter list containing the current parameters for this object.
Collection of types and exceptions used within the Belos solvers.
Belos::GmresPolyOp< ScalarType, MV, OP > gmres_poly_t
Teuchos::RCP< gmres_poly_t > poly_Op_
Teuchos::RCP< const Teuchos::ParameterList > validPL_
Cached default (valid) parameters.
ReturnType solve() override
This method performs possibly repeated calls to the underlying linear solver&#39;s iterate() routine unti...
virtual ~GmresPolySolMgr()
Destructor.
Defines the GMRES polynomial operator hybrid-GMRES iterative linear solver.
Belos::GmresPolyMv< ScalarType, MV > gmres_poly_mv_t
bool isLOADetected() const override
Return whether a loss of accuracy was detected by this solver during the most current solve...
static constexpr const char * label_default_
T & get(ParameterList &l, const std::string &name)
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
Teuchos::Array< Teuchos::RCP< Teuchos::Time > > getTimers() const
Return the timers for this object.
static RCP< Time > getNewCounter(const std::string &name)
static constexpr int maxDegree_default_
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Teuchos::RCP< Teuchos::ParameterList > outerParams_
Teuchos::RCP< Teuchos::ParameterList > params_
static const double polyTol
Relative residual tolerance for matrix polynomial construction.
Definition: BelosTypes.hpp:304
Traits class which defines basic operations on multivectors.
bool isParameter(const std::string &name) const
static constexpr std::ostream * outputStream_default_
ResetType
How to reset the solver.
Definition: BelosTypes.hpp:206
Alternative run-time polymorphic interface for operators.
GmresPolySolMgrLinearProblemFailure is thrown when the linear problem is not setup (i...
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Pure virtual base class which describes the basic interface for a solver manager. ...
bool is_null(const RCP< T > &p)
Teuchos::ScalarTraits< ScalarType > STS
static constexpr bool addRoots_default_
bool isSublist(const std::string &name) const
const LinearProblem< ScalarType, MV, OP > & getProblem() const override
Get current linear problem being solved for in this object.
A linear system to solve, and its associated information.
Class which describes the linear problem to be solved by the iterative solver.
void setProblem(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem) override
Set the linear problem that needs to be solved.
MultiVecTraits< ScalarType, MV > MVT
The GMRES polynomial can be created in conjunction with any standard preconditioner.
static constexpr const char * orthoType_default_
int getNumIters() const override
Get the iteration count for the most recent call to solve().
ReturnType
Whether the Belos solve converged for all linear systems.
Definition: BelosTypes.hpp:155
The Belos::SolverManager is a templated virtual base class that defines the basic interface that any ...
Belos&#39;s class for applying the GMRES polynomial operator that is used by the hybrid-GMRES linear solv...
void reset(const ResetType type) override
Reset the solver.
Teuchos::RCP< Teuchos::Time > timerPoly_
void validateParameters(ParameterList const &validParamList, int const depth=1000, EValidateUsed const validateUsed=VALIDATE_USED_ENABLED, EValidateDefaults const validateDefaults=VALIDATE_DEFAULTS_ENABLED) const
std::string description() const override
Method to return description of the hybrid block GMRES solver manager.
static constexpr bool randomRHS_default_
GmresPolySolMgrPolynomialFailure is thrown when their is a problem generating the GMRES polynomial fo...
virtual Teuchos::RCP< solver_base_type > create(const std::string &solverName, const Teuchos::RCP< Teuchos::ParameterList > &solverParams)
Create, configure, and return the specified solver.
GmresPolySolMgrPolynomialFailure(const std::string &what_arg)
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > &params) override
Set the parameters the solver manager should use to solve the linear problem.
bool isType(const std::string &name) const
Teuchos::ScalarTraits< ScalarType >::magnitudeType MagnitudeType
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
Get a parameter list containing the valid parameters for this object.
static constexpr const char * outerSolverType_default_
Interface for multivectors used by Belos&#39; linear solvers.
Teuchos::RCP< std::ostream > outputStream_
Parent class to all Belos exceptions.
Definition: BelosTypes.hpp:60
Default parameters common to most Belos solvers.
Definition: BelosTypes.hpp:291
GmresPolySolMgr()
Empty constructor for GmresPolySolMgr. This constructor takes no arguments and sets the default value...
static constexpr bool dampPoly_default_
Belos header file which uses auto-configuration information to include necessary C++ headers...
static constexpr int verbosity_default_
Teuchos::RCP< SolverManager< ScalarType, MV, OP > > clone() const override
clone for Inverted Injection (DII)
typename::Belos::Impl::SolverFactorySelector< SC, MV, OP >::type SolverFactory
GmresPolySolMgrLinearProblemFailure(const std::string &what_arg)
static constexpr const char * polyType_default_
Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > problem_
static std::string name()