13#ifndef BELOS_BLOCK_GCRODR_SOLMGR_HPP
14#define BELOS_BLOCK_GCRODR_SOLMGR_HPP
30#include "Teuchos_LAPACK.hpp"
31#include "Teuchos_as.hpp"
32#ifdef BELOS_TEUCHOS_TIME_MONITOR
33#include "Teuchos_TimeMonitor.hpp"
93template<
class ScalarType,
class MV,
class OP>
99 typedef Teuchos::ScalarTraits<ScalarType> SCT;
100 typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType MagnitudeType;
101 typedef Teuchos::ScalarTraits<MagnitudeType> MT;
102 typedef Teuchos::ScalarTraits<MagnitudeType> SMT;
104 typedef Teuchos::SerialDenseMatrix<int,ScalarType> SDM;
105 typedef Teuchos::SerialDenseVector<int,ScalarType> SDV;
144 const Teuchos::RCP<Teuchos::ParameterList> &pl);
191 TEUCHOS_TEST_FOR_EXCEPTION(problem.is_null(), std::invalid_argument,
192 "Belos::BlockGCRODRSolMgr::setProblem: The input LinearProblem cannot be null.");
195 if (! problem->isProblemSet()) {
196 const bool success = problem->setProblem();
197 TEUCHOS_TEST_FOR_EXCEPTION(success, std::runtime_error,
198 "Belos::BlockGCRODRSolMgr::setProblem: Calling the input LinearProblem's setProblem() method failed. This likely means that the "
199 "LinearProblem has a missing (null) matrix A, solution vector X, or right-hand side vector B. Please set these items in the LinearProblem and try again.");
206 void setParameters(
const Teuchos::RCP<Teuchos::ParameterList> ¶ms );
221 problem_->setProblem();
259 void initializeStateStorage();
277 int getHarmonicVecsKryl (
int m,
const SDM& HH, SDM& PP);
282 int getHarmonicVecsAugKryl (
int keff,
285 const Teuchos::RCP<const MV>& VV,
289 void sort (std::vector<MagnitudeType>& dlist,
int n, std::vector<int>& iperm);
292 Teuchos::LAPACK<int,ScalarType> lapack;
295 Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > problem_;
298 Teuchos::RCP<OutputManager<ScalarType> > printer_;
299 Teuchos::RCP<std::ostream> outputStream_;
302 Teuchos::RCP<StatusTest<ScalarType,MV,OP> > sTest_;
303 Teuchos::RCP<StatusTestMaxIters<ScalarType,MV,OP> > maxIterTest_;
304 Teuchos::RCP<StatusTest<ScalarType,MV,OP> > convTest_;
305 Teuchos::RCP<StatusTestGenResNorm<ScalarType,MV,OP> > expConvTest_, impConvTest_;
306 Teuchos::RCP<StatusTestOutput<ScalarType,MV,OP> > outputTest_;
309 ortho_factory_type orthoFactory_;
316 Teuchos::RCP<MatOrthoManager<ScalarType,MV,OP> > ortho_;
319 Teuchos::RCP<Teuchos::ParameterList> params_;
331 mutable Teuchos::RCP<const Teuchos::ParameterList> defaultParams_;
338 MagnitudeType convTol_, orthoKappa_, achievedTol_;
339 int blockSize_, maxRestarts_, maxIters_, numIters_;
340 int verbosity_, outputStyle_, outputFreq_;
341 bool adaptiveBlockSize_;
342 std::string orthoType_, recycleMethod_;
343 std::string impResScale_, expResScale_;
351 int numBlocks_, recycledBlocks_;
362 Teuchos::RCP<MV> U_, C_;
365 Teuchos::RCP<MV> U1_, C1_;
368 Teuchos::RCP<SDM > G_;
369 Teuchos::RCP<SDM > H_;
370 Teuchos::RCP<SDM > B_;
371 Teuchos::RCP<SDM > PP_;
372 Teuchos::RCP<SDM > HP_;
373 std::vector<ScalarType> tau_;
374 std::vector<ScalarType> work_;
375 Teuchos::RCP<SDM > F_;
376 std::vector<int> ipiv_;
379 Teuchos::RCP<Teuchos::Time> timerSolve_;
388 bool builtRecycleSpace_;
394 template<
class ScalarType,
class MV,
class OP>
395 const bool BlockGCRODRSolMgr<ScalarType,MV,OP>::adaptiveBlockSize_default_ =
true;
397 template<
class ScalarType,
class MV,
class OP>
398 const std::string BlockGCRODRSolMgr<ScalarType,MV,OP>::recycleMethod_default_ =
"harmvecs";
404 template<
class ScalarType,
class MV,
class OP>
410 template<
class ScalarType,
class MV,
class OP>
413 const Teuchos::RCP<Teuchos::ParameterList> &pl ) {
418 TEUCHOS_TEST_FOR_EXCEPTION(problem.is_null(), std::invalid_argument,
419 "Belos::BlockGCRODR constructor: The solver manager's constructor needs "
420 "the linear problem argument 'problem' to be nonnull.");
430 template<
class ScalarType,
class MV,
class OP>
431 void BlockGCRODRSolMgr<ScalarType,MV,OP>::init() {
435 loaDetected_ =
false;
436 builtRecycleSpace_ =
false;
507 template<
class ScalarType,
class MV,
class OP>
509 std::ostringstream oss;
510 oss <<
"Belos::BlockGCRODRSolMgr<" << SCT::name() <<
", ...>";
512 oss <<
"Ortho Type='"<<orthoType_ ;
513 oss <<
", Num Blocks=" <<numBlocks_;
514 oss <<
", Block Size=" <<blockSize_;
515 oss <<
", Num Recycle Blocks=" << recycledBlocks_;
516 oss <<
", Max Restarts=" << maxRestarts_;
521 template<
class ScalarType,
class MV,
class OP>
522 Teuchos::RCP<const Teuchos::ParameterList>
524 using Teuchos::ParameterList;
525 using Teuchos::parameterList;
527 using Teuchos::rcpFromRef;
529 if (defaultParams_.is_null()) {
530 RCP<ParameterList> pl = parameterList ();
532 const MagnitudeType convTol = SMT::squareroot (SCT::magnitude (SCT::eps()));
533 const int maxRestarts = 1000;
534 const int maxIters = 5000;
535 const int blockSize = 2;
536 const int numBlocks = 100;
537 const int numRecycledBlocks = 25;
540 const int outputFreq = 1;
541 RCP<std::ostream> outputStream = rcpFromRef (std::cout);
542 const std::string impResScale (
"Norm of Preconditioned Initial Residual");
543 const std::string expResScale (
"Norm of Initial Residual");
544 const std::string timerLabel (
"Belos");
545 const std::string orthoType (
"ICGS");
546 RCP<const ParameterList> orthoParams = orthoFactory_.getDefaultParameters (orthoType);
553 const MagnitudeType orthoKappa = -SMT::one();
556 pl->set (
"Convergence Tolerance", convTol,
557 "The tolerance that the solver needs to achieve in order for "
558 "the linear system(s) to be declared converged. The meaning "
559 "of this tolerance depends on the convergence test details.");
560 pl->set(
"Maximum Restarts", maxRestarts,
561 "The maximum number of restart cycles (not counting the first) "
562 "allowed for each set of right-hand sides solved.");
563 pl->set(
"Maximum Iterations", maxIters,
564 "The maximum number of iterations allowed for each set of "
565 "right-hand sides solved.");
566 pl->set(
"Block Size", blockSize,
567 "How many linear systems to solve at once.");
568 pl->set(
"Num Blocks", numBlocks,
569 "The maximum number of blocks allowed in the Krylov subspace "
570 "for each set of right-hand sides solved. This includes the "
571 "initial block. Total Krylov basis storage, not counting the "
572 "recycled subspace, is \"Num Blocks\" times \"Block Size\".");
573 pl->set(
"Num Recycled Blocks", numRecycledBlocks,
574 "The maximum number of vectors in the recycled subspace." );
575 pl->set(
"Verbosity", verbosity,
576 "What type(s) of solver information should be written "
577 "to the output stream.");
578 pl->set(
"Output Style", outputStyle,
579 "What style is used for the solver information to write "
580 "to the output stream.");
581 pl->set(
"Output Frequency", outputFreq,
582 "How often convergence information should be written "
583 "to the output stream.");
584 pl->set(
"Output Stream", outputStream,
585 "A reference-counted pointer to the output stream where all "
586 "solver output is sent.");
587 pl->set(
"Implicit Residual Scaling", impResScale,
588 "The type of scaling used in the implicit residual convergence test.");
589 pl->set(
"Explicit Residual Scaling", expResScale,
590 "The type of scaling used in the explicit residual convergence test.");
591 pl->set(
"Timer Label", timerLabel,
592 "The string to use as a prefix for the timer labels.");
593 pl->set(
"Orthogonalization", orthoType,
594 "The orthogonalization method to use. Valid options: " +
595 orthoFactory_.validNamesString());
596 pl->set (
"Orthogonalization Parameters", *orthoParams,
597 "Sublist of parameters specific to the orthogonalization method to use.");
598 pl->set(
"Orthogonalization Constant", orthoKappa,
599 "When using DGKS orthogonalization: the \"depTol\" constant, used "
600 "to determine whether another step of classical Gram-Schmidt is "
601 "necessary. Otherwise ignored. Nonpositive values are ignored.");
604 return defaultParams_;
607 template<
class ScalarType,
class MV,
class OP>
610 setParameters (
const Teuchos::RCP<Teuchos::ParameterList> ¶ms) {
611 using Teuchos::isParameterType;
612 using Teuchos::getParameter;
614 using Teuchos::ParameterList;
615 using Teuchos::parameterList;
618 using Teuchos::rcp_dynamic_cast;
619 using Teuchos::rcpFromRef;
620 using Teuchos::Exceptions::InvalidParameter;
621 using Teuchos::Exceptions::InvalidParameterName;
622 using Teuchos::Exceptions::InvalidParameterType;
627 if (params.is_null()) {
632 params_ = parameterList (*defaultParams);
644 params->validateParametersAndSetDefaults (*defaultParams);
660 const int maxRestarts = params_->get<
int> (
"Maximum Restarts");
661 TEUCHOS_TEST_FOR_EXCEPTION(maxRestarts <= 0, std::invalid_argument,
662 "Belos::BlockGCRODRSolMgr: The \"Maximum Restarts\" parameter "
663 "must be nonnegative, but you specified a negative value of "
664 << maxRestarts <<
".");
666 const int maxIters = params_->get<
int> (
"Maximum Iterations");
667 TEUCHOS_TEST_FOR_EXCEPTION(maxIters <= 0, std::invalid_argument,
668 "Belos::BlockGCRODRSolMgr: The \"Maximum Iterations\" parameter "
669 "must be positive, but you specified a nonpositive value of "
672 const int numBlocks = params_->get<
int> (
"Num Blocks");
673 TEUCHOS_TEST_FOR_EXCEPTION(numBlocks <= 0, std::invalid_argument,
674 "Belos::BlockGCRODRSolMgr: The \"Num Blocks\" parameter must be "
675 "positive, but you specified a nonpositive value of " << numBlocks
678 const int blockSize = params_->get<
int> (
"Block Size");
679 TEUCHOS_TEST_FOR_EXCEPTION(blockSize <= 0, std::invalid_argument,
680 "Belos::BlockGCRODRSolMgr: The \"Block Size\" parameter must be "
681 "positive, but you specified a nonpositive value of " << blockSize
684 const int recycledBlocks = params_->get<
int> (
"Num Recycled Blocks");
685 TEUCHOS_TEST_FOR_EXCEPTION(recycledBlocks <= 0, std::invalid_argument,
686 "Belos::BlockGCRODRSolMgr: The \"Num Recycled Blocks\" parameter must "
687 "be positive, but you specified a nonpositive value of "
688 << recycledBlocks <<
".");
689 TEUCHOS_TEST_FOR_EXCEPTION(recycledBlocks >= numBlocks,
690 std::invalid_argument,
"Belos::BlockGCRODRSolMgr: The \"Num Recycled "
691 "Blocks\" parameter must be less than the \"Num Blocks\" parameter, "
692 "but you specified \"Num Recycled Blocks\" = " << recycledBlocks
693 <<
" and \"Num Blocks\" = " << numBlocks <<
".");
697 maxRestarts_ = maxRestarts;
698 maxIters_ = maxIters;
699 numBlocks_ = numBlocks;
700 blockSize_ = blockSize;
701 recycledBlocks_ = recycledBlocks;
708 std::string tempLabel = params_->get<std::string> (
"Timer Label");
709 const bool labelChanged = (tempLabel != label_);
711#ifdef BELOS_TEUCHOS_TIME_MONITOR
712 std::string solveLabel = label_ +
": BlockGCRODRSolMgr total solve time";
713 if (timerSolve_.is_null()) {
715 timerSolve_ = Teuchos::TimeMonitor::getNewCounter (solveLabel);
716 }
else if (labelChanged) {
722 Teuchos::TimeMonitor::clearCounter (solveLabel);
723 timerSolve_ = Teuchos::TimeMonitor::getNewCounter (solveLabel);
731 if (params_->isParameter (
"Verbosity")) {
732 if (isParameterType<int> (*params_,
"Verbosity")) {
733 verbosity_ = params_->get<
int> (
"Verbosity");
736 verbosity_ = (int) getParameter<MsgType> (*params_,
"Verbosity");
743 if (params_->isParameter (
"Output Style")) {
744 if (isParameterType<int> (*params_,
"Output Style")) {
745 outputStyle_ = params_->get<
int> (
"Output Style");
748 outputStyle_ = (int) getParameter<OutputType> (*params_,
"Output Style");
776 if (params_->isParameter (
"Output Stream")) {
778 outputStream_ = getParameter<RCP<std::ostream> > (*params_,
"Output Stream");
780 catch (InvalidParameter&) {
781 outputStream_ = rcpFromRef (std::cout);
788 if (outputStream_.is_null()) {
789 outputStream_ = rcp (
new Teuchos::oblackholestream);
793 outputFreq_ = params_->get<
int> (
"Output Frequency");
796 if (! outputTest_.is_null()) {
797 outputTest_->setOutputFrequency (outputFreq_);
805 if (printer_.is_null()) {
808 printer_->setVerbosity (verbosity_);
809 printer_->setOStream (outputStream_);
820 if (params_->isParameter (
"Orthogonalization")) {
821 const std::string& tempOrthoType =
822 params_->get<std::string> (
"Orthogonalization");
824 if (! orthoFactory_.isValidName (tempOrthoType)) {
825 std::ostringstream os;
826 os <<
"Belos::BlockGCRODRSolMgr: Invalid orthogonalization name \""
827 << tempOrthoType <<
"\". The following are valid options "
828 <<
"for the \"Orthogonalization\" name parameter: ";
829 orthoFactory_.printValidNames (os);
830 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::invalid_argument, os.str());
832 if (tempOrthoType != orthoType_) {
834 orthoType_ = tempOrthoType;
851 RCP<ParameterList> orthoParams = sublist (params,
"Orthogonalization Parameters",
true);
852 TEUCHOS_TEST_FOR_EXCEPTION(orthoParams.is_null(), std::logic_error,
853 "Failed to get orthogonalization parameters. "
854 "Please report this bug to the Belos developers.");
880 ortho_ = orthoFactory_.makeMatOrthoManager (orthoType_, null, printer_,
881 label_, orthoParams);
893 if (params_->isParameter (
"Orthogonalization Constant")) {
894 const MagnitudeType orthoKappa =
895 params_->get<MagnitudeType> (
"Orthogonalization Constant");
896 if (orthoKappa > 0) {
897 orthoKappa_ = orthoKappa;
900 if (orthoType_ ==
"DGKS" && ! ortho_.is_null()) {
906 rcp_dynamic_cast<ortho_man_type>(ortho_)->
setDepTol (orthoKappa_);
916 convTol_ = params_->get<MagnitudeType> (
"Convergence Tolerance");
917 if (! impConvTest_.is_null()) {
918 impConvTest_->setTolerance (convTol_);
920 if (! expConvTest_.is_null()) {
921 expConvTest_->setTolerance (convTol_);
925 if (params_->isParameter (
"Implicit Residual Scaling")) {
926 std::string tempImpResScale =
927 getParameter<std::string> (*params_,
"Implicit Residual Scaling");
930 if (impResScale_ != tempImpResScale) {
932 impResScale_ = tempImpResScale;
934 if (! impConvTest_.is_null()) {
947 if (params_->isParameter(
"Explicit Residual Scaling")) {
948 std::string tempExpResScale =
949 getParameter<std::string> (*params_,
"Explicit Residual Scaling");
952 if (expResScale_ != tempExpResScale) {
954 expResScale_ = tempExpResScale;
956 if (! expConvTest_.is_null()) {
974 if (maxIterTest_.is_null()) {
977 maxIterTest_->setMaxIters (maxIters_);
982 if (impConvTest_.is_null()) {
983 impConvTest_ = rcp (
new StatusTestResNorm_t (convTol_));
989 if (expConvTest_.is_null()) {
990 expConvTest_ = rcp (
new StatusTestResNorm_t (convTol_));
991 expConvTest_->defineResForm (StatusTestResNorm_t::Explicit,
Belos::TwoNorm);
997 if (convTest_.is_null()) {
998 convTest_ = rcp (
new StatusTestCombo_t (StatusTestCombo_t::SEQ,
1006 sTest_ = rcp (
new StatusTestCombo_t (StatusTestCombo_t::OR,
1007 maxIterTest_, convTest_));
1011 outputTest_ = stoFactory.
create (printer_, sTest_, outputFreq_,
1015 std::string solverDesc =
"Block GCRODR ";
1016 outputTest_->setSolverDesc (solverDesc);
1023 template<
class ScalarType,
class MV,
class OP>
1025 BlockGCRODRSolMgr<ScalarType,MV,OP>::initializeStateStorage()
1028 ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
1031 Teuchos::RCP<const MV> rhsMV = problem_->getRHS();
1034 int KrylSpaDim = (numBlocks_ - 1) * blockSize_;
1037 int augSpaDim = KrylSpaDim + recycledBlocks_ + 1;
1040 int augSpaImgDim = KrylSpaDim + blockSize_ + recycledBlocks_+1;
1064 if (U_ == Teuchos::null) {
1065 U_ = MVT::Clone( *rhsMV, recycledBlocks_+1 );
1069 if (MVT::GetNumberVecs(*U_) < recycledBlocks_+1) {
1070 Teuchos::RCP<const MV> tmp = U_;
1071 U_ = MVT::Clone( *tmp, recycledBlocks_+1 );
1076 if (C_ == Teuchos::null) {
1077 C_ = MVT::Clone( *rhsMV, recycledBlocks_+1 );
1081 if (MVT::GetNumberVecs(*C_) < recycledBlocks_+1) {
1082 Teuchos::RCP<const MV> tmp = C_;
1083 C_ = MVT::Clone( *tmp, recycledBlocks_+1 );
1088 if (U1_ == Teuchos::null) {
1089 U1_ = MVT::Clone( *rhsMV, recycledBlocks_+1 );
1093 if (MVT::GetNumberVecs(*U1_) < recycledBlocks_+1) {
1094 Teuchos::RCP<const MV> tmp = U1_;
1095 U1_ = MVT::Clone( *tmp, recycledBlocks_+1 );
1100 if (C1_ == Teuchos::null) {
1101 C1_ = MVT::Clone( *rhsMV, recycledBlocks_+1 );
1105 if (MVT::GetNumberVecs(*U1_) < recycledBlocks_+1) {
1106 Teuchos::RCP<const MV> tmp = C1_;
1107 C1_ = MVT::Clone( *tmp, recycledBlocks_+1 );
1112 if (R_ == Teuchos::null){
1113 R_ = MVT::Clone( *rhsMV, blockSize_ );
1119 if (G_ == Teuchos::null){
1120 G_ = Teuchos::rcp(
new SDM( augSpaImgDim, augSpaDim ) );
1123 if ( (G_->numRows() != augSpaImgDim) || (G_->numCols() != augSpaDim) )
1125 G_->reshape( augSpaImgDim, augSpaDim );
1127 G_->putScalar(zero);
1131 if (H_ == Teuchos::null){
1132 H_ = Teuchos::rcp (
new SDM ( Teuchos::View, *G_, KrylSpaDim + blockSize_, KrylSpaDim, recycledBlocks_+1 ,recycledBlocks_+1 ) );
1136 if (F_ == Teuchos::null){
1137 F_ = Teuchos::rcp(
new SDM( recycledBlocks_+1, recycledBlocks_+1 ) );
1140 if ( (F_->numRows() != recycledBlocks_+1) || (F_->numCols() != recycledBlocks_+1) ){
1141 F_->reshape( recycledBlocks_+1, recycledBlocks_+1 );
1144 F_->putScalar(zero);
1147 if (PP_ == Teuchos::null){
1148 PP_ = Teuchos::rcp(
new SDM( augSpaImgDim, recycledBlocks_+1 ) );
1151 if ( (PP_->numRows() != augSpaImgDim) || (PP_->numCols() != recycledBlocks_+1) ){
1152 PP_->reshape( augSpaImgDim, recycledBlocks_+1 );
1157 if (HP_ == Teuchos::null)
1158 HP_ = Teuchos::rcp(
new SDM( augSpaImgDim, augSpaDim ) );
1160 if ( (HP_->numRows() != augSpaImgDim) || (HP_->numCols() != augSpaDim) ){
1161 HP_->reshape( augSpaImgDim, augSpaDim );
1166 tau_.resize(recycledBlocks_+1);
1169 work_.resize(recycledBlocks_+1);
1172 ipiv_.resize(recycledBlocks_+1);
1176template<
class ScalarType,
class MV,
class OP>
1179 ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
1180 ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
1182 int p = block_gmres_iter->getState().curDim;
1183 std::vector<int> index(keff);
1187 SDM HH(Teuchos::Copy, *H_, p+blockSize_, p);
1188 if(recycledBlocks_ >= p + blockSize_){
1192 for (
int ii=0; ii < p; ++ii) index[ii] = ii;
1193 Teuchos::RCP<MV> Utmp = MVT::CloneViewNonConst( *U_, index );
1194 MVT::SetBlock(*V_, index, *Utmp);
1199 Teuchos::RCP<SDM > PPtmp = Teuchos::rcp (
new SDM ( Teuchos::View, *PP_, p, recycledBlocks_+1 ) );
1200 if(recycleMethod_ ==
"harmvecs"){
1201 keff = getHarmonicVecsKryl(p, HH, *PPtmp);
1202 printer_->stream(
Debug) <<
"keff = " << keff << std::endl;
1205PPtmp = Teuchos::rcp (
new SDM ( Teuchos::View, *PP_, p, keff ) );
1208for (
int ii=0; ii<keff; ++ii) index[ii] = ii;
1209Teuchos::RCP<MV> Ctmp = MVT::CloneViewNonConst( *C_, index );
1210Teuchos::RCP<MV> Utmp = MVT::CloneViewNonConst( *U_, index );
1211Teuchos::RCP<MV> U1tmp = MVT::CloneViewNonConst( *U1_, index );
1213for (
int ii=0; ii < p; ++ii) index[ii] = ii;
1214Teuchos::RCP<const MV> Vtmp = MVT::CloneView( *V_, index );
1218MVT::MvTimesMatAddMv( one, *Vtmp, *PPtmp, zero, *U1tmp );
1221SDM HPtmp( Teuchos::View, *HP_, p+blockSize_, keff );
1222HPtmp.multiply( Teuchos::NO_TRANS, Teuchos::NO_TRANS, one, *H_, *PPtmp, zero );
1226lapack.GEQRF(HPtmp.numRows(),HPtmp.numCols(),HPtmp.values(),HPtmp.stride(),&tau_[0],&work_[0],lwork,&info);
1227TEUCHOS_TEST_FOR_EXCEPTION(info != 0,
BlockGCRODRSolMgrLAPACKFailure,
"Belos::BlockGCRODRSolMgr::solve(): LAPACK _GEQRF failed to compute a workspace size.");
1230 lwork =
static_cast<int> (Teuchos::ScalarTraits<ScalarType>::magnitude (work_[0]));
1232lapack.GEQRF(HPtmp.numRows(),HPtmp.numCols(),HPtmp.values(),HPtmp.stride(),&tau_[0],&work_[0],lwork,&info);
1233TEUCHOS_TEST_FOR_EXCEPTION(info != 0,
BlockGCRODRSolMgrLAPACKFailure,
"Belos::BlockGCRODRSolMgr::solve(): LAPACK _GEQRF failed to compute a QR factorization.");
1237SDM Rtmp( Teuchos::View, *F_, keff, keff );
1238for(
int ii=0;ii<keff;ii++) {
for(
int jj=ii;jj<keff;jj++) Rtmp(ii,jj) = HPtmp(ii,jj); }
1240lapack.UNGQR(HPtmp.numRows(),HPtmp.numCols(),HPtmp.numCols(),HPtmp.values(),HPtmp.stride(),&tau_[0],&work_[0],lwork,&info);
1241TEUCHOS_TEST_FOR_EXCEPTION(info != 0,
BlockGCRODRSolMgrLAPACKFailure,
"Belos::BlockGCRODRSolMgr::solve(): LAPACK _UNGQR failed to construct the Q factor.");
1245 index.resize( p+blockSize_ );
1246 for (
int ii=0; ii < (p+blockSize_); ++ii) { index[ii] = ii; }
1247 Vtmp = MVT::CloneView( *V_, index );
1248 MVT::MvTimesMatAddMv( one, *Vtmp, HPtmp, zero, *Ctmp );
1255ipiv_.resize(Rtmp.numRows());
1256lapack.GETRF(Rtmp.numRows(),Rtmp.numCols(),Rtmp.values(),Rtmp.stride(),&ipiv_[0],&info);
1257TEUCHOS_TEST_FOR_EXCEPTION(info != 0,
BlockGCRODRSolMgrLAPACKFailure,
"Belos::GCRODRSolMgr::solve(): LAPACK _GETRF failed to compute an LU factorization.");
1259lwork = Rtmp.numRows();
1261lapack.GETRI(Rtmp.numRows(),Rtmp.values(),Rtmp.stride(),&ipiv_[0],&work_[0],lwork,&info);
1264MVT::MvTimesMatAddMv( one, *U1tmp, Rtmp, zero, *Utmp );
1270template<
class ScalarType,
class MV,
class OP>
1272 const MagnitudeType one = Teuchos::ScalarTraits<ScalarType>::one();
1273 const ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
1275 std::vector<MagnitudeType> d(keff);
1276 std::vector<ScalarType> dscalar(keff);
1277 std::vector<int> index(numBlocks_+1);
1286 if(recycledBlocks_ >= keff + p){
1289 for (
int ii=0; ii < p; ++ii) { index[ii] = keff+ii; }
1290 Teuchos::RCP<MV> Utmp = MVT::CloneViewNonConst( *U_, index );
1291 for (
int ii=0; ii < p; ++ii) { index[ii] =ii; }
1292 MVT::SetBlock(*V_, index, *Utmp);
1299 for (
int ii=0; ii<keff; ++ii) { index[ii] = ii; }
1300 Teuchos::RCP<MV> Utmp = MVT::CloneViewNonConst( *U_, index );
1302 dscalar.resize(keff);
1303 MVT::MvNorm( *Utmp, d );
1304 for (
int i=0; i<keff; ++i) {
1306 dscalar[i] = (ScalarType)d[i];
1308 MVT::MvScale( *Utmp, dscalar );
1313 Teuchos::RCP<SDM> Gtmp = Teuchos::rcp(
new SDM( Teuchos::View, *G_, p+keff, p+keff-blockSize_ ) );
1316 for (
int i=0; i<keff; ++i)
1317 (*Gtmp)(i,i) = d[i];
1324 SDM PPtmp( Teuchos::View, *PP_, p+keff-blockSize_, recycledBlocks_+1 );
1325 keff_new = getHarmonicVecsAugKryl( keff, p-blockSize_, *Gtmp, oldState.V, PPtmp );
1332 Teuchos::RCP<MV> U1tmp;
1334 index.resize( keff );
1335 for (
int ii=0; ii<keff; ++ii) { index[ii] = ii; }
1336 Teuchos::RCP<const MV> Utmp = MVT::CloneView( *U_, index );
1337 index.resize( keff_new );
1338 for (
int ii=0; ii<keff_new; ++ii) { index[ii] = ii; }
1339 U1tmp = MVT::CloneViewNonConst( *U1_, index );
1340 SDM PPtmp( Teuchos::View, *PP_, keff, keff_new );
1341 MVT::MvTimesMatAddMv( one, *Utmp, PPtmp, zero, *U1tmp );
1346 index.resize(p-blockSize_);
1347 for (
int ii=0; ii < p-blockSize_; ii++) { index[ii] = ii; }
1348 Teuchos::RCP<const MV> Vtmp = MVT::CloneView( *V_, index );
1349 SDM PPtmp( Teuchos::View, *PP_, p-blockSize_, keff_new, keff );
1350 MVT::MvTimesMatAddMv( one, *Vtmp, PPtmp, one, *U1tmp );
1354 SDM HPtmp( Teuchos::View, *HP_, p+keff, keff_new );
1356 SDM PPtmp( Teuchos::View, *PP_, p-blockSize_+keff, keff_new );
1357 HPtmp.multiply(Teuchos::NO_TRANS,Teuchos::NO_TRANS,one,*Gtmp,PPtmp,zero);
1361 int info = 0, lwork = -1;
1362 tau_.resize(keff_new);
1363 lapack.GEQRF(HPtmp.numRows(),HPtmp.numCols(),HPtmp.values(),HPtmp.stride(),&tau_[0],&work_[0],lwork,&info);
1364 TEUCHOS_TEST_FOR_EXCEPTION(info != 0,
BlockGCRODRSolMgrLAPACKFailure,
"Belos::BlockGCRODRSolMgr::solve(): LAPACK _GEQRF failed to compute a workspace size.");
1366 lwork =
static_cast<int> (Teuchos::ScalarTraits<ScalarType>::real (work_[0]));
1367 work_.resize(lwork);
1368 lapack.GEQRF(HPtmp.numRows(),HPtmp.numCols(),HPtmp.values(),HPtmp.stride(),&tau_[0],&work_[0],lwork,&info);
1369 TEUCHOS_TEST_FOR_EXCEPTION(info != 0,
BlockGCRODRSolMgrLAPACKFailure,
"Belos::BlockGCRODRSolMgr::solve(): LAPACK _GEQRF failed to compute a QR factorization.");
1373 SDM Ftmp( Teuchos::View, *F_, keff_new, keff_new );
1374 for(
int i=0;i<keff_new;i++) {
for(
int j=i;j<keff_new;j++) Ftmp(i,j) = HPtmp(i,j); }
1376 lapack.UNGQR(HPtmp.numRows(),HPtmp.numCols(),HPtmp.numCols(),HPtmp.values(),HPtmp.stride(),&tau_[0],&work_[0],lwork,&info);
1377 TEUCHOS_TEST_FOR_EXCEPTION(info != 0,
BlockGCRODRSolMgrLAPACKFailure,
"Belos::BlockGCRODRSolMgr::solve(): LAPACK _UNGQR failed to construct the Q factor.");
1384 Teuchos::RCP<MV> C1tmp;
1387 for (
int i=0; i < keff; i++) { index[i] = i; }
1388 Teuchos::RCP<const MV> Ctmp = MVT::CloneView( *C_, index );
1389 index.resize(keff_new);
1390 for (
int i=0; i < keff_new; i++) { index[i] = i; }
1391 C1tmp = MVT::CloneViewNonConst( *C1_, index );
1392 SDM PPtmp( Teuchos::View, *HP_, keff, keff_new );
1393 MVT::MvTimesMatAddMv( one, *Ctmp, PPtmp, zero, *C1tmp );
1398 for (
int i=0; i < p; ++i) { index[i] = i; }
1399 Teuchos::RCP<const MV> Vtmp = MVT::CloneView( *V_, index );
1400 SDM PPtmp( Teuchos::View, *HP_, p, keff_new, keff, 0 );
1401 MVT::MvTimesMatAddMv( one, *Vtmp, PPtmp, one, *C1tmp );
1410 ipiv_.resize(Ftmp.numRows());
1411 lapack.GETRF(Ftmp.numRows(),Ftmp.numCols(),Ftmp.values(),Ftmp.stride(),&ipiv_[0],&info);
1412 TEUCHOS_TEST_FOR_EXCEPTION(info != 0,
BlockGCRODRSolMgrLAPACKFailure,
"Belos::BlockGCRODRSolMgr::solve(): LAPACK _GETRF failed to compute an LU factorization.");
1415 lwork = Ftmp.numRows();
1416 work_.resize(lwork);
1417 lapack.GETRI(Ftmp.numRows(),Ftmp.values(),Ftmp.stride(),&ipiv_[0],&work_[0],lwork,&info);
1418 TEUCHOS_TEST_FOR_EXCEPTION(info != 0,
BlockGCRODRSolMgrLAPACKFailure,
"Belos::BlockGCRODRSolMgr::solve(): LAPACK _GETRI failed to compute an LU factorization.");
1421 index.resize(keff_new);
1422 for (
int i=0; i < keff_new; i++) { index[i] = i; }
1423 Teuchos::RCP<MV> Utmp = MVT::CloneViewNonConst( *U_, index );
1424 MVT::MvTimesMatAddMv( one, *U1tmp, Ftmp, zero, *Utmp );
1429 if (keff != keff_new) {
1431 block_gcrodr_iter->setSize( keff, numBlocks_ );
1433 SDM b1( Teuchos::View, *G_, recycledBlocks_+2, 1, 0, recycledBlocks_ );
1439template<
class ScalarType,
class MV,
class OP>
1440int BlockGCRODRSolMgr<ScalarType,MV,OP>::getHarmonicVecsAugKryl(
int keff,
int m,
const SDM& GG,
const Teuchos::RCP<const MV>& VV, SDM& PP){
1442 int m2 = GG.numCols();
1443 bool xtraVec =
false;
1444 ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
1445 ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
1446 std::vector<int> index;
1449 std::vector<MagnitudeType> wr(m2), wi(m2);
1452 std::vector<MagnitudeType> w(m2);
1455 SDM vr(m2,m2,
false);
1458 std::vector<int> iperm(m2);
1461 builtRecycleSpace_ =
true;
1467 B.multiply(Teuchos::TRANS,Teuchos::NO_TRANS,one,GG,GG,zero);
1471 SDM A_tmp( keff+m+blockSize_, keff+m );
1476 for (
int i=0; i<keff; ++i) { index[i] = i; }
1477 Teuchos::RCP<const MV> Ctmp = MVT::CloneView( *C_, index );
1478 Teuchos::RCP<const MV> Utmp = MVT::CloneView( *U_, index );
1479 SDM A11( Teuchos::View, A_tmp, keff, keff );
1480 MVT::MvTransMv( one, *Ctmp, *Utmp, A11 );
1483 SDM A21( Teuchos::View, A_tmp, m+blockSize_, keff, keff );
1484 index.resize(m+blockSize_);
1485 for (i=0; i < m+blockSize_; i++) { index[i] = i; }
1486 Teuchos::RCP<const MV> Vp = MVT::CloneView( *VV, index );
1487 MVT::MvTransMv( one, *Vp, *Utmp, A21 );
1490 for( i=keff; i<keff+m; i++)
1494 SDM A( m2, A_tmp.numCols() );
1495 A.multiply( Teuchos::TRANS, Teuchos::NO_TRANS, one, GG, A_tmp, zero );
1503 char balanc=
'P', jobvl=
'N', jobvr=
'V', sense=
'N';
1504 int ld = A.numRows();
1506 int ldvl = ld, ldvr = ld;
1507 int info = 0,ilo = 0,ihi = 0;
1508 MagnitudeType abnrm = 0.0, bbnrm = 0.0;
1510 std::vector<ScalarType> beta(ld);
1511 std::vector<ScalarType> work(lwork);
1512 std::vector<MagnitudeType> rwork(lwork);
1513 std::vector<MagnitudeType> lscale(ld), rscale(ld);
1514 std::vector<MagnitudeType> rconde(ld), rcondv(ld);
1515 std::vector<int> iwork(ld+6);
1520 lapack.GGEVX(balanc, jobvl, jobvr, sense, ld, A.values(), ld, B.values(), ld, &wr[0], &wi[0],
1521 &beta[0], vl, ldvl, vr.values(), ldvr, &ilo, &ihi, &lscale[0], &rscale[0],
1522 &abnrm, &bbnrm, &rconde[0], &rcondv[0], &work[0], lwork, &rwork[0],
1523 &iwork[0], bwork, &info);
1524 TEUCHOS_TEST_FOR_EXCEPTION(info != 0,
BlockGCRODRSolMgrLAPACKFailure,
"Belos::BlockGCRODRSolMgr::solve(): LAPACK GGEVX failed to compute eigensolutions.");
1528 for( i=0; i<ld; i++ )
1529 w[i] = Teuchos::ScalarTraits<MagnitudeType>::squareroot( wr[i]*wr[i] + wi[i]*wi[i] ) / std::abs( beta[i] );
1531 this->sort(w,ld,iperm);
1533 bool scalarTypeIsComplex = Teuchos::ScalarTraits<ScalarType>::isComplex;
1536 for( i=0; i<recycledBlocks_; i++ )
1537 for( j=0; j<ld; j++ )
1538 PP(j,i) = vr(j,iperm[ld-recycledBlocks_+i]);
1540 if(scalarTypeIsComplex==
false) {
1543 if (wi[iperm[ld-recycledBlocks_]] != 0.0) {
1545 for ( i=ld-recycledBlocks_; i<ld; i++ )
1546 if (wi[iperm[i]] != 0.0) countImag++;
1548 if (countImag % 2) xtraVec =
true;
1552 if (wi[iperm[ld-recycledBlocks_]] > 0.0) {
1553 for( j=0; j<ld; j++ )
1554 PP(j,recycledBlocks_) = vr(j,iperm[ld-recycledBlocks_]+1);
1557 for( j=0; j<ld; j++ )
1558 PP(j,recycledBlocks_) = vr(j,iperm[ld-recycledBlocks_]-1);
1566 return recycledBlocks_+1;
1568 return recycledBlocks_;
1572template<
class ScalarType,
class MV,
class OP>
1573int BlockGCRODRSolMgr<ScalarType,MV,OP>::getHarmonicVecsKryl(
int m,
const SDM& HH, SDM& PP){
1574 bool xtraVec =
false;
1575 ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
1576 ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
1579 std::vector<MagnitudeType> wr(m), wi(m);
1585 std::vector<MagnitudeType> w(m);
1588 std::vector<int> iperm(m);
1594 builtRecycleSpace_ =
true;
1597 SDM HHt( HH, Teuchos::TRANS );
1598 Teuchos::RCP<SDM> harmRitzMatrix = Teuchos::rcp(
new SDM( m, blockSize_));
1601 for(
int i=0; i<=blockSize_-1; i++) (*harmRitzMatrix)[blockSize_-1-i][harmRitzMatrix->numRows()-1-i] = 1;
1604 lapack.GESV(m, blockSize_, HHt.values(), HHt.stride(), &iperm[0], harmRitzMatrix->values(), harmRitzMatrix->stride(), &info);
1606 TEUCHOS_TEST_FOR_EXCEPTION(info != 0,
BlockGCRODRSolMgrLAPACKFailure,
"Belos::BlockGCRODRSolMgr::solve(): LAPACK GESV failed to compute a solution.");
1610 Teuchos::SerialDenseMatrix<int, ScalarType> H_lbl(Teuchos::View, HH, blockSize_, blockSize_, (HH).numRows()-blockSize_, (HH).numCols()-blockSize_ );
1611 Teuchos::SerialDenseMatrix<int, ScalarType> H_lbl_t( H_lbl, Teuchos::TRANS );
1616 Teuchos::RCP<SDM> Htemp = Teuchos::null;
1617 Htemp = Teuchos::rcp(
new SDM(H_lbl_t.numRows(), H_lbl_t.numCols()));
1618 Htemp -> multiply(Teuchos::NO_TRANS, Teuchos::NO_TRANS, one, H_lbl_t, H_lbl, zero);
1619 H_lbl_t.assign(*Htemp);
1621 Htemp = Teuchos::rcp(
new SDM(harmRitzMatrix -> numRows(), harmRitzMatrix -> numCols()));
1622 Htemp -> multiply(Teuchos::NO_TRANS, Teuchos::NO_TRANS, one, *harmRitzMatrix, H_lbl_t, zero);
1623 harmRitzMatrix -> assign(*Htemp);
1626 int harmColIndex, HHColIndex;
1627 Htemp = Teuchos::rcp(
new SDM(Teuchos::Copy,HH,HH.numRows()-blockSize_,HH.numCols()));
1628 for(
int i = 0; i<blockSize_; i++){
1629 harmColIndex = harmRitzMatrix -> numCols() - i -1;
1631 for(
int j=0; j<m; j++) (*Htemp)[HHColIndex][j] += (*harmRitzMatrix)[harmColIndex][j];
1633 harmRitzMatrix = Htemp;
1641 std::vector<ScalarType> work(1);
1642 std::vector<MagnitudeType> rwork(2*m);
1648 lapack.GEEV(
'N',
'V', m, harmRitzMatrix->values(), harmRitzMatrix->stride(), &wr[0], &wi[0],
1649 vl, ldvl, vr.values(), vr.stride(), &work[0], lwork, &rwork[0], &info);
1651 lwork = std::abs (
static_cast<int> (Teuchos::ScalarTraits<ScalarType>::real (work[0])));
1652 work.resize( lwork );
1654 lapack.GEEV(
'N',
'V', m, harmRitzMatrix->values(), harmRitzMatrix->stride(), &wr[0], &wi[0],
1655 vl, ldvl, vr.values(), vr.stride(), &work[0], lwork, &rwork[0], &info);
1657 TEUCHOS_TEST_FOR_EXCEPTION(info != 0,
BlockGCRODRSolMgrLAPACKFailure,
"Belos::BlockGCRODRSolMgr::solve(): LAPACK GEEV failed to compute eigensolutions.");
1660 for(
int i=0; i<m; ++i ) w[i] = Teuchos::ScalarTraits<MagnitudeType>::squareroot( wr[i]*wr[i] + wi[i]*wi[i] );
1662 this->sort(w, m, iperm);
1664 bool scalarTypeIsComplex = Teuchos::ScalarTraits<ScalarType>::isComplex;
1667 for(
int i=0; i<recycledBlocks_; ++i )
1668 for(
int j=0; j<m; j++ )
1669 PP(j,i) = vr(j,iperm[i]);
1671 if(scalarTypeIsComplex==
false) {
1674 if (wi[iperm[recycledBlocks_-1]] != 0.0) {
1676 for (
int i=0; i<recycledBlocks_; ++i )
1677 if (wi[iperm[i]] != 0.0) countImag++;
1679 if (countImag % 2) xtraVec =
true;
1683 if (wi[iperm[recycledBlocks_-1]] > 0.0) {
1684 for(
int j=0; j<m; ++j )
1685 PP(j,recycledBlocks_) = vr(j,iperm[recycledBlocks_-1]+1);
1688 for(
int j=0; j<m; ++j )
1689 PP(j,recycledBlocks_) = vr(j,iperm[recycledBlocks_-1]-1);
1697 printer_->stream(
Debug) <<
"Recycled " << recycledBlocks_+1 <<
" vectors" << std::endl;
1698 return recycledBlocks_+1;
1701 printer_->stream(
Debug) <<
"Recycled " << recycledBlocks_ <<
" vectors" << std::endl;
1702 return recycledBlocks_;
1707template<
class ScalarType,
class MV,
class OP>
1708void BlockGCRODRSolMgr<ScalarType,MV,OP>::sort(std::vector<MagnitudeType>& dlist,
int n, std::vector<int>& iperm) {
1709 int l, r, j, i, flag;
1711 MagnitudeType dRR, dK;
1736 if (dlist[j] > dlist[j - 1]) j = j + 1;
1737 if (dlist[j - 1] > dK) {
1738 dlist[i - 1] = dlist[j - 1];
1739 iperm[i - 1] = iperm[j - 1];
1752 dlist[r] = dlist[0];
1753 iperm[r] = iperm[0];
1767template<
class ScalarType,
class MV,
class OP>
1771 using Teuchos::rcp_const_cast;
1775 ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
1776 ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
1777 std::vector<int> index(numBlocks_+1);
1794 int numCurrRHS = ( numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;
1798 std::vector<int> currIdx;
1800 if ( adaptiveBlockSize_ ) {
1801 blockSize_ = numCurrRHS;
1802 currIdx.resize( numCurrRHS );
1803 for (
int i=0; i<numCurrRHS; ++i)
1804 currIdx[i] = startPtr+i;
1807 currIdx.resize( blockSize_ );
1808 for (
int i=0; i<numCurrRHS; ++i)
1809 currIdx[i] = startPtr+i;
1810 for (
int i=numCurrRHS; i<blockSize_; ++i)
1815 problem_->setLSIndex( currIdx );
1821 loaDetected_ =
false;
1824 bool isConverged =
true;
1827 initializeStateStorage();
1830 Teuchos::ParameterList plist;
1832 while (numRHS2Solve > 0){
1842 printer_->stream(
Debug) <<
" Now solving RHS index " << currIdx[0] <<
" using recycled subspace of dimension " << keff << std::endl << std::endl;
1846 for (
int ii=0; ii<keff; ++ii) { index[ii] = ii; }
1849 problem_->apply( *Utmp, *Ctmp );
1855 SDM Ftmp( Teuchos::View, *F_, keff, keff );
1856 int rank = ortho_->normalize(*Ctmp, rcp(&Ftmp,
false));
1858 TEUCHOS_TEST_FOR_EXCEPTION(rank != keff,
BlockGCRODRSolMgrOrthoFailure,
"Belos::BlockGCRODRSolMgr::solve(): Failed to compute orthonormal basis for initial recycled subspace.");
1863 ipiv_.resize(Ftmp.numRows());
1864 lapack.GETRF(Ftmp.numRows(),Ftmp.numCols(),Ftmp.values(),Ftmp.stride(),&ipiv_[0],&info);
1865 TEUCHOS_TEST_FOR_EXCEPTION(info != 0,
BlockGCRODRSolMgrLAPACKFailure,
"Belos::BlockGCRODRSolMgr::solve(): LAPACK _GETRF failed to compute an LU factorization.");
1868 int lwork = Ftmp.numRows();
1869 work_.resize(lwork);
1870 lapack.GETRI(Ftmp.numRows(),Ftmp.values(),Ftmp.stride(),&ipiv_[0],&work_[0],lwork,&info);
1871 TEUCHOS_TEST_FOR_EXCEPTION(info != 0,
BlockGCRODRSolMgrLAPACKFailure,
"Belos::BlockGCRODRSolMgr::solve(): LAPACK _GETRI failed to invert triangular matrix.");
1879 for (
int ii=0; ii<keff; ++ii) { index[ii] = ii; }
1884 SDM Ctr(keff,blockSize_);
1885 problem_->computeCurrPrecResVec( &*R_ );
1889 RCP<MV> update =
MVT::Clone( *problem_->getCurrLHSVec(), blockSize_ );
1892 problem_->updateSolution( update,
true );
1903 if (V_ == Teuchos::null) {
1905 Teuchos::RCP<const MV> rhsMV = problem_->getRHS();
1906 V_ =
MVT::Clone( *rhsMV, (numBlocks_+1)*blockSize_ );
1911 Teuchos::RCP<const MV> tmp = V_;
1912 V_ =
MVT::Clone( *tmp, (numBlocks_+1)*blockSize_ );
1917 printer_->stream(
Debug) <<
" No recycled subspace available for RHS index " << std::endl << std::endl;
1919 Teuchos::ParameterList primeList;
1922 primeList.set(
"Num Blocks",numBlocks_-1);
1923 primeList.set(
"Block Size",blockSize_);
1924 primeList.set(
"Recycled Blocks",0);
1925 primeList.set(
"Keep Hessenberg",
true);
1926 primeList.set(
"Initialize Hessenberg",
true);
1929 if (blockSize_*
static_cast<ptrdiff_t
>(numBlocks_) > dim) {
1930 ptrdiff_t tmpNumBlocks = 0;
1931 if (blockSize_ == 1)
1932 tmpNumBlocks = dim / blockSize_;
1934 tmpNumBlocks = ( dim - blockSize_) / blockSize_;
1935 printer_->stream(
Warnings) <<
"Belos::BlockGCRODRSolMgr::solve(): Warning! Requested Krylov subspace dimension is larger than operator dimension!"
1936 << std::endl <<
"The maximum number of blocks allowed for the Krylov subspace will be adjusted to " << tmpNumBlocks << std::endl;
1937 primeList.set(
"Num Blocks",Teuchos::as<int>(tmpNumBlocks));
1941 primeList.set(
"Num Blocks",numBlocks_-1);
1944 Teuchos::RCP<BlockGmresIter<ScalarType,MV,OP> > block_gmres_iter;
1948 block_gmres_iter->setSize( blockSize_, numBlocks_-1 );
1951 Teuchos::RCP<MV> V_0;
1952 if (currIdx[blockSize_-1] == -1) {
1953 V_0 =
MVT::Clone( *(problem_->getInitPrecResVec()), blockSize_ );
1954 problem_->computeCurrPrecResVec( &*V_0 );
1957 V_0 =
MVT::CloneCopy( *(problem_->getInitPrecResVec()), currIdx );
1961 Teuchos::RCP<SDM > z_0 = Teuchos::rcp(
new SDM( blockSize_, blockSize_ ) );
1964 int rank = ortho_->normalize( *V_0, z_0 );
1973 block_gmres_iter->initializeGmres(newstate);
1975 bool primeConverged =
false;
1978 printer_->stream(
Debug) <<
" Preparing to Iterate!!!!" << std::endl << std::endl;
1979 block_gmres_iter->iterate();
1984 if ( convTest_->getStatus() ==
Passed ) {
1985 printer_->stream(
Debug) <<
"We converged during the prime the pump stage" << std::endl << std::endl;
1986 primeConverged = !(expConvTest_->getLOADetected());
1987 if ( expConvTest_->getLOADetected() ) {
1989 loaDetected_ =
true;
1990 printer_->stream(
Warnings) <<
"Belos::BlockGCRODRSolMgr::solve(): Warning! Solver has experienced a loss of accuracy!" << std::endl;
1996 else if( maxIterTest_->getStatus() ==
Passed ) {
1998 primeConverged =
false;
2004 printer_->stream(
Debug) <<
" We did not converge on priming cycle of Block GMRES. Therefore we recycle and restart. " << std::endl << std::endl;
2009 if (blockSize_ != 1) {
2010 printer_->stream(
Errors) <<
"Error! Caught std::exception in BlockGmresIter::iterate() at iteration "
2011 << block_gmres_iter->getNumIters() << std::endl
2012 << e.what() << std::endl;
2013 if (convTest_->getStatus() !=
Passed)
2014 primeConverged =
false;
2018 block_gmres_iter->updateLSQR( block_gmres_iter->getCurSubspaceDim() );
2020 sTest_->checkStatus( &*block_gmres_iter );
2021 if (convTest_->getStatus() !=
Passed)
2022 isConverged =
false;
2027 achievedTol_ = MT::one();
2028 Teuchos::RCP<MV> X = problem_->getLHS();
2030 printer_->stream(
Warnings) <<
"Belos::BlockGCRODRSolMgr::solve(): Warning! NaN has been detected!"
2034 catch (
const std::exception &e) {
2035 printer_->stream(
Errors) <<
"Error! Caught std::exception in BlockGmresIter::iterate() at iteration "
2036 << block_gmres_iter->getNumIters() << std::endl
2037 << e.what() << std::endl;
2045 RCP<MV> update = block_gmres_iter->getCurrentUpdate();
2046 problem_->updateSolution( update,
true );
2050 problem_->computeCurrPrecResVec( &*R_ );
2053 newstate = block_gmres_iter->getState();
2056 H_->assign(*(newstate.
H));
2065 V_ = rcp_const_cast<MV>(newstate.
V);
2066 newstate.
V.release();
2068 buildRecycleSpaceKryl(keff, block_gmres_iter);
2069 printer_->stream(
Debug) <<
"Generated recycled subspace using RHS index " << currIdx[0] <<
" of dimension " << keff << std::endl << std::endl;
2073 if (primeConverged) {
2094 problem_->setCurrLS();
2097 startPtr += numCurrRHS;
2098 numRHS2Solve -= numCurrRHS;
2099 if ( numRHS2Solve > 0 ) {
2100 numCurrRHS = ( numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;
2101 if ( adaptiveBlockSize_ ) {
2102 blockSize_ = numCurrRHS;
2103 currIdx.resize( numCurrRHS );
2104 for (
int i=0; i<numCurrRHS; ++i) currIdx[i] = startPtr+i;
2107 currIdx.resize( blockSize_ );
2108 for (
int i=0; i<numCurrRHS; ++i) currIdx[i] = startPtr+i;
2109 for (
int i=numCurrRHS; i<blockSize_; ++i) currIdx[i] = -1;
2112 problem_->setLSIndex( currIdx );
2115 currIdx.resize( numRHS2Solve );
2125 Teuchos::ParameterList blockgcrodrList;
2126 blockgcrodrList.set(
"Num Blocks",numBlocks_);
2127 blockgcrodrList.set(
"Block Size",blockSize_);
2128 blockgcrodrList.set(
"Recycled Blocks",keff);
2130 Teuchos::RCP<BlockGCRODRIter<ScalarType,MV,OP> > block_gcrodr_iter;
2134 index.resize( blockSize_ );
2135 for(
int ii = 0; ii < blockSize_; ii++) index[ii] = ii;
2142 for(
int i=0; i < keff; i++){ index[i] = i;};
2143 B_ = rcp(
new SDM(Teuchos::View, *G_, keff, numBlocks_*blockSize_, 0, keff));
2144 H_ = rcp(
new SDM(Teuchos::View, *G_, (numBlocks_-1)*blockSize_ + blockSize_, (numBlocks_-1)*blockSize_, keff ,keff ));
2151 newstate.
curDim = blockSize_;
2152 block_gcrodr_iter -> initialize(newstate);
2154 int numRestarts = 0;
2158 block_gcrodr_iter -> iterate();
2163 if( convTest_->getStatus() ==
Passed ) {
2171 else if(maxIterTest_->getStatus() ==
Passed ){
2173 isConverged =
false;
2180 else if (block_gcrodr_iter->getCurSubspaceDim() == block_gcrodr_iter->getMaxSubspaceDim()){
2185 Teuchos::RCP<MV> update = block_gcrodr_iter->getCurrentUpdate();
2186 problem_->updateSolution(update,
true);
2187 buildRecycleSpaceAugKryl(block_gcrodr_iter);
2189 printer_->stream(
Debug) <<
" Generated new recycled subspace using RHS index " << currIdx[0] <<
" of dimension " << keff << std::endl << std::endl;
2191 if(numRestarts >= maxRestarts_) {
2192 isConverged =
false;
2198 printer_ -> stream(
Debug) <<
" Performing restart number " << numRestarts <<
" of " << maxRestarts_ << std::endl << std::endl;
2201 problem_->computeCurrPrecResVec( &*R_ );
2202 index.resize( blockSize_ );
2203 for (
int ii=0; ii<blockSize_; ++ii) index[ii] = ii;
2209 index.resize( numBlocks_*blockSize_ );
2210 for (
int ii=0; ii<(numBlocks_*blockSize_); ++ii) index[ii] = ii;
2212 index.resize( keff );
2213 for (
int ii=0; ii<keff; ++ii) index[ii] = ii;
2216 B_ = rcp(
new SDM(Teuchos::View, *G_, keff, (numBlocks_-1)*blockSize_, 0, keff));
2217 H_ = rcp(
new SDM(Teuchos::View, *G_, numBlocks_*blockSize_, (numBlocks_-1)*blockSize_, keff ,keff ));
2218 restartState.
B = B_;
2219 restartState.
H = H_;
2220 restartState.
curDim = blockSize_;
2221 block_gcrodr_iter->initialize(restartState);
2230 TEUCHOS_TEST_FOR_EXCEPTION(
true,std::logic_error,
"Belos::BlockGCRODRSolMgr::solve(): Invalid return from BlockGCRODRIter::iterate().");
2236 block_gcrodr_iter->updateLSQR( block_gcrodr_iter->getCurSubspaceDim() );
2238 sTest_->checkStatus( &*block_gcrodr_iter );
2239 if (convTest_->getStatus() !=
Passed) isConverged =
false;
2242 catch(
const std::exception &e){
2243 printer_->stream(
Errors) <<
"Error! Caught exception in BlockGCRODRIter::iterate() at iteration "
2244 << block_gcrodr_iter->getNumIters() << std::endl
2245 << e.what() << std::endl;
2252 Teuchos::RCP<MV> update = block_gcrodr_iter->getCurrentUpdate();
2253 problem_->updateSolution( update,
true );
2272 problem_->setCurrLS();
2275 startPtr += numCurrRHS;
2276 numRHS2Solve -= numCurrRHS;
2277 if ( numRHS2Solve > 0 ) {
2278 numCurrRHS = ( numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;
2279 if ( adaptiveBlockSize_ ) {
2280 blockSize_ = numCurrRHS;
2281 currIdx.resize( numCurrRHS );
2282 for (
int i=0; i<numCurrRHS; ++i) currIdx[i] = startPtr+i;
2285 currIdx.resize( blockSize_ );
2286 for (
int i=0; i<numCurrRHS; ++i) currIdx[i] = startPtr+i;
2287 for (
int i=numCurrRHS; i<blockSize_; ++i) currIdx[i] = -1;
2290 problem_->setLSIndex( currIdx );
2293 currIdx.resize( numRHS2Solve );
2297 if (!builtRecycleSpace_) {
2298 buildRecycleSpaceAugKryl(block_gcrodr_iter);
2299 printer_->stream(
Debug) <<
" Generated new recycled subspace using RHS index " << currIdx[0] <<
" of dimension " << keff << std::endl << std::endl;
2307 #ifdef BELOS_TEUCHOS_TIME_MONITOR
2313 numIters_ = maxIterTest_->getNumIters();
2316 const std::vector<MagnitudeType>* pTestValues = NULL;
2317 pTestValues = impConvTest_->getTestValue();
2318 TEUCHOS_TEST_FOR_EXCEPTION(pTestValues == NULL, std::logic_error,
2319 "Belos::BlockGCRODRSolMgr::solve(): The implicit convergence test's "
2320 "getTestValue() method returned NULL. Please report this bug to the "
2321 "Belos developers.");
2322 TEUCHOS_TEST_FOR_EXCEPTION(pTestValues->size() < 1, std::logic_error,
2323 "Belos::BlockGCRODRSolMgr::solve(): The implicit convergence test's "
2324 "getTestValue() method returned a vector of length zero. Please report "
2325 "this bug to the Belos developers.");
2329 achievedTol_ = *std::max_element (pTestValues->begin(), pTestValues->end());
Belos concrete class for performing the block, flexible GMRES iteration.
Belos concrete class for performing the block GCRO-DR (block GMRES with recycling) iteration.
Belos concrete class for performing the block GMRES iteration.
Belos header file which uses auto-configuration information to include necessary C++ headers.
Pure virtual base class which augments the basic interface for a Gmres linear solver iteration.
Class which describes the linear problem to be solved by the iterative solver.
Class which manages the output and verbosity of the Belos solvers.
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)
BlockGCRODRIterOrthoFailure is thrown when the BlockGCRODRIter object is unable to compute independen...
Implementation of the Block GCRO-DR (Block Recycling GMRES) iteration.
Thrown when an LAPACK call fails.
BlockGCRODRSolMgrLAPACKFailure(const std::string &what_arg)
BlockGCRODRSolMgrLinearProblemFailure(const std::string &what_arg)
Thrown when the solution manager was unable to orthogonalize the basis vectors.
BlockGCRODRSolMgrOrthoFailure(const std::string &what_arg)
BlockGCRODRSolMgrRecyclingFailure(const std::string &what_arg)
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > ¶ms)
Set the parameters the solver should use to solve the linear problem.
BlockGCRODRSolMgr()
Default constructor.
Teuchos::RCP< const Teuchos::ParameterList > getCurrentParameters() const
Get a parameter list containing the current parameters for this object.
std::string description() const
A description of the Block GCRODR solver manager.
MagnitudeType achievedTol() const
Get the residual for the most recent call to solve().
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Get a parameter list containing the valid parameters for this object.
void reset(const ResetType type)
Performs a reset of the solver manager specified by the ResetType.
virtual ~BlockGCRODRSolMgr()
Destructor.
bool isLOADetected() const
Whether a loss of accuracy was detected during the most recent solve.
ReturnType solve()
Solve the current linear problem.
void setProblem(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem)
Set the linear problem to solve on the next call to solve().
int getNumIters() const
Get the iteration count for the most recent call to solve().
const LinearProblem< ScalarType, MV, OP > & getProblem() const
Get current linear problem being solved for in this object.
This class implements the block GMRES iteration, where a block Krylov subspace is constructed....
An implementation of the Belos::MatOrthoManager that performs orthogonalization using (potentially) m...
void setDepTol(const MagnitudeType dep_tol)
Set parameter for re-orthogonalization threshhold.
GmresIterationOrthoFailure is thrown when the GmresIteration object is unable to compute independent ...
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 ptrdiff_t GetGlobalLength(const MV &mv)
Return the number of rows in the given multivector mv.
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 MvInit(MV &mv, const ScalarType alpha=Teuchos::ScalarTraits< ScalarType >::zero())
Replace each element of the vectors in mv with alpha.
static void SetBlock(const MV &A, const std::vector< int > &index, MV &mv)
Copy the vectors in A to a set of vectors in mv indicated by the indices given in index.
static int GetNumberVecs(const MV &mv)
Obtain the number of vectors in mv.
static void Assign(const MV &A, MV &mv)
mv := A
Class which defines basic traits for the operator type.
Enumeration of all valid Belos (Mat)OrthoManager classes.
Belos's basic output manager for sending information of select verbosity levels to the appropriate ou...
SolverManager()
Empty constructor.
A class for extending the status testing capabilities of Belos via logical combinations.
Exception thrown to signal error in a status test during Belos::StatusTest::checkStatus().
An implementation of StatusTestResNorm using a family of residual norms.
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.
ScaleType convertStringToScaleType(const std::string &scaleType)
Convert the given string to its ScaleType enum value.
ReturnType
Whether the Belos solve converged for all linear systems.
const std::string BlockGCRODRSolMgr< ScalarType, MV, OP >::recycleMethod_default_
OutputType
Style of output used to display status test information.
ScaleType
The type of scaling to use on the residual norm value.
const bool BlockGCRODRSolMgr< ScalarType, MV, OP >::adaptiveBlockSize_default_
ResetType
How to reset the solver.
Structure to contain pointers to BlockGCRODRIter state variables.
Teuchos::RCP< MV > V
The current Krylov basis.
Teuchos::RCP< MV > U
The recycled subspace and its projection.
Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > B
The projection of the Krylov subspace against the recycled subspace *.
Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > H
The current Hessenberg matrix.
int curDim
The current dimension of the reduction.
Structure to contain pointers to GmresIteration state variables.
Teuchos::RCP< const MV > V
The current Krylov basis.
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > z
The current right-hand side of the least squares system RY = Z.
int curDim
The current dimension of the reduction.
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > H
The current Hessenberg matrix.