10 #ifndef BELOS_PSEUDO_BLOCK_GMRES_SOLMGR_HPP
11 #define BELOS_PSEUDO_BLOCK_GMRES_SOLMGR_HPP
28 #ifdef BELOS_TEUCHOS_TIME_MONITOR
90 template<
class ScalarType,
class MV,
class OP>
496 template<
class ScalarType,
class MV,
class OP>
498 outputStream_(Teuchos::rcpFromRef(std::cout)),
499 taggedTests_(Teuchos::null),
502 achievedTol_(Teuchos::ScalarTraits<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>::zero()),
503 maxRestarts_(maxRestarts_default_),
504 maxIters_(maxIters_default_),
506 blockSize_(blockSize_default_),
507 numBlocks_(numBlocks_default_),
508 verbosity_(verbosity_default_),
509 outputStyle_(outputStyle_default_),
510 outputFreq_(outputFreq_default_),
511 defQuorum_(defQuorum_default_),
512 showMaxResNormOnly_(showMaxResNormOnly_default_),
513 orthoType_(orthoType_default_),
514 impResScale_(impResScale_default_),
515 expResScale_(expResScale_default_),
517 label_(label_default_),
525 template<
class ScalarType,
class MV,
class OP>
530 outputStream_(Teuchos::rcpFromRef(std::cout)),
531 taggedTests_(Teuchos::null),
534 achievedTol_(Teuchos::ScalarTraits<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>::zero()),
535 maxRestarts_(maxRestarts_default_),
536 maxIters_(maxIters_default_),
538 blockSize_(blockSize_default_),
539 numBlocks_(numBlocks_default_),
540 verbosity_(verbosity_default_),
541 outputStyle_(outputStyle_default_),
542 outputFreq_(outputFreq_default_),
543 defQuorum_(defQuorum_default_),
544 showMaxResNormOnly_(showMaxResNormOnly_default_),
545 orthoType_(orthoType_default_),
546 impResScale_(impResScale_default_),
547 expResScale_(expResScale_default_),
549 label_(label_default_),
564 template<
class ScalarType,
class MV,
class OP>
570 using Teuchos::parameterList;
572 using Teuchos::rcp_dynamic_cast;
576 params_ = parameterList (*getValidParameters ());
588 maxRestarts_ = params->
get (
"Maximum Restarts", maxRestarts_default_);
591 params_->set (
"Maximum Restarts", maxRestarts_);
596 maxIters_ = params->
get (
"Maximum Iterations", maxIters_default_);
599 params_->set (
"Maximum Iterations", maxIters_);
600 if (! maxIterTest_.is_null ()) {
601 maxIterTest_->setMaxIters (maxIters_);
607 blockSize_ = params->
get (
"Block Size", blockSize_default_);
609 blockSize_ <= 0, std::invalid_argument,
610 "Belos::PseudoBlockGmresSolMgr::setParameters: "
611 "The \"Block Size\" parameter must be strictly positive, "
612 "but you specified a value of " << blockSize_ <<
".");
615 params_->set (
"Block Size", blockSize_);
620 numBlocks_ = params->
get (
"Num Blocks", numBlocks_default_);
622 numBlocks_ <= 0, std::invalid_argument,
623 "Belos::PseudoBlockGmresSolMgr::setParameters: "
624 "The \"Num Blocks\" parameter must be strictly positive, "
625 "but you specified a value of " << numBlocks_ <<
".");
628 params_->set (
"Num Blocks", numBlocks_);
633 const std::string tempLabel = params->
get (
"Timer Label", label_default_);
636 if (tempLabel != label_) {
638 params_->set (
"Timer Label", label_);
639 const std::string solveLabel =
640 label_ +
": PseudoBlockGmresSolMgr total solve time";
641 #ifdef BELOS_TEUCHOS_TIME_MONITOR
643 #endif // BELOS_TEUCHOS_TIME_MONITOR
645 ortho_->setLabel( label_ );
653 if (Teuchos::isParameterType<int> (*params,
"Verbosity")) {
654 verbosity_ = params->
get (
"Verbosity", verbosity_default_);
656 verbosity_ = (int) Teuchos::getParameter<Belos::MsgType> (*params,
"Verbosity");
660 params_->set (
"Verbosity", verbosity_);
661 if (! printer_.is_null ()) {
662 printer_->setVerbosity (verbosity_);
668 if (Teuchos::isParameterType<int> (*params,
"Output Style")) {
669 outputStyle_ = params->
get (
"Output Style", outputStyle_default_);
671 outputStyle_ = (int) Teuchos::getParameter<Belos::OutputType> (*params,
"Output Style");
675 params_->set (
"Output Style", outputStyle_);
676 if (! outputTest_.is_null ()) {
683 if (params->
isSublist (
"User Status Tests")) {
685 if ( userStatusTestsList.
numParams() > 0 ) {
686 std::string userCombo_string = params->
get<std::string>(
"User Status Tests Combo Type",
"SEQ");
688 setUserConvStatusTest( testFactory->buildStatusTests(userStatusTestsList), testFactory->stringToComboType(userCombo_string) );
689 taggedTests_ = testFactory->getTaggedTests();
696 outputStream_ = Teuchos::getParameter<Teuchos::RCP<std::ostream> > (*params,
"Output Stream");
699 params_->set(
"Output Stream", outputStream_);
700 if (! printer_.is_null ()) {
701 printer_->setOStream (outputStream_);
708 outputFreq_ = params->
get (
"Output Frequency", outputFreq_default_);
712 params_->set (
"Output Frequency", outputFreq_);
713 if (! outputTest_.is_null ()) {
714 outputTest_->setOutputFrequency (outputFreq_);
719 if (printer_.is_null ()) {
724 bool changedOrthoType =
false;
726 std::string tempOrthoType = params->
get (
"Orthogonalization", orthoType_default_);
727 if (tempOrthoType != orthoType_) {
728 orthoType_ = tempOrthoType;
729 changedOrthoType =
true;
732 params_->set(
"Orthogonalization", orthoType_);
735 if (params->
isParameter (
"Orthogonalization Constant")) {
737 orthoKappa_ = params->
get (
"Orthogonalization Constant",
741 orthoKappa_ = params->
get (
"Orthogonalization Constant",
746 params_->set (
"Orthogonalization Constant", orthoKappa_);
747 if (orthoType_ ==
"DGKS") {
748 if (orthoKappa_ > 0 && ! ortho_.is_null() && !changedOrthoType) {
750 rcp_dynamic_cast<ortho_type> (ortho_)->setDepTol (orthoKappa_);
756 if (ortho_.is_null() || changedOrthoType) {
759 if (orthoType_==
"DGKS" && orthoKappa_ > 0) {
761 paramsOrtho->
set (
"depTol", orthoKappa_ );
770 if (params->
isParameter (
"Convergence Tolerance")) {
772 convtol_ = params->
get (
"Convergence Tolerance",
780 params_->set (
"Convergence Tolerance", convtol_);
781 if (! impConvTest_.is_null ()) {
782 impConvTest_->setTolerance (convtol_);
784 if (! expConvTest_.is_null ()) {
785 expConvTest_->setTolerance (convtol_);
790 bool userDefinedResidualScalingUpdated =
false;
791 if (params->
isParameter (
"User Defined Residual Scaling")) {
794 tempResScaleFactor = params->
get (
"User Defined Residual Scaling",
798 tempResScaleFactor = params->
get (
"User Defined Residual Scaling",
803 if (resScaleFactor_ != tempResScaleFactor) {
804 resScaleFactor_ = tempResScaleFactor;
805 userDefinedResidualScalingUpdated =
true;
808 if(userDefinedResidualScalingUpdated)
810 if (! params->
isParameter (
"Implicit Residual Scaling") && ! impConvTest_.is_null ()) {
812 if(impResScale_ ==
"User Provided")
815 catch (std::exception& e) {
820 if (! params->
isParameter (
"Explicit Residual Scaling") && ! expConvTest_.is_null ()) {
822 if(expResScale_ ==
"User Provided")
825 catch (std::exception& e) {
834 if (params->
isParameter (
"Implicit Residual Scaling")) {
835 const std::string tempImpResScale =
836 Teuchos::getParameter<std::string> (*params,
"Implicit Residual Scaling");
839 if (impResScale_ != tempImpResScale) {
841 impResScale_ = tempImpResScale;
844 params_->set (
"Implicit Residual Scaling", impResScale_);
845 if (! impConvTest_.is_null ()) {
847 if(impResScale_ ==
"User Provided")
848 impConvTest_->defineScaleForm (impResScaleType,
Belos::TwoNorm, resScaleFactor_);
852 catch (std::exception& e) {
858 else if (userDefinedResidualScalingUpdated) {
861 if (! impConvTest_.is_null ()) {
863 if(impResScale_ ==
"User Provided")
864 impConvTest_->defineScaleForm (impResScaleType,
Belos::TwoNorm, resScaleFactor_);
866 catch (std::exception& e) {
874 if (params->
isParameter (
"Explicit Residual Scaling")) {
875 const std::string tempExpResScale =
876 Teuchos::getParameter<std::string> (*params,
"Explicit Residual Scaling");
879 if (expResScale_ != tempExpResScale) {
881 expResScale_ = tempExpResScale;
884 params_->set (
"Explicit Residual Scaling", expResScale_);
885 if (! expConvTest_.is_null ()) {
887 if(expResScale_ ==
"User Provided")
888 expConvTest_->defineScaleForm (expResScaleType,
Belos::TwoNorm, resScaleFactor_);
892 catch (std::exception& e) {
898 else if (userDefinedResidualScalingUpdated) {
901 if (! expConvTest_.is_null ()) {
903 if(expResScale_ ==
"User Provided")
904 expConvTest_->defineScaleForm (expResScaleType,
Belos::TwoNorm, resScaleFactor_);
906 catch (std::exception& e) {
915 if (params->
isParameter (
"Show Maximum Residual Norm Only")) {
916 showMaxResNormOnly_ =
917 Teuchos::getParameter<bool> (*params,
"Show Maximum Residual Norm Only");
920 params_->set (
"Show Maximum Residual Norm Only", showMaxResNormOnly_);
921 if (! impConvTest_.is_null ()) {
922 impConvTest_->setShowMaxResNormOnly (showMaxResNormOnly_);
924 if (! expConvTest_.is_null ()) {
925 expConvTest_->setShowMaxResNormOnly (showMaxResNormOnly_);
933 defQuorum_ = params->
get(
"Deflation Quorum", defQuorum_);
935 defQuorum_ > blockSize_, std::invalid_argument,
936 "Belos::PseudoBlockGmresSolMgr::setParameters: "
937 "The \"Deflation Quorum\" parameter (= " << defQuorum_ <<
") must not be "
938 "larger than \"Block Size\" (= " << blockSize_ <<
").");
939 params_->set (
"Deflation Quorum", defQuorum_);
940 if (! impConvTest_.is_null ()) {
941 impConvTest_->setQuorum (defQuorum_);
943 if (! expConvTest_.is_null ()) {
944 expConvTest_->setQuorum (defQuorum_);
950 std::string solveLabel = label_ +
": PseudoBlockGmresSolMgr total solve time";
951 #ifdef BELOS_TEUCHOS_TIME_MONITOR
961 template<
class ScalarType,
class MV,
class OP>
968 userConvStatusTest_ = userConvStatusTest;
969 comboType_ = comboType;
972 template<
class ScalarType,
class MV,
class OP>
978 debugStatusTest_ = debugStatusTest;
983 template<
class ScalarType,
class MV,
class OP>
996 "The relative residual tolerance that needs to be achieved by the\n"
997 "iterative solver in order for the linear system to be declared converged.");
998 pl->
set(
"Maximum Restarts", static_cast<int>(maxRestarts_default_),
999 "The maximum number of restarts allowed for each\n"
1000 "set of RHS solved.");
1001 pl->
set(
"Maximum Iterations", static_cast<int>(maxIters_default_),
1002 "The maximum number of block iterations allowed for each\n"
1003 "set of RHS solved.");
1004 pl->
set(
"Num Blocks", static_cast<int>(numBlocks_default_),
1005 "The maximum number of vectors allowed in the Krylov subspace\n"
1006 "for each set of RHS solved.");
1007 pl->
set(
"Block Size", static_cast<int>(blockSize_default_),
1008 "The number of RHS solved simultaneously.");
1009 pl->
set(
"Verbosity", static_cast<int>(verbosity_default_),
1010 "What type(s) of solver information should be outputted\n"
1011 "to the output stream.");
1012 pl->
set(
"Output Style", static_cast<int>(outputStyle_default_),
1013 "What style is used for the solver information outputted\n"
1014 "to the output stream.");
1015 pl->
set(
"Output Frequency", static_cast<int>(outputFreq_default_),
1016 "How often convergence information should be outputted\n"
1017 "to the output stream.");
1018 pl->
set(
"Deflation Quorum", static_cast<int>(defQuorum_default_),
1019 "The number of linear systems that need to converge before\n"
1020 "they are deflated. This number should be <= block size.");
1021 pl->
set(
"Output Stream", Teuchos::rcpFromRef(std::cout),
1022 "A reference-counted pointer to the output stream where all\n"
1023 "solver output is sent.");
1024 pl->
set(
"Show Maximum Residual Norm Only", static_cast<bool>(showMaxResNormOnly_default_),
1025 "When convergence information is printed, only show the maximum\n"
1026 "relative residual norm when the block size is greater than one.");
1027 pl->
set(
"Implicit Residual Scaling", static_cast<const char *>(impResScale_default_),
1028 "The type of scaling used in the implicit residual convergence test.");
1029 pl->
set(
"Explicit Residual Scaling", static_cast<const char *>(expResScale_default_),
1030 "The type of scaling used in the explicit residual convergence test.");
1031 pl->
set(
"Timer Label", static_cast<const char *>(label_default_),
1032 "The string to use as a prefix for the timer labels.");
1033 pl->
set(
"Orthogonalization", static_cast<const char *>(orthoType_default_),
1034 "The type of orthogonalization to use.");
1036 "The constant used by DGKS orthogonalization to determine\n"
1037 "whether another step of classical Gram-Schmidt is necessary.");
1038 pl->
sublist(
"User Status Tests");
1039 pl->
set(
"User Status Tests Combo Type",
"SEQ",
1040 "Type of logical combination operation of user-defined\n"
1041 "and/or solver-specific status tests.");
1048 template<
class ScalarType,
class MV,
class OP>
1068 Teuchos::rcp(
new StatusTestGenResNorm_t( convtol_, defQuorum_ ) );
1069 if(impResScale_ ==
"User Provided")
1074 tmpImpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
1075 impConvTest_ = tmpImpConvTest;
1079 Teuchos::rcp(
new StatusTestGenResNorm_t( convtol_, defQuorum_ ) );
1080 tmpExpConvTest->defineResForm( StatusTestGenResNorm_t::Explicit,
Belos::TwoNorm );
1081 if(expResScale_ ==
"User Provided")
1085 tmpExpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
1086 expConvTest_ = tmpExpConvTest;
1089 convTest_ =
Teuchos::rcp(
new StatusTestCombo_t( StatusTestCombo_t::SEQ, impConvTest_, expConvTest_ ) );
1096 Teuchos::rcp(
new StatusTestImpResNorm_t( convtol_, defQuorum_ ) );
1097 if(impResScale_ ==
"User Provided")
1101 tmpImpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
1102 impConvTest_ = tmpImpConvTest;
1105 expConvTest_ = impConvTest_;
1106 convTest_ = impConvTest_;
1109 if (
nonnull(userConvStatusTest_) ) {
1113 std::vector<Teuchos::RCP<StatusTest<ScalarType,MV,OP> > > tmpVec = tmpComboTest->getStatusTests();
1114 comboType_ = tmpComboTest->getComboType();
1115 const int numResTests =
static_cast<int>(tmpVec.size());
1116 auto newConvTest =
Teuchos::rcp(
new StatusTestCombo_t(comboType_));
1117 newConvTest->addStatusTest(convTest_);
1118 for (
int j = 0; j < numResTests; ++j) {
1119 newConvTest->addStatusTest(tmpVec[j]);
1121 convTest_ = newConvTest;
1126 new StatusTestCombo_t( comboType_, convTest_, userConvStatusTest_ ) );
1134 sTest_ =
Teuchos::rcp(
new StatusTestCombo_t( StatusTestCombo_t::OR, maxIterTest_, convTest_ ) );
1137 if (
nonnull(debugStatusTest_) ) {
1139 Teuchos::rcp_dynamic_cast<StatusTestCombo_t>(sTest_)->addStatusTest( debugStatusTest_ );
1148 std::string solverDesc =
" Pseudo Block Gmres ";
1149 outputTest_->setSolverDesc( solverDesc );
1160 template<
class ScalarType,
class MV,
class OP>
1166 if (!isSet_) { setParameters( params_ ); }
1169 "Belos::PseudoBlockGmresSolMgr::solve(): Linear problem is not ready, setProblem() has not been called.");
1172 if (!isSTSet_ || (!expResTest_ && !
Teuchos::is_null(problem_->getLeftPrec())) ) {
1174 "Belos::BlockGmresSolMgr::solve(): Linear problem and requested status tests are incompatible.");
1179 int numRHS2Solve = MVT::GetNumberVecs( *(problem_->getRHS()) );
1180 int numCurrRHS = ( numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;
1182 std::vector<int> currIdx( numCurrRHS );
1183 blockSize_ = numCurrRHS;
1184 for (
int i=0; i<numCurrRHS; ++i)
1185 { currIdx[i] = startPtr+i; }
1188 problem_->setLSIndex( currIdx );
1193 plist.
set(
"Num Blocks",numBlocks_);
1196 outputTest_->reset();
1197 loaDetected_ =
false;
1200 bool isConverged =
true;
1210 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1214 while ( numRHS2Solve > 0 ) {
1217 std::vector<int> convRHSIdx;
1218 std::vector<int> currRHSIdx( currIdx );
1219 currRHSIdx.resize(numCurrRHS);
1222 block_gmres_iter->setNumBlocks( numBlocks_ );
1225 block_gmres_iter->resetNumIters();
1228 outputTest_->resetNumCalls();
1234 std::vector<int> index(1);
1235 Teuchos::RCP<MV> tmpV, R_0 = MVT::CloneCopy( *(problem_->getInitPrecResVec()), currIdx );
1236 newState.
V.resize( blockSize_ );
1237 newState.
Z.resize( blockSize_ );
1238 for (
int i=0; i<blockSize_; ++i) {
1240 tmpV = MVT::CloneViewNonConst( *R_0, index );
1247 int rank = ortho_->normalize( *tmpV, tmpZ );
1249 "Belos::PseudoBlockGmresSolMgr::solve(): Failed to compute initial block of orthonormal vectors.");
1251 newState.
V[i] = tmpV;
1252 newState.
Z[i] = tmpZ;
1256 block_gmres_iter->initialize(newState);
1257 int numRestarts = 0;
1263 block_gmres_iter->iterate();
1270 if ( convTest_->getStatus() ==
Passed ) {
1272 if ( expConvTest_->getLOADetected() ) {
1282 loaDetected_ =
true;
1284 "Belos::PseudoBlockGmresSolMgr::solve(): Warning! Solver has experienced a loss of accuracy!" << std::endl;
1285 isConverged =
false;
1289 std::vector<int> convIdx = expConvTest_->convIndices();
1293 if (convIdx.size() == currRHSIdx.size())
1300 problem_->setCurrLS();
1304 int curDim = oldState.
curDim;
1309 std::vector<int> oldRHSIdx( currRHSIdx );
1310 std::vector<int> defRHSIdx;
1311 for (
unsigned int i=0; i<currRHSIdx.size(); ++i) {
1313 for (
unsigned int j=0; j<convIdx.size(); ++j) {
1314 if (currRHSIdx[i] == convIdx[j]) {
1320 defRHSIdx.push_back( i );
1323 defState.
V.push_back( Teuchos::rcp_const_cast<MV>( oldState.
V[i] ) );
1328 currRHSIdx[have] = currRHSIdx[i];
1332 defRHSIdx.resize(currRHSIdx.size()-have);
1333 currRHSIdx.resize(have);
1341 problem_->setLSIndex( convIdx );
1344 problem_->updateSolution( defUpdate,
true );
1348 problem_->setLSIndex( currRHSIdx );
1351 defState.
curDim = curDim;
1354 block_gmres_iter->initialize(defState);
1361 else if ( maxIterTest_->getStatus() ==
Passed ) {
1363 isConverged =
false;
1371 else if ( block_gmres_iter->getCurSubspaceDim() == block_gmres_iter->getMaxSubspaceDim() ) {
1373 if ( numRestarts >= maxRestarts_ ) {
1374 isConverged =
false;
1379 printer_->stream(
Debug) <<
" Performing restart number " << numRestarts <<
" of " << maxRestarts_ << std::endl << std::endl;
1383 problem_->updateSolution( update,
true );
1390 newstate.
V.resize(currRHSIdx.size());
1391 newstate.
Z.resize(currRHSIdx.size());
1395 R_0 = MVT::Clone( *(problem_->getInitPrecResVec()), currRHSIdx.size() );
1396 problem_->computeCurrPrecResVec( &*R_0 );
1397 for (
unsigned int i=0; i<currRHSIdx.size(); ++i) {
1400 tmpV = MVT::CloneViewNonConst( *R_0, index );
1407 int rank = ortho_->normalize( *tmpV, tmpZ );
1409 "Belos::PseudoBlockGmresSolMgr::solve(): Failed to compute initial block of orthonormal vectors after the restart.");
1411 newstate.
V[i] = tmpV;
1412 newstate.
Z[i] = tmpZ;
1417 block_gmres_iter->initialize(newstate);
1430 "Belos::PseudoBlockGmresSolMgr::solve(): Invalid return from PseudoBlockGmresIter::iterate().");
1436 block_gmres_iter->updateLSQR( block_gmres_iter->getCurSubspaceDim() );
1439 sTest_->checkStatus( &*block_gmres_iter );
1440 if (convTest_->getStatus() !=
Passed)
1441 isConverged =
false;
1446 achievedTol_ = MT::one();
1448 MVT::MvInit( *X, SCT::zero() );
1449 printer_->stream(
Warnings) <<
"Belos::PseudoBlockGmresSolMgr::solve(): Warning! NaN has been detected!"
1453 catch (
const std::exception &e) {
1454 printer_->stream(
Errors) <<
"Error! Caught std::exception in PseudoBlockGmresIter::iterate() at iteration "
1455 << block_gmres_iter->getNumIters() << std::endl
1456 << e.what() << std::endl;
1463 if (
nonnull(userConvStatusTest_)) {
1466 problem_->updateSolution( update,
true );
1468 else if (
nonnull(expConvTest_->getSolution())) {
1472 MVT::MvAddMv( 0.0, *newX, 1.0, *newX, *curX );
1477 problem_->updateSolution( update,
true );
1481 problem_->setCurrLS();
1484 startPtr += numCurrRHS;
1485 numRHS2Solve -= numCurrRHS;
1486 if ( numRHS2Solve > 0 ) {
1487 numCurrRHS = ( numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;
1489 blockSize_ = numCurrRHS;
1490 currIdx.resize( numCurrRHS );
1491 for (
int i=0; i<numCurrRHS; ++i)
1492 { currIdx[i] = startPtr+i; }
1495 if (defQuorum_ > blockSize_) {
1497 impConvTest_->setQuorum( blockSize_ );
1499 expConvTest_->setQuorum( blockSize_ );
1503 problem_->setLSIndex( currIdx );
1506 currIdx.resize( numRHS2Solve );
1517 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1526 numIters_ = maxIterTest_->getNumIters();
1539 const std::vector<MagnitudeType>* pTestValues = NULL;
1541 pTestValues = expConvTest_->getTestValue();
1542 if (pTestValues == NULL || pTestValues->size() < 1) {
1543 pTestValues = impConvTest_->getTestValue();
1548 pTestValues = impConvTest_->getTestValue();
1551 "Belos::PseudoBlockGmresSolMgr::solve(): The implicit convergence test's "
1552 "getTestValue() method returned NULL. Please report this bug to the "
1553 "Belos developers.");
1555 "Belos::PseudoBlockGmresSolMgr::solve(): The implicit convergence test's "
1556 "getTestValue() method returned a vector of length zero. Please report "
1557 "this bug to the Belos developers.");
1562 achievedTol_ = *std::max_element (pTestValues->begin(), pTestValues->end());
1565 if (!isConverged || loaDetected_) {
1572 template<
class ScalarType,
class MV,
class OP>
1575 std::ostringstream out;
1577 out <<
"\"Belos::PseudoBlockGmresSolMgr\": {";
1578 if (this->getObjectLabel () !=
"") {
1579 out <<
"Label: " << this->getObjectLabel () <<
", ";
1581 out <<
"Num Blocks: " << numBlocks_
1582 <<
", Maximum Iterations: " << maxIters_
1583 <<
", Maximum Restarts: " << maxRestarts_
1584 <<
", Convergence Tolerance: " << convtol_
1590 template<
class ScalarType,
class MV,
class OP>
1612 out <<
"\"Belos::PseudoBlockGmresSolMgr\":" << endl;
1614 out <<
"Template parameters:" << endl;
1617 out <<
"ScalarType: " << TypeNameTraits<ScalarType>::name () << endl
1618 <<
"MV: " << TypeNameTraits<MV>::name () << endl
1619 <<
"OP: " << TypeNameTraits<OP>::name () << endl;
1621 if (this->getObjectLabel () !=
"") {
1622 out <<
"Label: " << this->getObjectLabel () << endl;
1624 out <<
"Num Blocks: " << numBlocks_ << endl
1625 <<
"Maximum Iterations: " << maxIters_ << endl
1626 <<
"Maximum Restarts: " << maxRestarts_ << endl
1627 <<
"Convergence Tolerance: " << convtol_ << endl;
static const double orthoKappa
DGKS orthogonalization constant.
ScaleType convertStringToScaleType(const std::string &scaleType)
Convert the given string to its ScaleType enum value.
std::vector< Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > > H
The current Hessenberg matrix.
Collection of types and exceptions used within the Belos solvers.
ComboType
The test can be either the AND of all the component tests, or the OR of all the component tests...
virtual ~PseudoBlockGmresSolMgr()
Destructor.
Belos's basic output manager for sending information of select verbosity levels to the appropriate ou...
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > ¶ms) override
Set the parameters the solver manager should use to solve the linear problem.
Class which manages the output and verbosity of the Belos solvers.
Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > problem_
The current linear problem to solve.
bool is_null(const boost::shared_ptr< T > &p)
Teuchos::RCP< Teuchos::ParameterList > params_
PseudoBlockGmresSolMgr()
Empty constructor.
PseudoBlockGmresSolMgrLinearProblemFailure(const std::string &what_arg)
MagnitudeType resScaleFactor_
ScaleType
The type of scaling to use on the residual norm value.
Teuchos::RCP< StatusTestMaxIters< ScalarType, MV, OP > > maxIterTest_
Convergence test using the implicit residual norm(s), with an explicit residual norm(s) check for los...
Teuchos::RCP< StatusTest< ScalarType, MV, OP > > debugStatusTest_
T & get(ParameterList &l, const std::string &name)
PseudoBlockGmresSolMgrOrthoFailure(const std::string &what_arg)
bool is_null(const std::shared_ptr< T > &p)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
A factory class for generating StatusTestOutput objects.
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...
An abstract class of StatusTest for stopping criteria using residual norms.
std::vector< Teuchos::RCP< const MV > > V
The current Krylov basis.
An implementation of StatusTestResNorm using a family of residual norms.
RCP< ParameterList > sublist(const RCP< ParameterList > ¶mList, const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
static constexpr int outputFreq_default_
Ordinal numParams() const
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
Teuchos::RCP< StatusTestOutput< ScalarType, MV, OP > > outputTest_
static const double convTol
Default convergence tolerance.
A pure virtual class for defining the status tests for the Belos iterative solvers.
Structure to contain pointers to PseudoBlockGmresIter state variables.
Teuchos::Array< Teuchos::RCP< Teuchos::Time > > getTimers() const
Return the timers for this object.
static constexpr const char * label_default_
A factory class for generating StatusTestOutput objects.
Teuchos::RCP< std::ostream > outputStream_
Teuchos::RCP< StatusTest< ScalarType, MV, OP > > convTest_
Teuchos::RCP< MatOrthoManager< ScalarType, MV, OP > > ortho_
Traits class which defines basic operations on multivectors.
MagnitudeType achievedTol_
bool isParameter(const std::string &name) const
A Belos::StatusTest class for specifying a maximum number of iterations.
MultiVecTraits< ScalarType, MV > MVT
static constexpr int blockSize_default_
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. ...
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.
static constexpr int numBlocks_default_
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)
bool isSublist(const std::string &name) const
static constexpr int defQuorum_default_
static constexpr bool showMaxResNormOnly_default_
static constexpr int verbosity_default_
Teuchos::ScalarTraits< ScalarType >::magnitudeType MagnitudeType
PseudoBlockGmresSolMgrLinearProblemFailure is thrown when the linear problem is not setup (i...
std::vector< Teuchos::RCP< const Teuchos::SerialDenseVector< int, MagnitudeType > > > cs
A linear system to solve, and its associated information.
Teuchos::RCP< std::map< std::string, Teuchos::RCP< StatusTest< ScalarType, MV, OP > > > > taggedTests_
Class which describes the linear problem to be solved by the iterative solver.
void reset(const ResetType type) override
Performs a reset of the solver manager specified by the ResetType. This informs the solver manager th...
Interface to standard and "pseudoblock" GMRES.
int getNumIters() const override
Iteration count for the most recent call to solve().
Teuchos::RCP< Teuchos::Time > timerSolve_
ReturnType
Whether the Belos solve converged for all linear systems.
bool checkStatusTest()
Check current status tests against current linear problem.
The Belos::SolverManager is a templated virtual base class that defines the basic interface that any ...
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
A list of valid default parameters for this solver.
static constexpr const char * orthoType_default_
static constexpr int maxRestarts_default_
void validateParameters(ParameterList const &validParamList, int const depth=1000, EValidateUsed const validateUsed=VALIDATE_USED_ENABLED, EValidateDefaults const validateDefaults=VALIDATE_DEFAULTS_ENABLED) const
static const EVerbosityLevel verbLevel_default
static constexpr const char * impResScale_default_
bool isLOADetected() const override
Whether a "loss of accuracy" was detected during the last solve().
bool nonnull(const boost::shared_ptr< T > &p)
std::vector< Teuchos::RCP< const Teuchos::SerialDenseVector< int, ScalarType > > > Z
The current right-hand side of the least squares system RY = Z.
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.
void setProblem(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem) override
Set the linear problem to solve.
Belos concrete class for performing the pseudo-block GMRES iteration.
Teuchos::RCP< StatusTest< ScalarType, MV, OP > > sTest_
bool isType(const std::string &name) const
Teuchos::RCP< SolverManager< ScalarType, MV, OP > > clone() const override
clone for Inverted Injection (DII)
An implementation of the Belos::MatOrthoManager that performs orthogonalization using (potentially) m...
A class for extending the status testing capabilities of Belos via logical combinations.
std::string description() const override
Return a one-line description of this object.
const StatusTestResNorm< ScalarType, MV, OP > * getResidualStatusTest() const
Return the residual status test.
Class which defines basic traits for the operator type.
int curDim
The current dimension of the reduction.
Teuchos::ScalarTraits< MagnitudeType > MT
static constexpr const char * expResScale_default_
static constexpr int maxIters_default_
Parent class to all Belos exceptions.
ReturnType solve() override
This method performs possibly repeated calls to the underlying linear solver's iterate() routine unti...
Default parameters common to most Belos solvers.
std::vector< Teuchos::RCP< const Teuchos::SerialDenseVector< int, ScalarType > > > sn
The current Given's rotation coefficients.
static const double resScaleFactor
User-defined residual scaling factor.
static constexpr int outputStyle_default_
Belos header file which uses auto-configuration information to include necessary C++ headers...
Teuchos::RCP< StatusTestResNorm< ScalarType, MV, OP > > expConvTest_
OperatorTraits< ScalarType, MV, OP > OPT
This class implements the pseudo-block GMRES iteration, where a block Krylov subspace is constructed ...
Factory to build a set of status tests from a parameter list.
Teuchos::RCP< const Teuchos::ParameterList > getCurrentParameters() const override
The current parameters for this solver.
Teuchos::RCP< OutputManager< ScalarType > > printer_
PseudoBlockGmresIterOrthoFailure is thrown when the orthogonalization manager is unable to generate o...
const LinearProblem< ScalarType, MV, OP > & getProblem() const override
Return a reference to the linear problem being solved by this solver manager.
Teuchos::ScalarTraits< ScalarType > SCT
StatusTestCombo< ScalarType, MV, OP >::ComboType comboType_
Teuchos::RCP< StatusTestResNorm< ScalarType, MV, OP > > impConvTest_
PseudoBlockGmresSolMgrOrthoFailure is thrown when the orthogonalization manager is unable to generate...
virtual void setUserConvStatusTest(const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &userConvStatusTest, const typename StatusTestCombo< ScalarType, MV, OP >::ComboType &comboType=StatusTestCombo< ScalarType, MV, OP >::SEQ) override
Set a custom status test.
Teuchos::RCP< StatusTest< ScalarType, MV, OP > > userConvStatusTest_
MagnitudeType achievedTol() const override
Tolerance achieved by the last solve() invocation.
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 > ¶ms)
Return an instance of the specified MatOrthoManager subclass.
MagnitudeType orthoKappa_