42#ifndef BELOS_PCPG_SOLMGR_HPP
43#define BELOS_PCPG_SOLMGR_HPP
62#include "Teuchos_LAPACK.hpp"
63#ifdef BELOS_TEUCHOS_TIME_MONITOR
64# include "Teuchos_TimeMonitor.hpp"
66#if defined(HAVE_TEUCHOSCORE_CXX11)
67# include <type_traits>
69#include "Teuchos_TypeTraits.hpp"
138 template<
class ScalarType,
class MV,
class OP,
139 const bool supportsScalarType =
141 ! Teuchos::ScalarTraits<ScalarType>::isComplex>
144 Belos::Details::LapackSupportsScalar<ScalarType>::value &&
145 ! Teuchos::ScalarTraits<ScalarType>::isComplex>
149 ! Teuchos::ScalarTraits<ScalarType>::isComplex;
158 const Teuchos::RCP<Teuchos::ParameterList> &pl) :
164 Teuchos::RCP<SolverManager<ScalarType, MV, OP> >
clone ()
const override {
169 template<
class ScalarType,
class MV,
class OP>
175 typedef Teuchos::ScalarTraits<ScalarType>
SCT;
176 typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType
MagnitudeType;
177 typedef Teuchos::ScalarTraits<MagnitudeType>
MT;
227 const Teuchos::RCP<Teuchos::ParameterList> &pl );
233 virtual Teuchos::RCP<SolverManager<ScalarType, MV, OP> >
clone ()
const {
249 Teuchos::RCP<const Teuchos::ParameterList> getValidParameters()
const;
260 Teuchos::Array<Teuchos::RCP<Teuchos::Time> >
getTimers()
const {
291 void setParameters(
const Teuchos::RCP<Teuchos::ParameterList> ¶ms );
332 std::string description()
const;
340 int ARRQR(
int numVecs,
int numOrthVecs,
const Teuchos::SerialDenseMatrix<int,ScalarType>& D);
343 Teuchos::RCP<LinearProblem<ScalarType,MV,OP> >
problem_;
350 Teuchos::RCP<StatusTest<ScalarType,MV,OP> >
sTest_;
352 Teuchos::RCP<StatusTestGenResNorm<ScalarType,MV,OP> >
convTest_;
356 Teuchos::RCP<MatOrthoManager<ScalarType,MV,OP> >
ortho_;
409template<
class ScalarType,
class MV,
class OP>
430template<
class ScalarType,
class MV,
class OP>
433 const Teuchos::RCP<Teuchos::ParameterList> &pl ) :
452 TEUCHOS_TEST_FOR_EXCEPTION(
453 problem_.is_null (), std::invalid_argument,
454 "Belos::PCPGSolMgr two-argument constructor: "
455 "'problem' is null. You must supply a non-null Belos::LinearProblem "
456 "instance when calling this constructor.");
458 if (! pl.is_null ()) {
465template<
class ScalarType,
class MV,
class OP>
469 if (
params_ == Teuchos::null) {
477 if (params->isParameter(
"Maximum Iterations")) {
487 if (params->isParameter(
"Num Saved Blocks")) {
489 TEUCHOS_TEST_FOR_EXCEPTION(
savedBlocks_ <= 0, std::invalid_argument,
490 "Belos::PCPGSolMgr: \"Num Saved Blocks\" must be strictly positive.");
499 if (params->isParameter(
"Num Deflated Blocks")) {
502 "Belos::PCPGSolMgr: \"Num Deflated Blocks\" must be positive.");
505 "Belos::PCPGSolMgr: \"Num Deflated Blocks\" must be <= \"Num Saved Blocks\".");
513 if (params->isParameter(
"Timer Label")) {
514 std::string tempLabel = params->get(
"Timer Label",
label_default_);
517 if (tempLabel !=
label_) {
520 std::string solveLabel =
label_ +
": PCPGSolMgr total solve time";
521#ifdef BELOS_TEUCHOS_TIME_MONITOR
522 timerSolve_ = Teuchos::TimeMonitor::getNewCounter(solveLabel);
524 if (
ortho_ != Teuchos::null) {
531 if (params->isParameter(
"Verbosity")) {
532 if (Teuchos::isParameterType<int>(*params,
"Verbosity")) {
535 verbosity_ = (int)Teuchos::getParameter<Belos::MsgType>(*params,
"Verbosity");
545 if (params->isParameter(
"Output Style")) {
546 if (Teuchos::isParameterType<int>(*params,
"Output Style")) {
549 outputStyle_ = (int)Teuchos::getParameter<Belos::OutputType>(*params,
"Output Style");
558 if (params->isParameter(
"Output Stream")) {
559 outputStream_ = Teuchos::getParameter<Teuchos::RCP<std::ostream> >(*params,
"Output Stream");
569 if (params->isParameter(
"Output Frequency")) {
585 bool changedOrthoType =
false;
586 if (params->isParameter(
"Orthogonalization")) {
590 changedOrthoType =
true;
596 if (params->isParameter(
"Orthogonalization Constant")) {
597 if (params->isType<
MagnitudeType> (
"Orthogonalization Constant")) {
598 orthoKappa_ = params->get (
"Orthogonalization Constant",
602 orthoKappa_ = params->get (
"Orthogonalization Constant",
610 Teuchos::rcp_dynamic_cast<DGKSOrthoManager<ScalarType,MV,OP> >(
ortho_)->setDepTol(
orthoKappa_ );
616 if (
ortho_ == Teuchos::null || changedOrthoType) {
618 Teuchos::RCP<Teuchos::ParameterList> paramsOrtho;
631 if (params->isParameter(
"Convergence Tolerance")) {
632 if (params->isType<
MagnitudeType> (
"Convergence Tolerance")) {
633 convtol_ = params->get (
"Convergence Tolerance",
663 std::string solverDesc =
" PCPG ";
668 std::string solveLabel =
label_ +
": PCPGSolMgr total solve time";
669#ifdef BELOS_TEUCHOS_TIME_MONITOR
670 timerSolve_ = Teuchos::TimeMonitor::getNewCounter(solveLabel);
679template<
class ScalarType,
class MV,
class OP>
680Teuchos::RCP<const Teuchos::ParameterList>
683 static Teuchos::RCP<const Teuchos::ParameterList> validPL;
684 if (is_null(validPL)) {
685 Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
688 "The relative residual tolerance that needs to be achieved by the\n"
689 "iterative solver in order for the linear system to be declared converged.");
691 "The maximum number of iterations allowed for each\n"
692 "set of RHS solved.");
694 "The maximum number of vectors in the seed subspace." );
696 "The maximum number of vectors saved from old Krylov subspaces." );
698 "What type(s) of solver information should be outputted\n"
699 "to the output stream.");
701 "What style is used for the solver information outputted\n"
702 "to the output stream.");
704 "How often convergence information should be outputted\n"
705 "to the output stream.");
706 pl->set(
"Output Stream", Teuchos::rcpFromRef(std::cout),
707 "A reference-counted pointer to the output stream where all\n"
708 "solver output is sent.");
710 "The string to use as a prefix for the timer labels.");
712 "The type of orthogonalization to use: DGKS, ICGS, IMGS");
714 "The constant used by DGKS orthogonalization to determine\n"
715 "whether another step of classical Gram-Schmidt is necessary.");
723template<
class ScalarType,
class MV,
class OP>
729 Teuchos::LAPACK<int,ScalarType> lapack;
730 ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
731 ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
734 "Belos::PCPGSolMgr::solve(): Linear problem is not a valid object.");
737 "Belos::PCPGSolMgr::solve(): Linear problem is not ready, setProblem() has not been called.");
741 std::vector<int> currIdx(1);
750 bool isConverged =
true;
754 Teuchos::ParameterList plist;
756 plist.set(
"Block Size", 1);
757 plist.set(
"Keep Diagonal",
true);
758 plist.set(
"Initialize Diagonal",
true);
763 Teuchos::RCP<PCPGIter<ScalarType,MV,OP> > pcpg_iter;
769#ifdef BELOS_TEUCHOS_TIME_MONITOR
772 while ( numRHS2Solve > 0 ) {
778 if (
R_ == Teuchos::null)
787 if(
U_ != Teuchos::null ){
792 Teuchos::RCP<MV> cur_soln_vec =
problem_->getCurrLHSVec();
793 std::vector<MagnitudeType> rnorm0(1);
797 std::cout <<
"Solver Manager: dimU_ = " <<
dimU_ << std::endl;
798 Teuchos::SerialDenseMatrix<int,ScalarType> Z(
dimU_, 1 );
800 Teuchos::RCP<const MV> Uactive, Cactive;
801 std::vector<int> active_columns(
dimU_ );
802 for (
int i=0; i <
dimU_; ++i) active_columns[i] = i;
807 std::cout <<
" Solver Manager : check duality of seed basis " << std::endl;
808 Teuchos::SerialDenseMatrix<int,ScalarType> H(
dimU_,
dimU_ );
810 H.print( std::cout );
816 MVT::MvAddMv( one, *tempU, one, *cur_soln_vec, *cur_soln_vec );
819 std::vector<MagnitudeType> rnorm(1);
821 if( rnorm[0] < rnorm0[0] * .001 ){
824 MVT::MvAddMv( one, *tempU, one, *cur_soln_vec, *cur_soln_vec );
828 Uactive = Teuchos::null;
829 Cactive = Teuchos::null;
830 tempU = Teuchos::null;
841 if(
U_ != Teuchos::null ) pcpgState.
U =
U_;
842 if(
C_ != Teuchos::null ) pcpgState.
C =
C_;
844 pcpg_iter->initialize(pcpgState);
850 if( !
dimU_ )
printer_->stream(
Debug) <<
" No recycled subspace available for RHS index " << currIdx[0] << std::endl << std::endl;
851 pcpg_iter->resetNumIters();
854 std::cout <<
"Error: dimU_ = " <<
dimU_ <<
" > savedBlocks_ = " <<
savedBlocks_ << std::endl;
860 if( debug ) printf(
"********** Calling iterate...\n");
861 pcpg_iter->iterate();
891 TEUCHOS_TEST_FOR_EXCEPTION(
true,std::logic_error,
892 "Belos::PCPGSolMgr::solve(): Invalid return from PCPGIter::iterate().");
895 catch (
const std::exception &e) {
896 printer_->stream(
Errors) <<
"Error! Caught exception in PCPGIter::iterate() at iteration "
897 << pcpg_iter->getNumIters() << std::endl
898 << e.what() << std::endl;
904 Teuchos::RCP<MV> update = pcpg_iter->getCurrentUpdate();
905 problem_->updateSolution( update,
true );
916 std::cout <<
"SolverManager: dimU_ " <<
dimU_ <<
" prevUdim= " << q << std::endl;
919 std::cout <<
"SolverManager: Error deflatedBlocks = " <<
deflatedBlocks_ << std::endl;
932 std::cout <<
" rank decreased in ARRQR, something to do? " << std::endl;
947 bool Harmonic =
false;
949 Teuchos::RCP<MV> Uorth;
951 std::vector<int> active_cols(
dimU_ );
952 for (
int i=0; i <
dimU_; ++i) active_cols[i] = i;
962 Teuchos::SerialDenseMatrix<int,ScalarType> R(
dimU_,
dimU_);
963 rank =
ortho_->normalize(*Uorth, Teuchos::rcp(&R,
false));
964 Uorth = Teuchos::null;
970 "Belos::PCPGSolMgr::solve(): Failed to compute orthonormal basis for initial recycled subspace.");
974 Teuchos::SerialDenseMatrix<int,ScalarType> VT;
975 Teuchos::SerialDenseMatrix<int,ScalarType> Ur;
980 std::vector<ScalarType> work(lwork);
981 std::vector<ScalarType> Svec(
dimU_);
982 std::vector<ScalarType> rwork(lrwork);
983 lapack.GESVD(
'N',
'O',
984 R.numRows(),R.numCols(),R.values(), R.numRows(),
992 "Belos::PCPGSolMgr::solve(): LAPACK _GESVD failed to compute singular values.");
994 if( work[0] != 67. *
dimU_ )
995 std::cout <<
" SVD " <<
dimU_ <<
" lwork " << work[0] << std::endl;
996 for(
int i=0; i<
dimU_; i++)
997 std::cout << i <<
" " << Svec[i] << std::endl;
999 Teuchos::SerialDenseMatrix<int,ScalarType> wholeV( R, Teuchos::TRANS);
1001 int startRow = 0, startCol = 0;
1005 Teuchos::SerialDenseMatrix<int,ScalarType> V(Teuchos::Copy,
1011 std::vector<int> active_columns(
dimU_ );
1013 for (
int i=0; i <
dimU_; ++i) active_columns[i] = i;
1019 Ucopy = Teuchos::null;
1020 Uactive = Teuchos::null;
1024 Ccopy = Teuchos::null;
1025 Cactive = Teuchos::null;
1028 printer_->stream(
Debug) <<
" Generated recycled subspace using RHS index " << currIdx[0] <<
" of dimension " <<
dimU_ << std::endl << std::endl;
1035 if ( numRHS2Solve > 0 ) {
1042 currIdx.resize( numRHS2Solve );
1051#ifdef BELOS_TEUCHOS_TIME_MONITOR
1061 using Teuchos::rcp_dynamic_cast;
1064 const std::vector<MagnitudeType>* pTestValues =
1067 TEUCHOS_TEST_FOR_EXCEPTION(pTestValues == NULL, std::logic_error,
1068 "Belos::PCPGSolMgr::solve(): The convergence test's getTestValue() "
1069 "method returned NULL. Please report this bug to the Belos developers.");
1071 TEUCHOS_TEST_FOR_EXCEPTION(pTestValues->size() < 1, std::logic_error,
1072 "Belos::PCPGSolMgr::solve(): The convergence test's getTestValue() "
1073 "method returned a vector of length zero. Please report this bug to the "
1074 "Belos developers.");
1079 achievedTol_ = *std::max_element (pTestValues->begin(), pTestValues->end());
1094template<
class ScalarType,
class MV,
class OP>
1098 ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
1099 ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
1102 Teuchos::SerialDenseMatrix<int,ScalarType> alpha( 1, 1 );
1103 Teuchos::SerialDenseMatrix<int,ScalarType> gamma( 1, 1 );
1104 Teuchos::SerialDenseMatrix<int,ScalarType> anorm( 1, 1 );
1105 std::vector<int> curind(1);
1106 std::vector<int> ipiv(p - q);
1107 std::vector<ScalarType> Pivots(p);
1109 ScalarType rteps = 1.5e-8;
1112 for( i = q ; i < p ; i++ ){
1117 anorm(0,0) = one / Teuchos::ScalarTraits<ScalarType>::squareroot( D(i-q,i-q) ) ;
1123 for( i = q ; i < p ; i++ ){
1124 if( q < i && i < p-1 ){
1127 for( j = i+1 ; j < p ; j++ ){
1128 const int k = ipiv[j-q];
1129 if( Pivots[k] > Pivots[l] ){
1136 ipiv[imax-q] = ipiv[i-q];
1142 if( Pivots[k] > 1.5625e-2 ){
1143 anorm(0,0) = Pivots[k];
1150 anorm(0,0) = Teuchos::ScalarTraits<ScalarType>::squareroot( anorm(0,0) ) ;
1152 if( rteps <= anorm(0,0) && anorm(0,0) < 9.765625e-4){
1160 std::cout <<
"ARRQR: Bad case not implemented" << std::endl;
1162 if( anorm(0,0) < rteps ){
1163 std::cout <<
"ARRQR : deficient case not implemented " << std::endl;
1180 for( j = i+1 ; j < p ; j++ ){
1190 gamma(0,0) = ( Pivots[l] - alpha(0,0))*( Pivots[l] + alpha(0,0));
1191 if( gamma(0,0) > 0){
1192 Pivots[l] = Teuchos::ScalarTraits<ScalarType>::squareroot( gamma(0,0) );
1203template<
class ScalarType,
class MV,
class OP>
1206 std::ostringstream oss;
1207 oss <<
"Belos::PCPGSolMgr<...,"<<Teuchos::ScalarTraits<ScalarType>::name()<<
">";
Belos header file which uses auto-configuration information to include necessary C++ headers.
Class which describes the linear problem to be solved by the iterative solver.
Class which manages the output and verbosity of the Belos solvers.
Belos concrete class to iterate Preconditioned Conjugate Projected Gradients.
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)
Base class for Belos::SolverManager subclasses which normally can only compile with real ScalarType t...
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 MvTimesMatAddMv(const ScalarType alpha, const MV &A, const Teuchos::SerialDenseMatrix< int, ScalarType > &B, const ScalarType beta, MV &mv)
Update mv with .
static Teuchos::RCP< MV > CloneCopy(const MV &mv)
Creates a new MV and copies contents of mv into the new vector (deep copy).
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 MvAddMv(const ScalarType alpha, const MV &A, const ScalarType beta, const MV &B, MV &mv)
Replace mv with .
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.
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.
Belos's basic output manager for sending information of select verbosity levels to the appropriate ou...
This class implements the PCPG iteration, where a single-std::vector Krylov subspace is constructed....
PCPGSolMgrLAPACKFailure is thrown when a nonzero value is retuned from an LAPACK call.
PCPGSolMgrLAPACKFailure(const std::string &what_arg)
PCPGSolMgrLinearProblemFailure is thrown when the linear problem is not setup (i.e.
PCPGSolMgrLinearProblemFailure(const std::string &what_arg)
PCPGSolMgrOrthoFailure is thrown when the orthogonalization manager is unable to generate orthonormal...
PCPGSolMgrOrthoFailure(const std::string &what_arg)
const LinearProblem< ScalarType, MV, OP > & getProblem() const
Get current linear problem being solved for in this object.
Teuchos::RCP< const Teuchos::ParameterList > getCurrentParameters() const
Get a parameter list containing the current parameters for this object.
Teuchos::RCP< StatusTest< ScalarType, MV, OP > > sTest_
int maxIters_
Maximum iteration count (read from parameter list).
void setProblem(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem)
Set the linear problem that needs to be solved.
MagnitudeType orthoKappa_
Orthogonalization parameter (read from parameter list).
int ARRQR(int numVecs, int numOrthVecs, const Teuchos::SerialDenseMatrix< int, ScalarType > &D)
Teuchos::RCP< Teuchos::ParameterList > params_
MagnitudeType achievedTol_
Tolerance achieved by the last solve() invocation.
static constexpr int maxIters_default_
static constexpr const char * orthoType_default_
bool isLOADetected() const
Return whether a loss of accuracy was detected by this solver during the most current solve.
MultiVecTraits< ScalarType, MV > MVT
Teuchos::RCP< StatusTestMaxIters< ScalarType, MV, OP > > maxIterTest_
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > ¶ms)
Set the parameters the solver manager should use to solve the linear problem.
OperatorTraits< ScalarType, MV, OP > OPT
Teuchos::RCP< OutputManager< ScalarType > > printer_
static constexpr const char * label_default_
Teuchos::RCP< MatOrthoManager< ScalarType, MV, OP > > ortho_
static constexpr int outputFreq_default_
int getNumIters() const
Get the iteration count for the most recent call to solve().
Teuchos::RCP< StatusTestOutput< ScalarType, MV, OP > > outputTest_
virtual ~PCPGSolMgr()
Destructor.
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Get a parameter list containing the valid parameters for this object.
static constexpr int verbosity_default_
Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > problem_
int numIters_
Number of iterations taken by the last solve() invocation.
Teuchos::ScalarTraits< ScalarType >::magnitudeType MagnitudeType
PCPGSolMgr()
Empty constructor for PCPGSolMgr. This constructor takes no arguments and sets the default values for...
Teuchos::Array< Teuchos::RCP< Teuchos::Time > > getTimers() const
Return the timers for this object.
Teuchos::ScalarTraits< ScalarType > SCT
Teuchos::RCP< StatusTestGenResNorm< ScalarType, MV, OP > > convTest_
virtual Teuchos::RCP< SolverManager< ScalarType, MV, OP > > clone() const
clone for Inverted Injection (DII)
MagnitudeType convtol_
Convergence tolerance (read from parameter list).
MagnitudeType achievedTol() const
Tolerance achieved by the last solve() invocation.
void reset(const ResetType type)
Performs a reset of the solver manager specified by the ResetType. This informs the solver manager th...
Teuchos::RCP< std::ostream > outputStream_
Teuchos::RCP< Teuchos::Time > timerSolve_
static constexpr int deflatedBlocks_default_
static constexpr int savedBlocks_default_
Teuchos::ScalarTraits< MagnitudeType > MT
static constexpr int outputStyle_default_
Details::SolverManagerRequiresRealLapack< ScalarType, MV, OP, scalarTypeIsSupported > base_type
Teuchos::RCP< SolverManager< ScalarType, MV, OP > > clone() const override
clone for Inverted Injection (DII)
static const bool scalarTypeIsSupported
PCPGSolMgr(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem, const Teuchos::RCP< Teuchos::ParameterList > &pl)
A class for extending the status testing capabilities of Belos via logical combinations.
An implementation of StatusTestResNorm using a family of residual norms.
const std::vector< MagnitudeType > * getTestValue() const
Returns the test value, , computed in most recent call to CheckStatus.
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.
ReturnType
Whether the Belos solve converged for all linear systems.
ResetType
How to reset the solver.
Default parameters common to most Belos solvers.
static const double convTol
Default convergence tolerance.
static const double orthoKappa
DGKS orthogonalization constant.
Structure to contain pointers to PCPGIter state variables.
int prevUdim
Number of block columns in matrices C and U before current iteration.
Teuchos::RCP< MV > U
The recycled subspace.
Teuchos::RCP< MV > C
C = AU, U spans recycled subspace.
int curDim
The current dimension of the reduction.
Teuchos::RCP< MV > R
The current residual.
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > D
The current Hessenberg matrix.