10 #ifndef BELOS_BLOCK_CG_SOLMGR_HPP
11 #define BELOS_BLOCK_CG_SOLMGR_HPP
33 #ifdef BELOS_TEUCHOS_TIME_MONITOR
36 #if defined(HAVE_TEUCHOSCORE_CXX11)
37 # include <type_traits>
38 #endif // defined(HAVE_TEUCHOSCORE_CXX11)
75 template<
class ScalarType,
class MV,
class OP,
76 const bool lapackSupportsScalarType =
101 template<
class ScalarType,
class MV,
class OP>
213 return Teuchos::tuple(timerSolve_);
284 std::string description()
const override;
322 static constexpr
int maxIters_default_ = 1000;
323 static constexpr
bool adaptiveBlockSize_default_ =
true;
324 static constexpr
bool showMaxResNormOnly_default_ =
false;
325 static constexpr
bool useSingleReduction_default_ =
false;
326 static constexpr
int blockSize_default_ = 1;
329 static constexpr
int outputFreq_default_ = -1;
330 static constexpr
const char * resNorm_default_ =
"TwoNorm";
331 static constexpr
bool foldConvergenceDetectionIntoAllreduce_default_ =
false;
332 static constexpr
const char * resScale_default_ =
"Norm of Initial Residual";
333 static constexpr
const char * label_default_ =
"Belos";
334 static constexpr
const char * orthoType_default_ =
"ICGS";
335 static constexpr
bool assertPositiveDefiniteness_default_ =
true;
379 template<
class ScalarType,
class MV,
class OP>
381 outputStream_(Teuchos::rcpFromRef(std::cout)),
385 maxIters_(maxIters_default_),
387 blockSize_(blockSize_default_),
388 verbosity_(verbosity_default_),
389 outputStyle_(outputStyle_default_),
390 outputFreq_(outputFreq_default_),
391 adaptiveBlockSize_(adaptiveBlockSize_default_),
392 showMaxResNormOnly_(showMaxResNormOnly_default_),
393 useSingleReduction_(useSingleReduction_default_),
394 orthoType_(orthoType_default_),
395 resScale_(resScale_default_),
396 assertPositiveDefiniteness_(assertPositiveDefiniteness_default_),
397 foldConvergenceDetectionIntoAllreduce_(foldConvergenceDetectionIntoAllreduce_default_),
398 label_(label_default_),
404 template<
class ScalarType,
class MV,
class OP>
409 outputStream_(Teuchos::rcpFromRef(std::cout)),
413 maxIters_(maxIters_default_),
415 blockSize_(blockSize_default_),
416 verbosity_(verbosity_default_),
417 outputStyle_(outputStyle_default_),
418 outputFreq_(outputFreq_default_),
419 adaptiveBlockSize_(adaptiveBlockSize_default_),
420 showMaxResNormOnly_(showMaxResNormOnly_default_),
421 useSingleReduction_(useSingleReduction_default_),
422 orthoType_(orthoType_default_),
423 resScale_(resScale_default_),
424 assertPositiveDefiniteness_(assertPositiveDefiniteness_default_),
425 foldConvergenceDetectionIntoAllreduce_(foldConvergenceDetectionIntoAllreduce_default_),
426 label_(label_default_),
430 "BlockCGSolMgr's constructor requires a nonnull LinearProblem instance.");
440 template<
class ScalarType,
class MV,
class OP>
455 maxIters_ = params->
get(
"Maximum Iterations",maxIters_default_);
458 params_->set(
"Maximum Iterations", maxIters_);
460 maxIterTest_->setMaxIters( maxIters_ );
465 blockSize_ = params->
get(
"Block Size",blockSize_default_);
467 "Belos::BlockCGSolMgr: \"Block Size\" must be strictly positive.");
470 params_->set(
"Block Size", blockSize_);
475 adaptiveBlockSize_ = params->
get(
"Adaptive Block Size",adaptiveBlockSize_default_);
478 params_->set(
"Adaptive Block Size", adaptiveBlockSize_);
483 useSingleReduction_ = params->
get(
"Use Single Reduction", useSingleReduction_default_);
486 if (params->
isParameter(
"Fold Convergence Detection Into Allreduce")) {
487 foldConvergenceDetectionIntoAllreduce_ = params->
get(
"Fold Convergence Detection Into Allreduce",
488 foldConvergenceDetectionIntoAllreduce_default_);
493 std::string tempLabel = params->
get(
"Timer Label", label_default_);
496 if (tempLabel != label_) {
498 params_->set(
"Timer Label", label_);
499 std::string solveLabel = label_ +
": BlockCGSolMgr total solve time";
500 #ifdef BELOS_TEUCHOS_TIME_MONITOR
504 ortho_->setLabel( label_ );
511 if (Teuchos::isParameterType<int>(*params,
"Verbosity")) {
512 verbosity_ = params->
get(
"Verbosity", verbosity_default_);
514 verbosity_ = (int)Teuchos::getParameter<Belos::MsgType>(*params,
"Verbosity");
518 params_->set(
"Verbosity", verbosity_);
520 printer_->setVerbosity(verbosity_);
525 if (Teuchos::isParameterType<int>(*params,
"Output Style")) {
526 outputStyle_ = params->
get(
"Output Style", outputStyle_default_);
528 outputStyle_ = (int)Teuchos::getParameter<Belos::OutputType>(*params,
"Output Style");
532 params_->set(
"Output Style", outputStyle_);
538 outputStream_ = Teuchos::getParameter<Teuchos::RCP<std::ostream> >(*params,
"Output Stream");
541 params_->set(
"Output Stream", outputStream_);
543 printer_->setOStream( outputStream_ );
549 outputFreq_ = params->
get(
"Output Frequency", outputFreq_default_);
553 params_->set(
"Output Frequency", outputFreq_);
555 outputTest_->setOutputFrequency( outputFreq_ );
564 bool changedOrthoType =
false;
566 std::string tempOrthoType = params->
get(
"Orthogonalization",orthoType_default_);
567 if (tempOrthoType != orthoType_) {
568 orthoType_ = tempOrthoType;
569 changedOrthoType =
true;
572 params_->set(
"Orthogonalization", orthoType_);
575 if (params->
isParameter(
"Orthogonalization Constant")) {
577 orthoKappa_ = params->
get (
"Orthogonalization Constant",
581 orthoKappa_ = params->
get (
"Orthogonalization Constant",
586 params_->set(
"Orthogonalization Constant",orthoKappa_);
587 if (orthoType_==
"DGKS") {
588 if (orthoKappa_ > 0 && ortho_ !=
Teuchos::null && !changedOrthoType) {
598 if (orthoType_==
"DGKS" && orthoKappa_ > 0) {
600 paramsOrtho->
set (
"depTol", orthoKappa_ );
611 if (params->
isParameter(
"Convergence Tolerance")) {
613 convtol_ = params->
get (
"Convergence Tolerance",
621 params_->set(
"Convergence Tolerance", convtol_);
623 convTest_->setTolerance( convtol_ );
626 if (params->
isParameter(
"Show Maximum Residual Norm Only")) {
627 showMaxResNormOnly_ = Teuchos::getParameter<bool>(*params,
"Show Maximum Residual Norm Only");
630 params_->set(
"Show Maximum Residual Norm Only", showMaxResNormOnly_);
632 convTest_->setShowMaxResNormOnly( showMaxResNormOnly_ );
636 bool newResTest =
false;
638 std::string tempResScale = resScale_;
639 if (params->
isParameter (
"Implicit Residual Scaling")) {
640 tempResScale = params->
get<std::string> (
"Implicit Residual Scaling");
644 if (resScale_ != tempResScale) {
647 resScale_ = tempResScale;
650 params_->set (
"Implicit Residual Scaling", resScale_);
652 if (! convTest_.is_null ()) {
656 if (params->
isType<std::string> (
"Residual Norm")) {
660 convTest_->defineResForm(StatusTestResNorm_t::Implicit, normType);
663 catch (std::exception& e) {
678 if (convTest_.is_null () || newResTest) {
682 if (params->
isType<std::string> (
"Residual Norm")) {
687 convTest_ =
rcp (
new StatusTestResNorm_t (convtol_, 1, showMaxResNormOnly_));
688 convTest_->defineResForm(StatusTestResNorm_t::Implicit, normType);
692 if (sTest_.is_null () || newResTest)
693 sTest_ =
Teuchos::rcp(
new StatusTestCombo_t( StatusTestCombo_t::OR, maxIterTest_, convTest_ ) );
695 if (outputTest_.is_null () || newResTest) {
703 std::string solverDesc =
" Block CG ";
704 outputTest_->setSolverDesc( solverDesc );
709 if (params->
isParameter(
"Assert Positive Definiteness")) {
710 assertPositiveDefiniteness_ = Teuchos::getParameter<bool>(*params,
"Assert Positive Definiteness");
711 params_->set(
"Assert Positive Definiteness", assertPositiveDefiniteness_);
716 std::string solveLabel = label_ +
": BlockCGSolMgr total solve time";
717 #ifdef BELOS_TEUCHOS_TIME_MONITOR
727 template<
class ScalarType,
class MV,
class OP>
737 "The relative residual tolerance that needs to be achieved by the\n"
738 "iterative solver in order for the linear system to be declared converged.");
739 pl->
set(
"Maximum Iterations", static_cast<int>(maxIters_default_),
740 "The maximum number of block iterations allowed for each\n"
741 "set of RHS solved.");
742 pl->
set(
"Block Size", static_cast<int>(blockSize_default_),
743 "The number of vectors in each block.");
744 pl->
set(
"Adaptive Block Size", static_cast<bool>(adaptiveBlockSize_default_),
745 "Whether the solver manager should adapt to the block size\n"
746 "based on the number of RHS to solve.");
747 pl->
set(
"Verbosity", static_cast<int>(verbosity_default_),
748 "What type(s) of solver information should be outputted\n"
749 "to the output stream.");
750 pl->
set(
"Output Style", static_cast<int>(outputStyle_default_),
751 "What style is used for the solver information outputted\n"
752 "to the output stream.");
753 pl->
set(
"Output Frequency", static_cast<int>(outputFreq_default_),
754 "How often convergence information should be outputted\n"
755 "to the output stream.");
756 pl->
set(
"Output Stream", Teuchos::rcpFromRef(std::cout),
757 "A reference-counted pointer to the output stream where all\n"
758 "solver output is sent.");
759 pl->
set(
"Show Maximum Residual Norm Only", static_cast<bool>(showMaxResNormOnly_default_),
760 "When convergence information is printed, only show the maximum\n"
761 "relative residual norm when the block size is greater than one.");
762 pl->
set(
"Use Single Reduction", static_cast<bool>(useSingleReduction_default_),
763 "Use single reduction iteration when the block size is one.");
764 pl->
set(
"Implicit Residual Scaling", resScale_default_,
765 "The type of scaling used in the residual convergence test.");
766 pl->
set(
"Timer Label", static_cast<const char *>(label_default_),
767 "The string to use as a prefix for the timer labels.");
768 pl->
set(
"Orthogonalization", static_cast<const char *>(orthoType_default_),
769 "The type of orthogonalization to use: DGKS, ICGS, or IMGS.");
770 pl->
set(
"Assert Positive Definiteness",static_cast<bool>(assertPositiveDefiniteness_default_),
771 "Assert for positivity of p^H*A*p in CG iteration.");
773 "The constant used by DGKS orthogonalization to determine\n"
774 "whether another step of classical Gram-Schmidt is necessary.");
775 pl->
set(
"Residual Norm",static_cast<const char *>(resNorm_default_),
776 "Norm used for the convergence check on the residual.");
777 pl->
set(
"Fold Convergence Detection Into Allreduce",static_cast<bool>(foldConvergenceDetectionIntoAllreduce_default_),
778 "Merge the allreduce for convergence detection with the one for CG.\n"
779 "This saves one all-reduce, but incurs more computation.");
787 template<
class ScalarType,
class MV,
class OP>
791 using Teuchos::rcp_const_cast;
792 using Teuchos::rcp_dynamic_cast;
799 setParameters(Teuchos::parameterList(*getValidParameters()));
806 "Belos::BlockCGSolMgr::solve(): Linear problem is not ready, setProblem() "
807 "has not been called.");
811 int numRHS2Solve = MVT::GetNumberVecs( *(problem_->getRHS()) );
812 int numCurrRHS = ( numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;
814 std::vector<int> currIdx, currIdx2;
819 if ( adaptiveBlockSize_ ) {
820 blockSize_ = numCurrRHS;
821 currIdx.resize( numCurrRHS );
822 currIdx2.resize( numCurrRHS );
823 for (
int i=0; i<numCurrRHS; ++i)
824 { currIdx[i] = startPtr+i; currIdx2[i]=i; }
828 currIdx.resize( blockSize_ );
829 currIdx2.resize( blockSize_ );
830 for (
int i=0; i<numCurrRHS; ++i)
831 { currIdx[i] = startPtr+i; currIdx2[i]=i; }
832 for (
int i=numCurrRHS; i<blockSize_; ++i)
833 { currIdx[i] = -1; currIdx2[i] = i; }
837 problem_->setLSIndex( currIdx );
842 plist.
set(
"Block Size",blockSize_);
845 outputTest_->reset();
849 bool isConverged =
true;
854 plist.
set(
"Assert Positive Definiteness", assertPositiveDefiniteness_);
856 RCP<CGIteration<ScalarType,MV,OP> > block_cg_iter;
857 if (blockSize_ == 1) {
861 plist.
set(
"Fold Convergence Detection Into Allreduce",
862 foldConvergenceDetectionIntoAllreduce_);
863 if (useSingleReduction_) {
866 outputTest_, convTest_, plist));
871 outputTest_, convTest_, plist));
882 #ifdef BELOS_TEUCHOS_TIME_MONITOR
886 while ( numRHS2Solve > 0 ) {
889 std::vector<int> convRHSIdx;
890 std::vector<int> currRHSIdx( currIdx );
891 currRHSIdx.resize(numCurrRHS);
894 block_cg_iter->resetNumIters();
897 outputTest_->resetNumCalls();
900 RCP<MV> R_0 = MVT::CloneViewNonConst( *(rcp_const_cast<MV>(problem_->getInitResVec())), currIdx );
905 block_cg_iter->initializeCG(newstate);
911 block_cg_iter->iterate();
915 if (convTest_->getStatus() ==
Passed) {
920 std::vector<int> convIdx =
921 rcp_dynamic_cast<conv_test_type>(convTest_)->convIndices();
926 if (convIdx.size() == currRHSIdx.size())
931 problem_->setCurrLS();
936 std::vector<int> unconvIdx(currRHSIdx.size());
937 for (
unsigned int i=0; i<currRHSIdx.size(); ++i) {
939 for (
unsigned int j=0; j<convIdx.size(); ++j) {
940 if (currRHSIdx[i] == convIdx[j]) {
946 currIdx2[have] = currIdx2[i];
947 currRHSIdx[have++] = currRHSIdx[i];
952 currRHSIdx.resize(have);
953 currIdx2.resize(have);
956 problem_->setLSIndex( currRHSIdx );
959 std::vector<MagnitudeType> norms;
960 R_0 = MVT::CloneCopy( *(block_cg_iter->getNativeResiduals(&norms)),currIdx2 );
961 for (
int i=0; i<have; ++i) { currIdx2[i] = i; }
964 block_cg_iter->setBlockSize( have );
969 block_cg_iter->initializeCG(defstate);
975 else if (maxIterTest_->getStatus() ==
Passed) {
985 "Belos::BlockCGSolMgr::solve(): Neither the convergence test nor "
986 "the maximum iteration count test passed. Please report this bug "
987 "to the Belos developers.");
992 achievedTol_ = MT::one();
994 MVT::MvInit( *X, SCT::zero() );
995 printer_->stream(
Warnings) <<
"Belos::BlockCGSolMgr::solve(): Warning! NaN has been detected!"
999 catch (
const std::exception &e) {
1000 std::ostream& err = printer_->stream (
Errors);
1001 err <<
"Error! Caught std::exception in CGIteration::iterate() at "
1002 <<
"iteration " << block_cg_iter->getNumIters() << std::endl
1003 << e.what() << std::endl;
1010 problem_->setCurrLS();
1013 startPtr += numCurrRHS;
1014 numRHS2Solve -= numCurrRHS;
1015 if ( numRHS2Solve > 0 ) {
1016 numCurrRHS = ( numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;
1018 if ( adaptiveBlockSize_ ) {
1019 blockSize_ = numCurrRHS;
1020 currIdx.resize( numCurrRHS );
1021 currIdx2.resize( numCurrRHS );
1022 for (
int i=0; i<numCurrRHS; ++i)
1023 { currIdx[i] = startPtr+i; currIdx2[i] = i; }
1026 currIdx.resize( blockSize_ );
1027 currIdx2.resize( blockSize_ );
1028 for (
int i=0; i<numCurrRHS; ++i)
1029 { currIdx[i] = startPtr+i; currIdx2[i] = i; }
1030 for (
int i=numCurrRHS; i<blockSize_; ++i)
1031 { currIdx[i] = -1; currIdx2[i] = i; }
1034 problem_->setLSIndex( currIdx );
1037 block_cg_iter->setBlockSize( blockSize_ );
1040 currIdx.resize( numRHS2Solve );
1051 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1062 numIters_ = maxIterTest_->getNumIters();
1068 const std::vector<MagnitudeType>* pTestValues =
1069 rcp_dynamic_cast<conv_test_type>(convTest_)->getTestValue();
1072 "Belos::BlockCGSolMgr::solve(): The convergence test's getTestValue() "
1073 "method returned NULL. Please report this bug to the Belos developers.");
1076 "Belos::BlockCGSolMgr::solve(): The convergence test's getTestValue() "
1077 "method returned a vector of length zero. Please report this bug to the "
1078 "Belos developers.");
1083 achievedTol_ = *std::max_element (pTestValues->begin(), pTestValues->end());
1093 template<
class ScalarType,
class MV,
class OP>
1096 std::ostringstream oss;
1099 oss <<
"Ortho Type='"<<orthoType_<<
"\', Block Size=" << blockSize_;
static const double orthoKappa
DGKS orthogonalization constant.
ScaleType convertStringToScaleType(const std::string &scaleType)
Convert the given string to its ScaleType enum value.
Teuchos::RCP< const MV > R
The current residual.
bool isLOADetected() const override
Return whether a loss of accuracy was detected by this solver during the most current solve...
Collection of types and exceptions used within the Belos solvers.
Teuchos::Array< Teuchos::RCP< Teuchos::Time > > getTimers() const
Return the timers for this object.
Belos's basic output manager for sending information of select verbosity levels to the appropriate ou...
int maxIters_
Maximum iteration count (read from parameter list).
Class which manages the output and verbosity of the Belos solvers.
bool is_null(const boost::shared_ptr< T > &p)
bool isSet_
Whether or not the parameters have been set (via setParameters()).
Teuchos::RCP< Teuchos::ParameterList > params_
Current parameter list.
Teuchos::RCP< StatusTestOutput< ScalarType, MV, OP > > outputTest_
Output "status test" that controls all the other status tests.
ScaleType
The type of scaling to use on the residual norm value.
Teuchos::RCP< SolverManager< ScalarType, MV, OP > > clone() const override
clone for Inverted Injection (DII)
Teuchos::ScalarTraits< ScalarType > SCT
This class implements the preconditioned Conjugate Gradient (CG) iteration.
T & get(ParameterList &l, const std::string &name)
Structure to contain pointers to CGIteration state variables.
Belos concrete class for performing the conjugate-gradient (CG) iteration.
bool is_null(const std::shared_ptr< T > &p)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
A factory class for generating StatusTestOutput objects.
Belos concrete class for performing the block conjugate-gradient (CG) iteration.
An implementation of StatusTestResNorm using a family of residual norms.
void setProblem(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem) override
Set the linear problem that needs to be solved.
std::string label_
Prefix label for all the timers.
The Belos::BlockCGSolMgr provides a powerful and fully-featured solver manager over the CG and BlockC...
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.
static std::string name()
BlockCGSolMgrLinearProblemFailure(const std::string &what_arg)
Teuchos::RCP< const Teuchos::ParameterList > getCurrentParameters() const override
Get a parameter list containing the current parameters for this object.
BlockCGSolMgrLinearProblemFailure is thrown when the linear problem is not setup (i.e.
A factory class for generating StatusTestOutput objects.
Teuchos::RCP< MatOrthoManager< ScalarType, MV, OP > > ortho_
Orthogonalization manager.
BlockCGSolMgr(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem, const Teuchos::RCP< Teuchos::ParameterList > &pl)
Teuchos::RCP< StatusTestMaxIters< ScalarType, MV, OP > > maxIterTest_
Maximum iteration count stopping criterion.
Traits class which defines basic operations on multivectors.
Belos::StatusTest for logically combining several status tests.
Teuchos::RCP< StatusTest< ScalarType, MV, OP > > sTest_
Aggregate stopping criterion.
bool assertPositiveDefiniteness_
Details::SolverManagerRequiresLapack< ScalarType, MV, OP, requiresLapack > base_type
bool isParameter(const std::string &name) const
virtual ~BlockCGSolMgr()
Destructor.
OperatorTraits< ScalarType, MV, OP > OPT
MultiVecTraits< ScalarType, MV > MVT
A Belos::StatusTest class for specifying a maximum number of iterations.
ResetType
How to reset the solver.
void reset(const ResetType type) override
Performs a reset of the solver manager specified by the ResetType. This informs the solver manager th...
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. ...
Belos concrete class for performing a single-reduction conjugate-gradient (CG) iteration.
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)
NormType convertStringToNormType(const std::string &normType)
Convert the given string to its NormType enum value.
A linear system to solve, and its associated information.
Class which describes the linear problem to be solved by the iterative solver.
MagnitudeType achievedTol_
Tolerance achieved by the last solve() invocation.
MagnitudeType achievedTol() const override
Tolerance achieved by the last solve() invocation.
Type traits class that says whether Teuchos::LAPACK has a valid implementation for the given ScalarTy...
Teuchos::ScalarTraits< ScalarType >::magnitudeType MagnitudeType
ReturnType
Whether the Belos solve converged for all linear systems.
Teuchos::RCP< StatusTestGenResNorm< ScalarType, MV, OP > > convTest_
Convergence stopping criterion.
int getNumIters() const override
Get the iteration count for the most recent call to solve().
Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > problem_
The linear problem to solve.
void validateParameters(ParameterList const &validParamList, int const depth=1000, EValidateUsed const validateUsed=VALIDATE_USED_ENABLED, EValidateDefaults const validateDefaults=VALIDATE_DEFAULTS_ENABLED) const
NormType
The type of vector norm to compute.
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.
MagnitudeType convtol_
Convergence tolerance (read from parameter list).
bool isType(const std::string &name) const
Teuchos::RCP< Teuchos::Time > timerSolve_
Solve timer.
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.
const LinearProblem< ScalarType, MV, OP > & getProblem() const override
Return a reference to the linear problem being solved by this solver manager.
Class which defines basic traits for the operator type.
bool foldConvergenceDetectionIntoAllreduce_
Parent class to all Belos exceptions.
Default parameters common to most Belos solvers.
This class implements the preconditioned single-reduction Conjugate Gradient (CG) iteration...
Teuchos::RCP< OutputManager< ScalarType > > printer_
Output manager, that handles printing of different kinds of messages.
Base class for Belos::SolverManager subclasses which normally can only compile with ScalarType types ...
Belos header file which uses auto-configuration information to include necessary C++ headers...
int numIters_
Number of iterations taken by the last solve() invocation.
Teuchos::ScalarTraits< MagnitudeType > MT
static const bool requiresLapack
Teuchos::RCP< std::ostream > outputStream_
Output stream to which the output manager prints.
MagnitudeType orthoKappa_
Orthogonalization parameter (read from parameter list).
Stub implementation of BlockCGIter, for ScalarType types for which Teuchos::LAPACK does NOT have a va...
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.