42 #ifndef BELOS_BLOCK_CG_SOLMGR_HPP
43 #define BELOS_BLOCK_CG_SOLMGR_HPP
66 #include "Teuchos_BLAS.hpp"
67 #include "Teuchos_LAPACK.hpp"
68 #ifdef BELOS_TEUCHOS_TIME_MONITOR
69 # include "Teuchos_TimeMonitor.hpp"
71 #if defined(HAVE_TEUCHOSCORE_CXX11)
72 # include <type_traits>
73 #endif // defined(HAVE_TEUCHOSCORE_CXX11)
118 template<
class ScalarType,
class MV,
class OP,
119 const bool lapackSupportsScalarType =
124 static const bool requiresLapack =
134 const Teuchos::RCP<Teuchos::ParameterList>& pl) :
144 template<
class ScalarType,
class MV,
class OP>
171 typedef Teuchos::ScalarTraits<ScalarType> SCT;
172 typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType MagnitudeType;
173 typedef Teuchos::ScalarTraits<MagnitudeType> MT;
224 const Teuchos::RCP<Teuchos::ParameterList> &pl );
230 Teuchos::RCP<SolverManager<ScalarType, MV, OP> >
clone ()
const override {
244 Teuchos::RCP<const Teuchos::ParameterList> getValidParameters()
const override;
255 Teuchos::Array<Teuchos::RCP<Teuchos::Time> >
getTimers()
const {
256 return Teuchos::tuple(timerSolve_);
286 void setParameters(
const Teuchos::RCP<Teuchos::ParameterList> ¶ms )
override;
327 std::string description()
const override;
334 Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > problem_;
337 Teuchos::RCP<OutputManager<ScalarType> > printer_;
339 Teuchos::RCP<std::ostream> outputStream_;
345 Teuchos::RCP<StatusTest<ScalarType,MV,OP> > sTest_;
348 Teuchos::RCP<StatusTestMaxIters<ScalarType,MV,OP> > maxIterTest_;
351 Teuchos::RCP<StatusTestGenResNorm<ScalarType,MV,OP> > convTest_;
354 Teuchos::RCP<StatusTestOutput<ScalarType,MV,OP> > outputTest_;
357 Teuchos::RCP<MatOrthoManager<ScalarType,MV,OP> > ortho_;
360 Teuchos::RCP<Teuchos::ParameterList> params_;
365 static constexpr
int maxIters_default_ = 1000;
366 static constexpr
bool adaptiveBlockSize_default_ =
true;
367 static constexpr
bool showMaxResNormOnly_default_ =
false;
368 static constexpr
bool useSingleReduction_default_ =
false;
369 static constexpr
int blockSize_default_ = 1;
372 static constexpr
int outputFreq_default_ = -1;
373 static constexpr
const char * resScale_default_ =
"Norm of Initial Residual";
374 static constexpr
const char * label_default_ =
"Belos";
375 static constexpr
const char * orthoType_default_ =
"DGKS";
376 static constexpr std::ostream * outputStream_default_ = &std::cout;
383 MagnitudeType convtol_;
386 MagnitudeType orthoKappa_;
393 MagnitudeType achievedTol_;
402 int blockSize_, verbosity_, outputStyle_, outputFreq_;
403 bool adaptiveBlockSize_, showMaxResNormOnly_, useSingleReduction_;
404 std::string orthoType_, resScale_;
410 Teuchos::RCP<Teuchos::Time> timerSolve_;
418 template<
class ScalarType,
class MV,
class OP>
420 outputStream_(Teuchos::rcp(outputStream_default_,false)),
423 achievedTol_(Teuchos::ScalarTraits<MagnitudeType>::zero()),
424 maxIters_(maxIters_default_),
426 blockSize_(blockSize_default_),
427 verbosity_(verbosity_default_),
428 outputStyle_(outputStyle_default_),
429 outputFreq_(outputFreq_default_),
430 adaptiveBlockSize_(adaptiveBlockSize_default_),
431 showMaxResNormOnly_(showMaxResNormOnly_default_),
432 useSingleReduction_(useSingleReduction_default_),
433 orthoType_(orthoType_default_),
434 resScale_(resScale_default_),
435 label_(label_default_),
441 template<
class ScalarType,
class MV,
class OP>
444 const Teuchos::RCP<Teuchos::ParameterList> &pl) :
446 outputStream_(Teuchos::rcp(outputStream_default_,false)),
449 achievedTol_(Teuchos::ScalarTraits<MagnitudeType>::zero()),
450 maxIters_(maxIters_default_),
452 blockSize_(blockSize_default_),
453 verbosity_(verbosity_default_),
454 outputStyle_(outputStyle_default_),
455 outputFreq_(outputFreq_default_),
456 adaptiveBlockSize_(adaptiveBlockSize_default_),
457 showMaxResNormOnly_(showMaxResNormOnly_default_),
458 useSingleReduction_(useSingleReduction_default_),
459 orthoType_(orthoType_default_),
460 resScale_(resScale_default_),
461 label_(label_default_),
464 TEUCHOS_TEST_FOR_EXCEPTION(problem_.is_null(), std::invalid_argument,
465 "BlockCGSolMgr's constructor requires a nonnull LinearProblem instance.");
470 if (! pl.is_null()) {
475 template<
class ScalarType,
class MV,
class OP>
481 if (params_ == Teuchos::null) {
482 params_ = Teuchos::rcp(
new Teuchos::ParameterList(*getValidParameters()) );
485 params->validateParameters(*getValidParameters());
489 if (params->isParameter(
"Maximum Iterations")) {
490 maxIters_ = params->get(
"Maximum Iterations",maxIters_default_);
493 params_->set(
"Maximum Iterations", maxIters_);
494 if (maxIterTest_!=Teuchos::null)
495 maxIterTest_->setMaxIters( maxIters_ );
499 if (params->isParameter(
"Block Size")) {
500 blockSize_ = params->get(
"Block Size",blockSize_default_);
501 TEUCHOS_TEST_FOR_EXCEPTION(blockSize_ <= 0, std::invalid_argument,
502 "Belos::BlockCGSolMgr: \"Block Size\" must be strictly positive.");
505 params_->set(
"Block Size", blockSize_);
509 if (params->isParameter(
"Adaptive Block Size")) {
510 adaptiveBlockSize_ = params->get(
"Adaptive Block Size",adaptiveBlockSize_default_);
513 params_->set(
"Adaptive Block Size", adaptiveBlockSize_);
517 if (params->isParameter(
"Use Single Reduction")) {
518 useSingleReduction_ = params->get(
"Use Single Reduction", useSingleReduction_default_);
522 if (params->isParameter(
"Timer Label")) {
523 std::string tempLabel = params->get(
"Timer Label", label_default_);
526 if (tempLabel != label_) {
528 params_->set(
"Timer Label", label_);
529 std::string solveLabel = label_ +
": BlockCGSolMgr total solve time";
530 #ifdef BELOS_TEUCHOS_TIME_MONITOR
531 timerSolve_ = Teuchos::TimeMonitor::getNewCounter(solveLabel);
533 if (ortho_ != Teuchos::null) {
534 ortho_->setLabel( label_ );
540 if (params->isParameter(
"Orthogonalization")) {
541 std::string tempOrthoType = params->get(
"Orthogonalization",orthoType_default_);
542 TEUCHOS_TEST_FOR_EXCEPTION( tempOrthoType !=
"DGKS" && tempOrthoType !=
"ICGS" && tempOrthoType !=
"IMGS",
543 std::invalid_argument,
544 "Belos::BlockCGSolMgr: \"Orthogonalization\" must be either \"DGKS\", \"ICGS\", or \"IMGS\".");
545 if (tempOrthoType != orthoType_) {
546 orthoType_ = tempOrthoType;
547 params_->set(
"Orthogonalization", orthoType_);
549 if (orthoType_==
"DGKS") {
550 if (orthoKappa_ <= 0) {
558 else if (orthoType_==
"ICGS") {
561 else if (orthoType_==
"IMGS") {
568 if (params->isParameter(
"Orthogonalization Constant")) {
569 if (params->isType<MagnitudeType> (
"Orthogonalization Constant")) {
570 orthoKappa_ = params->get (
"Orthogonalization Constant",
574 orthoKappa_ = params->get (
"Orthogonalization Constant",
579 params_->set(
"Orthogonalization Constant",orthoKappa_);
580 if (orthoType_==
"DGKS") {
581 if (orthoKappa_ > 0 && ortho_ != Teuchos::null) {
588 if (params->isParameter(
"Verbosity")) {
589 if (Teuchos::isParameterType<int>(*params,
"Verbosity")) {
590 verbosity_ = params->get(
"Verbosity", verbosity_default_);
592 verbosity_ = (int)Teuchos::getParameter<Belos::MsgType>(*params,
"Verbosity");
596 params_->set(
"Verbosity", verbosity_);
597 if (printer_ != Teuchos::null)
598 printer_->setVerbosity(verbosity_);
602 if (params->isParameter(
"Output Style")) {
603 if (Teuchos::isParameterType<int>(*params,
"Output Style")) {
604 outputStyle_ = params->get(
"Output Style", outputStyle_default_);
606 outputStyle_ = (int)Teuchos::getParameter<Belos::OutputType>(*params,
"Output Style");
610 params_->set(
"Output Style", outputStyle_);
611 outputTest_ = Teuchos::null;
615 if (params->isParameter(
"Output Stream")) {
616 outputStream_ = Teuchos::getParameter<Teuchos::RCP<std::ostream> >(*params,
"Output Stream");
619 params_->set(
"Output Stream", outputStream_);
620 if (printer_ != Teuchos::null)
621 printer_->setOStream( outputStream_ );
626 if (params->isParameter(
"Output Frequency")) {
627 outputFreq_ = params->get(
"Output Frequency", outputFreq_default_);
631 params_->set(
"Output Frequency", outputFreq_);
632 if (outputTest_ != Teuchos::null)
633 outputTest_->setOutputFrequency( outputFreq_ );
637 if (printer_ == Teuchos::null) {
646 if (params->isParameter(
"Convergence Tolerance")) {
647 if (params->isType<MagnitudeType> (
"Convergence Tolerance")) {
648 convtol_ = params->get (
"Convergence Tolerance",
656 params_->set(
"Convergence Tolerance", convtol_);
657 if (convTest_ != Teuchos::null)
661 if (params->isParameter(
"Show Maximum Residual Norm Only")) {
662 showMaxResNormOnly_ = Teuchos::getParameter<bool>(*params,
"Show Maximum Residual Norm Only");
665 params_->set(
"Show Maximum Residual Norm Only", showMaxResNormOnly_);
666 if (convTest_ != Teuchos::null)
667 convTest_->setShowMaxResNormOnly( showMaxResNormOnly_ );
671 bool newResTest =
false;
673 std::string tempResScale = resScale_;
674 if (params->isParameter (
"Implicit Residual Scaling")) {
675 tempResScale = params->get<std::string> (
"Implicit Residual Scaling");
679 if (resScale_ != tempResScale) {
682 resScale_ = tempResScale;
685 params_->set (
"Implicit Residual Scaling", resScale_);
687 if (! convTest_.is_null ()) {
691 catch (std::exception& e) {
702 if (maxIterTest_ == Teuchos::null)
706 if (convTest_.is_null () || newResTest) {
707 convTest_ = rcp (
new StatusTestResNorm_t (convtol_, 1, showMaxResNormOnly_));
711 if (sTest_.is_null () || newResTest)
712 sTest_ = Teuchos::rcp(
new StatusTestCombo_t( StatusTestCombo_t::OR, maxIterTest_, convTest_ ) );
714 if (outputTest_.is_null () || newResTest) {
722 std::string solverDesc =
" Block CG ";
723 outputTest_->setSolverDesc( solverDesc );
728 if (ortho_ == Teuchos::null) {
729 params_->set(
"Orthogonalization", orthoType_);
730 if (orthoType_==
"DGKS") {
731 if (orthoKappa_ <= 0) {
739 else if (orthoType_==
"ICGS") {
742 else if (orthoType_==
"IMGS") {
746 TEUCHOS_TEST_FOR_EXCEPTION(orthoType_!=
"ICGS"&&orthoType_!=
"DGKS"&&orthoType_!=
"IMGS",std::logic_error,
747 "Belos::BlockCGSolMgr(): Invalid orthogonalization type.");
752 if (timerSolve_ == Teuchos::null) {
753 std::string solveLabel = label_ +
": BlockCGSolMgr total solve time";
754 #ifdef BELOS_TEUCHOS_TIME_MONITOR
755 timerSolve_ = Teuchos::TimeMonitor::getNewCounter(solveLabel);
764 template<
class ScalarType,
class MV,
class OP>
765 Teuchos::RCP<const Teuchos::ParameterList>
768 static Teuchos::RCP<const Teuchos::ParameterList> validPL;
771 if(is_null(validPL)) {
772 Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
774 "The relative residual tolerance that needs to be achieved by the\n"
775 "iterative solver in order for the linear system to be declared converged.");
776 pl->set(
"Maximum Iterations", static_cast<int>(maxIters_default_),
777 "The maximum number of block iterations allowed for each\n"
778 "set of RHS solved.");
779 pl->set(
"Block Size", static_cast<int>(blockSize_default_),
780 "The number of vectors in each block.");
781 pl->set(
"Adaptive Block Size", static_cast<bool>(adaptiveBlockSize_default_),
782 "Whether the solver manager should adapt to the block size\n"
783 "based on the number of RHS to solve.");
784 pl->set(
"Verbosity", static_cast<int>(verbosity_default_),
785 "What type(s) of solver information should be outputted\n"
786 "to the output stream.");
787 pl->set(
"Output Style", static_cast<int>(outputStyle_default_),
788 "What style is used for the solver information outputted\n"
789 "to the output stream.");
790 pl->set(
"Output Frequency", static_cast<int>(outputFreq_default_),
791 "How often convergence information should be outputted\n"
792 "to the output stream.");
793 pl->set(
"Output Stream", Teuchos::rcp(outputStream_default_,
false),
794 "A reference-counted pointer to the output stream where all\n"
795 "solver output is sent.");
796 pl->set(
"Show Maximum Residual Norm Only", static_cast<bool>(showMaxResNormOnly_default_),
797 "When convergence information is printed, only show the maximum\n"
798 "relative residual norm when the block size is greater than one.");
799 pl->set(
"Use Single Reduction", static_cast<bool>(useSingleReduction_default_),
800 "Use single reduction iteration when the block size is one.");
801 pl->set(
"Implicit Residual Scaling", resScale_default_,
802 "The type of scaling used in the residual convergence test.");
803 pl->set(
"Timer Label", static_cast<const char *>(label_default_),
804 "The string to use as a prefix for the timer labels.");
805 pl->set(
"Orthogonalization", static_cast<const char *>(orthoType_default_),
806 "The type of orthogonalization to use: DGKS, ICGS, or IMGS.");
808 "The constant used by DGKS orthogonalization to determine\n"
809 "whether another step of classical Gram-Schmidt is necessary.");
817 template<
class ScalarType,
class MV,
class OP>
821 using Teuchos::rcp_const_cast;
822 using Teuchos::rcp_dynamic_cast;
829 setParameters(Teuchos::parameterList(*getValidParameters()));
832 Teuchos::BLAS<int,ScalarType> blas;
833 Teuchos::LAPACK<int,ScalarType> lapack;
835 TEUCHOS_TEST_FOR_EXCEPTION( !problem_->isProblemSet(),
837 "Belos::BlockCGSolMgr::solve(): Linear problem is not ready, setProblem() "
838 "has not been called.");
842 int numRHS2Solve = MVT::GetNumberVecs( *(problem_->getRHS()) );
843 int numCurrRHS = ( numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;
845 std::vector<int> currIdx, currIdx2;
850 if ( adaptiveBlockSize_ ) {
851 blockSize_ = numCurrRHS;
852 currIdx.resize( numCurrRHS );
853 currIdx2.resize( numCurrRHS );
854 for (
int i=0; i<numCurrRHS; ++i)
855 { currIdx[i] = startPtr+i; currIdx2[i]=i; }
859 currIdx.resize( blockSize_ );
860 currIdx2.resize( blockSize_ );
861 for (
int i=0; i<numCurrRHS; ++i)
862 { currIdx[i] = startPtr+i; currIdx2[i]=i; }
863 for (
int i=numCurrRHS; i<blockSize_; ++i)
864 { currIdx[i] = -1; currIdx2[i] = i; }
868 problem_->setLSIndex( currIdx );
872 Teuchos::ParameterList plist;
873 plist.set(
"Block Size",blockSize_);
876 outputTest_->reset();
880 bool isConverged =
true;
885 RCP<CGIteration<ScalarType,MV,OP> > block_cg_iter;
886 if (blockSize_ == 1) {
890 if (useSingleReduction_) {
893 outputTest_, plist));
898 outputTest_, plist));
909 #ifdef BELOS_TEUCHOS_TIME_MONITOR
910 Teuchos::TimeMonitor slvtimer(*timerSolve_);
913 while ( numRHS2Solve > 0 ) {
916 std::vector<int> convRHSIdx;
917 std::vector<int> currRHSIdx( currIdx );
918 currRHSIdx.resize(numCurrRHS);
921 block_cg_iter->resetNumIters();
924 outputTest_->resetNumCalls();
927 RCP<MV> R_0 = MVT::CloneViewNonConst( *(rcp_const_cast<MV>(problem_->getInitResVec())), currIdx );
932 block_cg_iter->initializeCG(newstate);
938 block_cg_iter->iterate();
942 if (convTest_->getStatus() ==
Passed) {
947 std::vector<int> convIdx =
948 rcp_dynamic_cast<conv_test_type>(convTest_)->
convIndices();
953 if (convIdx.size() == currRHSIdx.size())
958 problem_->setCurrLS();
963 std::vector<int> unconvIdx(currRHSIdx.size());
964 for (
unsigned int i=0; i<currRHSIdx.size(); ++i) {
966 for (
unsigned int j=0; j<convIdx.size(); ++j) {
967 if (currRHSIdx[i] == convIdx[j]) {
973 currIdx2[have] = currIdx2[i];
974 currRHSIdx[have++] = currRHSIdx[i];
979 currRHSIdx.resize(have);
980 currIdx2.resize(have);
983 problem_->setLSIndex( currRHSIdx );
986 std::vector<MagnitudeType> norms;
987 R_0 = MVT::CloneCopy( *(block_cg_iter->getNativeResiduals(&norms)),currIdx2 );
988 for (
int i=0; i<have; ++i) { currIdx2[i] = i; }
991 block_cg_iter->setBlockSize( have );
996 block_cg_iter->initializeCG(defstate);
1002 else if (maxIterTest_->getStatus() ==
Passed) {
1003 isConverged =
false;
1011 TEUCHOS_TEST_FOR_EXCEPTION(
true,std::logic_error,
1012 "Belos::BlockCGSolMgr::solve(): Neither the convergence test nor "
1013 "the maximum iteration count test passed. Please report this bug "
1014 "to the Belos developers.");
1017 catch (
const std::exception &e) {
1018 std::ostream& err = printer_->stream (
Errors);
1019 err <<
"Error! Caught std::exception in CGIteration::iterate() at "
1020 <<
"iteration " << block_cg_iter->getNumIters() << std::endl
1021 << e.what() << std::endl;
1028 problem_->setCurrLS();
1031 startPtr += numCurrRHS;
1032 numRHS2Solve -= numCurrRHS;
1033 if ( numRHS2Solve > 0 ) {
1034 numCurrRHS = ( numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;
1036 if ( adaptiveBlockSize_ ) {
1037 blockSize_ = numCurrRHS;
1038 currIdx.resize( numCurrRHS );
1039 currIdx2.resize( numCurrRHS );
1040 for (
int i=0; i<numCurrRHS; ++i)
1041 { currIdx[i] = startPtr+i; currIdx2[i] = i; }
1044 currIdx.resize( blockSize_ );
1045 currIdx2.resize( blockSize_ );
1046 for (
int i=0; i<numCurrRHS; ++i)
1047 { currIdx[i] = startPtr+i; currIdx2[i] = i; }
1048 for (
int i=numCurrRHS; i<blockSize_; ++i)
1049 { currIdx[i] = -1; currIdx2[i] = i; }
1052 problem_->setLSIndex( currIdx );
1055 block_cg_iter->setBlockSize( blockSize_ );
1058 currIdx.resize( numRHS2Solve );
1069 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1075 Teuchos::TimeMonitor::summarize( printer_->stream(
TimingDetails) );
1080 numIters_ = maxIterTest_->getNumIters();
1086 const std::vector<MagnitudeType>* pTestValues =
1087 rcp_dynamic_cast<conv_test_type>(convTest_)->
getTestValue();
1089 TEUCHOS_TEST_FOR_EXCEPTION(pTestValues == NULL, std::logic_error,
1090 "Belos::BlockCGSolMgr::solve(): The convergence test's getTestValue() "
1091 "method returned NULL. Please report this bug to the Belos developers.");
1093 TEUCHOS_TEST_FOR_EXCEPTION(pTestValues->size() < 1, std::logic_error,
1094 "Belos::BlockCGSolMgr::solve(): The convergence test's getTestValue() "
1095 "method returned a vector of length zero. Please report this bug to the "
1096 "Belos developers.");
1101 achievedTol_ = *std::max_element (pTestValues->begin(), pTestValues->end());
1111 template<
class ScalarType,
class MV,
class OP>
1114 std::ostringstream oss;
1115 oss <<
"Belos::BlockCGSolMgr<...,"<<Teuchos::ScalarTraits<ScalarType>::name()<<
">";
1117 oss <<
"Ortho Type='"<<orthoType_<<
"\', Block Size=" << blockSize_;