42 #ifndef BELOS_TFQMR_SOLMGR_HPP
43 #define BELOS_TFQMR_SOLMGR_HPP
61 #ifdef BELOS_TEUCHOS_TIME_MONITOR
103 template<
class ScalarType,
class MV,
class OP>
302 template<
class ScalarType,
class MV,
class OP>
304 outputStream_(Teuchos::
rcp(outputStream_default_,false)),
307 achievedTol_(Teuchos::ScalarTraits<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>::zero()),
308 maxIters_(maxIters_default_),
310 verbosity_(verbosity_default_),
311 outputStyle_(outputStyle_default_),
312 outputFreq_(outputFreq_default_),
314 expResTest_(expResTest_default_),
315 impResScale_(impResScale_default_),
316 expResScale_(expResScale_default_),
317 label_(label_default_),
324 template<
class ScalarType,
class MV,
class OP>
329 outputStream_(Teuchos::
rcp(outputStream_default_,false)),
332 achievedTol_(Teuchos::ScalarTraits<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>::zero()),
333 maxIters_(maxIters_default_),
335 verbosity_(verbosity_default_),
336 outputStyle_(outputStyle_default_),
337 outputFreq_(outputFreq_default_),
339 expResTest_(expResTest_default_),
340 impResScale_(impResScale_default_),
341 expResScale_(expResScale_default_),
342 label_(label_default_),
354 template<
class ScalarType,
class MV,
class OP>
367 maxIters_ = params->
get(
"Maximum Iterations",maxIters_default_);
370 params_->set(
"Maximum Iterations", maxIters_);
372 maxIterTest_->setMaxIters( maxIters_ );
377 blockSize_ = params->
get(
"Block Size",1);
379 "Belos::TFQMRSolMgr: \"Block Size\" must be 1.");
382 params_->set(
"Block Size", blockSize_);
387 std::string tempLabel = params->
get(
"Timer Label", label_default_);
390 if (tempLabel != label_) {
392 params_->set(
"Timer Label", label_);
393 std::string solveLabel = label_ +
": TFQMRSolMgr total solve time";
394 #ifdef BELOS_TEUCHOS_TIME_MONITOR
402 if (Teuchos::isParameterType<int>(*params,
"Verbosity")) {
403 verbosity_ = params->
get(
"Verbosity", verbosity_default_);
405 verbosity_ = (int)Teuchos::getParameter<Belos::MsgType>(*params,
"Verbosity");
409 params_->set(
"Verbosity", verbosity_);
411 printer_->setVerbosity(verbosity_);
416 if (Teuchos::isParameterType<int>(*params,
"Output Style")) {
417 outputStyle_ = params->
get(
"Output Style", outputStyle_default_);
419 outputStyle_ = (int)Teuchos::getParameter<Belos::OutputType>(*params,
"Output Style");
423 params_->set(
"Output Style", outputStyle_);
429 outputStream_ = Teuchos::getParameter<Teuchos::RCP<std::ostream> >(*params,
"Output Stream");
432 params_->set(
"Output Stream", outputStream_);
434 printer_->setOStream( outputStream_ );
440 outputFreq_ = params->
get(
"Output Frequency", outputFreq_default_);
444 params_->set(
"Output Frequency", outputFreq_);
446 outputTest_->setOutputFrequency( outputFreq_ );
455 if (params->
isParameter(
"Convergence Tolerance")) {
457 convtol_ = params->
get (
"Convergence Tolerance",
465 params_->set(
"Convergence Tolerance", convtol_);
470 if (params->
isParameter(
"Implicit Tolerance Scale Factor")) {
472 impTolScale_ = params->
get (
"Implicit Tolerance Scale Factor",
477 impTolScale_ = params->
get (
"Implicit Tolerance Scale Factor",
482 params_->set(
"Implicit Tolerance Scale Factor", impTolScale_);
487 if (params->
isParameter(
"Implicit Residual Scaling")) {
488 std::string tempImpResScale = Teuchos::getParameter<std::string>( *params,
"Implicit Residual Scaling" );
491 if (impResScale_ != tempImpResScale) {
492 impResScale_ = tempImpResScale;
495 params_->set(
"Implicit Residual Scaling", impResScale_);
502 if (params->
isParameter(
"Explicit Residual Scaling")) {
503 std::string tempExpResScale = Teuchos::getParameter<std::string>( *params,
"Explicit Residual Scaling" );
506 if (expResScale_ != tempExpResScale) {
507 expResScale_ = tempExpResScale;
510 params_->set(
"Explicit Residual Scaling", expResScale_);
517 if (params->
isParameter(
"Explicit Residual Test")) {
518 expResTest_ = Teuchos::getParameter<bool>( *params,
"Explicit Residual Test" );
521 params_->set(
"Explicit Residual Test", expResTest_);
529 std::string solveLabel = label_ +
": TFQMRSolMgr total solve time";
530 #ifdef BELOS_TEUCHOS_TIME_MONITOR
541 template<
class ScalarType,
class MV,
class OP>
554 Teuchos::rcp(
new StatusTestGenResNorm_t( impTolScale_*convtol_ ) );
556 impConvTest_ = tmpImpConvTest;
561 tmpExpConvTest->defineResForm( StatusTestGenResNorm_t::Explicit,
Belos::TwoNorm );
563 expConvTest_ = tmpExpConvTest;
566 convTest_ =
Teuchos::rcp(
new StatusTestCombo_t( StatusTestCombo_t::SEQ, impConvTest_, expConvTest_ ) );
574 impConvTest_ = tmpImpConvTest;
577 expConvTest_ = impConvTest_;
578 convTest_ = impConvTest_;
580 sTest_ =
Teuchos::rcp(
new StatusTestCombo_t( StatusTestCombo_t::OR, maxIterTest_, convTest_ ) );
588 std::string solverDesc =
" TFQMR ";
589 outputTest_->setSolverDesc( solverDesc );
599 template<
class ScalarType,
class MV,
class OP>
612 "The relative residual tolerance that needs to be achieved by the\n"
613 "iterative solver in order for the linear system to be declared converged.");
615 "The scale factor used by the implicit residual test when explicit residual\n"
616 "testing is used. May enable faster convergence when TFQMR bound is too loose.");
617 pl->
set(
"Maximum Iterations", static_cast<int>(maxIters_default_),
618 "The maximum number of block iterations allowed for each\n"
619 "set of RHS solved.");
620 pl->
set(
"Verbosity", static_cast<int>(verbosity_default_),
621 "What type(s) of solver information should be outputted\n"
622 "to the output stream.");
623 pl->
set(
"Output Style", static_cast<int>(outputStyle_default_),
624 "What style is used for the solver information outputted\n"
625 "to the output stream.");
626 pl->
set(
"Output Frequency", static_cast<int>(outputFreq_default_),
627 "How often convergence information should be outputted\n"
628 "to the output stream.");
630 "A reference-counted pointer to the output stream where all\n"
631 "solver output is sent.");
632 pl->
set(
"Explicit Residual Test", static_cast<bool>(expResTest_default_),
633 "Whether the explicitly computed residual should be used in the convergence test.");
634 pl->
set(
"Implicit Residual Scaling", static_cast<const char *>(impResScale_default_),
635 "The type of scaling used in the implicit residual convergence test.");
636 pl->
set(
"Explicit Residual Scaling", static_cast<const char *>(expResScale_default_),
637 "The type of scaling used in the explicit residual convergence test.");
638 pl->
set(
"Timer Label", static_cast<const char *>(label_default_),
639 "The string to use as a prefix for the timer labels.");
647 template<
class ScalarType,
class MV,
class OP>
654 setParameters(Teuchos::parameterList(*getValidParameters()));
658 "Belos::TFQMRSolMgr::solve(): Linear problem is not a valid object.");
661 "Belos::TFQMRSolMgr::solve(): Linear problem is not ready, setProblem() has not been called.");
665 "Belos::TFQMRSolMgr::solve(): Linear problem and requested status tests are incompatible.");
670 int numRHS2Solve = MVT::GetNumberVecs( *(problem_->getRHS()) );
671 int numCurrRHS = blockSize_;
673 std::vector<int> currIdx, currIdx2;
676 currIdx.resize( blockSize_ );
677 currIdx2.resize( blockSize_ );
678 for (
int i=0; i<numCurrRHS; ++i)
679 { currIdx[i] = startPtr+i; currIdx2[i]=i; }
682 problem_->setLSIndex( currIdx );
687 plist.
set(
"Block Size",blockSize_);
690 outputTest_->reset();
693 bool isConverged =
true;
703 #ifdef BELOS_TEUCHOS_TIME_MONITOR
707 while ( numRHS2Solve > 0 ) {
710 std::vector<int> convRHSIdx;
711 std::vector<int> currRHSIdx( currIdx );
712 currRHSIdx.resize(numCurrRHS);
715 tfqmr_iter->resetNumIters();
718 outputTest_->resetNumCalls();
721 Teuchos::RCP<MV> R_0 = MVT::CloneViewNonConst( *(Teuchos::rcp_const_cast<MV>(problem_->getInitPrecResVec())), currIdx );
726 tfqmr_iter->initializeTFQMR(newstate);
732 tfqmr_iter->iterate();
739 if ( convTest_->getStatus() ==
Passed ) {
748 else if ( maxIterTest_->getStatus() ==
Passed ) {
763 "Belos::TFQMRSolMgr::solve(): Invalid return from TFQMRIter::iterate().");
766 catch (
const std::exception &e) {
767 printer_->stream(
Errors) <<
"Error! Caught std::exception in TFQMRIter::iterate() at iteration "
768 << tfqmr_iter->getNumIters() << std::endl
769 << e.what() << std::endl;
775 problem_->updateSolution( tfqmr_iter->getCurrentUpdate(), true );
778 problem_->setCurrLS();
781 startPtr += numCurrRHS;
782 numRHS2Solve -= numCurrRHS;
783 if ( numRHS2Solve > 0 ) {
784 numCurrRHS = blockSize_;
786 currIdx.resize( blockSize_ );
787 currIdx2.resize( blockSize_ );
788 for (
int i=0; i<numCurrRHS; ++i)
789 { currIdx[i] = startPtr+i; currIdx2[i] = i; }
791 problem_->setLSIndex( currIdx );
794 tfqmr_iter->setBlockSize( blockSize_ );
797 currIdx.resize( numRHS2Solve );
808 #ifdef BELOS_TEUCHOS_TIME_MONITOR
817 numIters_ = maxIterTest_->getNumIters();
830 const std::vector<MagnitudeType>* pTestValues = NULL;
832 pTestValues = expConvTest_->getTestValue();
833 if (pTestValues == NULL || pTestValues->size() < 1) {
834 pTestValues = impConvTest_->getTestValue();
839 pTestValues = impConvTest_->getTestValue();
842 "Belos::TFQMRSolMgr::solve(): The implicit convergence test's "
843 "getTestValue() method returned NULL. Please report this bug to the "
844 "Belos developers.");
846 "Belos::TMQMRSolMgr::solve(): The implicit convergence test's "
847 "getTestValue() method returned a vector of length zero. Please report "
848 "this bug to the Belos developers.");
853 achievedTol_ = *std::max_element (pTestValues->begin(), pTestValues->end());
863 template<
class ScalarType,
class MV,
class OP>
866 std::ostringstream oss;
ScaleType convertStringToScaleType(const std::string &scaleType)
Convert the given string to its ScaleType enum value.
Collection of types and exceptions used within the Belos solvers.
This class implements the preconditioned transpose-free QMR algorithm for solving non-Hermitian linea...
int getNumIters() const override
Get the iteration count for the most recent call to solve().
Belos's basic output manager for sending information of select verbosity levels to the appropriate ou...
Class which manages the output and verbosity of the Belos solvers.
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > ¶ms) override
Set the parameters the solver manager should use to solve the linear problem.
bool is_null(const boost::shared_ptr< T > &p)
static constexpr int outputStyle_default_
Teuchos::RCP< StatusTestOutput< ScalarType, MV, OP > > outputTest_
Teuchos::RCP< StatusTest< ScalarType, MV, OP > > convTest_
static constexpr bool expResTest_default_
Teuchos::Array< Teuchos::RCP< Teuchos::Time > > getTimers() const
Return the timers for this object.
bool isLOADetected() const override
Whether loss of accuracy was detected during the last solve() invocation.
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)
bool is_null(const std::shared_ptr< T > &p)
static constexpr int outputFreq_default_
MultiVecTraits< ScalarType, MV > MVT
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
A factory class for generating StatusTestOutput objects.
const LinearProblem< ScalarType, MV, OP > & getProblem() const override
Return a reference to the linear problem being solved by this solver manager.
Teuchos::RCP< Teuchos::Time > timerSolve_
Teuchos::RCP< SolverManager< ScalarType, MV, OP > > clone() const override
clone for Inverted Injection (DII)
static constexpr const char * label_default_
An implementation of StatusTestResNorm using a family of residual norms.
static constexpr int maxIters_default_
Belos concrete class for generating iterations with the preconditioned tranpose-free QMR (TFQMR) meth...
Teuchos::RCP< const MV > R
The current residual basis.
static const double convTol
Default convergence tolerance.
Belos::StatusTest class for specifying a maximum number of iterations.
Teuchos::RCP< StatusTestGenResNorm< ScalarType, MV, OP > > expConvTest_
Teuchos::RCP< OutputManager< ScalarType > > printer_
std::string description() const override
Method to return description of the TFQMR solver manager.
static std::string name()
Teuchos::ScalarTraits< ScalarType > SCT
A factory class for generating StatusTestOutput objects.
Teuchos::RCP< StatusTestGenResNorm< ScalarType, MV, OP > > impConvTest_
TFQMRSolMgrOrthoFailure(const std::string &what_arg)
ReturnType solve() override
This method performs possibly repeated calls to the underlying linear solver's iterate() routine unti...
Traits class which defines basic operations on multivectors.
Belos::StatusTest for logically combining several status tests.
static constexpr int verbosity_default_
bool isParameter(const std::string &name) const
MagnitudeType achievedTol() const override
Tolerance achieved by the last solve() invocation.
A Belos::StatusTest class for specifying a maximum number of iterations.
ResetType
How to reset the solver.
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. ...
static const double impTolScale
"Implicit Tolerance Scale Factor"
static void summarize(Ptr< const Comm< int > > comm, std::ostream &out=std::cout, const bool alwaysWriteLocal=false, const bool writeGlobalStats=true, const bool writeZeroTimers=true, const ECounterSetOp setOp=Intersection, const std::string &filter="", const bool ignoreZeroTimers=false)
Teuchos::RCP< StatusTestMaxIters< ScalarType, MV, OP > > maxIterTest_
Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > problem_
Teuchos::ScalarTraits< ScalarType >::magnitudeType MagnitudeType
A linear system to solve, and its associated information.
Structure to contain pointers to TFQMRIter state variables.
Class which describes the linear problem to be solved by the iterative solver.
TFQMRSolMgr()
Empty constructor for TFQMRSolMgr. This constructor takes no arguments and sets the default values fo...
static constexpr std::ostream * outputStream_default_
Teuchos::RCP< std::ostream > outputStream_
virtual ~TFQMRSolMgr()
Destructor.
ReturnType
Whether the Belos solve converged for all linear systems.
void reset(const ResetType type) override
Performs a reset of the solver manager specified by the ResetType. This informs the solver manager th...
The Belos::SolverManager is a templated virtual base class that defines the basic interface that any ...
The Belos::TFQMRSolMgr provides a powerful and fully-featured solver manager over the TFQMR linear so...
MagnitudeType impTolScale_
void validateParameters(ParameterList const &validParamList, int const depth=1000, EValidateUsed const validateUsed=VALIDATE_USED_ENABLED, EValidateDefaults const validateDefaults=VALIDATE_DEFAULTS_ENABLED) const
static constexpr const char * expResScale_default_
Teuchos::RCP< Teuchos::ParameterList > params_
Teuchos::RCP< StatusTestOutput< ScalarType, MV, OP > > create(const Teuchos::RCP< OutputManager< ScalarType > > &printer, Teuchos::RCP< StatusTest< ScalarType, MV, OP > > test, int mod, int printStates)
Create the StatusTestOutput object specified by the outputStyle.
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
Get a parameter list containing the valid parameters for this object.
Belos::StatusTestResNorm for specifying general residual norm stopping criteria.
bool isType(const std::string &name) const
A class for extending the status testing capabilities of Belos via logical combinations.
MagnitudeType achievedTol_
TFQMRSolMgrLinearProblemFailure(const std::string &what_arg)
static constexpr const char * impResScale_default_
Class which defines basic traits for the operator type.
TFQMRSolMgrOrthoFailure is thrown when the orthogonalization manager is unable to generate orthonorma...
Parent class to all Belos exceptions.
Default parameters common to most Belos solvers.
TFQMRSolMgrLinearProblemFailure is thrown when the linear problem is not setup (i.e.
Teuchos::RCP< StatusTest< ScalarType, MV, OP > > sTest_
Belos header file which uses auto-configuration information to include necessary C++ headers...
OperatorTraits< ScalarType, MV, OP > OPT
void setProblem(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem) override
Set the linear problem that needs to be solved.
Teuchos::RCP< const Teuchos::ParameterList > getCurrentParameters() const override
Get a parameter list containing the current parameters for this object.
Teuchos::ScalarTraits< MagnitudeType > MT