10#ifndef BELOS_BLOCK_GMRES_SOLMGR_HPP
11#define BELOS_BLOCK_GMRES_SOLMGR_HPP
34#ifdef BELOS_TEUCHOS_TIME_MONITOR
35#include "Teuchos_TimeMonitor.hpp"
95template<
class ScalarType,
class MV,
class OP>
101 typedef Teuchos::ScalarTraits<ScalarType> SCT;
102 typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType MagnitudeType;
103 typedef Teuchos::ScalarTraits<MagnitudeType> MT;
138 const Teuchos::RCP<Teuchos::ParameterList> &
pl );
144 Teuchos::RCP<SolverManager<ScalarType, MV, OP> >
clone ()
const override {
171 Teuchos::Array<Teuchos::RCP<Teuchos::Time> >
getTimers()
const {
172 return Teuchos::tuple(timerSolve_);
259 const Teuchos::EVerbosityLevel
verbLevel =
260 Teuchos::Describable::verbLevel_default)
const override;
270 bool checkStatusTest();
273 Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > problem_;
276 Teuchos::RCP<OutputManager<ScalarType> > printer_;
277 Teuchos::RCP<std::ostream> outputStream_;
280 Teuchos::RCP<StatusTest<ScalarType,MV,OP> > debugStatusTest_;
281 Teuchos::RCP<StatusTest<ScalarType,MV,OP> > sTest_;
282 Teuchos::RCP<StatusTestMaxIters<ScalarType,MV,OP> > maxIterTest_;
283 Teuchos::RCP<StatusTest<ScalarType,MV,OP> > convTest_;
284 Teuchos::RCP<StatusTestResNorm<ScalarType,MV,OP> > expConvTest_, impConvTest_;
285 Teuchos::RCP<StatusTestOutput<ScalarType,MV,OP> > outputTest_;
288 Teuchos::RCP<MatOrthoManager<ScalarType,MV,OP> > ortho_;
291 Teuchos::RCP<Teuchos::ParameterList> params_;
294 static constexpr int maxRestarts_default_ = 20;
295 static constexpr int maxIters_default_ = 1000;
296 static constexpr bool adaptiveBlockSize_default_ =
true;
297 static constexpr bool showMaxResNormOnly_default_ =
false;
298 static constexpr bool flexibleGmres_default_ =
false;
299 static constexpr bool keepHessenberg_default_ =
false;
300 static constexpr bool expResTest_default_ =
false;
301 static constexpr int blockSize_default_ = 1;
302 static constexpr int numBlocks_default_ = 300;
305 static constexpr int outputFreq_default_ = -1;
306 static constexpr const char * impResScale_default_ =
"Norm of Preconditioned Initial Residual";
307 static constexpr const char * expResScale_default_ =
"Norm of Initial Residual";
308 static constexpr const char * label_default_ =
"Belos";
309 static constexpr const char * orthoType_default_ =
"ICGS";
312 MagnitudeType convtol_, orthoKappa_, achievedTol_;
313 int maxRestarts_, maxIters_, numIters_;
314 int blockSize_, numBlocks_, verbosity_, outputStyle_, outputFreq_;
315 bool adaptiveBlockSize_, showMaxResNormOnly_, isFlexible_, keepHessenberg_, expResTest_;
316 std::string orthoType_;
317 std::string impResScale_, expResScale_;
321 Teuchos::RCP<Teuchos::Time> timerSolve_;
324 Teuchos::RCP<GmresIteration<ScalarType,MV,OP> > block_gmres_iter_;
327 bool isSet_, isSTSet_;
328 bool needsIterRebuild_;
334template<
class ScalarType,
class MV,
class OP>
340 maxRestarts_(maxRestarts_default_),
341 maxIters_(maxIters_default_),
343 blockSize_(blockSize_default_),
344 numBlocks_(numBlocks_default_),
345 verbosity_(verbosity_default_),
346 outputStyle_(outputStyle_default_),
347 outputFreq_(outputFreq_default_),
348 adaptiveBlockSize_(adaptiveBlockSize_default_),
349 showMaxResNormOnly_(showMaxResNormOnly_default_),
350 isFlexible_(flexibleGmres_default_),
351 keepHessenberg_(keepHessenberg_default_),
352 expResTest_(expResTest_default_),
353 orthoType_(orthoType_default_),
354 impResScale_(impResScale_default_),
355 expResScale_(expResScale_default_),
356 label_(label_default_),
359 needsIterRebuild_(
true),
365template<
class ScalarType,
class MV,
class OP>
368 const Teuchos::RCP<Teuchos::ParameterList> &
pl) :
374 maxRestarts_(maxRestarts_default_),
375 maxIters_(maxIters_default_),
377 blockSize_(blockSize_default_),
378 numBlocks_(numBlocks_default_),
379 verbosity_(verbosity_default_),
380 outputStyle_(outputStyle_default_),
381 outputFreq_(outputFreq_default_),
382 adaptiveBlockSize_(adaptiveBlockSize_default_),
383 showMaxResNormOnly_(showMaxResNormOnly_default_),
384 isFlexible_(flexibleGmres_default_),
385 keepHessenberg_(keepHessenberg_default_),
386 expResTest_(expResTest_default_),
387 orthoType_(orthoType_default_),
388 impResScale_(impResScale_default_),
389 expResScale_(expResScale_default_),
390 label_(label_default_),
393 needsIterRebuild_(
true),
407template<
class ScalarType,
class MV,
class OP>
408Teuchos::RCP<const Teuchos::ParameterList>
411 static Teuchos::RCP<const Teuchos::ParameterList>
validPL;
413 Teuchos::RCP<Teuchos::ParameterList>
pl = Teuchos::parameterList();
418 "The relative residual tolerance that needs to be achieved by the\n"
419 "iterative solver in order for the linear system to be declared converged." );
420 pl->set(
"Maximum Restarts",
static_cast<int>(maxRestarts_default_),
421 "The maximum number of restarts allowed for each\n"
422 "set of RHS solved.");
423 pl->set(
"Maximum Iterations",
static_cast<int>(maxIters_default_),
424 "The maximum number of block iterations allowed for each\n"
425 "set of RHS solved.");
426 pl->set(
"Num Blocks",
static_cast<int>(numBlocks_default_),
427 "The maximum number of blocks allowed in the Krylov subspace\n"
428 "for each set of RHS solved.");
429 pl->set(
"Block Size",
static_cast<int>(blockSize_default_),
430 "The number of vectors in each block. This number times the\n"
431 "number of blocks is the total Krylov subspace dimension.");
432 pl->set(
"Adaptive Block Size",
static_cast<bool>(adaptiveBlockSize_default_),
433 "Whether the solver manager should adapt the block size\n"
434 "based on the number of RHS to solve.");
435 pl->set(
"Verbosity",
static_cast<int>(verbosity_default_),
436 "What type(s) of solver information should be outputted\n"
437 "to the output stream.");
438 pl->set(
"Output Style",
static_cast<int>(outputStyle_default_),
439 "What style is used for the solver information outputted\n"
440 "to the output stream.");
441 pl->set(
"Output Frequency",
static_cast<int>(outputFreq_default_),
442 "How often convergence information should be outputted\n"
443 "to the output stream.");
444 pl->set(
"Output Stream", Teuchos::rcpFromRef(std::cout),
445 "A reference-counted pointer to the output stream where all\n"
446 "solver output is sent.");
447 pl->set(
"Show Maximum Residual Norm Only",
static_cast<bool>(showMaxResNormOnly_default_),
448 "When convergence information is printed, only show the maximum\n"
449 "relative residual norm when the block size is greater than one.");
450 pl->set(
"Flexible Gmres",
static_cast<bool>(flexibleGmres_default_),
451 "Whether the solver manager should use the flexible variant\n"
453 pl->set(
"Keep Hessenberg",
static_cast<bool>(keepHessenberg_default_),
454 "Whether the raw upper Hessenberg matrix should be stored separately\n"
455 "from the QR-factored least squares system. Useful for harmonic Ritz\n"
456 "pair computation from the GMRES iteration state.");
457 pl->set(
"Explicit Residual Test",
static_cast<bool>(expResTest_default_),
458 "Whether the explicitly computed residual should be used in the convergence test.");
459 pl->set(
"Implicit Residual Scaling",
static_cast<const char *
>(impResScale_default_),
460 "The type of scaling used in the implicit residual convergence test.");
461 pl->set(
"Explicit Residual Scaling",
static_cast<const char *
>(expResScale_default_),
462 "The type of scaling used in the explicit residual convergence test.");
463 pl->set(
"Timer Label",
static_cast<const char *
>(label_default_),
464 "The string to use as a prefix for the timer labels.");
465 pl->set(
"Orthogonalization",
static_cast<const char *
>(orthoType_default_),
466 "The type of orthogonalization to use: DGKS, ICGS, or IMGS.");
468 "The constant used by DGKS orthogonalization to determine\n"
469 "whether another step of classical Gram-Schmidt is necessary.");
476template<
class ScalarType,
class MV,
class OP>
481 if (params_ == Teuchos::null) {
482 params_ = Teuchos::rcp(
new Teuchos::ParameterList(*getValidParameters()) );
485 params->validateParameters(*getValidParameters());
489 if (
params->isParameter(
"Maximum Restarts")) {
490 maxRestarts_ =
params->get(
"Maximum Restarts",maxRestarts_default_);
493 params_->set(
"Maximum Restarts", maxRestarts_);
497 if (
params->isParameter(
"Maximum Iterations")) {
498 maxIters_ =
params->get(
"Maximum Iterations",maxIters_default_);
501 params_->set(
"Maximum Iterations", maxIters_);
502 if (maxIterTest_!=Teuchos::null)
503 maxIterTest_->setMaxIters( maxIters_ );
507 if (
params->isParameter(
"Block Size")) {
508 blockSize_ =
params->get(
"Block Size",blockSize_default_);
510 "Belos::BlockGmresSolMgr: \"Block Size\" must be strictly positive.");
513 params_->set(
"Block Size", blockSize_);
517 if (
params->isParameter(
"Adaptive Block Size")) {
518 adaptiveBlockSize_ =
params->get(
"Adaptive Block Size",adaptiveBlockSize_default_);
521 params_->set(
"Adaptive Block Size", adaptiveBlockSize_);
525 if (
params->isParameter(
"Num Blocks")) {
526 numBlocks_ =
params->get(
"Num Blocks",numBlocks_default_);
528 "Belos::BlockGmresSolMgr: \"Num Blocks\" must be strictly positive.");
531 params_->set(
"Num Blocks", numBlocks_);
535 if (
params->isParameter(
"Timer Label")) {
541 params_->set(
"Timer Label", label_);
542 std::string
solveLabel = label_ +
": BlockGmresSolMgr total solve time";
543#ifdef BELOS_TEUCHOS_TIME_MONITOR
544 timerSolve_ = Teuchos::TimeMonitor::getNewCounter(
solveLabel);
546 if (ortho_ != Teuchos::null) {
547 ortho_->setLabel( label_ );
553 if (
params->isParameter(
"Keep Hessenberg")) {
554 keepHessenberg_ = Teuchos::getParameter<bool>(*
params,
"Keep Hessenberg");
555 params_->set(
"Keep Hessenberg", keepHessenberg_);
559 if (
params->isParameter(
"Flexible Gmres")) {
560 isFlexible_ = Teuchos::getParameter<bool>(*
params,
"Flexible Gmres");
561 params_->set(
"Flexible Gmres", isFlexible_);
562 if (isFlexible_ && expResTest_) {
569 if (
params->isParameter(
"Verbosity")) {
570 if (Teuchos::isParameterType<int>(*
params,
"Verbosity")) {
571 verbosity_ =
params->get(
"Verbosity", verbosity_default_);
573 verbosity_ = (
int)Teuchos::getParameter<Belos::MsgType>(*
params,
"Verbosity");
577 params_->set(
"Verbosity", verbosity_);
578 if (printer_ != Teuchos::null)
579 printer_->setVerbosity(verbosity_);
583 if (
params->isParameter(
"Output Style")) {
584 if (Teuchos::isParameterType<int>(*
params,
"Output Style")) {
585 outputStyle_ =
params->get(
"Output Style", outputStyle_default_);
587 outputStyle_ = (
int)Teuchos::getParameter<Belos::OutputType>(*
params,
"Output Style");
591 params_->set(
"Output Style", outputStyle_);
592 if (outputTest_ != Teuchos::null) {
598 if (
params->isParameter(
"Output Stream")) {
599 outputStream_ = Teuchos::getParameter<Teuchos::RCP<std::ostream> >(*
params,
"Output Stream");
602 params_->set(
"Output Stream", outputStream_);
603 if (printer_ != Teuchos::null)
604 printer_->setOStream( outputStream_ );
609 if (
params->isParameter(
"Output Frequency")) {
610 outputFreq_ =
params->get(
"Output Frequency", outputFreq_default_);
614 params_->set(
"Output Frequency", outputFreq_);
615 if (outputTest_ != Teuchos::null)
616 outputTest_->setOutputFrequency( outputFreq_ );
620 if (printer_ == Teuchos::null) {
626 if (
params->isParameter(
"Orthogonalization")) {
633 params_->set(
"Orthogonalization", orthoType_);
636 if (
params->isParameter(
"Orthogonalization Constant")) {
637 if (
params->isType<MagnitudeType> (
"Orthogonalization Constant")) {
638 orthoKappa_ =
params->get (
"Orthogonalization Constant",
642 orthoKappa_ =
params->get (
"Orthogonalization Constant",
647 params_->set(
"Orthogonalization Constant",orthoKappa_);
648 if (orthoType_==
"DGKS") {
650 Teuchos::rcp_dynamic_cast<DGKSOrthoManager<ScalarType,MV,OP> >(ortho_)->setDepTol( orthoKappa_ );
659 if (orthoType_==
"DGKS" && orthoKappa_ > 0) {
660 paramsOrtho = Teuchos::rcp(
new Teuchos::ParameterList());
664 ortho_ =
factory.makeMatOrthoManager (orthoType_, Teuchos::null, printer_, label_,
paramsOrtho);
668 if (
params->isParameter(
"Convergence Tolerance")) {
669 if (
params->isType<MagnitudeType> (
"Convergence Tolerance")) {
670 convtol_ =
params->get (
"Convergence Tolerance",
678 params_->set(
"Convergence Tolerance", convtol_);
679 if (impConvTest_ != Teuchos::null)
680 impConvTest_->setTolerance( convtol_ );
681 if (expConvTest_ != Teuchos::null)
682 expConvTest_->setTolerance( convtol_ );
686 if (
params->isParameter(
"Implicit Residual Scaling")) {
687 std::string
tempImpResScale = Teuchos::getParameter<std::string>( *
params,
"Implicit Residual Scaling" );
695 params_->set(
"Implicit Residual Scaling", impResScale_);
696 if (impConvTest_ != Teuchos::null) {
700 catch (std::exception&
e) {
708 if (
params->isParameter(
"Explicit Residual Scaling")) {
709 std::string
tempExpResScale = Teuchos::getParameter<std::string>( *
params,
"Explicit Residual Scaling" );
717 params_->set(
"Explicit Residual Scaling", expResScale_);
718 if (expConvTest_ != Teuchos::null) {
722 catch (std::exception&
e) {
730 if (
params->isParameter(
"Explicit Residual Test")) {
731 expResTest_ = Teuchos::getParameter<bool>( *
params,
"Explicit Residual Test" );
734 params_->set(
"Explicit Residual Test", expResTest_);
735 if (expConvTest_ == Teuchos::null) {
740 if (
params->isParameter(
"Show Maximum Residual Norm Only")) {
741 showMaxResNormOnly_ = Teuchos::getParameter<bool>(*
params,
"Show Maximum Residual Norm Only");
744 params_->set(
"Show Maximum Residual Norm Only", showMaxResNormOnly_);
745 if (impConvTest_ != Teuchos::null)
746 impConvTest_->setShowMaxResNormOnly( showMaxResNormOnly_ );
747 if (expConvTest_ != Teuchos::null)
748 expConvTest_->setShowMaxResNormOnly( showMaxResNormOnly_ );
753 if (timerSolve_ == Teuchos::null) {
754 std::string
solveLabel = label_ +
": BlockGmresSolMgr total solve time";
755#ifdef BELOS_TEUCHOS_TIME_MONITOR
756 timerSolve_ = Teuchos::TimeMonitor::getNewCounter(
solveLabel);
762 needsIterRebuild_ =
true;
767template<
class ScalarType,
class MV,
class OP>
780 if (isFlexible_ && Teuchos::is_null(problem_->getRightPrec())) {
782 params_->set(
"Flexible Gmres", isFlexible_);
786 "Belos::BlockGmresSolMgr::solve(): Linear problem has a left preconditioner, not a right preconditioner, which is required for flexible GMRES.");
791 if (!Teuchos::is_null(problem_->getLeftPrec()) && !isFlexible_) {
813 convTest_ = Teuchos::rcp(
new StatusTestCombo_t( StatusTestCombo_t::SEQ, impConvTest_, expConvTest_ ) );
836 expConvTest_ = impConvTest_;
837 convTest_ = impConvTest_;
841 sTest_ = Teuchos::rcp(
new StatusTestCombo_t( StatusTestCombo_t::OR, maxIterTest_, convTest_ ) );
844 if (
nonnull(debugStatusTest_) ) {
846 Teuchos::rcp_dynamic_cast<StatusTestCombo_t>(sTest_)->addStatusTest( debugStatusTest_ );
866template<
class ScalarType,
class MV,
class OP>
876template<
class ScalarType,
class MV,
class OP>
883 setParameters(Teuchos::parameterList(*getValidParameters()));
887 "Belos::BlockGmresSolMgr::solve(): Linear problem is not a valid object.");
890 "Belos::BlockGmresSolMgr::solve(): Linear problem is not ready, setProblem() has not been called.");
892 if (!isSTSet_ || (!expResTest_ && !Teuchos::is_null(problem_->getLeftPrec())) ) {
894 "Belos::BlockGmresSolMgr::solve(): Linear problem and requested status tests are incompatible.");
899 int numRHS2Solve = MVT::GetNumberVecs( *(problem_->getRHS()) );
905 if ( adaptiveBlockSize_ ) {
921 problem_->setLSIndex(
currIdx );
925 Teuchos::ParameterList
plist;
926 plist.set(
"Block Size",blockSize_);
927 plist.set(
"Keep Hessenberg",keepHessenberg_);
929 ptrdiff_t dim = MVT::GetGlobalLength( *(problem_->getRHS()) );
930 if (blockSize_*
static_cast<ptrdiff_t>(numBlocks_) >
dim) {
937 "Belos::BlockGmresSolMgr::solve(): Warning! Requested Krylov subspace dimension is larger than operator dimension!"
938 << std::endl <<
" The maximum number of blocks allowed for the Krylov subspace will be adjusted to " <<
tmpNumBlocks << std::endl;
942 plist.set(
"Num Blocks",numBlocks_);
945 outputTest_->reset();
946 loaDetected_ =
false;
954 if (needsIterRebuild_) {
959 needsIterRebuild_ =
false;
961 Teuchos::RCP<GmresIteration<ScalarType,MV,OP> > &
block_gmres_iter = block_gmres_iter_;
965#ifdef BELOS_TEUCHOS_TIME_MONITOR
966 Teuchos::TimeMonitor
slvtimer(*timerSolve_);
972 if (blockSize_*numBlocks_ >
dim) {
987 outputTest_->resetNumCalls();
990 Teuchos::RCP<MV>
V_0;
993 if (
currIdx[blockSize_-1] == -1) {
994 V_0 = MVT::Clone( *(problem_->getInitResVec()), blockSize_ );
995 problem_->computeCurrResVec( &*
V_0 );
998 V_0 = MVT::CloneCopy( *(problem_->getInitResVec()),
currIdx );
1003 if (
currIdx[blockSize_-1] == -1) {
1004 V_0 = MVT::Clone( *(problem_->getInitPrecResVec()), blockSize_ );
1005 problem_->computeCurrPrecResVec( &*
V_0 );
1008 V_0 = MVT::CloneCopy( *(problem_->getInitPrecResVec()),
currIdx );
1013 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> >
z_0 =
1014 Teuchos::rcp(
new Teuchos::SerialDenseMatrix<int,ScalarType>( blockSize_, blockSize_ ) );
1019 "Belos::BlockGmresSolMgr::solve(): Failed to compute initial block of orthonormal vectors.");
1039 if ( convTest_->getStatus() ==
Passed ) {
1040 if ( expConvTest_->getLOADetected() ) {
1042 loaDetected_ =
true;
1044 "Belos::BlockGmresSolMgr::solve(): Warning! Solver has experienced a loss of accuracy!" << std::endl;
1054 else if ( maxIterTest_->getStatus() ==
Passed ) {
1072 printer_->stream(
Debug) <<
" Performing restart number " <<
numRestarts <<
" of " << maxRestarts_ << std::endl << std::endl;
1078 Teuchos::RCP<MV>
curX = problem_->getCurrLHSVec();
1082 problem_->updateSolution(
update,
true );
1091 problem_->computeCurrResVec( &*
V_0 );
1093 problem_->computeCurrPrecResVec( &*
V_0 );
1096 z_0 = Teuchos::rcp(
new Teuchos::SerialDenseMatrix<int,ScalarType>( blockSize_, blockSize_ ) );
1101 "Belos::BlockGmresSolMgr::solve(): Failed to compute initial block of orthonormal vectors after restart.");
1120 "Belos::BlockGmresSolMgr::solve(): Invalid return from BlockGmresIter::iterate().");
1125 if (blockSize_ != 1) {
1126 printer_->stream(
Errors) <<
"Error! Caught std::exception in BlockGmresIter::iterate() at iteration "
1128 <<
e.what() << std::endl;
1129 if (convTest_->getStatus() !=
Passed)
1139 if (convTest_->getStatus() !=
Passed)
1146 achievedTol_ = MT::one();
1147 Teuchos::RCP<MV>
X = problem_->getLHS();
1148 MVT::MvInit( *
X, SCT::zero() );
1149 printer_->stream(
Warnings) <<
"Belos::BlockGmresSolMgr::solve(): Warning! NaN has been detected!"
1153 catch (
const std::exception &
e) {
1154 printer_->stream(
Errors) <<
"Error! Caught std::exception in BlockGmresIter::iterate() at iteration "
1156 <<
e.what() << std::endl;
1166 Teuchos::RCP<MV>
curX = problem_->getCurrLHSVec();
1168 if (
update != Teuchos::null)
1173 if ( !Teuchos::is_null(expConvTest_->getSolution()) ) {
1174 Teuchos::RCP<MV>
newX = expConvTest_->getSolution();
1175 Teuchos::RCP<MV>
curX = problem_->getCurrLHSVec();
1180 problem_->updateSolution(
update,
true );
1185 problem_->setCurrLS();
1193 if ( adaptiveBlockSize_ ) {
1207 problem_->setLSIndex(
currIdx );
1221#ifdef BELOS_TEUCHOS_TIME_MONITOR
1226 Teuchos::TimeMonitor::summarize( printer_->stream(
TimingDetails) );
1230 numIters_ = maxIterTest_->getNumIters();
1256 "Belos::BlockGmresSolMgr::solve(): The implicit convergence test's "
1257 "getTestValue() method returned NULL. Please report this bug to the "
1258 "Belos developers.");
1260 "Belos::BlockGmresSolMgr::solve(): The implicit convergence test's "
1261 "getTestValue() method returned a vector of length zero. Please report "
1262 "this bug to the Belos developers.");
1277template<
class ScalarType,
class MV,
class OP>
1280 std::ostringstream
out;
1281 out <<
"\"Belos::BlockGmresSolMgr\": {";
1285 out <<
"Flexible: " << (isFlexible_ ?
"true" :
"false")
1286 <<
", Num Blocks: " << numBlocks_
1287 <<
", Maximum Iterations: " << maxIters_
1288 <<
", Maximum Restarts: " << maxRestarts_
1289 <<
", Convergence Tolerance: " << convtol_
1295template<
class ScalarType,
class MV,
class OP>
1299 const Teuchos::EVerbosityLevel
verbLevel)
const
1301 using Teuchos::TypeNameTraits;
1302 using Teuchos::VERB_DEFAULT;
1303 using Teuchos::VERB_NONE;
1304 using Teuchos::VERB_LOW;
1311 const Teuchos::EVerbosityLevel
vl =
1317 out <<
"\"Belos::BlockGmresSolMgr\":" << endl;
1319 out <<
"Template parameters:" << endl;
1329 out <<
"Flexible: " << (isFlexible_ ?
"true" :
"false") << endl
1330 <<
"Num Blocks: " << numBlocks_ << endl
1331 <<
"Maximum Iterations: " << maxIters_ << endl
1332 <<
"Maximum Restarts: " << maxRestarts_ << endl
1333 <<
"Convergence Tolerance: " << convtol_ << endl;
Belos concrete class for performing the block, flexible GMRES iteration.
Belos concrete class for performing the block GMRES iteration.
Belos header file which uses auto-configuration information to include necessary C++ headers.
Pure virtual base class which augments the basic interface for a Gmres linear solver iteration.
Class which describes the linear problem to be solved by the iterative solver.
Class which manages the output and verbosity of the Belos solvers.
Pure virtual base class which describes the basic interface for a solver manager.
Belos::StatusTest for logically combining several status tests.
Belos::StatusTestResNorm for specifying general residual norm stopping criteria.
Belos::StatusTest for specifying an implicit residual norm stopping criteria that checks for loss of ...
Belos::StatusTest class for specifying a maximum number of iterations.
Virtual base class for StatusTest that printing status tests.
A factory class for generating StatusTestOutput objects.
Collection of types and exceptions used within the Belos solvers.
Parent class to all Belos exceptions.
Interface to Block GMRES and Flexible GMRES.
void setProblem(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem) override
Set the linear problem that needs to be solved.
int getNumIters() const override
Get the iteration count for the most recent call to solve().
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const override
Print the object with the given verbosity level to a FancyOStream.
void setDebugStatusTest(const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &debugStatusTest) override
Set a debug status test that will be checked at the same time as the top-level status test.
ReturnType solve() override
This method performs possibly repeated calls to the underlying linear solver's iterate() routine unti...
Teuchos::Array< Teuchos::RCP< Teuchos::Time > > getTimers() const
Return the timers for this object.
std::string description() const override
Return a one-line description of this object.
MagnitudeType achievedTol() const override
Tolerance achieved by the last solve() invocation.
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > ¶ms) override
Set the parameters the solver manager should use to solve the linear problem.
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
Get a parameter list containing the valid parameters for this object.
Teuchos::RCP< SolverManager< ScalarType, MV, OP > > clone() const override
clone for Inverted Injection (DII)
BlockGmresSolMgr()
Empty constructor for BlockGmresSolMgr. This constructor takes no arguments and sets the default valu...
virtual ~BlockGmresSolMgr()
Destructor.
Teuchos::RCP< const Teuchos::ParameterList > getCurrentParameters() const override
Get a parameter list containing the current parameters for this object.
void reset(const ResetType type) override
Performs a reset of the solver manager specified by the ResetType. This informs the solver manager th...
bool isLOADetected() const override
Return whether a loss of accuracy was detected by this solver during the most current solve.
const LinearProblem< ScalarType, MV, OP > & getProblem() const override
Get current linear problem being solved for in this object.
BlockGmresSolMgrLinearProblemFailure is thrown when the linear problem is not setup (i....
BlockGmresSolMgrLinearProblemFailure(const std::string &what_arg)
BlockGmresSolMgrOrthoFailure is thrown when the orthogonalization manager is unable to generate ortho...
BlockGmresSolMgrOrthoFailure(const std::string &what_arg)
GmresIterationOrthoFailure is thrown when the GmresIteration object is unable to compute independent ...
Alternative run-time polymorphic interface for operators.
Operator()
Default constructor (does nothing).
The Belos::SolverManager is a templated virtual base class that defines the basic interface that any ...
ScaleType convertStringToScaleType(const std::string &scaleType)
Convert the given string to its ScaleType enum value.
ReturnType
Whether the Belos solve converged for all linear systems.
ScaleType
The type of scaling to use on the residual norm value.
ResetType
How to reset the solver.
Default parameters common to most Belos solvers.
static const double convTol
Default convergence tolerance.
static const double orthoKappa
DGKS orthogonalization constant.