42 #ifndef BELOS_PCPG_SOLMGR_HPP
43 #define BELOS_PCPG_SOLMGR_HPP
66 #ifdef BELOS_TEUCHOS_TIME_MONITOR
69 #if defined(HAVE_TEUCHOSCORE_CXX11)
70 # include <type_traits>
71 #endif // defined(HAVE_TEUCHOSCORE_CXX11)
151 template<
class ScalarType,
class MV,
class OP,
152 const bool supportsScalarType =
157 Belos::Details::LapackSupportsScalar<ScalarType>::value &&
158 ! Teuchos::ScalarTraits<ScalarType>::isComplex>
182 template<
class ScalarType,
class MV,
class OP>
274 return Teuchos::tuple(timerSolve_);
345 std::string description()
const;
375 static constexpr
int maxIters_default_ = 1000;
376 static constexpr
int deflatedBlocks_default_ = 2;
377 static constexpr
int savedBlocks_default_ = 16;
380 static constexpr
int outputFreq_default_ = -1;
381 static constexpr
const char * label_default_ =
"Belos";
382 static constexpr
const char * orthoType_default_ =
"DGKS";
383 static constexpr std::ostream * outputStream_default_ = &std::cout;
404 int deflatedBlocks_, savedBlocks_,
verbosity_, outputStyle_, outputFreq_;
423 template<
class ScalarType,
class MV,
class OP>
425 outputStream_(Teuchos::
rcp(outputStream_default_,false)),
430 maxIters_(maxIters_default_),
431 deflatedBlocks_(deflatedBlocks_default_),
432 savedBlocks_(savedBlocks_default_),
433 verbosity_(verbosity_default_),
434 outputStyle_(outputStyle_default_),
435 outputFreq_(outputFreq_default_),
436 orthoType_(orthoType_default_),
438 label_(label_default_),
444 template<
class ScalarType,
class MV,
class OP>
449 outputStream_(Teuchos::
rcp(outputStream_default_,false)),
455 maxIters_(maxIters_default_),
456 deflatedBlocks_(deflatedBlocks_default_),
457 savedBlocks_(savedBlocks_default_),
458 verbosity_(verbosity_default_),
459 outputStyle_(outputStyle_default_),
460 outputFreq_(outputFreq_default_),
461 orthoType_(orthoType_default_),
463 label_(label_default_),
467 problem_.is_null (), std::invalid_argument,
468 "Belos::PCPGSolMgr two-argument constructor: "
469 "'problem' is null. You must supply a non-null Belos::LinearProblem "
470 "instance when calling this constructor.");
479 template<
class ScalarType,
class MV,
class OP>
492 maxIters_ = params->
get(
"Maximum Iterations",maxIters_default_);
495 params_->set(
"Maximum Iterations", maxIters_);
497 maxIterTest_->setMaxIters( maxIters_ );
502 savedBlocks_ = params->
get(
"Num Saved Blocks",savedBlocks_default_);
504 "Belos::PCPGSolMgr: \"Num Saved Blocks\" must be strictly positive.");
511 params_->set(
"Num Saved Blocks", savedBlocks_);
514 deflatedBlocks_ = params->
get(
"Num Deflated Blocks",deflatedBlocks_default_);
516 "Belos::PCPGSolMgr: \"Num Deflated Blocks\" must be positive.");
519 "Belos::PCPGSolMgr: \"Num Deflated Blocks\" must be <= \"Num Saved Blocks\".");
523 params_->set(
"Num Deflated Blocks", static_cast<int>(deflatedBlocks_));
528 std::string tempLabel = params->
get(
"Timer Label", label_default_);
531 if (tempLabel != label_) {
533 params_->set(
"Timer Label", label_);
534 std::string solveLabel = label_ +
": PCPGSolMgr total solve time";
535 #ifdef BELOS_TEUCHOS_TIME_MONITOR
539 ortho_->setLabel( label_ );
546 std::string tempOrthoType = params->
get(
"Orthogonalization",orthoType_default_);
548 std::invalid_argument,
549 "Belos::PCPGSolMgr: \"Orthogonalization\" must be either \"DGKS\", \"ICGS\", or \"IMGS\".");
550 if (tempOrthoType != orthoType_) {
551 orthoType_ = tempOrthoType;
552 params_->set(
"Orthogonalization", orthoType_);
554 if (orthoType_==
"DGKS") {
555 if (orthoKappa_ <= 0) {
563 else if (orthoType_==
"ICGS") {
566 else if (orthoType_==
"IMGS") {
573 if (params->
isParameter(
"Orthogonalization Constant")) {
575 orthoKappa_ = params->
get (
"Orthogonalization Constant",
579 orthoKappa_ = params->
get (
"Orthogonalization Constant",
584 params_->set(
"Orthogonalization Constant",orthoKappa_);
585 if (orthoType_==
"DGKS") {
594 if (Teuchos::isParameterType<int>(*params,
"Verbosity")) {
595 verbosity_ = params->
get(
"Verbosity", verbosity_default_);
597 verbosity_ = (int)Teuchos::getParameter<Belos::MsgType>(*params,
"Verbosity");
601 params_->set(
"Verbosity", verbosity_);
603 printer_->setVerbosity(verbosity_);
608 if (Teuchos::isParameterType<int>(*params,
"Output Style")) {
609 outputStyle_ = params->
get(
"Output Style", outputStyle_default_);
611 outputStyle_ = (int)Teuchos::getParameter<Belos::OutputType>(*params,
"Output Style");
615 params_->set(
"Output Style", outputStyle_);
621 outputStream_ = Teuchos::getParameter<Teuchos::RCP<std::ostream> >(*params,
"Output Stream");
624 params_->set(
"Output Stream", outputStream_);
626 printer_->setOStream( outputStream_ );
632 outputFreq_ = params->
get(
"Output Frequency", outputFreq_default_);
636 params_->set(
"Output Frequency", outputFreq_);
638 outputTest_->setOutputFrequency( outputFreq_ );
651 if (params->
isParameter(
"Convergence Tolerance")) {
653 convtol_ = params->
get (
"Convergence Tolerance",
661 params_->set(
"Convergence Tolerance", convtol_);
663 convTest_->setTolerance( convtol_ );
673 convTest_ =
Teuchos::rcp(
new StatusTestResNorm_t( convtol_, 1 ) );
675 sTest_ =
Teuchos::rcp(
new StatusTestCombo_t( StatusTestCombo_t::OR, maxIterTest_, convTest_ ) );
683 std::string solverDesc =
" PCPG ";
684 outputTest_->setSolverDesc( solverDesc );
689 params_->set(
"Orthogonalization", orthoType_);
690 if (orthoType_==
"DGKS") {
691 if (orthoKappa_ <= 0) {
699 else if (orthoType_==
"ICGS") {
702 else if (orthoType_==
"IMGS") {
707 "Belos::PCPGSolMgr(): Invalid orthogonalization type.");
713 std::string solveLabel = label_ +
": PCPGSolMgr total solve time";
714 #ifdef BELOS_TEUCHOS_TIME_MONITOR
724 template<
class ScalarType,
class MV,
class OP>
733 "The relative residual tolerance that needs to be achieved by the\n"
734 "iterative solver in order for the linear system to be declared converged.");
735 pl->
set(
"Maximum Iterations", static_cast<int>(maxIters_default_),
736 "The maximum number of iterations allowed for each\n"
737 "set of RHS solved.");
738 pl->
set(
"Num Deflated Blocks", static_cast<int>(deflatedBlocks_default_),
739 "The maximum number of vectors in the seed subspace." );
740 pl->
set(
"Num Saved Blocks", static_cast<int>(savedBlocks_default_),
741 "The maximum number of vectors saved from old Krylov subspaces." );
742 pl->
set(
"Verbosity", static_cast<int>(verbosity_default_),
743 "What type(s) of solver information should be outputted\n"
744 "to the output stream.");
745 pl->
set(
"Output Style", static_cast<int>(outputStyle_default_),
746 "What style is used for the solver information outputted\n"
747 "to the output stream.");
748 pl->
set(
"Output Frequency", static_cast<int>(outputFreq_default_),
749 "How often convergence information should be outputted\n"
750 "to the output stream.");
752 "A reference-counted pointer to the output stream where all\n"
753 "solver output is sent.");
754 pl->
set(
"Timer Label", static_cast<const char *>(label_default_),
755 "The string to use as a prefix for the timer labels.");
756 pl->
set(
"Orthogonalization", static_cast<const char *>(orthoType_default_),
757 "The type of orthogonalization to use: DGKS, ICGS, IMGS");
759 "The constant used by DGKS orthogonalization to determine\n"
760 "whether another step of classical Gram-Schmidt is necessary.");
768 template<
class ScalarType,
class MV,
class OP>
772 if (!isSet_) { setParameters( params_ ); }
780 "Belos::PCPGSolMgr::solve(): Linear problem is not a valid object.");
783 "Belos::PCPGSolMgr::solve(): Linear problem is not ready, setProblem() has not been called.");
786 int numRHS2Solve = MVT::GetNumberVecs( *(problem_->getRHS()) );
787 std::vector<int> currIdx(1);
793 problem_->setLSIndex( currIdx );
796 bool isConverged =
true;
801 plist.
set(
"Saved Blocks", savedBlocks_);
802 plist.
set(
"Block Size", 1);
803 plist.
set(
"Keep Diagonal",
true);
804 plist.
set(
"Initialize Diagonal",
true);
815 #ifdef BELOS_TEUCHOS_TIME_MONITOR
818 while ( numRHS2Solve > 0 ) {
821 outputTest_->reset();
825 R_ = MVT::Clone( *(problem_->getRHS()), 1 );
827 problem_->computeCurrResVec( &*R_ );
839 std::vector<MagnitudeType> rnorm0(1);
840 MVT::MvNorm( *R_, rnorm0 );
843 std::cout <<
"Solver Manager: dimU_ = " << dimU_ << std::endl;
847 std::vector<int> active_columns( dimU_ );
848 for (
int i=0; i < dimU_; ++i) active_columns[i] = i;
849 Uactive = MVT::CloneView(*U_, active_columns);
850 Cactive = MVT::CloneView(*C_, active_columns);
853 std::cout <<
" Solver Manager : check duality of seed basis " << std::endl;
855 MVT::MvTransMv( one, *Uactive, *Cactive, H );
856 H.
print( std::cout );
859 MVT::MvTransMv( one, *Uactive, *R_, Z );
861 MVT::MvTimesMatAddMv( one, *Uactive, Z, zero, *tempU );
862 MVT::MvAddMv( one, *tempU, one, *cur_soln_vec, *cur_soln_vec );
863 MVT::MvTimesMatAddMv( one, *Cactive, Z, zero, *tempU );
864 MVT::MvAddMv( -one, *tempU, one, *R_, *R_ );
865 std::vector<MagnitudeType> rnorm(1);
866 MVT::MvNorm( *R_, rnorm );
867 if( rnorm[0] < rnorm0[0] * .001 ){
868 MVT::MvTransMv( one, *Uactive, *R_, Z );
869 MVT::MvTimesMatAddMv( one, *Uactive, Z, zero, *tempU );
870 MVT::MvAddMv( one, *tempU, one, *cur_soln_vec, *cur_soln_vec );
871 MVT::MvTimesMatAddMv( one, *Cactive, Z, zero, *tempU );
872 MVT::MvAddMv( -one, *tempU, one, *R_, *R_ );
889 if( dimU_ > 0 ) pcpgState.
curDim = dimU_;
890 pcpg_iter->initialize(pcpgState);
896 if( !dimU_ ) printer_->stream(
Debug) <<
" No recycled subspace available for RHS index " << currIdx[0] << std::endl << std::endl;
897 pcpg_iter->resetNumIters();
899 if( dimU_ > savedBlocks_ )
900 std::cout <<
"Error: dimU_ = " << dimU_ <<
" > savedBlocks_ = " << savedBlocks_ << std::endl;
906 if( debug ) printf(
"********** Calling iterate...\n");
907 pcpg_iter->iterate();
914 if ( convTest_->getStatus() ==
Passed ) {
923 else if ( maxIterTest_->getStatus() ==
Passed ) {
938 "Belos::PCPGSolMgr::solve(): Invalid return from PCPGIter::iterate().");
944 sTest_->checkStatus( &*pcpg_iter );
945 if (convTest_->getStatus() !=
Passed)
949 catch (
const std::exception &e) {
950 printer_->stream(
Errors) <<
"Error! Caught exception in PCPGIter::iterate() at iteration "
951 << pcpg_iter->getNumIters() << std::endl
952 << e.what() << std::endl;
959 problem_->updateSolution( update,
true );
962 problem_->setCurrLS();
970 std::cout <<
"SolverManager: dimU_ " << dimU_ <<
" prevUdim= " << q << std::endl;
972 if( q > deflatedBlocks_ )
973 std::cout <<
"SolverManager: Error deflatedBlocks = " << deflatedBlocks_ << std::endl;
984 rank = ARRQR(dimU_,q, *oldState.
D );
986 std::cout <<
" rank decreased in ARRQR, something to do? " << std::endl;
992 if( dimU_ > deflatedBlocks_ ){
994 if( !deflatedBlocks_ ){
997 dimU_ = deflatedBlocks_;
1001 bool Harmonic =
false;
1005 std::vector<int> active_cols( dimU_ );
1006 for (
int i=0; i < dimU_; ++i) active_cols[i] = i;
1009 Uorth = MVT::CloneCopy(*C_, active_cols);
1012 Uorth = MVT::CloneCopy(*U_, active_cols);
1017 rank = ortho_->normalize(*Uorth,
Teuchos::rcp(&R,
false));
1024 "Belos::PCPGSolMgr::solve(): Failed to compute orthonormal basis for initial recycled subspace.");
1030 int lwork = 5*dimU_;
1033 if( problem_->isHermitian() ) lrwork = dimU_;
1034 std::vector<ScalarType> work(lwork);
1035 std::vector<ScalarType> Svec(dimU_);
1036 std::vector<ScalarType> rwork(lrwork);
1037 lapack.
GESVD(
'N',
'O',
1046 "Belos::PCPGSolMgr::solve(): LAPACK _GESVD failed to compute singular values.");
1048 if( work[0] != 67. * dimU_ )
1049 std::cout <<
" SVD " << dimU_ <<
" lwork " << work[0] << std::endl;
1050 for(
int i=0; i< dimU_; i++)
1051 std::cout << i <<
" " << Svec[i] << std::endl;
1055 int startRow = 0, startCol = 0;
1057 startCol = dimU_ - deflatedBlocks_;
1065 std::vector<int> active_columns( dimU_ );
1066 std::vector<int> def_cols( deflatedBlocks_ );
1067 for (
int i=0; i < dimU_; ++i) active_columns[i] = i;
1068 for (
int i=0; i < deflatedBlocks_; ++i) def_cols[i] = i;
1072 MVT::MvTimesMatAddMv( one, *Ucopy, V, zero, *Uactive );
1077 MVT::MvTimesMatAddMv( one, *Ccopy, V, zero, *Cactive );
1080 dimU_ = deflatedBlocks_;
1082 printer_->stream(
Debug) <<
" Generated recycled subspace using RHS index " << currIdx[0] <<
" of dimension " << dimU_ << std::endl << std::endl;
1085 problem_->setCurrLS();
1089 if ( numRHS2Solve > 0 ) {
1093 problem_->setLSIndex( currIdx );
1096 currIdx.resize( numRHS2Solve );
1105 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1115 using Teuchos::rcp_dynamic_cast;
1118 const std::vector<MagnitudeType>* pTestValues =
1119 rcp_dynamic_cast<conv_test_type>(convTest_)->getTestValue();
1122 "Belos::PCPGSolMgr::solve(): The convergence test's getTestValue() "
1123 "method returned NULL. Please report this bug to the Belos developers.");
1126 "Belos::PCPGSolMgr::solve(): The convergence test's getTestValue() "
1127 "method returned a vector of length zero. Please report this bug to the "
1128 "Belos developers.");
1133 achievedTol_ = *std::max_element (pTestValues->begin(), pTestValues->end());
1137 numIters_ = maxIterTest_->getNumIters();
1148 template<
class ScalarType,
class MV,
class OP>
1159 std::vector<int> curind(1);
1160 std::vector<int> ipiv(p - q);
1161 std::vector<ScalarType> Pivots(p);
1163 ScalarType rteps = 1.5e-8;
1166 for( i = q ; i < p ; i++ ){
1169 RCP<MV> P = MVT::CloneViewNonConst(*U_,curind);
1170 RCP<MV> AP = MVT::CloneViewNonConst(*C_,curind);
1172 MVT::MvAddMv( anorm(0,0), *P, zero, *AP, *P );
1173 MVT::MvAddMv( zero, *P, anorm(0,0), *AP, *AP );
1177 for( i = q ; i < p ; i++ ){
1178 if( q < i && i < p-1 ){
1181 for( j = i+1 ; j < p ; j++ ){
1182 const int k = ipiv[j-q];
1183 if( Pivots[k] > Pivots[l] ){
1190 ipiv[imax-q] = ipiv[i-q];
1196 if( Pivots[k] > 1.5625e-2 ){
1197 anorm(0,0) = Pivots[k];
1201 RCP<const MV> P = MVT::CloneView(*U_,curind);
1202 RCP<const MV> AP = MVT::CloneView(*C_,curind);
1203 MVT::MvTransMv( one, *P, *AP, anorm );
1206 if( rteps <= anorm(0,0) && anorm(0,0) < 9.765625e-4){
1214 std::cout <<
"ARRQR: Bad case not implemented" << std::endl;
1216 if( anorm(0,0) < rteps ){
1217 std::cout <<
"ARRQR : deficient case not implemented " << std::endl;
1225 RCP<MV> P = MVT::CloneViewNonConst(*U_,curind);
1226 RCP<MV> AP = MVT::CloneViewNonConst(*C_,curind);
1227 MVT::MvAddMv( anorm(0,0), *P, zero, *AP, *P );
1228 MVT::MvAddMv( zero, *P, anorm(0,0), *AP, *AP );
1232 P = MVT::CloneViewNonConst(*U_,curind);
1233 AP = MVT::CloneViewNonConst(*C_,curind);
1234 for( j = i+1 ; j < p ; j++ ){
1237 RCP<MV> Q = MVT::CloneViewNonConst(*U_,curind);
1238 MVT::MvTransMv( one, *Q, *AP, alpha);
1239 MVT::MvAddMv( -alpha(0,0), *P, one, *Q, *Q );
1241 RCP<MV> AQ = MVT::CloneViewNonConst(*C_,curind);
1242 MVT::MvAddMv( -alpha(0,0), *AP, one, *AQ, *AQ );
1244 gamma(0,0) = ( Pivots[l] - alpha(0,0))*( Pivots[l] + alpha(0,0));
1245 if( gamma(0,0) > 0){
1257 template<
class ScalarType,
class MV,
class OP>
1260 std::ostringstream oss;
1263 oss <<
"Ortho Type='"<<orthoType_;
static const double orthoKappa
DGKS orthogonalization constant.
ScalarType * values() const
Collection of types and exceptions used within the Belos solvers.
int prevUdim
Number of block columns in matrices C and U before current iteration.
bool isLOADetected() const
Return whether a loss of accuracy was detected by this solver during the most current solve...
Belos's basic output manager for sending information of select verbosity levels to the appropriate ou...
PCPGSolMgr(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem, const Teuchos::RCP< Teuchos::ParameterList > &pl)
Class which manages the output and verbosity of the Belos solvers.
static const bool scalarTypeIsSupported
bool is_null(const boost::shared_ptr< T > &p)
virtual void print(std::ostream &os) const
Teuchos::RCP< StatusTestOutput< ScalarType, MV, OP > > outputTest_
Belos concrete class to iterate Preconditioned Conjugate Projected Gradients.
int numIters_
Number of iterations taken by the last solve() invocation.
PCPGSolMgrOrthoFailure is thrown when the orthogonalization manager is unable to generate orthonormal...
Teuchos::RCP< MV > R
The current residual.
T & get(ParameterList &l, const std::string &name)
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
bool is_null(const std::shared_ptr< T > &p)
Teuchos::RCP< MV > C
C = AU, U spans recycled subspace.
Base class for Belos::SolverManager subclasses which normally can only compile with real ScalarType t...
PCPGSolMgrRecyclingFailure is thrown when any problem occurs in using/creating the recycling subspace...
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
A factory class for generating StatusTestOutput objects.
int curDim
The current dimension of the reduction.
An implementation of StatusTestResNorm using a family of residual norms.
PCPGSolMgrOrthoFailure(const std::string &what_arg)
int maxIters_
Maximum iteration count (read from parameter list).
const LinearProblem< ScalarType, MV, OP > & getProblem() const
Get current linear problem being solved for in this object.
MagnitudeType orthoKappa_
Orthogonalization parameter (read from parameter list).
Structure to contain pointers to PCPGIter state variables.
Teuchos::RCP< Teuchos::Time > timerSolve_
PCPGSolMgrLinearProblemFailure is thrown when the linear problem is not setup (i.e.
static const double convTol
Default convergence tolerance.
Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > problem_
Belos::StatusTest class for specifying a maximum number of iterations.
Teuchos::RCP< std::ostream > outputStream_
static std::string name()
A factory class for generating StatusTestOutput objects.
Iterated Modified Gram-Schmidt (IMGS) implementation of the Belos::OrthoManager class.
Traits class which defines basic operations on multivectors.
Belos::StatusTest for logically combining several status tests.
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > D
The current Hessenberg matrix.
bool isParameter(const std::string &name) const
Teuchos::RCP< const Teuchos::ParameterList > getCurrentParameters() const
Get a parameter list containing the current parameters for this object.
OperatorTraits< ScalarType, MV, OP > OPT
Classical Gram-Schmidt (with DGKS correction) implementation of the Belos::OrthoManager class...
A Belos::StatusTest class for specifying a maximum number of iterations.
MagnitudeType achievedTol() const
Tolerance achieved by the last solve() invocation.
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. ...
MagnitudeType convtol_
Convergence tolerance (read from parameter list).
void reset(const ResetType type)
Performs a reset of the solver manager specified by the ResetType. This informs the solver manager th...
Teuchos::RCP< StatusTest< ScalarType, MV, OP > > sTest_
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)
MagnitudeType achievedTol_
Tolerance achieved by the last solve() invocation.
A linear system to solve, and its associated information.
MultiVecTraits< ScalarType, MV > MVT
void setProblem(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem)
Set the linear problem that needs to be solved.
Class which describes the linear problem to be solved by the iterative solver.
PCPGSolMgrLAPACKFailure is thrown when a nonzero value is retuned from an LAPACK call.
virtual Teuchos::RCP< SolverManager< ScalarType, MV, OP > > clone() const
clone for Inverted Injection (DII)
PCPGSolMgrLAPACKFailure(const std::string &what_arg)
int getNumIters() const
Get the iteration count for the most recent call to solve().
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.
Iterated Classical Gram-Schmidt (ICGS) implementation of the Belos::OrthoManager class.
void validateParameters(ParameterList const &validParamList, int const depth=1000, EValidateUsed const validateUsed=VALIDATE_USED_ENABLED, EValidateDefaults const validateDefaults=VALIDATE_DEFAULTS_ENABLED) const
An implementation of the Belos::MatOrthoManager that performs orthogonalization using multiple steps ...
Teuchos::RCP< OutputManager< ScalarType > > printer_
Teuchos::RCP< MV > U
The recycled subspace.
Teuchos::Array< Teuchos::RCP< Teuchos::Time > > getTimers() const
Return the timers for this object.
An implementation of the Belos::MatOrthoManager that performs orthogonalization using multiple steps ...
OrdinalType numCols() const
Teuchos::ScalarTraits< ScalarType >::magnitudeType MagnitudeType
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.
PCPGIterOrthoFailure is thrown when the PCPGIter object is unable to compute independent direction ve...
void GESVD(const char &JOBU, const char &JOBVT, const OrdinalType &m, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, MagnitudeType *S, ScalarType *U, const OrdinalType &ldu, ScalarType *V, const OrdinalType &ldv, ScalarType *WORK, const OrdinalType &lwork, MagnitudeType *RWORK, OrdinalType *info) const
Teuchos::ScalarTraits< ScalarType > SCT
Belos::StatusTestResNorm for specifying general residual norm stopping criteria.
PCPG iterative linear solver.
bool isType(const std::string &name) const
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.
Teuchos::RCP< StatusTestGenResNorm< ScalarType, MV, OP > > convTest_
Teuchos::RCP< StatusTestMaxIters< ScalarType, MV, OP > > maxIterTest_
virtual ~PCPGSolMgr()
Destructor.
Class which defines basic traits for the operator type.
Parent class to all Belos exceptions.
Default parameters common to most Belos solvers.
Teuchos::RCP< MatOrthoManager< ScalarType, MV, OP > > ortho_
Details::SolverManagerRequiresRealLapack< ScalarType, MV, OP, scalarTypeIsSupported > base_type
Teuchos::RCP< SolverManager< ScalarType, MV, OP > > clone() const override
clone for Inverted Injection (DII)
PCPGSolMgrRecyclingFailure(const std::string &what_arg)
This class implements the PCPG iteration, where a single-std::vector Krylov subspace is constructed...
Belos header file which uses auto-configuration information to include necessary C++ headers...
Teuchos::RCP< Teuchos::ParameterList > params_
PCPGSolMgrLinearProblemFailure(const std::string &what_arg)
OrdinalType numRows() const
Teuchos::ScalarTraits< MagnitudeType > MT