42#ifndef BELOS_GCRODR_SOLMGR_HPP
43#define BELOS_GCRODR_SOLMGR_HPP
62#include "Teuchos_BLAS.hpp"
63#include "Teuchos_LAPACK.hpp"
64#include "Teuchos_as.hpp"
65#ifdef BELOS_TEUCHOS_TIME_MONITOR
66# include "Teuchos_TimeMonitor.hpp"
68#if defined(HAVE_TEUCHOSCORE_CXX11)
69# include <type_traits>
154 template<
class ScalarType,
class MV,
class OP,
155 const bool lapackSupportsScalarType =
170 const Teuchos::RCP<Teuchos::ParameterList>& pl) :
179 template<
class ScalarType,
class MV,
class OP>
184#if defined(HAVE_TEUCHOSCORE_CXX11)
185# if defined(HAVE_TEUCHOS_COMPLEX)
186 #if defined(HAVE_TEUCHOS_LONG_DOUBLE)
187 static_assert (std::is_same<ScalarType, std::complex<float> >::value ||
188 std::is_same<ScalarType, std::complex<double> >::value ||
189 std::is_same<ScalarType, float>::value ||
190 std::is_same<ScalarType, double>::value ||
191 std::is_same<ScalarType, long double>::value,
192 "Belos::GCRODRSolMgr: ScalarType must be one of the four "
193 "types (S,D,C,Z) supported by LAPACK or long double (largely not impl'd).");
195 static_assert (std::is_same<ScalarType, std::complex<float> >::value ||
196 std::is_same<ScalarType, std::complex<double> >::value ||
197 std::is_same<ScalarType, float>::value ||
198 std::is_same<ScalarType, double>::value,
199 "Belos::GCRODRSolMgr: ScalarType must be one of the four "
200 "types (S,D,C,Z) supported by LAPACK.");
203 #if defined(HAVE_TEUCHOS_LONG_DOUBLE)
204 static_assert (std::is_same<ScalarType, float>::value ||
205 std::is_same<ScalarType, double>::value ||
206 std::is_same<ScalarType, long double>::value,
207 "Belos::GCRODRSolMgr: ScalarType must be float, double or long double. "
208 "Complex arithmetic support is currently disabled. To "
209 "enable it, set Teuchos_ENABLE_COMPLEX=ON.");
211 static_assert (std::is_same<ScalarType, float>::value ||
212 std::is_same<ScalarType, double>::value,
213 "Belos::GCRODRSolMgr: ScalarType must be float or double. "
214 "Complex arithmetic support is currently disabled. To "
215 "enable it, set Teuchos_ENABLE_COMPLEX=ON.");
223 typedef Teuchos::ScalarTraits<ScalarType>
SCT;
224 typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType
MagnitudeType;
225 typedef Teuchos::ScalarTraits<MagnitudeType>
MT;
292 const Teuchos::RCP<Teuchos::ParameterList> &pl);
298 Teuchos::RCP<SolverManager<ScalarType, MV, OP> >
clone ()
const override {
314 Teuchos::RCP<const Teuchos::ParameterList> getValidParameters()
const override;
327 Teuchos::Array<Teuchos::RCP<Teuchos::Time> >
getTimers()
const {
360 void setParameters(
const Teuchos::RCP<Teuchos::ParameterList> ¶ms )
override;
374 throw "Could not set problem.";
418 std::string description()
const override;
428 void initializeStateStorage();
437 int getHarmonicVecs1(
int m,
438 const Teuchos::SerialDenseMatrix<int,ScalarType>& HH,
439 Teuchos::SerialDenseMatrix<int,ScalarType>& PP);
446 int getHarmonicVecs2(
int keff,
int m,
447 const Teuchos::SerialDenseMatrix<int,ScalarType>& HH,
448 const Teuchos::RCP<const MV>& VV,
449 Teuchos::SerialDenseMatrix<int,ScalarType>& PP);
452 void sort(std::vector<MagnitudeType>& dlist,
int n, std::vector<int>& iperm);
458 Teuchos::RCP<LinearProblem<ScalarType,MV,OP> >
problem_;
465 Teuchos::RCP<StatusTest<ScalarType,MV,OP> >
sTest_;
474 Teuchos::RCP<MatOrthoManager<ScalarType,MV,OP> >
ortho_;
523 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> >
H2_;
524 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> >
H_;
525 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> >
B_;
526 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> >
PP_;
527 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> >
HP_;
530 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> >
R_;
547template<
class ScalarType,
class MV,
class OP>
557template<
class ScalarType,
class MV,
class OP>
560 const Teuchos::RCP<Teuchos::ParameterList>& pl):
568 TEUCHOS_TEST_FOR_EXCEPTION(
569 problem == Teuchos::null, std::invalid_argument,
570 "Belos::GCRODRSolMgr constructor: The solver manager's "
571 "constructor needs the linear problem argument 'problem' "
580 if (! pl.is_null ()) {
586template<
class ScalarType,
class MV,
class OP>
619template<
class ScalarType,
class MV,
class OP>
622setParameters (
const Teuchos::RCP<Teuchos::ParameterList> ¶ms)
624 using Teuchos::isParameterType;
625 using Teuchos::getParameter;
627 using Teuchos::ParameterList;
628 using Teuchos::parameterList;
631 using Teuchos::rcp_dynamic_cast;
632 using Teuchos::rcpFromRef;
633 using Teuchos::Exceptions::InvalidParameter;
634 using Teuchos::Exceptions::InvalidParameterName;
635 using Teuchos::Exceptions::InvalidParameterType;
658 params_ = parameterList (*defaultParams);
672 params_ = parameterList (*params);
707 params_->validateParametersAndSetDefaults (*defaultParams);
711 if (params->isParameter (
"Maximum Restarts")) {
719 if (params->isParameter (
"Maximum Iterations")) {
729 if (params->isParameter (
"Num Blocks")) {
731 TEUCHOS_TEST_FOR_EXCEPTION(
numBlocks_ <= 0, std::invalid_argument,
732 "Belos::GCRODRSolMgr: The \"Num Blocks\" parameter must "
733 "be strictly positive, but you specified a value of "
740 if (params->isParameter (
"Num Recycled Blocks")) {
743 TEUCHOS_TEST_FOR_EXCEPTION(
recycledBlocks_ <= 0, std::invalid_argument,
744 "Belos::GCRODRSolMgr: The \"Num Recycled Blocks\" "
745 "parameter must be strictly positive, but you specified "
748 "Belos::GCRODRSolMgr: The \"Num Recycled Blocks\" "
749 "parameter must be less than the \"Num Blocks\" "
750 "parameter, but you specified \"Num Recycled Blocks\" "
760 if (params->isParameter (
"Timer Label")) {
761 std::string tempLabel = params->get (
"Timer Label",
label_default_);
764 if (tempLabel !=
label_) {
767 std::string solveLabel =
label_ +
": GCRODRSolMgr total solve time";
768#ifdef BELOS_TEUCHOS_TIME_MONITOR
769 timerSolve_ = Teuchos::TimeMonitor::getNewCounter (solveLabel);
771 if (
ortho_ != Teuchos::null) {
778 if (params->isParameter (
"Verbosity")) {
779 if (isParameterType<int> (*params,
"Verbosity")) {
782 verbosity_ = (int) getParameter<Belos::MsgType> (*params,
"Verbosity");
793 if (params->isParameter (
"Output Style")) {
794 if (isParameterType<int> (*params,
"Output Style")) {
797 outputStyle_ = (int) getParameter<OutputType> (*params,
"Output Style");
819 if (params->isParameter (
"Output Stream")) {
821 outputStream_ = getParameter<RCP<std::ostream> > (*params,
"Output Stream");
822 }
catch (InvalidParameter&) {
844 if (params->isParameter (
"Output Frequency")) {
869 bool changedOrthoType =
false;
870 if (params->isParameter (
"Orthogonalization")) {
871 const std::string& tempOrthoType =
875 std::ostringstream os;
876 os <<
"Belos::GCRODRSolMgr: Invalid orthogonalization name \""
877 << tempOrthoType <<
"\". The following are valid options "
878 <<
"for the \"Orthogonalization\" name parameter: ";
880 throw std::invalid_argument (os.str());
883 changedOrthoType =
true;
902 RCP<ParameterList> orthoParams;
905 using Teuchos::sublist;
907 const std::string paramName (
"Orthogonalization Parameters");
910 orthoParams = sublist (
params_, paramName,
true);
911 }
catch (InvalidParameter&) {
918 orthoParams = sublist (
params_, paramName,
true);
921 TEUCHOS_TEST_FOR_EXCEPTION(orthoParams.is_null(), std::logic_error,
922 "Failed to get orthogonalization parameters. "
923 "Please report this bug to the Belos developers.");
928 if (
ortho_.is_null() || changedOrthoType) {
942 typedef Teuchos::ParameterListAcceptor PLA;
943 RCP<PLA> pla = rcp_dynamic_cast<PLA> (
ortho_);
951 pla->setParameterList (orthoParams);
963 if (params->isParameter (
"Orthogonalization Constant")) {
965 if (params->isType<
MagnitudeType> (
"Orthogonalization Constant")) {
966 orthoKappa = params->get (
"Orthogonalization Constant", orthoKappa);
972 if (orthoKappa > 0) {
994 if (params->isParameter(
"Convergence Tolerance")) {
995 if (params->isType<
MagnitudeType> (
"Convergence Tolerance")) {
996 convTol_ = params->get (
"Convergence Tolerance",
1012 if (params->isParameter (
"Implicit Residual Scaling")) {
1013 std::string tempImpResScale =
1014 getParameter<std::string> (*params,
"Implicit Residual Scaling");
1045 if (params->isParameter(
"Explicit Residual Scaling")) {
1046 std::string tempExpResScale =
1047 getParameter<std::string> (*params,
"Explicit Residual Scaling");
1096 convTest_ = rcp (
new StatusTestCombo_t (StatusTestCombo_t::SEQ,
1104 sTest_ = rcp (
new StatusTestCombo_t (StatusTestCombo_t::OR,
1114 std::string solverDesc =
" GCRODR ";
1119 std::string solveLabel =
label_ +
": GCRODRSolMgr total solve time";
1120#ifdef BELOS_TEUCHOS_TIME_MONITOR
1121 timerSolve_ = Teuchos::TimeMonitor::getNewCounter(solveLabel);
1130template<
class ScalarType,
class MV,
class OP>
1131Teuchos::RCP<const Teuchos::ParameterList>
1134 using Teuchos::ParameterList;
1135 using Teuchos::parameterList;
1138 static RCP<const ParameterList> validPL;
1139 if (is_null(validPL)) {
1140 RCP<ParameterList> pl = parameterList ();
1144 "The relative residual tolerance that needs to be achieved by the\n"
1145 "iterative solver in order for the linear system to be declared converged.");
1147 "The maximum number of cycles allowed for each\n"
1148 "set of RHS solved.");
1150 "The maximum number of iterations allowed for each\n"
1151 "set of RHS solved.");
1156 "Block Size Parameter -- currently must be 1 for GCRODR");
1158 "The maximum number of vectors allowed in the Krylov subspace\n"
1159 "for each set of RHS solved.");
1161 "The maximum number of vectors in the recycled subspace." );
1163 "What type(s) of solver information should be outputted\n"
1164 "to the output stream.");
1166 "What style is used for the solver information outputted\n"
1167 "to the output stream.");
1169 "How often convergence information should be outputted\n"
1170 "to the output stream.");
1171 pl->set(
"Output Stream", Teuchos::rcpFromRef(std::cout),
1172 "A reference-counted pointer to the output stream where all\n"
1173 "solver output is sent.");
1175 "The type of scaling used in the implicit residual convergence test.");
1177 "The type of scaling used in the explicit residual convergence test.");
1178 pl->set(
"Timer Label",
static_cast<const char *
>(
label_default_),
1179 "The string to use as a prefix for the timer labels.");
1183 "The type of orthogonalization to use. Valid options: " +
1185 RCP<const ParameterList> orthoParams =
1187 pl->set (
"Orthogonalization Parameters", *orthoParams,
1188 "Parameters specific to the type of orthogonalization used.");
1191 "When using DGKS orthogonalization: the \"depTol\" constant, used "
1192 "to determine whether another step of classical Gram-Schmidt is "
1193 "necessary. Otherwise ignored.");
1200template<
class ScalarType,
class MV,
class OP>
1203 ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
1206 Teuchos::RCP<const MV> rhsMV =
problem_->getRHS();
1207 if (rhsMV == Teuchos::null) {
1215 "Belos::GCRODRSolMgr::initializeStateStorage(): Cannot generate a Krylov basis with dimension larger the operator!");
1218 if (
U_ == Teuchos::null) {
1224 Teuchos::RCP<const MV> tmp =
U_;
1230 if (
C_ == Teuchos::null) {
1236 Teuchos::RCP<const MV> tmp =
C_;
1242 if (
V_ == Teuchos::null) {
1248 Teuchos::RCP<const MV> tmp =
V_;
1254 if (
U1_ == Teuchos::null) {
1260 Teuchos::RCP<const MV> tmp =
U1_;
1266 if (
C1_ == Teuchos::null) {
1272 Teuchos::RCP<const MV> tmp =
C1_;
1278 if (
r_ == Teuchos::null)
1291 if (
H2_ == Teuchos::null)
1297 H2_->putScalar(zero);
1300 if (
R_ == Teuchos::null)
1306 R_->putScalar(zero);
1309 if (
PP_ == Teuchos::null)
1317 if (
HP_ == Teuchos::null)
1329template<
class ScalarType,
class MV,
class OP>
1339 ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
1340 ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
1349 std::vector<int> currIdx(1);
1357 if (
static_cast<ptrdiff_t
>(
numBlocks_) > dim) {
1360 "Warning! Requested Krylov subspace dimension is larger than operator dimension!" << std::endl <<
1361 " The maximum number of blocks allowed for the Krylov subspace will be adjusted to " <<
numBlocks_ << std::endl;
1366 bool isConverged =
true;
1373 Teuchos::ParameterList plist;
1381 RCP<GCRODRIter<ScalarType,MV,OP> > gcrodr_iter;
1384 int prime_iterations = 0;
1388#ifdef BELOS_TEUCHOS_TIME_MONITOR
1392 while ( numRHS2Solve > 0 ) {
1406 "Belos::GCRODRSolMgr::solve(): Requested size of recycled subspace is not consistent with the current recycle subspace.");
1408 printer_->stream(
Debug) <<
" Now solving RHS index " << currIdx[0] <<
" using recycled subspace of dimension " <<
keff << std::endl << std::endl;
1411 for (
int ii=0; ii<
keff; ++ii) { index[ii] = ii; }
1420 Teuchos::SerialDenseMatrix<int,ScalarType> Rtmp( Teuchos::View, *
R_,
keff,
keff );
1421 int rank =
ortho_->normalize(*Ctmp, rcp(&Rtmp,
false));
1423 TEUCHOS_TEST_FOR_EXCEPTION(rank !=
keff,
GCRODRSolMgrOrthoFailure,
"Belos::GCRODRSolMgr::solve(): Failed to compute orthonormal basis for initial recycled subspace.");
1428 ipiv_.resize(Rtmp.numRows());
1429 lapack.GETRF(Rtmp.numRows(),Rtmp.numCols(),Rtmp.values(),Rtmp.stride(),&
ipiv_[0],&info);
1430 TEUCHOS_TEST_FOR_EXCEPTION(info != 0,
GCRODRSolMgrLAPACKFailure,
"Belos::GCRODRSolMgr::solve(): LAPACK _GETRF failed to compute an LU factorization.");
1433 int lwork = Rtmp.numRows();
1434 work_.resize(lwork);
1435 lapack.GETRI(Rtmp.numRows(),Rtmp.values(),Rtmp.stride(),&
ipiv_[0],&
work_[0],lwork,&info);
1436 TEUCHOS_TEST_FOR_EXCEPTION(info != 0,
GCRODRSolMgrLAPACKFailure,
"Belos::GCRODRSolMgr::solve(): LAPACK _GETRI failed to invert triangular matrix.");
1444 for (
int ii=0; ii<
keff; ++ii) { index[ii] = ii; }
1449 Teuchos::SerialDenseMatrix<int,ScalarType> Ctr(
keff,1);
1457 problem_->updateSolution( update,
true );
1463 prime_iterations = 0;
1469 printer_->stream(
Debug) <<
" No recycled subspace available for RHS index " << currIdx[0] << std::endl << std::endl;
1471 Teuchos::ParameterList primeList;
1475 primeList.set(
"Recycled Blocks",0);
1478 RCP<GCRODRIter<ScalarType,MV,OP> > gcrodr_prime_iter;
1483 index.resize( 1 ); index[0] = 0;
1490 for (
int ii=0; ii<(
numBlocks_+1); ++ii) { index[ii] = ii; }
1492 newstate.
U = Teuchos::null;
1493 newstate.
C = Teuchos::null;
1495 newstate.
B = Teuchos::null;
1497 gcrodr_prime_iter->initialize(newstate);
1500 bool primeConverged =
false;
1502 gcrodr_prime_iter->iterate();
1507 primeConverged =
true;
1512 gcrodr_prime_iter->updateLSQR( gcrodr_prime_iter->getCurSubspaceDim() );
1515 sTest_->checkStatus( &*gcrodr_prime_iter );
1517 primeConverged =
true;
1519 catch (
const std::exception &e) {
1520 printer_->stream(
Errors) <<
"Error! Caught exception in GCRODRIter::iterate() at iteration "
1521 << gcrodr_prime_iter->getNumIters() << std::endl
1522 << e.what() << std::endl;
1526 prime_iterations = gcrodr_prime_iter->getNumIters();
1529 RCP<MV> update = gcrodr_prime_iter->getCurrentUpdate();
1530 problem_->updateSolution( update,
true );
1533 newstate = gcrodr_prime_iter->getState();
1543 RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > PPtmp = rcp (
new Teuchos::SerialDenseMatrix<int,ScalarType> ( Teuchos::View, *
PP_, p,
recycledBlocks_+1 ) );
1547 PPtmp = rcp (
new Teuchos::SerialDenseMatrix<int,ScalarType> ( Teuchos::View, *
PP_, p,
keff ) );
1550 for (
int ii=0; ii<
keff; ++ii) { index[ii] = ii; }
1555 for (
int ii=0; ii < p; ++ii) { index[ii] = ii; }
1568 Teuchos::SerialDenseMatrix<int,ScalarType> HPtmp( Teuchos::View, *
HP_, p+1,
keff );
1569 HPtmp.multiply( Teuchos::NO_TRANS, Teuchos::NO_TRANS, one, Htmp, *PPtmp, zero );
1576 lapack.GEQRF (HPtmp.numRows (), HPtmp.numCols (), HPtmp.values (),
1577 HPtmp.stride (), &
tau_[0], &
work_[0], lwork, &info);
1578 TEUCHOS_TEST_FOR_EXCEPTION(
1580 " LAPACK's _GEQRF failed to compute a workspace size.");
1588 lwork = std::abs (
static_cast<int> (Teuchos::ScalarTraits<ScalarType>::real (
work_[0])));
1589 work_.resize (lwork);
1590 lapack.GEQRF (HPtmp.numRows (), HPtmp.numCols (), HPtmp.values (),
1591 HPtmp.stride (), &
tau_[0], &
work_[0], lwork, &info);
1592 TEUCHOS_TEST_FOR_EXCEPTION(
1594 " LAPACK's _GEQRF failed to compute a QR factorization.");
1598 Teuchos::SerialDenseMatrix<int,ScalarType> Rtmp( Teuchos::View, *
R_,
keff,
keff );
1599 for (
int ii = 0; ii <
keff; ++ii) {
1600 for (
int jj = ii; jj <
keff; ++jj) {
1601 Rtmp(ii,jj) = HPtmp(ii,jj);
1608 lapack.UNGQR (HPtmp.numRows (), HPtmp.numCols (), HPtmp.numCols (),
1609 HPtmp.values (), HPtmp.stride (), &
tau_[0], &
work_[0],
1611 TEUCHOS_TEST_FOR_EXCEPTION(
1613 "LAPACK's _UNGQR failed to construct the Q factor.");
1618 index.resize (p + 1);
1619 for (
int ii = 0; ii < (p+1); ++ii) {
1630 ipiv_.resize(Rtmp.numRows());
1631 lapack.GETRF(Rtmp.numRows(),Rtmp.numCols(),Rtmp.values(),Rtmp.stride(),&
ipiv_[0],&info);
1632 TEUCHOS_TEST_FOR_EXCEPTION(
1634 "LAPACK's _GETRF failed to compute an LU factorization.");
1643 lwork = Rtmp.numRows();
1644 work_.resize(lwork);
1645 lapack.GETRI(Rtmp.numRows(),Rtmp.values(),Rtmp.stride(),&
ipiv_[0],&
work_[0],lwork,&info);
1646 TEUCHOS_TEST_FOR_EXCEPTION(
1648 "LAPACK's _GETRI failed to invert triangular matrix.");
1654 <<
" Generated recycled subspace using RHS index " << currIdx[0]
1655 <<
" of dimension " <<
keff << std::endl << std::endl;
1660 if (primeConverged) {
1666 if (numRHS2Solve > 0) {
1671 currIdx.resize (numRHS2Solve);
1684 gcrodr_iter->resetNumIters(prime_iterations);
1691 index.resize( 1 ); index[0] = 0;
1698 for (
int ii=0; ii<(
numBlocks_+1); ++ii) { index[ii] = ii; }
1700 index.resize(
keff );
1701 for (
int ii=0; ii<
keff; ++ii) { index[ii] = ii; }
1704 newstate.
B = rcp(
new Teuchos::SerialDenseMatrix<int,ScalarType>( Teuchos::View, *
H2_,
keff,
numBlocks_, 0,
keff ) );
1707 gcrodr_iter->initialize(newstate);
1710 int numRestarts = 0;
1715 gcrodr_iter->iterate();
1733 isConverged =
false;
1741 else if ( gcrodr_iter->getCurSubspaceDim() == gcrodr_iter->getMaxSubspaceDim() ) {
1746 RCP<MV> update = gcrodr_iter->getCurrentUpdate();
1747 problem_->updateSolution( update,
true );
1752 <<
" Generated new recycled subspace using RHS index "
1753 << currIdx[0] <<
" of dimension " <<
keff << std::endl
1758 isConverged =
false;
1764 <<
" Performing restart number " << numRestarts <<
" of "
1769 index.resize( 1 ); index[0] = 0;
1776 for (
int ii=0; ii<(
numBlocks_+1); ++ii) { index[ii] = ii; }
1778 index.resize(
keff );
1779 for (
int ii=0; ii<
keff; ++ii) { index[ii] = ii; }
1782 restartState.
B = rcp(
new Teuchos::SerialDenseMatrix<int,ScalarType>( Teuchos::View, *
H2_,
keff,
numBlocks_, 0,
keff ) );
1785 gcrodr_iter->initialize(restartState);
1798 TEUCHOS_TEST_FOR_EXCEPTION(
1799 true, std::logic_error,
"Belos::GCRODRSolMgr::solve: "
1800 "Invalid return from GCRODRIter::iterate().");
1805 gcrodr_iter->updateLSQR( gcrodr_iter->getCurSubspaceDim() );
1808 sTest_->checkStatus( &*gcrodr_iter );
1810 isConverged =
false;
1813 catch (
const std::exception& e) {
1815 <<
"Error! Caught exception in GCRODRIter::iterate() at iteration "
1816 << gcrodr_iter->getNumIters() << std::endl << e.what() << std::endl;
1823 RCP<MV> update = gcrodr_iter->getCurrentUpdate();
1824 problem_->updateSolution( update,
true );
1835 <<
" Generated new recycled subspace using RHS index " << currIdx[0]
1836 <<
" of dimension " <<
keff << std::endl << std::endl;
1841 if (numRHS2Solve > 0) {
1846 currIdx.resize (numRHS2Solve);
1854#ifdef BELOS_TEUCHOS_TIME_MONITOR
1875 const std::vector<MagnitudeType>* pTestValues =
expConvTest_->getTestValue();
1876 if (pTestValues == NULL || pTestValues->size() < 1) {
1879 TEUCHOS_TEST_FOR_EXCEPTION(pTestValues == NULL, std::logic_error,
1880 "Belos::GCRODRSolMgr::solve(): The implicit convergence test's getTestValue() "
1881 "method returned NULL. Please report this bug to the Belos developers.");
1882 TEUCHOS_TEST_FOR_EXCEPTION(pTestValues->size() < 1, std::logic_error,
1883 "Belos::GCRODRSolMgr::solve(): The implicit convergence test's getTestValue() "
1884 "method returned a vector of length zero. Please report this bug to the "
1885 "Belos developers.");
1890 achievedTol_ = *std::max_element (pTestValues->begin(), pTestValues->end());
1897template<
class ScalarType,
class MV,
class OP>
1900 MagnitudeType one = Teuchos::ScalarTraits<MagnitudeType>::one();
1901 ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
1903 std::vector<MagnitudeType> d(
keff);
1904 std::vector<ScalarType> dscalar(
keff);
1917 for (
int ii=0; ii<
keff; ++ii) { index[ii] = ii; }
1920 dscalar.resize(
keff);
1922 for (
int i=0; i<
keff; ++i) {
1924 dscalar[i] = (ScalarType)d[i];
1930 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > H2tmp = Teuchos::rcp(
new Teuchos::SerialDenseMatrix<int,ScalarType>( Teuchos::View, *
H2_, p+
keff+1, p+
keff ) );
1933 for (
int i=0; i<
keff; ++i) {
1934 (*H2tmp)(i,i) = d[i];
1949 Teuchos::RCP<MV> U1tmp;
1951 index.resize(
keff );
1952 for (
int ii=0; ii<
keff; ++ii) { index[ii] = ii; }
1954 index.resize( keff_new );
1955 for (
int ii=0; ii<keff_new; ++ii) { index[ii] = ii; }
1957 Teuchos::SerialDenseMatrix<int,ScalarType> PPtmp( Teuchos::View, *
PP_,
keff, keff_new );
1964 for (
int ii=0; ii < p; ii++) { index[ii] = ii; }
1966 Teuchos::SerialDenseMatrix<int,ScalarType> PPtmp( Teuchos::View, *
PP_, p, keff_new,
keff );
1971 Teuchos::SerialDenseMatrix<int,ScalarType> HPtmp( Teuchos::View, *
HP_, p+
keff+1, keff_new );
1973 Teuchos::SerialDenseMatrix<int,ScalarType> PPtmp( Teuchos::View, *
PP_, p+
keff, keff_new );
1974 HPtmp.multiply(Teuchos::NO_TRANS,Teuchos::NO_TRANS,one,*H2tmp,PPtmp,zero);
1978 int info = 0, lwork = -1;
1979 tau_.resize (keff_new);
1980 lapack.GEQRF (HPtmp.numRows (), HPtmp.numCols (), HPtmp.values (),
1981 HPtmp.stride (), &
tau_[0], &
work_[0], lwork, &info);
1982 TEUCHOS_TEST_FOR_EXCEPTION(
1984 "LAPACK's _GEQRF failed to compute a workspace size.");
1990 lwork = std::abs (
static_cast<int> (Teuchos::ScalarTraits<ScalarType>::real (
work_[0])));
1991 work_.resize (lwork);
1992 lapack.GEQRF (HPtmp.numRows (), HPtmp.numCols (), HPtmp.values (),
1993 HPtmp.stride (), &
tau_[0], &
work_[0], lwork, &info);
1994 TEUCHOS_TEST_FOR_EXCEPTION(
1996 "LAPACK's _GEQRF failed to compute a QR factorization.");
2000 Teuchos::SerialDenseMatrix<int,ScalarType> Rtmp( Teuchos::View, *
R_, keff_new, keff_new );
2001 for(
int i=0;i<keff_new;i++) {
for(
int j=i;j<keff_new;j++) Rtmp(i,j) = HPtmp(i,j); }
2007 lapack.UNGQR (HPtmp.numRows (), HPtmp.numCols (), HPtmp.numCols (),
2008 HPtmp.values (), HPtmp.stride (), &
tau_[0], &
work_[0],
2010 TEUCHOS_TEST_FOR_EXCEPTION(
2012 "LAPACK's _UNGQR failed to construct the Q factor.");
2019 Teuchos::RCP<MV> C1tmp;
2022 for (
int i=0; i <
keff; i++) { index[i] = i; }
2024 index.resize(keff_new);
2025 for (
int i=0; i < keff_new; i++) { index[i] = i; }
2027 Teuchos::SerialDenseMatrix<int,ScalarType> PPtmp( Teuchos::View, *
HP_,
keff, keff_new );
2032 index.resize( p+1 );
2033 for (
int i=0; i < p+1; ++i) { index[i] = i; }
2035 Teuchos::SerialDenseMatrix<int,ScalarType> PPtmp( Teuchos::View, *
HP_, p+1, keff_new,
keff, 0 );
2045 ipiv_.resize(Rtmp.numRows());
2046 lapack.GETRF(Rtmp.numRows(),Rtmp.numCols(),Rtmp.values(),Rtmp.stride(),&
ipiv_[0],&info);
2047 TEUCHOS_TEST_FOR_EXCEPTION(info != 0,
GCRODRSolMgrLAPACKFailure,
"Belos::GCRODRSolMgr::solve(): LAPACK _GETRF failed to compute an LU factorization.");
2050 lwork = Rtmp.numRows();
2051 work_.resize(lwork);
2052 lapack.GETRI(Rtmp.numRows(),Rtmp.values(),Rtmp.stride(),&
ipiv_[0],&
work_[0],lwork,&info);
2053 TEUCHOS_TEST_FOR_EXCEPTION(info != 0,
GCRODRSolMgrLAPACKFailure,
"Belos::GCRODRSolMgr::solve(): LAPACK _GETRI failed to compute an LU factorization.");
2056 index.resize(keff_new);
2057 for (
int i=0; i < keff_new; i++) { index[i] = i; }
2063 if (
keff != keff_new) {
2075template<
class ScalarType,
class MV,
class OP>
2077 const Teuchos::SerialDenseMatrix<int,ScalarType>& HH,
2078 Teuchos::SerialDenseMatrix<int,ScalarType>& PP) {
2080 bool xtraVec =
false;
2081 ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
2084 std::vector<MagnitudeType> wr(m), wi(m);
2087 Teuchos::SerialDenseMatrix<int,ScalarType> vr(m,m,
false);
2090 std::vector<MagnitudeType> w(m);
2093 std::vector<int> iperm(m);
2102 Teuchos::SerialDenseMatrix<int, ScalarType> HHt( HH, Teuchos::TRANS );
2103 Teuchos::SerialDenseVector<int, ScalarType> e_m( m );
2105 lapack.GESV(m, 1, HHt.values(), HHt.stride(), &iperm[0], e_m.values(), e_m.stride(), &info);
2106 TEUCHOS_TEST_FOR_EXCEPTION(info != 0,
GCRODRSolMgrLAPACKFailure,
"Belos::GCRODRSolMgr::solve(): LAPACK GESV failed to compute a solution.");
2109 ScalarType d = HH(m, m-1) * HH(m, m-1);
2110 Teuchos::SerialDenseMatrix<int, ScalarType> harmHH( Teuchos::Copy, HH, m, m );
2111 for( i=0; i<m; ++i )
2112 harmHH(i, m-1) += d * e_m[i];
2121 std::vector<ScalarType> work(1);
2122 std::vector<MagnitudeType> rwork(2*m);
2125 lapack.GEEV(
'N',
'V', m, harmHH.values(), harmHH.stride(), &wr[0], &wi[0],
2126 vl, ldvl, vr.values(), vr.stride(), &work[0], lwork, &rwork[0], &info);
2128 lwork = std::abs (
static_cast<int> (Teuchos::ScalarTraits<ScalarType>::real (work[0])));
2129 work.resize( lwork );
2131 lapack.GEEV(
'N',
'V', m, harmHH.values(), harmHH.stride(), &wr[0], &wi[0],
2132 vl, ldvl, vr.values(), vr.stride(), &work[0], lwork, &rwork[0], &info);
2133 TEUCHOS_TEST_FOR_EXCEPTION(info != 0,
GCRODRSolMgrLAPACKFailure,
"Belos::GCRODRSolMgr::solve(): LAPACK GEEV failed to compute eigensolutions.");
2136 for( i=0; i<m; ++i )
2137 w[i] = Teuchos::ScalarTraits<MagnitudeType>::squareroot( wr[i]*wr[i] + wi[i]*wi[i] );
2140 this->
sort(w, m, iperm);
2142 const bool scalarTypeIsComplex = Teuchos::ScalarTraits<ScalarType>::isComplex;
2146 for( j=0; j<m; j++ ) {
2147 PP(j,i) = vr(j,iperm[i]);
2151 if(!scalarTypeIsComplex) {
2157 if (wi[iperm[i]] != 0.0)
2167 for( j=0; j<m; ++j ) {
2172 for( j=0; j<m; ++j ) {
2191template<
class ScalarType,
class MV,
class OP>
2193 const Teuchos::SerialDenseMatrix<int,ScalarType>& HH,
2194 const Teuchos::RCP<const MV>& VV,
2195 Teuchos::SerialDenseMatrix<int,ScalarType>& PP) {
2197 int m2 = HH.numCols();
2198 bool xtraVec =
false;
2199 ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
2200 ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
2201 std::vector<int> index;
2204 std::vector<MagnitudeType> wr(m2), wi(m2);
2207 std::vector<MagnitudeType> w(m2);
2210 Teuchos::SerialDenseMatrix<int,ScalarType> vr(m2,m2,
false);
2213 std::vector<int> iperm(m2);
2221 Teuchos::SerialDenseMatrix<int,ScalarType> B(m2,m2,
false);
2222 B.multiply(Teuchos::TRANS,Teuchos::NO_TRANS,one,HH,HH,zero);
2226 Teuchos::SerialDenseMatrix<int,ScalarType> A_tmp( keffloc+m+1, keffloc+m );
2229 index.resize(keffloc);
2230 for (i=0; i<keffloc; ++i) { index[i] = i; }
2233 Teuchos::SerialDenseMatrix<int,ScalarType> A11( Teuchos::View, A_tmp, keffloc, keffloc );
2237 Teuchos::SerialDenseMatrix<int,ScalarType> A21( Teuchos::View, A_tmp, m+1, keffloc, keffloc );
2239 for (i=0; i < m+1; i++) { index[i] = i; }
2244 for( i=keffloc; i<keffloc+m; i++ ) {
2249 Teuchos::SerialDenseMatrix<int,ScalarType> A( m2, A_tmp.numCols() );
2250 A.multiply( Teuchos::TRANS, Teuchos::NO_TRANS, one, HH, A_tmp, zero );
2258 char balanc=
'P', jobvl=
'N', jobvr=
'V', sense=
'N';
2259 int ld = A.numRows();
2261 int ldvl = ld, ldvr = ld;
2262 int info = 0,ilo = 0,ihi = 0;
2265 std::vector<ScalarType> beta(ld);
2266 std::vector<ScalarType> work(lwork);
2267 std::vector<MagnitudeType> rwork(lwork);
2268 std::vector<MagnitudeType> lscale(ld), rscale(ld);
2269 std::vector<MagnitudeType> rconde(ld), rcondv(ld);
2270 std::vector<int> iwork(ld+6);
2275 lapack.GGEVX(balanc, jobvl, jobvr, sense, ld, A.values(), ld, B.values(), ld, &wr[0], &wi[0],
2276 &beta[0], vl, ldvl, vr.values(), ldvr, &ilo, &ihi, &lscale[0], &rscale[0],
2277 &abnrm, &bbnrm, &rconde[0], &rcondv[0], &work[0], lwork, &rwork[0],
2278 &iwork[0], bwork, &info);
2279 TEUCHOS_TEST_FOR_EXCEPTION(info != 0,
GCRODRSolMgrLAPACKFailure,
"Belos::GCRODRSolMgr::solve(): LAPACK GGEVX failed to compute eigensolutions.");
2283 for( i=0; i<ld; i++ ) {
2284 w[i] = Teuchos::ScalarTraits<MagnitudeType>::squareroot (wr[i]*wr[i] + wi[i]*wi[i]) /
2285 Teuchos::ScalarTraits<ScalarType>::magnitude (beta[i]);
2289 this->
sort(w,ld,iperm);
2291 const bool scalarTypeIsComplex = Teuchos::ScalarTraits<ScalarType>::isComplex;
2295 for( j=0; j<ld; j++ ) {
2300 if(!scalarTypeIsComplex) {
2306 if (wi[iperm[i]] != 0.0)
2316 for( j=0; j<ld; j++ ) {
2321 for( j=0; j<ld; j++ ) {
2341template<
class ScalarType,
class MV,
class OP>
2343 int l, r, j, i, flag;
2372 if (dlist[j] > dlist[j - 1]) j = j + 1;
2374 if (dlist[j - 1] > dK) {
2375 dlist[i - 1] = dlist[j - 1];
2376 iperm[i - 1] = iperm[j - 1];
2390 dlist[r] = dlist[0];
2391 iperm[r] = iperm[0];
2406template<
class ScalarType,
class MV,
class OP>
2408 std::ostringstream out;
2409 out <<
"Belos::GCRODRSolMgr<...,"<<Teuchos::ScalarTraits<ScalarType>::name()<<
">";
2411 out <<
"Ortho Type: \"" <<
orthoType_ <<
"\"";
Belos concrete class for performing the block, flexible GMRES iteration.
Belos header file which uses auto-configuration information to include necessary C++ headers.
Belos concrete class for performing the GCRO-DR 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 class for specifying a maximum number of iterations.
A factory class for generating StatusTestOutput objects.
Collection of types and exceptions used within the Belos solvers.
BelosError(const std::string &what_arg)
An implementation of the Belos::MatOrthoManager that performs orthogonalization using (potentially) m...
void setDepTol(const MagnitudeType dep_tol)
Set parameter for re-orthogonalization threshhold.
Base class for Belos::SolverManager subclasses which normally can only compile with ScalarType types ...
GCRODRIterOrthoFailure is thrown when the GCRODRIter object is unable to compute independent directio...
This class implements the GCRODR iteration, where a single-stdvector Krylov subspace is constructed....
GCRODRSolMgrLAPACKFailure is thrown when a nonzero value is retuned from an LAPACK call.
GCRODRSolMgrLAPACKFailure(const std::string &what_arg)
GCRODRSolMgrLinearProblemFailure is thrown when the linear problem is not setup (i....
GCRODRSolMgrLinearProblemFailure(const std::string &what_arg)
GCRODRSolMgrOrthoFailure is thrown when the orthogonalization manager is unable to generate orthonorm...
GCRODRSolMgrOrthoFailure(const std::string &what_arg)
GCRODRSolMgrRecyclingFailure is thrown when any problem occurs in using/creating the recycling subspa...
GCRODRSolMgrRecyclingFailure(const std::string &what_arg)
std::vector< ScalarType > work_
Teuchos::RCP< StatusTest< ScalarType, MV, OP > > convTest_
const LinearProblem< ScalarType, MV, OP > & getProblem() const override
Get current linear problem being solved for in 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...
static constexpr int verbosity_default_
virtual ~GCRODRSolMgr()
Destructor.
Teuchos::RCP< SolverManager< ScalarType, MV, OP > > clone() const override
clone for Inverted Injection (DII)
MagnitudeType orthoKappa_
Teuchos::Array< Teuchos::RCP< Teuchos::Time > > getTimers() const
Return the timers for this object.
static constexpr int maxRestarts_default_
std::vector< ScalarType > tau_
Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > H_
Teuchos::RCP< StatusTestGenResNorm< ScalarType, MV, OP > > impConvTest_
MagnitudeType achievedTol_
Teuchos::RCP< OutputManager< ScalarType > > printer_
OperatorTraits< ScalarType, MV, OP > OPT
Teuchos::RCP< StatusTestOutput< ScalarType, MV, OP > > outputTest_
Teuchos::RCP< std::ostream > outputStream_
static constexpr double orthoKappa_default_
Teuchos::LAPACK< int, ScalarType > lapack
OrthoManagerFactory< ScalarType, MV, OP > ortho_factory_type
static constexpr int blockSize_default_
Teuchos::RCP< MatOrthoManager< ScalarType, MV, OP > > ortho_
Orthogonalization manager.
static constexpr const char * orthoType_default_
Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > problem_
Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > PP_
void initializeStateStorage()
static constexpr int outputFreq_default_
int getNumIters() const override
Get the iteration count for the most recent call to solve().
Teuchos::ScalarTraits< ScalarType >::magnitudeType MagnitudeType
static constexpr int recycledBlocks_default_
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
Get a parameter list containing the valid parameters for this object.
Teuchos::ScalarTraits< MagnitudeType > MT
void setProblem(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem) override
Set the linear problem that needs to be solved.
Teuchos::RCP< Teuchos::Time > timerSolve_
void sort(std::vector< MagnitudeType > &dlist, int n, std::vector< int > &iperm)
MultiVecTraits< ScalarType, MV > MVT
static constexpr const char * label_default_
Teuchos::RCP< StatusTestGenResNorm< ScalarType, MV, OP > > expConvTest_
Teuchos::RCP< Teuchos::ParameterList > params_
Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > H2_
int getHarmonicVecs1(int m, const Teuchos::SerialDenseMatrix< int, ScalarType > &HH, Teuchos::SerialDenseMatrix< int, ScalarType > &PP)
GCRODRSolMgr()
Empty constructor for GCRODRSolMgr. This constructor takes no arguments and sets the default values f...
static constexpr const char * expResScale_default_
int getHarmonicVecs2(int keff, int m, const Teuchos::SerialDenseMatrix< int, ScalarType > &HH, const Teuchos::RCP< const MV > &VV, Teuchos::SerialDenseMatrix< int, ScalarType > &PP)
Teuchos::RCP< const Teuchos::ParameterList > getCurrentParameters() const override
Get a parameter list containing the current parameters for this object.
static constexpr int maxIters_default_
MagnitudeType achievedTol() const override
Tolerance achieved by the last solve() invocation.
static constexpr const char * impResScale_default_
static constexpr int numBlocks_default_
Teuchos::ScalarTraits< ScalarType > SCT
static constexpr int outputStyle_default_
void buildRecycleSpace2(Teuchos::RCP< GCRODRIter< ScalarType, MV, OP > > gcrodr_iter)
bool isLOADetected() const override
Return whether a loss of accuracy was detected by this solver during the most current solve.
Teuchos::RCP< StatusTestMaxIters< ScalarType, MV, OP > > maxIterTest_
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > ¶ms) override
Set the parameters the solver manager should use to solve the linear problem.
Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > R_
Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > B_
Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > HP_
Teuchos::RCP< StatusTest< ScalarType, MV, OP > > sTest_
static const bool requiresLapack
Details::SolverManagerRequiresLapack< ScalarType, MV, OP, requiresLapack > base_type
GCRODRSolMgr(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem, const Teuchos::RCP< Teuchos::ParameterList > &pl)
Traits class which defines basic operations on multivectors.
static void MvTransMv(const ScalarType alpha, const MV &A, const MV &mv, Teuchos::SerialDenseMatrix< int, ScalarType > &B)
Compute a dense matrix B through the matrix-matrix multiply .
static void MvScale(MV &mv, const ScalarType alpha)
Scale each element of the vectors in mv with alpha.
static void MvTimesMatAddMv(const ScalarType alpha, const MV &A, const Teuchos::SerialDenseMatrix< int, ScalarType > &B, const ScalarType beta, MV &mv)
Update mv with .
static ptrdiff_t GetGlobalLength(const MV &mv)
Return the number of rows in the given multivector mv.
static void MvNorm(const MV &mv, std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > &normvec, NormType type=TwoNorm)
Compute the norm of each individual vector of mv. Upon return, normvec[i] holds the value of ,...
static Teuchos::RCP< MV > Clone(const MV &mv, const int numvecs)
Creates a new empty MV containing numvecs columns.
static Teuchos::RCP< const MV > CloneView(const MV &mv, const std::vector< int > &index)
Creates a new const MV that shares the selected contents of mv (shallow copy).
static Teuchos::RCP< MV > CloneViewNonConst(MV &mv, const std::vector< int > &index)
Creates a new MV that shares the selected contents of mv (shallow copy).
static void MvInit(MV &mv, const ScalarType alpha=Teuchos::ScalarTraits< ScalarType >::zero())
Replace each element of the vectors in mv with alpha.
static void SetBlock(const MV &A, const std::vector< int > &index, MV &mv)
Copy the vectors in A to a set of vectors in mv indicated by the indices given in index.
static int GetNumberVecs(const MV &mv)
Obtain the number of vectors in mv.
Class which defines basic traits for the operator type.
Enumeration of all valid Belos (Mat)OrthoManager classes.
std::ostream & printValidNames(std::ostream &out) const
Print all recognized MatOrthoManager names to the given ostream.
bool isValidName(const std::string &name) const
Whether this factory recognizes the MatOrthoManager with the given name.
std::string validNamesString() const
List (as a string) of recognized MatOrthoManager names.
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.
Teuchos::RCP< const Teuchos::ParameterList > getDefaultParameters(const std::string &name) const
Default parameters for the given MatOrthoManager subclass.
Belos's basic output manager for sending information of select verbosity levels to the appropriate ou...
A class for extending the status testing capabilities of Belos via logical combinations.
Exception thrown to signal error in a status test during Belos::StatusTest::checkStatus().
An implementation of StatusTestResNorm using a family of residual norms.
A Belos::StatusTest class for specifying a maximum number of iterations.
A factory class for generating StatusTestOutput objects.
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.
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.
static const double convTol
Default convergence tolerance.
Structure to contain pointers to GCRODRIter state variables.
Teuchos::RCP< MV > U
The recycled subspace and its projection.
Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > B
The projection of the Krylov subspace against the recycled subspace.
Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > H
The current Hessenberg matrix.
Teuchos::RCP< MV > V
The current Krylov basis.
int curDim
The current dimension of the reduction.