10#ifndef BELOS_SOLVERMANAGER_HPP
11#define BELOS_SOLVERMANAGER_HPP
22#include "Teuchos_ParameterList.hpp"
23#include "Teuchos_RCP.hpp"
24#include "Teuchos_Describable.hpp"
34template <
class ScalarType,
class MV,
class OP>
38template<
class ScalarType,
class MV,
class OP>
55 virtual Teuchos::RCP<SolverManager<ScalarType, MV, OP> >
clone ()
const = 0;
81 virtual typename Teuchos::ScalarTraits<ScalarType>::magnitudeType
achievedTol()
const {
121 <<
" overridden for the class" << this->description() <<
" yet!");
130 <<
" overridden for the class" << this->description() <<
" yet!");
190 const bool isComplex = Teuchos::ScalarTraits<ScalarType>::isComplex>
194 template<
class ScalarType,
class MV,
class OP>
209 template<
class ScalarType,
class MV,
class OP>
221 "This solver is not implemented for complex ScalarType." );
225 "This solver is not implemented for complex ScalarType." );
229 "This solver is not implemented for complex ScalarType." );
233 "This solver is not implemented for complex ScalarType." );
237 "This solver is not implemented for complex ScalarType." );
241 "This solver is not implemented for complex ScalarType." );
245 "This solver is not implemented for complex ScalarType." );
249 "This solver is not implemented for complex ScalarType." );
253 "This solver is not implemented for complex ScalarType." );
261 template<
class ScalarType>
279#ifdef HAVE_TEUCHOS_LONG_DOUBLE
283 const static bool value =
true;
288#ifdef HAVE_TEUCHOS_COMPLEX
290 class LapackSupportsScalar<std::
complex<float> > {
292 const static bool value =
true;
296 class LapackSupportsScalar<std::complex<double> > {
298 const static bool value =
true;
306 template<
class ScalarType,
309 const bool lapackSupportsScalarType =
317 template<
class ScalarType,
class MV,
class OP>
331 template<
class ScalarType,
class MV,
class OP>
337 (
true, std::logic_error,
"This solver is not implemented for ScalarType"
338 " types for which Teuchos::LAPACK does not have a valid implementation. "
339 "ScalarType = " << Teuchos::TypeNameTraits<ScalarType>::name () <<
".");
345 (
true, std::logic_error,
"This solver is not implemented for ScalarType"
346 " types for which Teuchos::LAPACK does not have a valid implementation. "
347 "ScalarType = " << Teuchos::TypeNameTraits<ScalarType>::name () <<
".");
351 (
true, std::logic_error,
"This solver is not implemented for ScalarType"
352 " types for which Teuchos::LAPACK does not have a valid implementation. "
353 "ScalarType = " << Teuchos::TypeNameTraits<ScalarType>::name () <<
".");
357 (
true, std::logic_error,
"This solver is not implemented for ScalarType"
358 " types for which Teuchos::LAPACK does not have a valid implementation. "
359 "ScalarType = " << Teuchos::TypeNameTraits<ScalarType>::name () <<
".");
363 (
true, std::logic_error,
"This solver is not implemented for ScalarType"
364 " types for which Teuchos::LAPACK does not have a valid implementation. "
365 "ScalarType = " << Teuchos::TypeNameTraits<ScalarType>::name () <<
".");
369 (
true, std::logic_error,
"This solver is not implemented for ScalarType"
370 " types for which Teuchos::LAPACK does not have a valid implementation. "
371 "ScalarType = " << Teuchos::TypeNameTraits<ScalarType>::name () <<
".");
375 (
true, std::logic_error,
"This solver is not implemented for ScalarType"
376 " types for which Teuchos::LAPACK does not have a valid implementation. "
377 "ScalarType = " << Teuchos::TypeNameTraits<ScalarType>::name () <<
".");
381 (
true, std::logic_error,
"This solver is not implemented for ScalarType"
382 " types for which Teuchos::LAPACK does not have a valid implementation. "
383 "ScalarType = " << Teuchos::TypeNameTraits<ScalarType>::name () <<
".");
387 (
true, std::logic_error,
"This solver is not implemented for ScalarType"
388 " types for which Teuchos::LAPACK does not have a valid implementation. "
389 "ScalarType = " << Teuchos::TypeNameTraits<ScalarType>::name () <<
".");
393 (
true, std::logic_error,
"This solver is not implemented for ScalarType"
394 " types for which Teuchos::LAPACK does not have a valid implementation. "
395 "ScalarType = " << Teuchos::TypeNameTraits<ScalarType>::name () <<
".");
408 ! Teuchos::ScalarTraits<ScalarType>::isComplex>
418 template<
class ScalarType,
class MV,
class OP>
433 template<
class ScalarType,
class MV,
class OP>
445 (
true, std::logic_error,
"This solver is not implemented for complex "
446 "ScalarType types, or for ScalarType types for which Teuchos::LAPACK "
447 "does not have a valid implementation."
448 "ScalarType = " << Teuchos::TypeNameTraits<ScalarType>::name () <<
".");
452 (
true, std::logic_error,
"This solver is not implemented for complex "
453 "ScalarType types, or for ScalarType types for which Teuchos::LAPACK "
454 "does not have a valid implementation."
455 "ScalarType = " << Teuchos::TypeNameTraits<ScalarType>::name () <<
".");
459 (
true, std::logic_error,
"This solver is not implemented for complex "
460 "ScalarType types, or for ScalarType types for which Teuchos::LAPACK "
461 "does not have a valid implementation."
462 "ScalarType = " << Teuchos::TypeNameTraits<ScalarType>::name () <<
".");
466 (
true, std::logic_error,
"This solver is not implemented for complex "
467 "ScalarType types, or for ScalarType types for which Teuchos::LAPACK "
468 "does not have a valid implementation."
469 "ScalarType = " << Teuchos::TypeNameTraits<ScalarType>::name () <<
".");
473 (
true, std::logic_error,
"This solver is not implemented for complex "
474 "ScalarType types, or for ScalarType types for which Teuchos::LAPACK "
475 "does not have a valid implementation."
476 "ScalarType = " << Teuchos::TypeNameTraits<ScalarType>::name () <<
".");
481 (
true, std::logic_error,
"This solver is not implemented for complex "
482 "ScalarType types, or for ScalarType types for which Teuchos::LAPACK "
483 "does not have a valid implementation."
484 "ScalarType = " << Teuchos::TypeNameTraits<ScalarType>::name () <<
".");
488 (
true, std::logic_error,
"This solver is not implemented for complex "
489 "ScalarType types, or for ScalarType types for which Teuchos::LAPACK "
490 "does not have a valid implementation."
491 "ScalarType = " << Teuchos::TypeNameTraits<ScalarType>::name () <<
".");
495 (
true, std::logic_error,
"This solver is not implemented for complex "
496 "ScalarType types, or for ScalarType types for which Teuchos::LAPACK "
497 "does not have a valid implementation."
498 "ScalarType = " << Teuchos::TypeNameTraits<ScalarType>::name () <<
".");
502 (
true, std::logic_error,
"This solver is not implemented for complex "
503 "ScalarType types, or for ScalarType types for which Teuchos::LAPACK "
504 "does not have a valid implementation."
505 "ScalarType = " << Teuchos::TypeNameTraits<ScalarType>::name () <<
".");
Belos header file which uses auto-configuration information to include necessary C++ headers.
Class which describes the linear problem to be solved by the iterative solver.
Belos::StatusTest for logically combining several status tests.
Collection of types and exceptions used within the Belos solvers.
Type traits class that says whether Teuchos::LAPACK has a valid implementation for the given ScalarTy...
virtual ~RealSolverManager()
virtual ~RealSolverManager()
virtual Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Return the valid parameters for this solver manager.
virtual const LinearProblem< ScalarType, MV, OP > & getProblem() const
Return a reference to the linear problem being solved by this solver manager.
virtual Teuchos::RCP< const Teuchos::ParameterList > getCurrentParameters() const
Return the current parameters being used for this solver manager.
virtual void setParameters(const Teuchos::RCP< Teuchos::ParameterList > ¶ms)
Set the parameters to use when solving the linear problem.
virtual void reset(const ResetType type)
Reset the solver manager.
virtual ReturnType solve()
Iterate until the status test tells us to stop.
virtual int getNumIters() const
Get the iteration count for the most recent call to solve().
virtual void setProblem(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem)
Set the linear problem that needs to be solved.
virtual bool isLOADetected() const
Returns whether a loss of accuracy was detected in the solver.
Base class for Belos::SolverManager subclasses which normally can only compile for real ScalarType.
virtual const LinearProblem< ScalarType, MV, OP > & getProblem() const
Return a reference to the linear problem being solved by this solver manager.
virtual Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Return the valid parameters for this solver manager.
virtual void setProblem(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem)
Set the linear problem that needs to be solved.
virtual ~SolverManagerRequiresLapack()
virtual ReturnType solve()
Iterate until the status test tells us to stop.
SolverManagerRequiresLapack()
virtual void reset(const ResetType type)
Reset the solver manager.
virtual bool isLOADetected() const
Returns whether a loss of accuracy was detected in the solver.
virtual int getNumIters() const
Get the iteration count for the most recent call to solve().
virtual void setParameters(const Teuchos::RCP< Teuchos::ParameterList > ¶ms)
Set the parameters to use when solving the linear problem.
virtual Teuchos::RCP< const Teuchos::ParameterList > getCurrentParameters() const
Return the current parameters being used for this solver manager.
SolverManagerRequiresLapack()
virtual ~SolverManagerRequiresLapack()
Base class for Belos::SolverManager subclasses which normally can only compile with ScalarType types ...
virtual int getNumIters() const
Get the iteration count for the most recent call to solve().
virtual Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Return the valid parameters for this solver manager.
virtual ~SolverManagerRequiresRealLapack()
virtual void setParameters(const Teuchos::RCP< Teuchos::ParameterList > ¶ms)
Set the parameters to use when solving the linear problem.
virtual bool isLOADetected() const
Returns whether a loss of accuracy was detected in the solver.
SolverManagerRequiresRealLapack()
virtual const LinearProblem< ScalarType, MV, OP > & getProblem() const
Return a reference to the linear problem being solved by this solver manager.
virtual void setProblem(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &)
Set the linear problem that needs to be solved.
virtual ReturnType solve()
Iterate until the status test tells us to stop.
virtual void reset(const ResetType type)
Reset the solver manager.
virtual Teuchos::RCP< const Teuchos::ParameterList > getCurrentParameters() const
Return the current parameters being used for this solver manager.
virtual ~SolverManagerRequiresRealLapack()
SolverManagerRequiresRealLapack()
Base class for Belos::SolverManager subclasses which normally can only compile with real ScalarType t...
Alternative run-time polymorphic interface for operators.
The Belos::SolverManager is a templated virtual base class that defines the basic interface that any ...
virtual void setUserConvStatusTest(const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &, const typename StatusTestCombo< ScalarType, MV, OP >::ComboType &=StatusTestCombo< ScalarType, MV, OP >::SEQ)
Set user-defined convergence status test.
virtual ReturnType solve()=0
Iterate until the status test tells us to stop.
virtual Teuchos::RCP< const Teuchos::ParameterList > getCurrentParameters() const =0
Return the current parameters being used for this solver manager.
virtual void reset(const ResetType type)=0
Reset the solver manager.
SolverManager()
Empty constructor.
virtual int getNumIters() const =0
Get the iteration count for the most recent call to solve().
virtual void setParameters(const Teuchos::RCP< Teuchos::ParameterList > ¶ms)=0
Set the parameters to use when solving the linear problem.
virtual Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const =0
Return the valid parameters for this solver manager.
virtual void setProblem(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem)=0
Set the linear problem that needs to be solved.
virtual Teuchos::RCP< SolverManager< ScalarType, MV, OP > > clone() const =0
clone the solver manager.
virtual bool isLOADetected() const =0
Returns whether a loss of accuracy was detected in the solver.
virtual ~SolverManager()
Destructor.
virtual const LinearProblem< ScalarType, MV, OP > & getProblem() const =0
Return a reference to the linear problem being solved by this solver manager.
virtual Teuchos::ScalarTraits< ScalarType >::magnitudeType achievedTol() const
Tolerance achieved by the last solve() invocation.
virtual void setDebugStatusTest(const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &)
Set user-defined debug status test.
A class for extending the status testing capabilities of Belos via logical combinations.
ComboType
The test can be either the AND of all the component tests, or the OR of all the component tests,...
ReturnType
Whether the Belos solve converged for all linear systems.
ResetType
How to reset the solver.