10 #ifndef BELOS_GCRODR_SOLMGR_HPP
11 #define BELOS_GCRODR_SOLMGR_HPP
33 #ifdef BELOS_TEUCHOS_TIME_MONITOR
35 #endif // BELOS_TEUCHOS_TIME_MONITOR
36 #if defined(HAVE_TEUCHOSCORE_CXX11)
37 # include <type_traits>
38 #endif // defined(HAVE_TEUCHOSCORE_CXX11)
125 template<
class ScalarType,
class MV,
class OP,
126 const bool lapackSupportsScalarType =
150 template<
class ScalarType,
class MV,
class OP>
155 #if defined(HAVE_TEUCHOSCORE_CXX11)
156 # if defined(HAVE_TEUCHOS_COMPLEX)
157 #if defined(HAVE_TEUCHOS_LONG_DOUBLE)
158 static_assert (std::is_same<ScalarType, std::complex<float> >::value ||
159 std::is_same<ScalarType, std::complex<double> >::value ||
160 std::is_same<ScalarType, Kokkos::complex<double> >::value ||
161 std::is_same<ScalarType, float>::value ||
162 std::is_same<ScalarType, double>::value ||
163 std::is_same<ScalarType, long double>::value,
164 "Belos::GCRODRSolMgr: ScalarType must be one of the four "
165 "types (S,D,C,Z) supported by LAPACK or long double (largely not impl'd).");
167 static_assert (std::is_same<ScalarType, std::complex<float> >::value ||
168 std::is_same<ScalarType, std::complex<double> >::value ||
169 std::is_same<ScalarType, Kokkos::complex<double> >::value ||
170 std::is_same<ScalarType, float>::value ||
171 std::is_same<ScalarType, double>::value,
172 "Belos::GCRODRSolMgr: ScalarType must be one of the four "
173 "types (S,D,C,Z) supported by LAPACK.");
176 #if defined(HAVE_TEUCHOS_LONG_DOUBLE)
177 static_assert (std::is_same<ScalarType, float>::value ||
178 std::is_same<ScalarType, double>::value ||
179 std::is_same<ScalarType, long double>::value,
180 "Belos::GCRODRSolMgr: ScalarType must be float, double or long double. "
181 "Complex arithmetic support is currently disabled. To "
182 "enable it, set Teuchos_ENABLE_COMPLEX=ON.");
184 static_assert (std::is_same<ScalarType, float>::value ||
185 std::is_same<ScalarType, double>::value,
186 "Belos::GCRODRSolMgr: ScalarType must be float or double. "
187 "Complex arithmetic support is currently disabled. To "
188 "enable it, set Teuchos_ENABLE_COMPLEX=ON.");
190 # endif // defined(HAVE_TEUCHOS_COMPLEX)
191 #endif // defined(HAVE_TEUCHOSCORE_CXX11)
301 return Teuchos::tuple(timerSolve_);
345 bool set = problem_->setProblem();
347 throw "Could not set problem.";
391 std::string description()
const override;
401 void initializeStateStorage();
410 int getHarmonicVecs1(
int m,
419 int getHarmonicVecs2(
int keff,
int m,
425 void sort(std::vector<MagnitudeType>& dlist,
int n, std::vector<int>& iperm);
453 static constexpr
double orthoKappa_default_ = 0.0;
454 static constexpr
int maxRestarts_default_ = 100;
455 static constexpr
int maxIters_default_ = 1000;
456 static constexpr
int numBlocks_default_ = 50;
457 static constexpr
int blockSize_default_ = 1;
458 static constexpr
int recycledBlocks_default_ = 5;
461 static constexpr
int outputFreq_default_ = -1;
462 static constexpr
const char * impResScale_default_ =
"Norm of Preconditioned Initial Residual";
463 static constexpr
const char * expResScale_default_ =
"Norm of Initial Residual";
464 static constexpr
const char * label_default_ =
"Belos";
465 static constexpr
const char * orthoType_default_ =
"ICGS";
520 template<
class ScalarType,
class MV,
class OP>
530 template<
class ScalarType,
class MV,
class OP>
543 "Belos::GCRODRSolMgr constructor: The solver manager's "
544 "constructor needs the linear problem argument 'problem' "
559 template<
class ScalarType,
class MV,
class OP>
561 outputStream_ = Teuchos::rcpFromRef(std::cout);
563 orthoKappa_ = orthoKappa_default_;
564 maxRestarts_ = maxRestarts_default_;
565 maxIters_ = maxIters_default_;
566 numBlocks_ = numBlocks_default_;
567 recycledBlocks_ = recycledBlocks_default_;
568 verbosity_ = verbosity_default_;
569 outputStyle_ = outputStyle_default_;
570 outputFreq_ = outputFreq_default_;
571 orthoType_ = orthoType_default_;
572 impResScale_ = impResScale_default_;
573 expResScale_ = expResScale_default_;
574 label_ = label_default_;
576 builtRecycleSpace_ =
false;
592 template<
class ScalarType,
class MV,
class OP>
597 using Teuchos::isParameterType;
598 using Teuchos::getParameter;
601 using Teuchos::parameterList;
604 using Teuchos::rcp_dynamic_cast;
605 using Teuchos::rcpFromRef;
612 RCP<const ParameterList> defaultParams = getValidParameters();
630 if (params_.is_null()) {
631 params_ = parameterList (*defaultParams);
639 if (params_ != params) {
645 params_ = parameterList (*params);
680 params_->validateParametersAndSetDefaults (*defaultParams);
685 maxRestarts_ = params->
get(
"Maximum Restarts", maxRestarts_default_);
688 params_->set (
"Maximum Restarts", maxRestarts_);
693 maxIters_ = params->
get (
"Maximum Iterations", maxIters_default_);
696 params_->set (
"Maximum Iterations", maxIters_);
697 if (! maxIterTest_.is_null())
698 maxIterTest_->setMaxIters (maxIters_);
703 numBlocks_ = params->
get (
"Num Blocks", numBlocks_default_);
705 "Belos::GCRODRSolMgr: The \"Num Blocks\" parameter must "
706 "be strictly positive, but you specified a value of "
707 << numBlocks_ <<
".");
709 params_->set (
"Num Blocks", numBlocks_);
714 recycledBlocks_ = params->
get (
"Num Recycled Blocks",
715 recycledBlocks_default_);
717 "Belos::GCRODRSolMgr: The \"Num Recycled Blocks\" "
718 "parameter must be strictly positive, but you specified "
719 "a value of " << recycledBlocks_ <<
".");
721 "Belos::GCRODRSolMgr: The \"Num Recycled Blocks\" "
722 "parameter must be less than the \"Num Blocks\" "
723 "parameter, but you specified \"Num Recycled Blocks\" "
724 "= " << recycledBlocks_ <<
" and \"Num Blocks\" = "
725 << numBlocks_ <<
".");
727 params_->set(
"Num Recycled Blocks", recycledBlocks_);
734 std::string tempLabel = params->
get (
"Timer Label", label_default_);
737 if (tempLabel != label_) {
739 params_->set (
"Timer Label", label_);
740 std::string solveLabel = label_ +
": GCRODRSolMgr total solve time";
741 #ifdef BELOS_TEUCHOS_TIME_MONITOR
745 ortho_->setLabel( label_ );
752 if (isParameterType<int> (*params,
"Verbosity")) {
753 verbosity_ = params->
get (
"Verbosity", verbosity_default_);
755 verbosity_ = (int) getParameter<Belos::MsgType> (*params,
"Verbosity");
758 params_->set (
"Verbosity", verbosity_);
761 if (! printer_.is_null())
762 printer_->setVerbosity (verbosity_);
767 if (isParameterType<int> (*params,
"Output Style")) {
768 outputStyle_ = params->
get (
"Output Style", outputStyle_default_);
770 outputStyle_ = (int) getParameter<OutputType> (*params,
"Output Style");
774 params_->set (
"Output Style", outputStyle_);
794 outputStream_ = getParameter<RCP<std::ostream> > (*params,
"Output Stream");
795 }
catch (InvalidParameter&) {
796 outputStream_ = rcpFromRef (std::cout);
803 if (outputStream_.is_null()) {
807 params_->set (
"Output Stream", outputStream_);
810 if (! printer_.is_null()) {
811 printer_->setOStream (outputStream_);
818 outputFreq_ = params->
get (
"Output Frequency", outputFreq_default_);
822 params_->set(
"Output Frequency", outputFreq_);
823 if (! outputTest_.is_null())
824 outputTest_->setOutputFrequency (outputFreq_);
831 if (printer_.is_null()) {
842 bool changedOrthoType =
false;
844 const std::string& tempOrthoType =
845 params->
get (
"Orthogonalization", orthoType_default_);
848 std::ostringstream os;
849 os <<
"Belos::GCRODRSolMgr: Invalid orthogonalization name \""
850 << tempOrthoType <<
"\". The following are valid options "
851 <<
"for the \"Orthogonalization\" name parameter: ";
853 throw std::invalid_argument (os.str());
855 if (tempOrthoType != orthoType_) {
856 changedOrthoType =
true;
857 orthoType_ = tempOrthoType;
859 params_->set (
"Orthogonalization", orthoType_);
875 RCP<ParameterList> orthoParams;
878 using Teuchos::sublist;
880 const std::string paramName (
"Orthogonalization Parameters");
883 orthoParams = sublist (params_, paramName,
true);
884 }
catch (InvalidParameter&) {
891 orthoParams = sublist (params_, paramName,
true);
895 "Failed to get orthogonalization parameters. "
896 "Please report this bug to the Belos developers.");
901 if (ortho_.is_null() || changedOrthoType) {
907 label_, orthoParams);
916 RCP<PLA> pla = rcp_dynamic_cast<PLA> (ortho_);
922 label_, orthoParams);
924 pla->setParameterList (orthoParams);
936 if (params->
isParameter (
"Orthogonalization Constant")) {
939 orthoKappa = params->
get (
"Orthogonalization Constant", orthoKappa);
942 orthoKappa = params->
get (
"Orthogonalization Constant", orthoKappa_default_);
945 if (orthoKappa > 0) {
946 orthoKappa_ = orthoKappa;
948 params_->set(
"Orthogonalization Constant", orthoKappa_);
950 if (orthoType_ ==
"DGKS" && ! ortho_.is_null()) {
957 rcp_dynamic_cast<ortho_man_type>(ortho_)->setDepTol (orthoKappa_);
967 if (params->
isParameter(
"Convergence Tolerance")) {
969 convTol_ = params->
get (
"Convergence Tolerance",
977 params_->set (
"Convergence Tolerance", convTol_);
978 if (! impConvTest_.is_null())
979 impConvTest_->setTolerance (convTol_);
980 if (! expConvTest_.is_null())
981 expConvTest_->setTolerance (convTol_);
985 if (params->
isParameter (
"Implicit Residual Scaling")) {
986 std::string tempImpResScale =
987 getParameter<std::string> (*params,
"Implicit Residual Scaling");
990 if (impResScale_ != tempImpResScale) {
992 impResScale_ = tempImpResScale;
995 params_->set(
"Implicit Residual Scaling", impResScale_);
1005 if (! impConvTest_.is_null()) {
1011 impConvTest_ = null;
1018 if (params->
isParameter(
"Explicit Residual Scaling")) {
1019 std::string tempExpResScale =
1020 getParameter<std::string> (*params,
"Explicit Residual Scaling");
1023 if (expResScale_ != tempExpResScale) {
1025 expResScale_ = tempExpResScale;
1028 params_->set(
"Explicit Residual Scaling", expResScale_);
1031 if (! expConvTest_.is_null()) {
1037 expConvTest_ = null;
1048 if (maxIterTest_.is_null())
1053 if (impConvTest_.is_null()) {
1054 impConvTest_ =
rcp (
new StatusTestResNorm_t (convTol_));
1060 if (expConvTest_.is_null()) {
1061 expConvTest_ =
rcp (
new StatusTestResNorm_t (convTol_));
1062 expConvTest_->defineResForm (StatusTestResNorm_t::Explicit,
Belos::TwoNorm);
1068 if (convTest_.is_null()) {
1069 convTest_ =
rcp (
new StatusTestCombo_t (StatusTestCombo_t::SEQ,
1077 sTest_ =
rcp (
new StatusTestCombo_t (StatusTestCombo_t::OR,
1083 outputTest_ = stoFactory.
create (printer_, sTest_, outputFreq_,
1087 std::string solverDesc =
" GCRODR ";
1088 outputTest_->setSolverDesc( solverDesc );
1091 if (timerSolve_.is_null()) {
1092 std::string solveLabel = label_ +
": GCRODRSolMgr total solve time";
1093 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1103 template<
class ScalarType,
class MV,
class OP>
1108 using Teuchos::parameterList;
1111 static RCP<const ParameterList> validPL;
1113 RCP<ParameterList> pl = parameterList ();
1117 "The relative residual tolerance that needs to be achieved by the\n"
1118 "iterative solver in order for the linear system to be declared converged.");
1119 pl->set(
"Maximum Restarts", static_cast<int>(maxRestarts_default_),
1120 "The maximum number of cycles allowed for each\n"
1121 "set of RHS solved.");
1122 pl->set(
"Maximum Iterations", static_cast<int>(maxIters_default_),
1123 "The maximum number of iterations allowed for each\n"
1124 "set of RHS solved.");
1128 pl->set(
"Block Size", static_cast<int>(blockSize_default_),
1129 "Block Size Parameter -- currently must be 1 for GCRODR");
1130 pl->set(
"Num Blocks", static_cast<int>(numBlocks_default_),
1131 "The maximum number of vectors allowed in the Krylov subspace\n"
1132 "for each set of RHS solved.");
1133 pl->set(
"Num Recycled Blocks", static_cast<int>(recycledBlocks_default_),
1134 "The maximum number of vectors in the recycled subspace." );
1135 pl->set(
"Verbosity", static_cast<int>(verbosity_default_),
1136 "What type(s) of solver information should be outputted\n"
1137 "to the output stream.");
1138 pl->set(
"Output Style", static_cast<int>(outputStyle_default_),
1139 "What style is used for the solver information outputted\n"
1140 "to the output stream.");
1141 pl->set(
"Output Frequency", static_cast<int>(outputFreq_default_),
1142 "How often convergence information should be outputted\n"
1143 "to the output stream.");
1144 pl->set(
"Output Stream", Teuchos::rcpFromRef(std::cout),
1145 "A reference-counted pointer to the output stream where all\n"
1146 "solver output is sent.");
1147 pl->set(
"Implicit Residual Scaling", static_cast<const char *>(impResScale_default_),
1148 "The type of scaling used in the implicit residual convergence test.");
1149 pl->set(
"Explicit Residual Scaling", static_cast<const char *>(expResScale_default_),
1150 "The type of scaling used in the explicit residual convergence test.");
1151 pl->set(
"Timer Label", static_cast<const char *>(label_default_),
1152 "The string to use as a prefix for the timer labels.");
1155 pl->set(
"Orthogonalization", static_cast<const char *>(orthoType_default_),
1156 "The type of orthogonalization to use. Valid options: " +
1158 RCP<const ParameterList> orthoParams =
1160 pl->
set (
"Orthogonalization Parameters", *orthoParams,
1161 "Parameters specific to the type of orthogonalization used.");
1163 pl->
set(
"Orthogonalization Constant",static_cast<MagnitudeType>(orthoKappa_default_),
1164 "When using DGKS orthogonalization: the \"depTol\" constant, used "
1165 "to determine whether another step of classical Gram-Schmidt is "
1166 "necessary. Otherwise ignored.");
1173 template<
class ScalarType,
class MV,
class OP>
1188 "Belos::GCRODRSolMgr::initializeStateStorage(): Cannot generate a Krylov basis with dimension larger the operator!");
1192 U_ = MVT::Clone( *rhsMV, recycledBlocks_+1 );
1196 if (MVT::GetNumberVecs(*U_) < recycledBlocks_+1) {
1198 U_ = MVT::Clone( *tmp, recycledBlocks_+1 );
1204 C_ = MVT::Clone( *rhsMV, recycledBlocks_+1 );
1208 if (MVT::GetNumberVecs(*C_) < recycledBlocks_+1) {
1210 C_ = MVT::Clone( *tmp, recycledBlocks_+1 );
1216 V_ = MVT::Clone( *rhsMV, numBlocks_+1 );
1220 if (MVT::GetNumberVecs(*V_) < numBlocks_+1) {
1222 V_ = MVT::Clone( *tmp, numBlocks_+1 );
1228 U1_ = MVT::Clone( *rhsMV, recycledBlocks_+1 );
1232 if (MVT::GetNumberVecs(*U1_) < recycledBlocks_+1) {
1234 U1_ = MVT::Clone( *tmp, recycledBlocks_+1 );
1240 C1_ = MVT::Clone( *rhsMV, recycledBlocks_+1 );
1244 if (MVT::GetNumberVecs(*C1_) < recycledBlocks_+1) {
1246 C1_ = MVT::Clone( *tmp, recycledBlocks_+1 );
1252 r_ = MVT::Clone( *rhsMV, 1 );
1255 tau_.resize(recycledBlocks_+1);
1258 work_.resize(recycledBlocks_+1);
1261 ipiv_.resize(recycledBlocks_+1);
1267 if ( (H2_->numRows() != numBlocks_+recycledBlocks_+2) || (H2_->numCols() != numBlocks_+recycledBlocks_+1) )
1268 H2_->reshape( numBlocks_+recycledBlocks_+2, numBlocks_+recycledBlocks_+1 );
1270 H2_->putScalar(zero);
1276 if ( (R_->numRows() != recycledBlocks_+1) || (R_->numCols() != recycledBlocks_+1) )
1277 R_->reshape( recycledBlocks_+1, recycledBlocks_+1 );
1279 R_->putScalar(zero);
1285 if ( (PP_->numRows() != numBlocks_+recycledBlocks_+2) || (PP_->numCols() != recycledBlocks_+1) )
1286 PP_->reshape( numBlocks_+recycledBlocks_+2, recycledBlocks_+1 );
1293 if ( (HP_->numRows() != numBlocks_+recycledBlocks_+2) || (HP_->numCols() != numBlocks_+recycledBlocks_+1) )
1294 HP_->reshape( numBlocks_+recycledBlocks_+2, numBlocks_+recycledBlocks_+1 );
1302 template<
class ScalarType,
class MV,
class OP>
1310 if (!isSet_) { setParameters( params_ ); }
1314 std::vector<int> index(numBlocks_+1);
1321 int numRHS2Solve = MVT::GetNumberVecs( *(problem_->getRHS()) );
1322 std::vector<int> currIdx(1);
1326 problem_->setLSIndex( currIdx );
1329 ptrdiff_t dim = MVT::GetGlobalLength( *(problem_->getRHS()) );
1330 if (static_cast<ptrdiff_t>(numBlocks_) > dim) {
1331 numBlocks_ = Teuchos::as<int>(dim);
1333 "Warning! Requested Krylov subspace dimension is larger than operator dimension!" << std::endl <<
1334 " The maximum number of blocks allowed for the Krylov subspace will be adjusted to " << numBlocks_ << std::endl;
1335 params_->set(
"Num Blocks", numBlocks_);
1339 bool isConverged =
true;
1342 initializeStateStorage();
1348 plist.
set(
"Num Blocks",numBlocks_);
1349 plist.
set(
"Recycled Blocks",recycledBlocks_);
1354 RCP<GCRODRIter<ScalarType,MV,OP> > gcrodr_iter;
1357 int prime_iterations = 0;
1361 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1365 while ( numRHS2Solve > 0 ) {
1368 builtRecycleSpace_ =
false;
1371 outputTest_->reset();
1379 "Belos::GCRODRSolMgr::solve(): Requested size of recycled subspace is not consistent with the current recycle subspace.");
1381 printer_->stream(
Debug) <<
" Now solving RHS index " << currIdx[0] <<
" using recycled subspace of dimension " << keff << std::endl << std::endl;
1384 for (
int ii=0; ii<keff; ++ii) { index[ii] = ii; }
1385 RCP<const MV> Utmp = MVT::CloneView( *U_, index );
1386 RCP<MV> Ctmp = MVT::CloneViewNonConst( *C_, index );
1387 problem_->apply( *Utmp, *Ctmp );
1389 RCP<MV> U1tmp = MVT::CloneViewNonConst( *U1_, index );
1394 int rank = ortho_->normalize(*Ctmp,
rcp(&Rtmp,
false));
1407 work_.resize(lwork);
1412 MVT::MvTimesMatAddMv( one, *Utmp, Rtmp, zero, *U1tmp );
1417 for (
int ii=0; ii<keff; ++ii) { index[ii] = ii; }
1418 Ctmp = MVT::CloneViewNonConst( *C_, index );
1419 Utmp = MVT::CloneView( *U_, index );
1423 problem_->computeCurrPrecResVec( &*r_ );
1424 MVT::MvTransMv( one, *Ctmp, *r_, Ctr );
1427 RCP<MV> update = MVT::Clone( *problem_->getCurrLHSVec(), 1 );
1428 MVT::MvInit( *update, 0.0 );
1429 MVT::MvTimesMatAddMv( one, *Utmp, Ctr, one, *update );
1430 problem_->updateSolution( update,
true );
1433 MVT::MvTimesMatAddMv( -one, *Ctmp, Ctr, one, *r_ );
1436 prime_iterations = 0;
1442 printer_->stream(
Debug) <<
" No recycled subspace available for RHS index " << currIdx[0] << std::endl << std::endl;
1447 primeList.
set(
"Num Blocks",numBlocks_);
1448 primeList.
set(
"Recycled Blocks",0);
1451 RCP<GCRODRIter<ScalarType,MV,OP> > gcrodr_prime_iter;
1455 problem_->computeCurrPrecResVec( &*r_ );
1456 index.resize( 1 ); index[0] = 0;
1457 RCP<MV> v0 = MVT::CloneViewNonConst( *V_, index );
1458 MVT::SetBlock(*r_,index,*v0);
1462 index.resize( numBlocks_+1 );
1463 for (
int ii=0; ii<(numBlocks_+1); ++ii) { index[ii] = ii; }
1464 newstate.
V = MVT::CloneViewNonConst( *V_, index );
1470 gcrodr_prime_iter->initialize(newstate);
1473 bool primeConverged =
false;
1475 gcrodr_prime_iter->iterate();
1478 if ( convTest_->getStatus() ==
Passed ) {
1480 primeConverged =
true;
1485 gcrodr_prime_iter->updateLSQR( gcrodr_prime_iter->getCurSubspaceDim() );
1488 sTest_->checkStatus( &*gcrodr_prime_iter );
1489 if (convTest_->getStatus() ==
Passed)
1490 primeConverged =
true;
1494 achievedTol_ = MT::one();
1496 MVT::MvInit( *X, SCT::zero() );
1497 printer_->stream(
Warnings) <<
"Belos::GCRODRSolMgr::solve(): Warning! NaN has been detected!"
1501 catch (
const std::exception &e) {
1502 printer_->stream(
Errors) <<
"Error! Caught exception in GCRODRIter::iterate() at iteration "
1503 << gcrodr_prime_iter->getNumIters() << std::endl
1504 << e.what() << std::endl;
1508 prime_iterations = gcrodr_prime_iter->getNumIters();
1511 RCP<MV> update = gcrodr_prime_iter->getCurrentUpdate();
1512 problem_->updateSolution( update,
true );
1515 newstate = gcrodr_prime_iter->getState();
1523 if (recycledBlocks_ < p+1) {
1527 keff = getHarmonicVecs1( p, *newstate.
H, *PPtmp );
1532 for (
int ii=0; ii<keff; ++ii) { index[ii] = ii; }
1533 RCP<MV> Ctmp = MVT::CloneViewNonConst( *C_, index );
1534 RCP<MV> Utmp = MVT::CloneViewNonConst( *U_, index );
1535 RCP<MV> U1tmp = MVT::CloneViewNonConst( *U1_, index );
1537 for (
int ii=0; ii < p; ++ii) { index[ii] = ii; }
1538 RCP<const MV> Vtmp = MVT::CloneView( *V_, index );
1542 MVT::MvTimesMatAddMv( one, *Vtmp, *PPtmp, zero, *U1tmp );
1559 HPtmp.
stride (), &tau_[0], &work_[0], lwork, &info);
1562 " LAPACK's _GEQRF failed to compute a workspace size.");
1571 work_.resize (lwork);
1573 HPtmp.
stride (), &tau_[0], &work_[0], lwork, &info);
1576 " LAPACK's _GEQRF failed to compute a QR factorization.");
1581 for (
int ii = 0; ii < keff; ++ii) {
1582 for (
int jj = ii; jj < keff; ++jj) {
1583 Rtmp(ii,jj) = HPtmp(ii,jj);
1595 "LAPACK's _UNGQR failed to construct the Q factor.");
1600 index.resize (p + 1);
1601 for (
int ii = 0; ii < (p+1); ++ii) {
1604 Vtmp = MVT::CloneView( *V_, index );
1605 MVT::MvTimesMatAddMv( one, *Vtmp, HPtmp, zero, *Ctmp );
1616 "LAPACK's _GETRF failed to compute an LU factorization.");
1626 work_.resize(lwork);
1630 "LAPACK's _GETRI failed to invert triangular matrix.");
1633 MVT::MvTimesMatAddMv( one, *U1tmp, Rtmp, zero, *Utmp );
1635 printer_->stream(
Debug)
1636 <<
" Generated recycled subspace using RHS index " << currIdx[0]
1637 <<
" of dimension " << keff << std::endl << std::endl;
1642 if (primeConverged) {
1644 problem_->setCurrLS();
1648 if (numRHS2Solve > 0) {
1650 problem_->setLSIndex (currIdx);
1653 currIdx.resize (numRHS2Solve);
1663 gcrodr_iter->setSize( keff, numBlocks_ );
1666 gcrodr_iter->resetNumIters(prime_iterations);
1669 outputTest_->resetNumCalls();
1672 problem_->computeCurrPrecResVec( &*r_ );
1673 index.resize( 1 ); index[0] = 0;
1674 RCP<MV> v0 = MVT::CloneViewNonConst( *V_, index );
1675 MVT::SetBlock(*r_,index,*v0);
1679 index.resize( numBlocks_+1 );
1680 for (
int ii=0; ii<(numBlocks_+1); ++ii) { index[ii] = ii; }
1681 newstate.
V = MVT::CloneViewNonConst( *V_, index );
1682 index.resize( keff );
1683 for (
int ii=0; ii<keff; ++ii) { index[ii] = ii; }
1684 newstate.
C = MVT::CloneViewNonConst( *C_, index );
1685 newstate.
U = MVT::CloneViewNonConst( *U_, index );
1689 gcrodr_iter->initialize(newstate);
1692 int numRestarts = 0;
1697 gcrodr_iter->iterate();
1704 if ( convTest_->getStatus() ==
Passed ) {
1713 else if ( maxIterTest_->getStatus() ==
Passed ) {
1715 isConverged =
false;
1723 else if ( gcrodr_iter->getCurSubspaceDim() == gcrodr_iter->getMaxSubspaceDim() ) {
1728 RCP<MV> update = gcrodr_iter->getCurrentUpdate();
1729 problem_->updateSolution( update,
true );
1731 buildRecycleSpace2(gcrodr_iter);
1733 printer_->stream(
Debug)
1734 <<
" Generated new recycled subspace using RHS index "
1735 << currIdx[0] <<
" of dimension " << keff << std::endl
1739 if (numRestarts >= maxRestarts_) {
1740 isConverged =
false;
1745 printer_->stream(
Debug)
1746 <<
" Performing restart number " << numRestarts <<
" of "
1747 << maxRestarts_ << std::endl << std::endl;
1750 problem_->computeCurrPrecResVec( &*r_ );
1751 index.resize( 1 ); index[0] = 0;
1752 RCP<MV> v00 = MVT::CloneViewNonConst( *V_, index );
1753 MVT::SetBlock(*r_,index,*v00);
1757 index.resize( numBlocks_+1 );
1758 for (
int ii=0; ii<(numBlocks_+1); ++ii) { index[ii] = ii; }
1759 restartState.
V = MVT::CloneViewNonConst( *V_, index );
1760 index.resize( keff );
1761 for (
int ii=0; ii<keff; ++ii) { index[ii] = ii; }
1762 restartState.
U = MVT::CloneViewNonConst( *U_, index );
1763 restartState.
C = MVT::CloneViewNonConst( *C_, index );
1767 gcrodr_iter->initialize(restartState);
1781 true, std::logic_error,
"Belos::GCRODRSolMgr::solve: "
1782 "Invalid return from GCRODRIter::iterate().");
1787 gcrodr_iter->updateLSQR( gcrodr_iter->getCurSubspaceDim() );
1790 sTest_->checkStatus( &*gcrodr_iter );
1791 if (convTest_->getStatus() !=
Passed)
1792 isConverged =
false;
1795 catch (
const std::exception& e) {
1797 <<
"Error! Caught exception in GCRODRIter::iterate() at iteration "
1798 << gcrodr_iter->getNumIters() << std::endl << e.what() << std::endl;
1805 RCP<MV> update = gcrodr_iter->getCurrentUpdate();
1806 problem_->updateSolution( update,
true );
1809 problem_->setCurrLS();
1814 if (!builtRecycleSpace_) {
1815 buildRecycleSpace2(gcrodr_iter);
1816 printer_->stream(
Debug)
1817 <<
" Generated new recycled subspace using RHS index " << currIdx[0]
1818 <<
" of dimension " << keff << std::endl << std::endl;
1823 if (numRHS2Solve > 0) {
1825 problem_->setLSIndex (currIdx);
1828 currIdx.resize (numRHS2Solve);
1836 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1842 #endif // BELOS_TEUCHOS_TIME_MONITOR
1845 numIters_ = maxIterTest_->getNumIters ();
1857 const std::vector<MagnitudeType>* pTestValues = expConvTest_->getTestValue();
1858 if (pTestValues == NULL || pTestValues->size() < 1) {
1859 pTestValues = impConvTest_->getTestValue();
1862 "Belos::GCRODRSolMgr::solve(): The implicit convergence test's getTestValue() "
1863 "method returned NULL. Please report this bug to the Belos developers.");
1865 "Belos::GCRODRSolMgr::solve(): The implicit convergence test's getTestValue() "
1866 "method returned a vector of length zero. Please report this bug to the "
1867 "Belos developers.");
1872 achievedTol_ = *std::max_element (pTestValues->begin(), pTestValues->end());
1879 template<
class ScalarType,
class MV,
class OP>
1885 std::vector<MagnitudeType> d(keff);
1886 std::vector<ScalarType> dscalar(keff);
1887 std::vector<int> index(numBlocks_+1);
1899 for (
int ii=0; ii<keff; ++ii) { index[ii] = ii; }
1902 dscalar.resize(keff);
1903 MVT::MvNorm( *Utmp, d );
1904 for (
int i=0; i<keff; ++i) {
1906 dscalar[i] = (ScalarType)d[i];
1908 MVT::MvScale( *Utmp, dscalar );
1915 for (
int i=0; i<keff; ++i) {
1916 (*H2tmp)(i,i) = d[i];
1924 keff_new = getHarmonicVecs2( keff, p, *H2tmp, oldState.
V, PPtmp );
1933 index.resize( keff );
1934 for (
int ii=0; ii<keff; ++ii) { index[ii] = ii; }
1936 index.resize( keff_new );
1937 for (
int ii=0; ii<keff_new; ++ii) { index[ii] = ii; }
1938 U1tmp = MVT::CloneViewNonConst( *U1_, index );
1940 MVT::MvTimesMatAddMv( one, *Utmp, PPtmp, zero, *U1tmp );
1946 for (
int ii=0; ii < p; ii++) { index[ii] = ii; }
1949 MVT::MvTimesMatAddMv( one, *Vtmp, PPtmp, one, *U1tmp );
1960 int info = 0, lwork = -1;
1961 tau_.resize (keff_new);
1963 HPtmp.
stride (), &tau_[0], &work_[0], lwork, &info);
1966 "LAPACK's _GEQRF failed to compute a workspace size.");
1973 work_.resize (lwork);
1975 HPtmp.
stride (), &tau_[0], &work_[0], lwork, &info);
1978 "LAPACK's _GEQRF failed to compute a QR factorization.");
1983 for(
int i=0;i<keff_new;i++) {
for(
int j=i;j<keff_new;j++) Rtmp(i,j) = HPtmp(i,j); }
1994 "LAPACK's _UNGQR failed to construct the Q factor.");
2004 for (
int i=0; i < keff; i++) { index[i] = i; }
2006 index.resize(keff_new);
2007 for (
int i=0; i < keff_new; i++) { index[i] = i; }
2008 C1tmp = MVT::CloneViewNonConst( *C1_, index );
2010 MVT::MvTimesMatAddMv( one, *Ctmp, PPtmp, zero, *C1tmp );
2014 index.resize( p+1 );
2015 for (
int i=0; i < p+1; ++i) { index[i] = i; }
2018 MVT::MvTimesMatAddMv( one, *Vtmp, PPtmp, one, *C1tmp );
2033 work_.resize(lwork);
2038 index.resize(keff_new);
2039 for (
int i=0; i < keff_new; i++) { index[i] = i; }
2041 MVT::MvTimesMatAddMv( one, *U1tmp, Rtmp, zero, *Utmp );
2045 if (keff != keff_new) {
2047 gcrodr_iter->setSize( keff, numBlocks_ );
2057 template<
class ScalarType,
class MV,
class OP>
2062 bool xtraVec =
false;
2066 std::vector<MagnitudeType> wr(m), wi(m);
2072 std::vector<MagnitudeType> w(m);
2075 std::vector<int> iperm(m);
2081 builtRecycleSpace_ =
true;
2091 ScalarType d = HH(m, m-1) * HH(m, m-1);
2093 for( i=0; i<m; ++i )
2094 harmHH(i, m-1) += d * e_m[i];
2103 std::vector<ScalarType> work(1);
2104 std::vector<MagnitudeType> rwork(2*m);
2107 lapack.GEEV(
'N',
'V', m, harmHH.
values(), harmHH.
stride(), &wr[0], &wi[0],
2108 vl, ldvl, vr.
values(), vr.
stride(), &work[0], lwork, &rwork[0], &info);
2111 work.resize( lwork );
2113 lapack.GEEV(
'N',
'V', m, harmHH.
values(), harmHH.
stride(), &wr[0], &wi[0],
2114 vl, ldvl, vr.
values(), vr.
stride(), &work[0], lwork, &rwork[0], &info);
2118 for( i=0; i<m; ++i )
2122 this->sort(w, m, iperm);
2127 for( i=0; i<recycledBlocks_; ++i ) {
2128 for( j=0; j<m; j++ ) {
2129 PP(j,i) = vr(j,iperm[i]);
2133 if(!scalarTypeIsComplex) {
2136 if (wi[iperm[recycledBlocks_-1]] != 0.0) {
2138 for ( i=0; i<recycledBlocks_; ++i ) {
2139 if (wi[iperm[i]] != 0.0)
2148 if (wi[iperm[recycledBlocks_-1]] > 0.0) {
2149 for( j=0; j<m; ++j ) {
2150 PP(j,recycledBlocks_) = vr(j,iperm[recycledBlocks_-1]+1);
2154 for( j=0; j<m; ++j ) {
2155 PP(j,recycledBlocks_) = vr(j,iperm[recycledBlocks_-1]-1);
2164 return recycledBlocks_+1;
2167 return recycledBlocks_;
2173 template<
class ScalarType,
class MV,
class OP>
2180 bool xtraVec =
false;
2183 std::vector<int> index;
2186 std::vector<MagnitudeType> wr(m2), wi(m2);
2189 std::vector<MagnitudeType> w(m2);
2195 std::vector<int> iperm(m2);
2198 builtRecycleSpace_ =
true;
2211 index.resize(keffloc);
2212 for (i=0; i<keffloc; ++i) { index[i] = i; }
2216 MVT::MvTransMv( one, *Ctmp, *Utmp, A11 );
2221 for (i=0; i < m+1; i++) { index[i] = i; }
2223 MVT::MvTransMv( one, *Vp, *Utmp, A21 );
2226 for( i=keffloc; i<keffloc+m; i++ ) {
2240 char balanc=
'P', jobvl=
'N', jobvr=
'V', sense=
'N';
2241 int ld =
A.numRows();
2243 int ldvl = ld, ldvr = ld;
2244 int info = 0,ilo = 0,ihi = 0;
2247 std::vector<ScalarType> beta(ld);
2248 std::vector<ScalarType> work(lwork);
2249 std::vector<MagnitudeType> rwork(lwork);
2250 std::vector<MagnitudeType> lscale(ld), rscale(ld);
2251 std::vector<MagnitudeType> rconde(ld), rcondv(ld);
2252 std::vector<int> iwork(ld+6);
2257 lapack.GGEVX(balanc, jobvl, jobvr, sense, ld,
A.values(), ld, B.
values(), ld, &wr[0], &wi[0],
2258 &beta[0], vl, ldvl, vr.
values(), ldvr, &ilo, &ihi, &lscale[0], &rscale[0],
2259 &abnrm, &bbnrm, &rconde[0], &rcondv[0], &work[0], lwork, &rwork[0],
2260 &iwork[0], bwork, &info);
2265 for( i=0; i<ld; i++ ) {
2271 this->sort(w,ld,iperm);
2276 for( i=0; i<recycledBlocks_; i++ ) {
2277 for( j=0; j<ld; j++ ) {
2278 PP(j,i) = vr(j,iperm[ld-recycledBlocks_+i]);
2282 if(!scalarTypeIsComplex) {
2285 if (wi[iperm[ld-recycledBlocks_]] != 0.0) {
2287 for ( i=ld-recycledBlocks_; i<ld; i++ ) {
2288 if (wi[iperm[i]] != 0.0)
2297 if (wi[iperm[ld-recycledBlocks_]] > 0.0) {
2298 for( j=0; j<ld; j++ ) {
2299 PP(j,recycledBlocks_) = vr(j,iperm[ld-recycledBlocks_]+1);
2303 for( j=0; j<ld; j++ ) {
2304 PP(j,recycledBlocks_) = vr(j,iperm[ld-recycledBlocks_]-1);
2313 return recycledBlocks_+1;
2316 return recycledBlocks_;
2323 template<
class ScalarType,
class MV,
class OP>
2325 int l, r, j, i, flag;
2354 if (dlist[j] > dlist[j - 1]) j = j + 1;
2356 if (dlist[j - 1] > dK) {
2357 dlist[i - 1] = dlist[j - 1];
2358 iperm[i - 1] = iperm[j - 1];
2372 dlist[r] = dlist[0];
2373 iperm[r] = iperm[0];
2388 template<
class ScalarType,
class MV,
class OP>
2390 std::ostringstream out;
2393 out <<
"Ortho Type: \"" << orthoType_ <<
"\"";
2394 out <<
", Num Blocks: " <<numBlocks_;
2395 out <<
", Num Recycle Blocks: " << recycledBlocks_;
2396 out <<
", Max Restarts: " << maxRestarts_;
ScaleType convertStringToScaleType(const std::string &scaleType)
Convert the given string to its ScaleType enum value.
ScalarType * values() const
Teuchos::RCP< StatusTestMaxIters< ScalarType, MV, OP > > maxIterTest_
Collection of types and exceptions used within the Belos solvers.
Teuchos::ScalarTraits< ScalarType >::magnitudeType MagnitudeType
Teuchos::ScalarTraits< MagnitudeType > MT
Partial specialization for ScalarType types for which Teuchos::LAPACK has a valid implementation...
Belos's basic output manager for sending information of select verbosity levels to the appropriate ou...
Teuchos::RCP< OutputManager< ScalarType > > printer_
GCRODRSolMgrLAPACKFailure(const std::string &what_arg)
Class which manages the output and verbosity of the Belos solvers.
bool is_null(const boost::shared_ptr< T > &p)
static const bool requiresLapack
ScaleType
The type of scaling to use on the residual norm value.
Exception thrown to signal error in a status test during Belos::StatusTest::checkStatus().
T & get(ParameterList &l, const std::string &name)
GCRODRSolMgrLAPACKFailure is thrown when a nonzero value is retuned from an LAPACK call...
Teuchos::RCP< std::ostream > outputStream_
bool is_null(const std::shared_ptr< T > &p)
MultiVecTraits< ScalarType, MV > MVT
int multiply(ETransp transa, ETransp transb, ScalarType alpha, const SerialDenseMatrix< OrdinalType, ScalarType > &A, const SerialDenseMatrix< OrdinalType, ScalarType > &B, ScalarType beta)
Teuchos::ScalarTraits< ScalarType > SCT
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
A factory class for generating StatusTestOutput objects.
Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > H_
Teuchos::RCP< const Teuchos::ParameterList > getDefaultParameters(const std::string &name) const
Default parameters for the given MatOrthoManager subclass.
An implementation of StatusTestResNorm using a family of residual norms.
GCRODRSolMgrLinearProblemFailure is thrown when the linear problem is not setup (i.e.
GCRODRSolMgrRecyclingFailure(const std::string &what_arg)
This class implements the GCRODR iteration, where a single-std::vector Krylov subspace is constructed...
std::vector< ScalarType > tau_
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
static const double convTol
Default convergence tolerance.
Belos::StatusTest class for specifying a maximum number of iterations.
Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > PP_
int curDim
The current dimension of the reduction.
static std::string name()
Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > HP_
A factory class for generating StatusTestOutput objects.
std::string validNamesString() const
List (as a string) of recognized MatOrthoManager names.
Traits class which defines basic operations on multivectors.
Belos::StatusTest for logically combining several status tests.
bool isParameter(const std::string &name) const
A Belos::StatusTest class for specifying a maximum number of iterations.
GCRODRSolMgrOrthoFailure(const std::string &what_arg)
Teuchos::RCP< SolverManager< ScalarType, MV, OP > > clone() const override
clone for Inverted Injection (DII)
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. ...
Teuchos::RCP< MV > V
The current Krylov basis.
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 is_null(const RCP< T > &p)
bool isValidName(const std::string &name) const
Whether this factory recognizes the MatOrthoManager with the given name.
Teuchos::RCP< StatusTest< ScalarType, MV, OP > > convTest_
Teuchos::RCP< MV > U
The recycled subspace and its projection.
std::ostream & printValidNames(std::ostream &out) const
Print all recognized MatOrthoManager names to the given ostream.
std::vector< ScalarType > work_
int putScalar(const ScalarType value=Teuchos::ScalarTraits< ScalarType >::zero())
Teuchos::RCP< StatusTestOutput< ScalarType, MV, OP > > outputTest_
A linear system to solve, and its associated information.
GCRODRSolMgrOrthoFailure is thrown when the orthogonalization manager is unable to generate orthonorm...
Class which describes the linear problem to be solved by the iterative solver.
Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > H2_
GCRODRIterOrthoFailure is thrown when the GCRODRIter object is unable to compute independent directio...
Implementation of the GCRODR (Recycling GMRES) iterative linear solver.
Type traits class that says whether Teuchos::LAPACK has a valid implementation for the given ScalarTy...
ReturnType
Whether the Belos solve converged for all linear systems.
Teuchos::RCP< Teuchos::Time > timerSolve_
Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > B_
Belos concrete class for performing the GCRO-DR iteration.
Structure to contain pointers to GCRODRIter state variables.
MagnitudeType achievedTol() const override
Tolerance achieved by the last solve() invocation.
Teuchos::Array< Teuchos::RCP< Teuchos::Time > > getTimers() const
Return the timers 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...
void setProblem(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem) override
Set the linear problem that needs to be solved.
GCRODRSolMgrRecyclingFailure is thrown when any problem occurs in using/creating the recycling subspa...
Teuchos::RCP< StatusTestGenResNorm< ScalarType, MV, OP > > impConvTest_
Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > problem_
static magnitudeType magnitude(T a)
virtual ~GCRODRSolMgr()
Destructor.
OrdinalType numCols() const
GCRODRSolMgrLinearProblemFailure(const std::string &what_arg)
GCRODRSolMgr(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem, const Teuchos::RCP< Teuchos::ParameterList > &pl)
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.
Belos::StatusTestResNorm for specifying general residual norm stopping criteria.
Teuchos::RCP< const Teuchos::ParameterList > getCurrentParameters() const override
Get a parameter list containing the current parameters for this object.
Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > B
The projection of the Krylov subspace against the recycled subspace.
OrthoManagerFactory< ScalarType, MV, OP > ortho_factory_type
bool isType(const std::string &name) const
An implementation of the Belos::MatOrthoManager that performs orthogonalization using (potentially) m...
const LinearProblem< ScalarType, MV, OP > & getProblem() const override
Get current linear problem being solved for in this object.
A class for extending the status testing capabilities of Belos via logical combinations.
int getNumIters() const override
Get the iteration count for the most recent call to solve().
Details::SolverManagerRequiresLapack< ScalarType, MV, OP, requiresLapack > base_type
Class which defines basic traits for the operator type.
OperatorTraits< ScalarType, MV, OP > OPT
Teuchos::LAPACK< int, ScalarType > lapack
Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > H
The current Hessenberg matrix.
Parent class to all Belos exceptions.
MagnitudeType orthoKappa_
Teuchos::RCP< StatusTest< ScalarType, MV, OP > > sTest_
Base class for Belos::SolverManager subclasses which normally can only compile with ScalarType types ...
Teuchos::RCP< MatOrthoManager< ScalarType, MV, OP > > ortho_
Orthogonalization manager.
Belos header file which uses auto-configuration information to include necessary C++ headers...
bool isLOADetected() const override
Return whether a loss of accuracy was detected by this solver during the most current solve...
Teuchos::RCP< Teuchos::ParameterList > params_
Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > R_
OrdinalType stride() const
OrdinalType numRows() const
Belos concrete class for performing the block, flexible GMRES iteration.
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.