42 #ifndef BELOS_GMRES_POLY_SOLMGR_HPP 43 #define BELOS_GMRES_POLY_SOLMGR_HPP 67 #include "Teuchos_BLAS.hpp" 68 #include "Teuchos_as.hpp" 69 #ifdef BELOS_TEUCHOS_TIME_MONITOR 70 #include "Teuchos_TimeMonitor.hpp" 154 template<
class ScalarType,
class MV,
class OP>
159 typedef Teuchos::ScalarTraits<ScalarType> STS;
160 typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType MagnitudeType;
161 typedef Teuchos::ScalarTraits<MagnitudeType> MT;
189 const Teuchos::RCP<Teuchos::ParameterList> &pl );
206 Teuchos::RCP<const Teuchos::ParameterList> getValidParameters()
const;
217 Teuchos::Array<Teuchos::RCP<Teuchos::Time> >
getTimers()
const {
218 return Teuchos::tuple(timerSolve_, timerPoly_);
240 void setParameters(
const Teuchos::RCP<Teuchos::ParameterList> ¶ms );
256 problem_->setProblem ();
257 isPolyBuilt_ =
false;
290 std::string description()
const;
297 bool checkStatusTest();
303 Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > problem_;
306 Teuchos::RCP<OutputManager<ScalarType> > printer_;
307 Teuchos::RCP<std::ostream> outputStream_;
310 Teuchos::RCP<StatusTest<ScalarType,MV,OP> > sTest_;
311 Teuchos::RCP<StatusTestMaxIters<ScalarType,MV,OP> > maxIterTest_;
312 Teuchos::RCP<StatusTest<ScalarType,MV,OP> > convTest_;
313 Teuchos::RCP<StatusTestResNorm<ScalarType,MV,OP> > expConvTest_, impConvTest_;
314 Teuchos::RCP<StatusTestOutput<ScalarType,MV,OP> > outputTest_;
317 Teuchos::RCP<MatOrthoManager<ScalarType,MV,OP> > ortho_;
320 Teuchos::RCP<Teuchos::ParameterList> params_;
323 static const MagnitudeType polytol_default_;
324 static const MagnitudeType convtol_default_;
325 static const MagnitudeType orthoKappa_default_;
326 static const int maxDegree_default_;
327 static const int maxRestarts_default_;
328 static const int maxIters_default_;
329 static const bool strictConvTol_default_;
330 static const bool showMaxResNormOnly_default_;
331 static const int blockSize_default_;
332 static const int numBlocks_default_;
333 static const int verbosity_default_;
334 static const int outputStyle_default_;
335 static const int outputFreq_default_;
336 static const std::string impResScale_default_;
337 static const std::string expResScale_default_;
338 static const std::string label_default_;
339 static const std::string orthoType_default_;
340 static const Teuchos::RCP<std::ostream> outputStream_default_;
343 MagnitudeType polytol_, convtol_, orthoKappa_;
344 int maxDegree_, maxRestarts_, maxIters_, numIters_;
345 int blockSize_, numBlocks_, verbosity_, outputStyle_, outputFreq_;
346 bool strictConvTol_, showMaxResNormOnly_;
347 std::string orthoType_;
348 std::string impResScale_, expResScale_;
352 Teuchos::RCP<Teuchos::SerialDenseMatrix<int, ScalarType> > poly_H_, poly_y_;
353 Teuchos::RCP<Teuchos::SerialDenseVector<int, ScalarType> > poly_r0_;
354 Teuchos::RCP<Belos::GmresPolyOp<ScalarType, MV, OP> > poly_Op_;
358 Teuchos::RCP<Teuchos::Time> timerSolve_, timerPoly_;
362 bool isSet_, isSTSet_, expResTest_;
366 mutable Teuchos::RCP<const Teuchos::ParameterList> validPL_;
371 template<
class ScalarType,
class MV,
class OP>
372 const typename GmresPolySolMgr<ScalarType,MV,OP>::MagnitudeType
375 template<
class ScalarType,
class MV,
class OP>
376 const typename GmresPolySolMgr<ScalarType,MV,OP>::MagnitudeType
379 template<
class ScalarType,
class MV,
class OP>
380 const typename GmresPolySolMgr<ScalarType,MV,OP>::MagnitudeType
382 -Teuchos::ScalarTraits<MagnitudeType>::one();
384 template<
class ScalarType,
class MV,
class OP>
387 template<
class ScalarType,
class MV,
class OP>
390 template<
class ScalarType,
class MV,
class OP>
393 template<
class ScalarType,
class MV,
class OP>
396 template<
class ScalarType,
class MV,
class OP>
399 template<
class ScalarType,
class MV,
class OP>
402 template<
class ScalarType,
class MV,
class OP>
405 template<
class ScalarType,
class MV,
class OP>
408 template<
class ScalarType,
class MV,
class OP>
411 template<
class ScalarType,
class MV,
class OP>
414 template<
class ScalarType,
class MV,
class OP>
417 template<
class ScalarType,
class MV,
class OP>
420 template<
class ScalarType,
class MV,
class OP>
423 template<
class ScalarType,
class MV,
class OP>
426 template<
class ScalarType,
class MV,
class OP>
427 const Teuchos::RCP<std::ostream>
431 template<
class ScalarType,
class MV,
class OP>
433 outputStream_ (outputStream_default_),
434 polytol_ (polytol_default_),
435 convtol_ (convtol_default_),
436 orthoKappa_ (orthoKappa_default_),
437 maxDegree_ (maxDegree_default_),
438 maxRestarts_ (maxRestarts_default_),
439 maxIters_ (maxIters_default_),
441 blockSize_ (blockSize_default_),
442 numBlocks_ (numBlocks_default_),
443 verbosity_ (verbosity_default_),
444 outputStyle_ (outputStyle_default_),
445 outputFreq_ (outputFreq_default_),
446 strictConvTol_ (strictConvTol_default_),
447 showMaxResNormOnly_ (showMaxResNormOnly_default_),
448 orthoType_ (orthoType_default_),
449 impResScale_ (impResScale_default_),
450 expResScale_ (expResScale_default_),
452 label_ (label_default_),
453 isPolyBuilt_ (false),
461 template<
class ScalarType,
class MV,
class OP>
464 const Teuchos::RCP<Teuchos::ParameterList> &pl) :
466 outputStream_ (outputStream_default_),
467 polytol_ (polytol_default_),
468 convtol_ (convtol_default_),
469 orthoKappa_ (orthoKappa_default_),
470 maxDegree_ (maxDegree_default_),
471 maxRestarts_ (maxRestarts_default_),
472 maxIters_ (maxIters_default_),
474 blockSize_ (blockSize_default_),
475 numBlocks_ (numBlocks_default_),
476 verbosity_ (verbosity_default_),
477 outputStyle_ (outputStyle_default_),
478 outputFreq_ (outputFreq_default_),
479 strictConvTol_ (strictConvTol_default_),
480 showMaxResNormOnly_ (showMaxResNormOnly_default_),
481 orthoType_ (orthoType_default_),
482 impResScale_ (impResScale_default_),
483 expResScale_ (expResScale_default_),
485 label_ (label_default_),
486 isPolyBuilt_ (false),
492 TEUCHOS_TEST_FOR_EXCEPTION(
493 problem_.is_null (), std::invalid_argument,
494 "Belos::GmresPolySolMgr: The given linear problem is null. " 495 "Please call this constructor with a nonnull LinearProblem argument, " 496 "or call the constructor that does not take a LinearProblem.");
500 if (! pl.is_null ()) {
506 template<
class ScalarType,
class MV,
class OP>
507 Teuchos::RCP<const Teuchos::ParameterList>
510 if (validPL_.is_null ()) {
511 Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList ();
512 pl->set(
"Polynomial Tolerance", polytol_default_,
513 "The relative residual tolerance that used to construct the GMRES polynomial.");
514 pl->set(
"Maximum Degree", maxDegree_default_,
515 "The maximum degree allowed for any GMRES polynomial.");
516 pl->set(
"Convergence Tolerance", convtol_default_,
517 "The relative residual tolerance that needs to be achieved by the\n" 518 "iterative solver in order for the linear system to be declared converged." );
519 pl->set(
"Maximum Restarts", maxRestarts_default_,
520 "The maximum number of restarts allowed for each\n" 521 "set of RHS solved.");
522 pl->set(
"Maximum Iterations", maxIters_default_,
523 "The maximum number of block iterations allowed for each\n" 524 "set of RHS solved.");
525 pl->set(
"Num Blocks", numBlocks_default_,
526 "The maximum number of blocks allowed in the Krylov subspace\n" 527 "for each set of RHS solved.");
528 pl->set(
"Block Size", blockSize_default_,
529 "The number of vectors in each block. This number times the\n" 530 "number of blocks is the total Krylov subspace dimension.");
531 pl->set(
"Verbosity", verbosity_default_,
532 "What type(s) of solver information should be outputted\n" 533 "to the output stream.");
534 pl->set(
"Output Style", outputStyle_default_,
535 "What style is used for the solver information outputted\n" 536 "to the output stream.");
537 pl->set(
"Output Frequency", outputFreq_default_,
538 "How often convergence information should be outputted\n" 539 "to the output stream.");
540 pl->set(
"Output Stream", outputStream_default_,
541 "A reference-counted pointer to the output stream where all\n" 542 "solver output is sent.");
543 pl->set(
"Strict Convergence", strictConvTol_default_,
544 "After polynomial is applied, whether solver should try to achieve\n" 545 "the relative residual tolerance.");
546 pl->set(
"Show Maximum Residual Norm Only", showMaxResNormOnly_default_,
547 "When convergence information is printed, only show the maximum\n" 548 "relative residual norm when the block size is greater than one.");
549 pl->set(
"Implicit Residual Scaling", impResScale_default_,
550 "The type of scaling used in the implicit residual convergence test.");
551 pl->set(
"Explicit Residual Scaling", expResScale_default_,
552 "The type of scaling used in the explicit residual convergence test.");
553 pl->set(
"Timer Label", label_default_,
554 "The string to use as a prefix for the timer labels.");
555 pl->set(
"Orthogonalization", orthoType_default_,
556 "The type of orthogonalization to use: DGKS, ICGS, or IMGS.");
557 pl->set(
"Orthogonalization Constant",orthoKappa_default_,
558 "The constant used by DGKS orthogonalization to determine\n" 559 "whether another step of classical Gram-Schmidt is necessary.");
566 template<
class ScalarType,
class MV,
class OP>
571 if (params_.is_null ()) {
579 if (params->isParameter(
"Maximum Degree")) {
580 maxDegree_ = params->get(
"Maximum Degree",maxDegree_default_);
583 params_->set(
"Maximum Degree", maxDegree_);
587 if (params->isParameter(
"Maximum Restarts")) {
588 maxRestarts_ = params->get(
"Maximum Restarts",maxRestarts_default_);
591 params_->set(
"Maximum Restarts", maxRestarts_);
595 if (params->isParameter(
"Maximum Iterations")) {
596 maxIters_ = params->get(
"Maximum Iterations",maxIters_default_);
599 params_->set(
"Maximum Iterations", maxIters_);
600 if (maxIterTest_!=Teuchos::null)
601 maxIterTest_->setMaxIters( maxIters_ );
605 if (params->isParameter(
"Block Size")) {
606 blockSize_ = params->get(
"Block Size",blockSize_default_);
607 TEUCHOS_TEST_FOR_EXCEPTION(blockSize_ <= 0, std::invalid_argument,
608 "Belos::GmresPolySolMgr: \"Block Size\" must be strictly positive.");
611 params_->set(
"Block Size", blockSize_);
615 if (params->isParameter(
"Num Blocks")) {
616 numBlocks_ = params->get(
"Num Blocks",numBlocks_default_);
617 TEUCHOS_TEST_FOR_EXCEPTION(numBlocks_ <= 0, std::invalid_argument,
618 "Belos::GmresPolySolMgr: \"Num Blocks\" must be strictly positive.");
621 params_->set(
"Num Blocks", numBlocks_);
625 if (params->isParameter(
"Timer Label")) {
626 std::string tempLabel = params->get(
"Timer Label", label_default_);
629 if (tempLabel != label_) {
631 params_->set(
"Timer Label", label_);
632 std::string solveLabel = label_ +
": GmresPolySolMgr total solve time";
633 #ifdef BELOS_TEUCHOS_TIME_MONITOR 634 timerSolve_ = Teuchos::TimeMonitor::getNewCounter(solveLabel);
636 std::string polyLabel = label_ +
": GmresPolySolMgr polynomial creation time";
637 #ifdef BELOS_TEUCHOS_TIME_MONITOR 638 timerPoly_ = Teuchos::TimeMonitor::getNewCounter(polyLabel);
640 if (ortho_ != Teuchos::null) {
641 ortho_->setLabel( label_ );
647 if (params->isParameter(
"Orthogonalization")) {
648 std::string tempOrthoType = params->get(
"Orthogonalization",orthoType_default_);
649 TEUCHOS_TEST_FOR_EXCEPTION( tempOrthoType !=
"DGKS" && tempOrthoType !=
"ICGS" && tempOrthoType !=
"IMGS",
650 std::invalid_argument,
651 "Belos::GmresPolySolMgr: \"Orthogonalization\" must be either \"DGKS\", \"ICGS\", or \"IMGS\".");
652 if (tempOrthoType != orthoType_) {
653 params_->set(
"Orthogonalization", orthoType_);
654 orthoType_ = tempOrthoType;
656 if (orthoType_==
"DGKS") {
657 if (orthoKappa_ <= 0) {
665 else if (orthoType_==
"ICGS") {
668 else if (orthoType_==
"IMGS") {
675 if (params->isParameter(
"Orthogonalization Constant")) {
676 orthoKappa_ = params->get(
"Orthogonalization Constant",orthoKappa_default_);
679 params_->set(
"Orthogonalization Constant",orthoKappa_);
680 if (orthoType_==
"DGKS") {
681 if (orthoKappa_ > 0 && ortho_ != Teuchos::null) {
688 if (params->isParameter(
"Verbosity")) {
689 if (Teuchos::isParameterType<int>(*params,
"Verbosity")) {
690 verbosity_ = params->get(
"Verbosity", verbosity_default_);
692 verbosity_ = (int)Teuchos::getParameter<Belos::MsgType>(*params,
"Verbosity");
696 params_->set(
"Verbosity", verbosity_);
697 if (printer_ != Teuchos::null)
698 printer_->setVerbosity(verbosity_);
702 if (params->isParameter(
"Output Style")) {
703 if (Teuchos::isParameterType<int>(*params,
"Output Style")) {
704 outputStyle_ = params->get(
"Output Style", outputStyle_default_);
706 outputStyle_ = (int)Teuchos::getParameter<Belos::OutputType>(*params,
"Output Style");
710 params_->set(
"Output Style", outputStyle_);
711 if (outputTest_ != Teuchos::null) {
717 if (params->isParameter(
"Output Stream")) {
718 outputStream_ = Teuchos::getParameter<Teuchos::RCP<std::ostream> >(*params,
"Output Stream");
721 params_->set(
"Output Stream", outputStream_);
722 if (printer_ != Teuchos::null)
723 printer_->setOStream( outputStream_ );
728 if (params->isParameter(
"Output Frequency")) {
729 outputFreq_ = params->get(
"Output Frequency", outputFreq_default_);
733 params_->set(
"Output Frequency", outputFreq_);
734 if (outputTest_ != Teuchos::null)
735 outputTest_->setOutputFrequency( outputFreq_ );
739 if (printer_ == Teuchos::null) {
748 if (params->isParameter(
"Polynomial Tolerance")) {
749 polytol_ = params->get(
"Polynomial Tolerance",polytol_default_);
752 params_->set(
"Polynomial Tolerance", polytol_);
756 if (params->isParameter(
"Convergence Tolerance")) {
757 convtol_ = params->get(
"Convergence Tolerance",convtol_default_);
760 params_->set(
"Convergence Tolerance", convtol_);
761 if (impConvTest_ != Teuchos::null)
762 impConvTest_->setTolerance( convtol_ );
763 if (expConvTest_ != Teuchos::null)
764 expConvTest_->setTolerance( convtol_ );
768 if (params->isParameter(
"Strict Convergence")) {
769 strictConvTol_ = params->get(
"Strict Convergence",strictConvTol_default_);
772 params_->set(
"Strict Convergence", strictConvTol_);
776 if (params->isParameter(
"Implicit Residual Scaling")) {
777 std::string tempImpResScale = Teuchos::getParameter<std::string>( *params,
"Implicit Residual Scaling" );
780 if (impResScale_ != tempImpResScale) {
782 impResScale_ = tempImpResScale;
785 params_->set(
"Implicit Residual Scaling", impResScale_);
786 if (impConvTest_ != Teuchos::null) {
790 catch (std::exception& e) {
798 if (params->isParameter(
"Explicit Residual Scaling")) {
799 std::string tempExpResScale = Teuchos::getParameter<std::string>( *params,
"Explicit Residual Scaling" );
802 if (expResScale_ != tempExpResScale) {
804 expResScale_ = tempExpResScale;
807 params_->set(
"Explicit Residual Scaling", expResScale_);
808 if (expConvTest_ != Teuchos::null) {
812 catch (std::exception& e) {
821 if (params->isParameter(
"Show Maximum Residual Norm Only")) {
822 showMaxResNormOnly_ = Teuchos::getParameter<bool>(*params,
"Show Maximum Residual Norm Only");
825 params_->set(
"Show Maximum Residual Norm Only", showMaxResNormOnly_);
826 if (impConvTest_ != Teuchos::null)
827 impConvTest_->setShowMaxResNormOnly( showMaxResNormOnly_ );
828 if (expConvTest_ != Teuchos::null)
829 expConvTest_->setShowMaxResNormOnly( showMaxResNormOnly_ );
833 if (ortho_ == Teuchos::null) {
834 params_->set(
"Orthogonalization", orthoType_);
835 if (orthoType_==
"DGKS") {
836 if (orthoKappa_ <= 0) {
844 else if (orthoType_==
"ICGS") {
847 else if (orthoType_==
"IMGS") {
851 TEUCHOS_TEST_FOR_EXCEPTION(orthoType_!=
"ICGS"&&orthoType_!=
"DGKS"&&orthoType_!=
"IMGS",std::logic_error,
852 "Belos::GmresPolySolMgr(): Invalid orthogonalization type.");
857 if (timerSolve_ == Teuchos::null) {
858 std::string solveLabel = label_ +
": GmresPolySolMgr total solve time";
859 #ifdef BELOS_TEUCHOS_TIME_MONITOR 860 timerSolve_ = Teuchos::TimeMonitor::getNewCounter(solveLabel);
864 if (timerPoly_ == Teuchos::null) {
865 std::string polyLabel = label_ +
": GmresPolySolMgr polynomial creation time";
866 #ifdef BELOS_TEUCHOS_TIME_MONITOR 867 timerPoly_ = Teuchos::TimeMonitor::getNewCounter(polyLabel);
876 template<
class ScalarType,
class MV,
class OP>
888 if (!Teuchos::is_null(problem_->getLeftPrec())) {
895 Teuchos::RCP<StatusTestGenResNorm_t> tmpImpConvTest =
896 Teuchos::rcp(
new StatusTestGenResNorm_t( convtol_ ) );
898 tmpImpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
899 impConvTest_ = tmpImpConvTest;
902 Teuchos::RCP<StatusTestGenResNorm_t> tmpExpConvTest =
903 Teuchos::rcp(
new StatusTestGenResNorm_t( convtol_ ) );
904 tmpExpConvTest->defineResForm( StatusTestGenResNorm_t::Explicit,
Belos::TwoNorm );
906 tmpExpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
907 expConvTest_ = tmpExpConvTest;
910 convTest_ = Teuchos::rcp(
new StatusTestCombo_t( StatusTestCombo_t::SEQ, impConvTest_, expConvTest_ ) );
916 Teuchos::RCP<StatusTestImpResNorm_t> tmpImpConvTest =
917 Teuchos::rcp(
new StatusTestImpResNorm_t( convtol_ ) );
919 tmpImpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
920 impConvTest_ = tmpImpConvTest;
923 expConvTest_ = impConvTest_;
924 convTest_ = impConvTest_;
927 sTest_ = Teuchos::rcp(
new StatusTestCombo_t( StatusTestCombo_t::OR, maxIterTest_, convTest_ ) );
935 std::string solverDesc =
" Gmres Polynomial ";
936 outputTest_->setSolverDesc( solverDesc );
945 template<
class ScalarType,
class MV,
class OP>
949 Teuchos::RCP<MV> newX =
MVT::Clone( *(problem_->getLHS()), 1 );
950 Teuchos::RCP<MV> newB =
MVT::Clone( *(problem_->getRHS()), 1 );
953 Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > newProblem =
955 newProblem->setLeftPrec( problem_->getLeftPrec() );
956 newProblem->setRightPrec( problem_->getRightPrec() );
957 newProblem->setLabel(
"Belos GMRES Poly Generation");
958 newProblem->setProblem();
959 std::vector<int> idx(1,0);
960 newProblem->setLSIndex( idx );
963 Teuchos::ParameterList polyList;
966 polyList.set(
"Num Blocks",maxDegree_);
967 polyList.set(
"Block Size",1);
968 polyList.set(
"Keep Hessenberg",
true);
971 Teuchos::RCP<StatusTestMaxIters<ScalarType,MV,OP> > maxItrTst =
975 Teuchos::RCP<StatusTestGenResNorm<ScalarType,MV,OP> > convTst =
980 Teuchos::RCP<StatusTestCombo<ScalarType,MV,OP> > polyTest =
984 Teuchos::RCP<BlockGmresIter<ScalarType,MV,OP> > gmres_iter;
988 Teuchos::RCP<MV> V_0 =
MVT::Clone( *(newProblem->getRHS()), 1 );
989 newProblem->computeCurrPrecResVec( &*V_0 );
992 poly_r0_ = Teuchos::rcp(
new Teuchos::SerialDenseVector<int,ScalarType>(1) );
995 int rank = ortho_->normalize( *V_0, poly_r0_ );
997 "Belos::GmresPolySolMgr::generatePoly(): Failed to compute initial block of orthonormal vectors for polynomial generation.");
1002 newstate.
z = poly_r0_;
1004 gmres_iter->initializeGmres(newstate);
1007 bool polyConverged =
false;
1009 gmres_iter->iterate();
1012 if ( convTst->getStatus() ==
Passed ) {
1014 polyConverged =
true;
1019 gmres_iter->updateLSQR( gmres_iter->getCurSubspaceDim() );
1022 polyTest->checkStatus( &*gmres_iter );
1023 if (convTst->getStatus() ==
Passed)
1024 polyConverged =
true;
1026 catch (std::exception e) {
1027 printer_->stream(
Errors) <<
"Error! Caught exception in BlockGmresIter::iterate() at iteration " 1028 << gmres_iter->getNumIters() << std::endl
1029 << e.what() << std::endl;
1036 (void) polyConverged;
1039 Teuchos::RCP<MV> currX = gmres_iter->getCurrentUpdate();
1045 poly_dim_ = gmresState.
curDim;
1046 if (poly_dim_ == 0) {
1052 poly_y_ = Teuchos::rcp(
new Teuchos::SerialDenseMatrix<int,ScalarType>( Teuchos::Copy, *gmresState.
z, poly_dim_, 1 ) );
1053 poly_H_ = Teuchos::rcp(
new Teuchos::SerialDenseMatrix<int,ScalarType>( *gmresState.
H ) );
1057 const ScalarType one = STS::one ();
1058 Teuchos::BLAS<int,ScalarType> blas;
1059 blas.TRSM( Teuchos::LEFT_SIDE, Teuchos::UPPER_TRI, Teuchos::NO_TRANS,
1060 Teuchos::NON_UNIT_DIAG, poly_dim_, 1, one,
1061 gmresState.
R->values(), gmresState.
R->stride(),
1062 poly_y_->values(), poly_y_->stride() );
1066 poly_Op_ = Teuchos::rcp(
1073 template<
class ScalarType,
class MV,
class OP>
1078 typedef Teuchos::SerialDenseMatrix<int, ScalarType> SDM;
1088 TEUCHOS_TEST_FOR_EXCEPTION(
1090 "Belos::GmresPolySolMgr::solve: The linear problem has not been set yet, " 1091 "or was set to null. Please call setProblem() with a nonnull input before " 1092 "calling solve().");
1094 TEUCHOS_TEST_FOR_EXCEPTION(
1096 "Belos::GmresPolySolMgr::solve: The linear problem is not ready. Please " 1097 "call setProblem() on the LinearProblem object before calling solve().");
1099 if (! isSTSet_ || (! expResTest_ && ! problem_->getLeftPrec ().is_null ())) {
1104 const bool check = checkStatusTest ();
1105 TEUCHOS_TEST_FOR_EXCEPTION(
1107 "Belos::GmresPolySolMgr::solve: Linear problem and requested status " 1108 "tests are incompatible.");
1113 if (! isPolyBuilt_) {
1114 #ifdef BELOS_TEUCHOS_TIME_MONITOR 1115 Teuchos::TimeMonitor slvtimer(*timerPoly_);
1117 isPolyBuilt_ = generatePoly();
1119 "Belos::GmresPolySolMgr::generatePoly(): Failed to generate a non-trivial polynomial, reduce polynomial tolerance.");
1121 "Belos::GmresPolySolMgr::generatePoly(): Failed to generate polynomial that satisfied requirements.");
1125 bool isConverged =
true;
1129 #ifdef BELOS_TEUCHOS_TIME_MONITOR 1130 Teuchos::TimeMonitor slvtimer(*timerSolve_);
1134 poly_Op_->Apply( *problem_->getRHS(), *problem_->getLHS() );
1137 problem_->setProblem ();
1140 if (strictConvTol_) {
1145 int numCurrRHS = ( numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;
1153 std::vector<int> currIdx (blockSize_);
1154 for (
int i = 0; i < numCurrRHS; ++i) {
1155 currIdx[i] = startPtr+i;
1157 for (
int i = numCurrRHS; i < blockSize_; ++i) {
1162 problem_->setLSIndex (currIdx);
1166 Teuchos::ParameterList plist;
1167 plist.set (
"Block Size", blockSize_);
1170 if (blockSize_*static_cast<ptrdiff_t>(numBlocks_) > dim) {
1171 ptrdiff_t tmpNumBlocks = 0;
1172 if (blockSize_ == 1) {
1173 tmpNumBlocks = dim / blockSize_;
1175 tmpNumBlocks = ( dim - blockSize_) / blockSize_;
1178 <<
"Warning! Requested Krylov subspace dimension is larger than " 1179 <<
"operator dimension!" << std::endl <<
"The maximum number of " 1180 <<
"blocks allowed for the Krylov subspace will be adjusted to " 1181 << tmpNumBlocks << std::endl;
1182 plist.set (
"Num Blocks", Teuchos::asSafe<int> (tmpNumBlocks));
1185 plist.set (
"Num Blocks", numBlocks_);
1189 outputTest_->reset ();
1190 loaDetected_ =
false;
1200 RCP<GmresIteration<ScalarType,MV,OP> > block_gmres_iter =
1202 outputTest_, ortho_, plist));
1205 while (numRHS2Solve > 0) {
1208 if (blockSize_*numBlocks_ > dim) {
1209 int tmpNumBlocks = 0;
1210 if (blockSize_ == 1) {
1212 tmpNumBlocks = dim / blockSize_;
1215 tmpNumBlocks = (dim - blockSize_) / blockSize_;
1217 block_gmres_iter->setSize (blockSize_, tmpNumBlocks);
1220 block_gmres_iter->setSize (blockSize_, numBlocks_);
1224 block_gmres_iter->resetNumIters ();
1227 outputTest_->resetNumCalls ();
1230 RCP<MV> V_0 =
MVT::CloneCopy (*(problem_->getInitPrecResVec ()), currIdx);
1233 RCP<SDM> z_0 = rcp (
new SDM (blockSize_, blockSize_));
1236 int rank = ortho_->normalize (*V_0, z_0);
1237 TEUCHOS_TEST_FOR_EXCEPTION(
1239 "Belos::GmresPolySolMgr::solve: Failed to compute initial block of " 1240 "orthonormal vectors.");
1247 block_gmres_iter->initializeGmres(newstate);
1248 int numRestarts = 0;
1252 block_gmres_iter->iterate();
1255 if ( convTest_->getStatus() ==
Passed ) {
1256 if ( expConvTest_->getLOADetected() ) {
1258 loaDetected_ =
true;
1259 isConverged =
false;
1265 else if ( maxIterTest_->getStatus() ==
Passed ) {
1267 isConverged =
false;
1271 else if (block_gmres_iter->getCurSubspaceDim () ==
1272 block_gmres_iter->getMaxSubspaceDim ()) {
1273 if (numRestarts >= maxRestarts_) {
1274 isConverged =
false;
1279 printer_->stream(
Debug)
1280 <<
" Performing restart number " << numRestarts <<
" of " 1281 << maxRestarts_ << std::endl << std::endl;
1284 RCP<MV> update = block_gmres_iter->getCurrentUpdate();
1285 problem_->updateSolution (update,
true);
1295 RCP<MV> VV_0 =
MVT::Clone (*(oldState.
V), blockSize_);
1296 problem_->computeCurrPrecResVec (&*VV_0);
1302 RCP<SDM> zz_0 = rcp (
new SDM (blockSize_, blockSize_));
1305 const int theRank = ortho_->normalize( *VV_0, zz_0 );
1306 TEUCHOS_TEST_FOR_EXCEPTION(
1308 "Belos::GmresPolySolMgr::solve: Failed to compute initial " 1309 "block of orthonormal vectors after restart.");
1313 theNewState.
V = VV_0;
1314 theNewState.
z = zz_0;
1316 block_gmres_iter->initializeGmres (theNewState);
1324 TEUCHOS_TEST_FOR_EXCEPTION(
1325 true, std::logic_error,
1326 "Belos::GmresPolySolMgr::solve: Invalid return from " 1327 "BlockGmresIter::iterate(). Please report this bug " 1328 "to the Belos developers.");
1333 if (blockSize_ != 1) {
1335 <<
"Error! Caught std::exception in BlockGmresIter::iterate() " 1336 <<
"at iteration " << block_gmres_iter->getNumIters()
1337 << std::endl << e.what() << std::endl;
1338 if (convTest_->getStatus() !=
Passed) {
1339 isConverged =
false;
1346 block_gmres_iter->updateLSQR (block_gmres_iter->getCurSubspaceDim ());
1350 sTest_->checkStatus (&*block_gmres_iter);
1351 if (convTest_->getStatus() !=
Passed) {
1352 isConverged =
false;
1357 catch (
const std::exception &e) {
1359 <<
"Error! Caught std::exception in BlockGmresIter::iterate() " 1360 <<
"at iteration " << block_gmres_iter->getNumIters() << std::endl
1361 << e.what() << std::endl;
1369 if (! Teuchos::is_null (expConvTest_->getSolution ()) ) {
1370 RCP<MV> newX = expConvTest_->getSolution ();
1371 RCP<MV> curX = problem_->getCurrLHSVec ();
1372 MVT::MvAddMv (STS::zero (), *newX, STS::one (), *newX, *curX);
1375 RCP<MV> update = block_gmres_iter->getCurrentUpdate ();
1376 problem_->updateSolution (update,
true);
1380 problem_->setCurrLS ();
1383 startPtr += numCurrRHS;
1384 numRHS2Solve -= numCurrRHS;
1385 if (numRHS2Solve > 0) {
1386 numCurrRHS = (numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;
1388 currIdx.resize (blockSize_);
1389 for (
int i=0; i<numCurrRHS; ++i) {
1390 currIdx[i] = startPtr+i;
1392 for (
int i=numCurrRHS; i<blockSize_; ++i) {
1397 problem_->setLSIndex( currIdx );
1400 currIdx.resize( numRHS2Solve );
1412 #ifdef BELOS_TEUCHOS_TIME_MONITOR 1417 Teuchos::TimeMonitor::summarize( printer_->stream(TimingDetails) );
1420 if (!isConverged || loaDetected_) {
1427 template<
class ScalarType,
class MV,
class OP>
1430 std::ostringstream out;
1432 out <<
"\"Belos::GmresPolySolMgr\": {" 1433 <<
"ScalarType: " << Teuchos::TypeNameTraits<ScalarType>::name ()
1434 <<
", Ortho Type: " << orthoType_
1435 <<
", Block Size: " << blockSize_
1436 <<
", Num Blocks: " << numBlocks_
1437 <<
", Max Restarts: " << maxRestarts_;
1444 #endif // BELOS_GMRES_POLY_SOLMGR_HPP ScaleType convertStringToScaleType(const std::string &scaleType)
Convert the given string to its ScaleType enum value.
Collection of types and exceptions used within the Belos solvers.
Belos's basic output manager for sending information of select verbosity levels to the appropriate ou...
Class which manages the output and verbosity of the Belos solvers.
virtual ~GmresPolySolMgr()
Destructor.
Defines the GMRES polynomial operator hybrid-GMRES iterative linear solver.
ScaleType
The type of scaling to use on the residual norm value.
static void MvRandom(MV &mv)
Replace the vectors in mv with random vectors.
This class implements the block GMRES iteration, where a block Krylov subspace is constructed...
Convergence test using the implicit residual norm(s), with an explicit residual norm(s) check for los...
Teuchos::RCP< const MV > V
The current Krylov basis.
Teuchos::RCP< const Teuchos::ParameterList > getCurrentParameters() const
Get a parameter list containing the current parameters for this object.
A factory class for generating StatusTestOutput objects.
void setProblem(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem)
Set the linear problem that needs to be solved.
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > z
The current right-hand side of the least squares system RY = Z.
An implementation of StatusTestResNorm using a family of residual norms.
Structure to contain pointers to GmresIteration state variables.
Belos::StatusTest class for specifying a maximum number of iterations.
static void MvInit(MV &mv, const ScalarType alpha=Teuchos::ScalarTraits< ScalarType >::zero())
Replace each element of the vectors in mv with alpha.
static int GetNumberVecs(const MV &mv)
Obtain the number of vectors in mv.
A factory class for generating StatusTestOutput objects.
Iterated Modified Gram-Schmidt (IMGS) implementation of the Belos::OrthoManager class.
void reset(const ResetType type)
Reset the solver.
static void MvAddMv(const ScalarType alpha, const MV &A, const ScalarType beta, const MV &B, MV &mv)
Replace mv with .
Traits class which defines basic operations on multivectors.
Belos::StatusTest for logically combining several status tests.
std::string description() const
Method to return description of the hybrid block GMRES solver manager.
Classical Gram-Schmidt (with DGKS correction) implementation of the Belos::OrthoManager class...
A Belos::StatusTest class for specifying a maximum number of iterations.
static Teuchos::RCP< MV > Clone(const MV &mv, const int numvecs)
Creates a new empty MV containing numvecs columns.
ResetType
How to reset the solver.
GmresPolySolMgrLinearProblemFailure is thrown when the linear problem is not setup (i...
Pure virtual base class which describes the basic interface for a solver manager. ...
int curDim
The current dimension of the reduction.
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Get a parameter list containing the valid parameters for this object.
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > R
The current upper-triangular matrix from the QR reduction of H.
A linear system to solve, and its associated information.
Teuchos::Array< Teuchos::RCP< Teuchos::Time > > getTimers() const
Return the timers for this object.
Class which describes the linear problem to be solved by the iterative solver.
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > ¶ms)
Set the parameters the solver manager should use to solve the linear problem.
const LinearProblem< ScalarType, MV, OP > & getProblem() const
Get current linear problem being solved for in this object.
Hybrid block GMRES iterative linear solver.
static Teuchos::RCP< MV > CloneCopy(const MV &mv)
Creates a new MV and copies contents of mv into the new vector (deep copy).
ReturnType
Whether the Belos solve converged for all linear systems.
The Belos::SolverManager is a templated virtual base class that defines the basic interface that any ...
ReturnType solve()
This method performs possibly repeated calls to the underlying linear solver's iterate() routine unti...
Belos's class for applying the GMRES polynomial operator that is used by the hybrid-GMRES linear solv...
Belos concrete class for performing the block GMRES iteration.
Iterated Classical Gram-Schmidt (ICGS) implementation of the Belos::OrthoManager class.
An implementation of the Belos::MatOrthoManager that performs orthogonalization using multiple steps ...
An implementation of the Belos::MatOrthoManager that performs orthogonalization using multiple steps ...
static ptrdiff_t GetGlobalLength(const MV &mv)
Return the number of rows in the given multivector mv.
GmresPolySolMgrPolynomialFailure is thrown when their is a problem generating the GMRES polynomial fo...
GmresPolySolMgrPolynomialFailure(const std::string &what_arg)
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.
Belos::StatusTest for specifying an implicit residual norm stopping criteria that checks for loss of ...
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.
Class which defines basic traits for the operator type.
Parent class to all Belos exceptions.
GmresPolySolMgr()
Empty constructor for GmresPolySolMgr. This constructor takes no arguments and sets the default value...
Pure virtual base class which augments the basic interface for a Gmres linear solver iteration...
GmresIterationOrthoFailure is thrown when the GmresIteration object is unable to compute independent ...
int getNumIters() const
Get the iteration count for the most recent call to solve().
Belos header file which uses auto-configuration information to include necessary C++ headers...
bool isLOADetected() const
Return whether a loss of accuracy was detected by this solver during the most current solve...
GmresPolySolMgrOrthoFailure is thrown when the orthogonalization manager is unable to generate orthon...
GmresPolySolMgrOrthoFailure(const std::string &what_arg)
GmresPolySolMgrLinearProblemFailure(const std::string &what_arg)
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > H
The current Hessenberg matrix.