42 #ifndef BELOS_BLOCK_GMRES_SOLMGR_HPP 
   43 #define BELOS_BLOCK_GMRES_SOLMGR_HPP 
   69 #ifdef BELOS_TEUCHOS_TIME_MONITOR 
  124 template<
class ScalarType, 
class MV, 
class OP>
 
  201     return Teuchos::tuple(timerSolve_);
 
  299   bool checkStatusTest();
 
  323   static constexpr 
int maxRestarts_default_ = 20;
 
  324   static constexpr 
int maxIters_default_ = 1000;
 
  325   static constexpr 
bool adaptiveBlockSize_default_ = 
true;
 
  326   static constexpr 
bool showMaxResNormOnly_default_ = 
false;
 
  327   static constexpr 
bool flexibleGmres_default_ = 
false;
 
  328   static constexpr 
bool expResTest_default_ = 
false;
 
  329   static constexpr 
int blockSize_default_ = 1;
 
  330   static constexpr 
int numBlocks_default_ = 300;
 
  333   static constexpr 
int outputFreq_default_ = -1;
 
  334   static constexpr 
const char * impResScale_default_ = 
"Norm of Preconditioned Initial Residual";
 
  335   static constexpr 
const char * expResScale_default_ = 
"Norm of Initial Residual";
 
  336   static constexpr 
const char * label_default_ = 
"Belos";
 
  337   static constexpr 
const char * orthoType_default_ = 
"DGKS";
 
  338   static constexpr std::ostream * outputStream_default_ = &std::cout;
 
  341   MagnitudeType convtol_, orthoKappa_, achievedTol_;
 
  342   int maxRestarts_, maxIters_, numIters_;
 
  343   int blockSize_, numBlocks_, verbosity_, outputStyle_, outputFreq_;
 
  344   bool adaptiveBlockSize_, showMaxResNormOnly_, isFlexible_, expResTest_;
 
  345   std::string orthoType_;
 
  346   std::string impResScale_, expResScale_;
 
  353   bool isSet_, isSTSet_;
 
  359 template<
class ScalarType, 
class MV, 
class OP>
 
  361   outputStream_(Teuchos::
rcp(outputStream_default_,false)),
 
  364   achievedTol_(Teuchos::ScalarTraits<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>::zero()),
 
  365   maxRestarts_(maxRestarts_default_),
 
  366   maxIters_(maxIters_default_),
 
  368   blockSize_(blockSize_default_),
 
  369   numBlocks_(numBlocks_default_),
 
  370   verbosity_(verbosity_default_),
 
  371   outputStyle_(outputStyle_default_),
 
  372   outputFreq_(outputFreq_default_),
 
  373   adaptiveBlockSize_(adaptiveBlockSize_default_),
 
  374   showMaxResNormOnly_(showMaxResNormOnly_default_),
 
  375   isFlexible_(flexibleGmres_default_),
 
  376   expResTest_(expResTest_default_),
 
  377   orthoType_(orthoType_default_),
 
  378   impResScale_(impResScale_default_),
 
  379   expResScale_(expResScale_default_),
 
  380   label_(label_default_),
 
  388 template<
class ScalarType, 
class MV, 
class OP>
 
  393   outputStream_(Teuchos::
rcp(outputStream_default_,false)),
 
  396   achievedTol_(Teuchos::ScalarTraits<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>::zero()),
 
  397   maxRestarts_(maxRestarts_default_),
 
  398   maxIters_(maxIters_default_),
 
  400   blockSize_(blockSize_default_),
 
  401   numBlocks_(numBlocks_default_),
 
  402   verbosity_(verbosity_default_),
 
  403   outputStyle_(outputStyle_default_),
 
  404   outputFreq_(outputFreq_default_),
 
  405   adaptiveBlockSize_(adaptiveBlockSize_default_),
 
  406   showMaxResNormOnly_(showMaxResNormOnly_default_),
 
  407   isFlexible_(flexibleGmres_default_),
 
  408   expResTest_(expResTest_default_),
 
  409   orthoType_(orthoType_default_),
 
  410   impResScale_(impResScale_default_),
 
  411   expResScale_(expResScale_default_),
 
  412   label_(label_default_),
 
  428 template<
class ScalarType, 
class MV, 
class OP>
 
  439       "The relative residual tolerance that needs to be achieved by the\n" 
  440       "iterative solver in order for the linear system to be declared converged." );
 
  441     pl->
set(
"Maximum Restarts", static_cast<int>(maxRestarts_default_),
 
  442       "The maximum number of restarts allowed for each\n" 
  443       "set of RHS solved.");
 
  444     pl->
set(
"Maximum Iterations", static_cast<int>(maxIters_default_),
 
  445       "The maximum number of block iterations allowed for each\n" 
  446       "set of RHS solved.");
 
  447     pl->
set(
"Num Blocks", static_cast<int>(numBlocks_default_),
 
  448       "The maximum number of blocks allowed in the Krylov subspace\n" 
  449       "for each set of RHS solved.");
 
  450     pl->
set(
"Block Size", static_cast<int>(blockSize_default_),
 
  451       "The number of vectors in each block.  This number times the\n" 
  452       "number of blocks is the total Krylov subspace dimension.");
 
  453     pl->
set(
"Adaptive Block Size", static_cast<bool>(adaptiveBlockSize_default_),
 
  454       "Whether the solver manager should adapt the block size\n" 
  455       "based on the number of RHS to solve.");
 
  456     pl->
set(
"Verbosity", static_cast<int>(verbosity_default_),
 
  457       "What type(s) of solver information should be outputted\n" 
  458       "to the output stream.");
 
  459     pl->
set(
"Output Style", static_cast<int>(outputStyle_default_),
 
  460       "What style is used for the solver information outputted\n" 
  461       "to the output stream.");
 
  462     pl->
set(
"Output Frequency", static_cast<int>(outputFreq_default_),
 
  463       "How often convergence information should be outputted\n" 
  464       "to the output stream.");
 
  466       "A reference-counted pointer to the output stream where all\n" 
  467       "solver output is sent.");
 
  468     pl->
set(
"Show Maximum Residual Norm Only", static_cast<bool>(showMaxResNormOnly_default_),
 
  469       "When convergence information is printed, only show the maximum\n" 
  470       "relative residual norm when the block size is greater than one.");
 
  471     pl->
set(
"Flexible Gmres", static_cast<bool>(flexibleGmres_default_),
 
  472       "Whether the solver manager should use the flexible variant\n" 
  474     pl->
set(
"Explicit Residual Test", static_cast<bool>(expResTest_default_),
 
  475       "Whether the explicitly computed residual should be used in the convergence test.");
 
  476     pl->
set(
"Implicit Residual Scaling", static_cast<const char *>(impResScale_default_),
 
  477       "The type of scaling used in the implicit residual convergence test.");
 
  478     pl->
set(
"Explicit Residual Scaling", static_cast<const char *>(expResScale_default_),
 
  479       "The type of scaling used in the explicit residual convergence test.");
 
  480     pl->
set(
"Timer Label", static_cast<const char *>(label_default_),
 
  481       "The string to use as a prefix for the timer labels.");
 
  482     pl->
set(
"Orthogonalization", static_cast<const char *>(orthoType_default_),
 
  483       "The type of orthogonalization to use: DGKS, ICGS, or IMGS.");
 
  485       "The constant used by DGKS orthogonalization to determine\n" 
  486       "whether another step of classical Gram-Schmidt is necessary.");
 
  493 template<
class ScalarType, 
class MV, 
class OP>
 
  498   if (params_ == Teuchos::null) {
 
  507     maxRestarts_ = params->
get(
"Maximum Restarts",maxRestarts_default_);
 
  510     params_->set(
"Maximum Restarts", maxRestarts_);
 
  515     maxIters_ = params->
get(
"Maximum Iterations",maxIters_default_);
 
  518     params_->set(
"Maximum Iterations", maxIters_);
 
  519     if (maxIterTest_!=Teuchos::null)
 
  520       maxIterTest_->setMaxIters( maxIters_ );
 
  525     blockSize_ = params->
get(
"Block Size",blockSize_default_);
 
  527       "Belos::BlockGmresSolMgr: \"Block Size\" must be strictly positive.");
 
  530     params_->set(
"Block Size", blockSize_);
 
  535     adaptiveBlockSize_ = params->
get(
"Adaptive Block Size",adaptiveBlockSize_default_);
 
  538     params_->set(
"Adaptive Block Size", adaptiveBlockSize_);
 
  543     numBlocks_ = params->
get(
"Num Blocks",numBlocks_default_);
 
  545       "Belos::BlockGmresSolMgr: \"Num Blocks\" must be strictly positive.");
 
  548     params_->set(
"Num Blocks", numBlocks_);
 
  553     std::string tempLabel = params->
get(
"Timer Label", label_default_);
 
  556     if (tempLabel != label_) {
 
  558       params_->set(
"Timer Label", label_);
 
  559       std::string solveLabel = label_ + 
": BlockGmresSolMgr total solve time";
 
  560 #ifdef BELOS_TEUCHOS_TIME_MONITOR 
  563       if (ortho_ != Teuchos::null) {
 
  564         ortho_->setLabel( label_ );
 
  571     isFlexible_ = Teuchos::getParameter<bool>(*params,
"Flexible Gmres");
 
  572     params_->set(
"Flexible Gmres", isFlexible_);
 
  573     if (isFlexible_ && expResTest_) {
 
  582     std::string tempOrthoType = params->
get(
"Orthogonalization",orthoType_default_);
 
  584                         std::invalid_argument,
 
  585                         "Belos::BlockGmresSolMgr: \"Orthogonalization\" must be either \"DGKS\", \"ICGS\", or \"IMGS\".");
 
  586     if (tempOrthoType != orthoType_) {
 
  587       orthoType_ = tempOrthoType;
 
  588       params_->set(
"Orthogonalization", orthoType_);
 
  590       if (orthoType_==
"DGKS") {
 
  591         if (orthoKappa_ <= 0) {
 
  599       else if (orthoType_==
"ICGS") {
 
  602       else if (orthoType_==
"IMGS") {
 
  609   if (params->
isParameter(
"Orthogonalization Constant")) {
 
  610     if (params->
isType<MagnitudeType> (
"Orthogonalization Constant")) {
 
  611       orthoKappa_ = params->
get (
"Orthogonalization Constant",
 
  615       orthoKappa_ = params->
get (
"Orthogonalization Constant",
 
  620     params_->set(
"Orthogonalization Constant",orthoKappa_);
 
  621     if (orthoType_==
"DGKS") {
 
  622       if (orthoKappa_ > 0 && ortho_ != Teuchos::null) {
 
  630     if (Teuchos::isParameterType<int>(*params,
"Verbosity")) {
 
  631       verbosity_ = params->
get(
"Verbosity", verbosity_default_);
 
  633       verbosity_ = (int)Teuchos::getParameter<Belos::MsgType>(*params,
"Verbosity");
 
  637     params_->set(
"Verbosity", verbosity_);
 
  638     if (printer_ != Teuchos::null)
 
  639       printer_->setVerbosity(verbosity_);
 
  644     if (Teuchos::isParameterType<int>(*params,
"Output Style")) {
 
  645       outputStyle_ = params->
get(
"Output Style", outputStyle_default_);
 
  647       outputStyle_ = (int)Teuchos::getParameter<Belos::OutputType>(*params,
"Output Style");
 
  651     params_->set(
"Output Style", outputStyle_);
 
  652     if (outputTest_ != Teuchos::null) {
 
  659     outputStream_ = Teuchos::getParameter<Teuchos::RCP<std::ostream> >(*params,
"Output Stream");
 
  662     params_->set(
"Output Stream", outputStream_);
 
  663     if (printer_ != Teuchos::null)
 
  664       printer_->setOStream( outputStream_ );
 
  670       outputFreq_ = params->
get(
"Output Frequency", outputFreq_default_);
 
  674     params_->set(
"Output Frequency", outputFreq_);
 
  675     if (outputTest_ != Teuchos::null)
 
  676       outputTest_->setOutputFrequency( outputFreq_ );
 
  680   if (printer_ == Teuchos::null) {
 
  685   if (params->
isParameter(
"Convergence Tolerance")) {
 
  686     if (params->
isType<MagnitudeType> (
"Convergence Tolerance")) {
 
  687       convtol_ = params->
get (
"Convergence Tolerance",
 
  695     params_->set(
"Convergence Tolerance", convtol_);
 
  696     if (impConvTest_ != Teuchos::null)
 
  697       impConvTest_->setTolerance( convtol_ );
 
  698     if (expConvTest_ != Teuchos::null)
 
  699       expConvTest_->setTolerance( convtol_ );
 
  703   if (params->
isParameter(
"Implicit Residual Scaling")) {
 
  704     std::string tempImpResScale = Teuchos::getParameter<std::string>( *params, 
"Implicit Residual Scaling" );
 
  707     if (impResScale_ != tempImpResScale) {
 
  709       impResScale_ = tempImpResScale;
 
  712       params_->set(
"Implicit Residual Scaling", impResScale_);
 
  713       if (impConvTest_ != Teuchos::null) {
 
  717         catch (std::exception& e) {
 
  725   if (params->
isParameter(
"Explicit Residual Scaling")) {
 
  726     std::string tempExpResScale = Teuchos::getParameter<std::string>( *params, 
"Explicit Residual Scaling" );
 
  729     if (expResScale_ != tempExpResScale) {
 
  731       expResScale_ = tempExpResScale;
 
  734       params_->set(
"Explicit Residual Scaling", expResScale_);
 
  735       if (expConvTest_ != Teuchos::null) {
 
  739         catch (std::exception& e) {
 
  747   if (params->
isParameter(
"Explicit Residual Test")) {
 
  748     expResTest_ = Teuchos::getParameter<bool>( *params,
"Explicit Residual Test" );
 
  751     params_->set(
"Explicit Residual Test", expResTest_);
 
  752     if (expConvTest_ == Teuchos::null) {
 
  758   if (params->
isParameter(
"Show Maximum Residual Norm Only")) {
 
  759     showMaxResNormOnly_ = Teuchos::getParameter<bool>(*params,
"Show Maximum Residual Norm Only");
 
  762     params_->set(
"Show Maximum Residual Norm Only", showMaxResNormOnly_);
 
  763     if (impConvTest_ != Teuchos::null)
 
  764       impConvTest_->setShowMaxResNormOnly( showMaxResNormOnly_ );
 
  765     if (expConvTest_ != Teuchos::null)
 
  766       expConvTest_->setShowMaxResNormOnly( showMaxResNormOnly_ );
 
  770   if (ortho_ == Teuchos::null) {
 
  771     params_->set(
"Orthogonalization", orthoType_);
 
  772     if (orthoType_==
"DGKS") {
 
  773       if (orthoKappa_ <= 0) {
 
  781     else if (orthoType_==
"ICGS") {
 
  784     else if (orthoType_==
"IMGS") {
 
  789         "Belos::BlockGmresSolMgr(): Invalid orthogonalization type.");
 
  794   if (timerSolve_ == Teuchos::null) {
 
  795     std::string solveLabel = label_ + 
": BlockGmresSolMgr total solve time";
 
  796 #ifdef BELOS_TEUCHOS_TIME_MONITOR 
  806 template<
class ScalarType, 
class MV, 
class OP>
 
  821     params_->set(
"Flexible Gmres", isFlexible_);
 
  825       "Belos::BlockGmresSolMgr::solve(): Linear problem has a left preconditioner, not a right preconditioner, which is required for flexible GMRES.");
 
  840     tmpImpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
 
  841     impConvTest_ = tmpImpConvTest;
 
  846     tmpExpConvTest->defineResForm( StatusTestGenResNorm_t::Explicit, 
Belos::TwoNorm );
 
  848     tmpExpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
 
  849     expConvTest_ = tmpExpConvTest;
 
  852     convTest_ = 
Teuchos::rcp( 
new StatusTestCombo_t( StatusTestCombo_t::SEQ, impConvTest_, expConvTest_ ) );
 
  861       tmpImpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
 
  862       impConvTest_ = tmpImpConvTest;
 
  870       tmpImpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
 
  871       impConvTest_ = tmpImpConvTest;
 
  875     expConvTest_ = impConvTest_;
 
  876     convTest_ = impConvTest_;
 
  880   sTest_ = 
Teuchos::rcp( 
new StatusTestCombo_t( StatusTestCombo_t::OR, maxIterTest_, convTest_ ) );
 
  883   if (
nonnull(debugStatusTest_) ) {
 
  885     Teuchos::rcp_dynamic_cast<StatusTestCombo_t>(sTest_)->addStatusTest( debugStatusTest_ );
 
  890   StatusTestOutputFactory<ScalarType,MV,OP> stoFactory( outputStyle_ );
 
  894   std::string solverDesc = 
" Block Gmres ";
 
  896     solverDesc = 
"Flexible" + solverDesc;
 
  897   outputTest_->setSolverDesc( solverDesc );
 
  905 template<
class ScalarType, 
class MV, 
class OP>
 
  910   debugStatusTest_ = debugStatusTest;
 
  915 template<
class ScalarType, 
class MV, 
class OP>
 
  922     setParameters(Teuchos::parameterList(*getValidParameters()));
 
  928     "Belos::BlockGmresSolMgr::solve(): Linear problem is not a valid object.");
 
  931     "Belos::BlockGmresSolMgr::solve(): Linear problem is not ready, setProblem() has not been called.");
 
  933   if (!isSTSet_ || (!expResTest_ && !
Teuchos::is_null(problem_->getLeftPrec())) ) {
 
  935       "Belos::BlockGmresSolMgr::solve(): Linear problem and requested status tests are incompatible.");
 
  940   int numRHS2Solve = MVT::GetNumberVecs( *(problem_->getRHS()) );
 
  941   int numCurrRHS = ( numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;
 
  943   std::vector<int> currIdx;
 
  946   if ( adaptiveBlockSize_ ) {
 
  947     blockSize_ = numCurrRHS;
 
  948     currIdx.resize( numCurrRHS  );
 
  949     for (
int i=0; i<numCurrRHS; ++i)
 
  950     { currIdx[i] = startPtr+i; }
 
  954     currIdx.resize( blockSize_ );
 
  955     for (
int i=0; i<numCurrRHS; ++i)
 
  956     { currIdx[i] = startPtr+i; }
 
  957     for (
int i=numCurrRHS; i<blockSize_; ++i)
 
  962   problem_->setLSIndex( currIdx );
 
  967   plist.
set(
"Block Size",blockSize_);
 
  969   ptrdiff_t dim = MVT::GetGlobalLength( *(problem_->getRHS()) );
 
  970   if (blockSize_*static_cast<ptrdiff_t>(numBlocks_) > dim) {
 
  971     int tmpNumBlocks = 0;
 
  973       tmpNumBlocks = dim / blockSize_;  
 
  975       tmpNumBlocks = ( dim - blockSize_) / blockSize_;  
 
  977       "Belos::BlockGmresSolMgr::solve():  Warning! Requested Krylov subspace dimension is larger than operator dimension!" 
  978       << std::endl << 
" The maximum number of blocks allowed for the Krylov subspace will be adjusted to " << tmpNumBlocks << std::endl;
 
  979     plist.
set(
"Num Blocks",tmpNumBlocks);
 
  982     plist.
set(
"Num Blocks",numBlocks_);
 
  985   outputTest_->reset();
 
  986   loaDetected_ = 
false;
 
  989   bool isConverged = 
true;
 
 1003 #ifdef BELOS_TEUCHOS_TIME_MONITOR 
 1007     while ( numRHS2Solve > 0 ) {
 
 1010       if (blockSize_*numBlocks_ > dim) {
 
 1011         int tmpNumBlocks = 0;
 
 1012         if (blockSize_ == 1)
 
 1013           tmpNumBlocks = dim / blockSize_;  
 
 1015           tmpNumBlocks = ( dim - blockSize_) / blockSize_;  
 
 1016         block_gmres_iter->setSize( blockSize_, tmpNumBlocks );
 
 1019         block_gmres_iter->setSize( blockSize_, numBlocks_ );
 
 1022       block_gmres_iter->resetNumIters();
 
 1025       outputTest_->resetNumCalls();
 
 1031         if (currIdx[blockSize_-1] == -1) {
 
 1032           V_0 = MVT::Clone( *(problem_->getInitResVec()), blockSize_ );
 
 1033           problem_->computeCurrResVec( &*V_0 );
 
 1036           V_0 = MVT::CloneCopy( *(problem_->getInitResVec()), currIdx );
 
 1041         if (currIdx[blockSize_-1] == -1) {
 
 1042           V_0 = MVT::Clone( *(problem_->getInitPrecResVec()), blockSize_ );
 
 1043           problem_->computeCurrPrecResVec( &*V_0 );
 
 1046           V_0 = MVT::CloneCopy( *(problem_->getInitPrecResVec()), currIdx );
 
 1055       int rank = ortho_->normalize( *V_0, z_0 );
 
 1057         "Belos::BlockGmresSolMgr::solve(): Failed to compute initial block of orthonormal vectors.");
 
 1064       block_gmres_iter->initializeGmres(newstate);
 
 1065       int numRestarts = 0;
 
 1070           block_gmres_iter->iterate();
 
 1077           if ( convTest_->getStatus() == 
Passed ) {
 
 1078             if ( expConvTest_->getLOADetected() ) {
 
 1080               loaDetected_ = 
true;
 
 1082                 "Belos::BlockGmresSolMgr::solve(): Warning! Solver has experienced a loss of accuracy!" << std::endl;
 
 1083               isConverged = 
false;
 
 1092           else if ( maxIterTest_->getStatus() == 
Passed ) {
 
 1094             isConverged = 
false;
 
 1102           else if ( block_gmres_iter->getCurSubspaceDim() == block_gmres_iter->getMaxSubspaceDim() ) {
 
 1104             if ( numRestarts >= maxRestarts_ ) {
 
 1105               isConverged = 
false;
 
 1110             printer_->stream(
Debug) << 
" Performing restart number " << numRestarts << 
" of " << maxRestarts_ << std::endl << std::endl;
 
 1117               MVT::MvAddMv( 1.0, *curX, 1.0, *update, *curX );
 
 1120               problem_->updateSolution( update, 
true );
 
 1127             V_0  = MVT::Clone( *(oldState.
V), blockSize_ );
 
 1129               problem_->computeCurrResVec( &*V_0 );
 
 1131               problem_->computeCurrPrecResVec( &*V_0 );
 
 1137             rank = ortho_->normalize( *V_0, z_0 );
 
 1139               "Belos::BlockGmresSolMgr::solve(): Failed to compute initial block of orthonormal vectors after restart.");
 
 1145             block_gmres_iter->initializeGmres(newstate);
 
 1158               "Belos::BlockGmresSolMgr::solve(): Invalid return from BlockGmresIter::iterate().");
 
 1163           if (blockSize_ != 1) {
 
 1164             printer_->stream(
Errors) << 
"Error! Caught std::exception in BlockGmresIter::iterate() at iteration " 
 1165                                      << block_gmres_iter->getNumIters() << std::endl
 
 1166                                      << e.what() << std::endl;
 
 1167             if (convTest_->getStatus() != 
Passed)
 
 1168               isConverged = 
false;
 
 1173             block_gmres_iter->updateLSQR( block_gmres_iter->getCurSubspaceDim() );
 
 1176             sTest_->checkStatus( &*block_gmres_iter );
 
 1177             if (convTest_->getStatus() != 
Passed)
 
 1178               isConverged = 
false;
 
 1182         catch (
const std::exception &e) {
 
 1183           printer_->stream(
Errors) << 
"Error! Caught std::exception in BlockGmresIter::iterate() at iteration " 
 1184                                    << block_gmres_iter->getNumIters() << std::endl
 
 1185                                    << e.what() << std::endl;
 
 1197         if (update != Teuchos::null)
 
 1198           MVT::MvAddMv( 1.0, *curX, 1.0, *update, *curX );
 
 1205           MVT::MvAddMv( 0.0, *newX, 1.0, *newX, *curX );
 
 1209           problem_->updateSolution( update, 
true );
 
 1214       problem_->setCurrLS();
 
 1217       startPtr += numCurrRHS;
 
 1218       numRHS2Solve -= numCurrRHS;
 
 1219       if ( numRHS2Solve > 0 ) {
 
 1220         numCurrRHS = ( numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;
 
 1222         if ( adaptiveBlockSize_ ) {
 
 1223           blockSize_ = numCurrRHS;
 
 1224           currIdx.resize( numCurrRHS  );
 
 1225           for (
int i=0; i<numCurrRHS; ++i)
 
 1226           { currIdx[i] = startPtr+i; }
 
 1229           currIdx.resize( blockSize_ );
 
 1230           for (
int i=0; i<numCurrRHS; ++i)
 
 1231           { currIdx[i] = startPtr+i; }
 
 1232           for (
int i=numCurrRHS; i<blockSize_; ++i)
 
 1233           { currIdx[i] = -1; }
 
 1236         problem_->setLSIndex( currIdx );
 
 1239         currIdx.resize( numRHS2Solve );
 
 1250 #ifdef BELOS_TEUCHOS_TIME_MONITOR 
 1259   numIters_ = maxIterTest_->getNumIters();
 
 1273     const std::vector<MagnitudeType>* pTestValues = NULL;
 
 1275       pTestValues = expConvTest_->getTestValue();
 
 1276       if (pTestValues == NULL || pTestValues->size() < 1) {
 
 1277         pTestValues = impConvTest_->getTestValue();
 
 1282       pTestValues = impConvTest_->getTestValue();
 
 1285       "Belos::BlockGmresSolMgr::solve(): The implicit convergence test's " 
 1286       "getTestValue() method returned NULL.  Please report this bug to the " 
 1287       "Belos developers.");
 
 1289       "Belos::BlockGmresSolMgr::solve(): The implicit convergence test's " 
 1290       "getTestValue() method returned a vector of length zero.  Please report " 
 1291       "this bug to the Belos developers.");
 
 1296     achievedTol_ = *std::max_element (pTestValues->begin(), pTestValues->end());
 
 1299   if (!isConverged || loaDetected_) {
 
 1306 template<
class ScalarType, 
class MV, 
class OP>
 
 1309   std::ostringstream out;
 
 1310   out << 
"\"Belos::BlockGmresSolMgr\": {";
 
 1311   if (this->getObjectLabel () != 
"") {
 
 1312     out << 
"Label: " << this->getObjectLabel () << 
", ";
 
 1314   out << 
"Flexible: " << (isFlexible_ ? 
"true" : 
"false")
 
 1315       << 
", Num Blocks: " << numBlocks_
 
 1316       << 
", Maximum Iterations: " << maxIters_
 
 1317       << 
", Maximum Restarts: " << maxRestarts_
 
 1318       << 
", Convergence Tolerance: " << convtol_
 
 1324 template<
class ScalarType, 
class MV, 
class OP>
 
 1346     out << 
"\"Belos::BlockGmresSolMgr\":" << endl;
 
 1348     out << 
"Template parameters:" << endl;
 
 1351       out << 
"ScalarType: " << TypeNameTraits<ScalarType>::name () << endl
 
 1352           << 
"MV: " << TypeNameTraits<MV>::name () << endl
 
 1353           << 
"OP: " << TypeNameTraits<OP>::name () << endl;
 
 1355     if (this->getObjectLabel () != 
"") {
 
 1356       out << 
"Label: " << this->getObjectLabel () << endl;
 
 1358     out << 
"Flexible: " << (isFlexible_ ? 
"true" : 
"false") << endl
 
 1359         << 
"Num Blocks: " << numBlocks_ << endl
 
 1360         << 
"Maximum Iterations: " << maxIters_ << endl
 
 1361         << 
"Maximum Restarts: " << maxRestarts_ << endl
 
 1362         << 
"Convergence Tolerance: " << convtol_ << endl;
 
static const double orthoKappa
DGKS orthogonalization constant. 
 
ScaleType convertStringToScaleType(const std::string &scaleType)
Convert the given string to its ScaleType enum value. 
 
Teuchos::Array< Teuchos::RCP< Teuchos::Time > > getTimers() const 
Return the timers for this object. 
 
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. 
 
bool is_null(const boost::shared_ptr< T > &p)
 
void setDebugStatusTest(const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &debugStatusTest) override
Set a debug status test that will be checked at the same time as the top-level status test...
 
ScaleType
The type of scaling to use on the residual norm value. 
 
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. 
 
T & get(const std::string &name, T def_value)
 
Virtual base class for StatusTest that printing status tests. 
 
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
 
bool is_null(const std::shared_ptr< T > &p)
 
MagnitudeType achievedTol() const override
Tolerance achieved by the last solve() invocation. 
 
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
 
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
Get a parameter list containing the valid parameters for this object. 
 
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. 
 
BlockGmresSolMgrLinearProblemFailure is thrown when the linear problem is not setup (i...
 
Structure to contain pointers to GmresIteration state variables. 
 
static const double convTol
Default convergence tolerance. 
 
Belos::StatusTest class for specifying a maximum number of iterations. 
 
Interface to Block GMRES and Flexible GMRES. 
 
A pure virtual class for defining the status tests for the Belos iterative solvers. 
 
virtual ~BlockGmresSolMgr()
Destructor. 
 
void setProblem(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem) override
Set the linear problem that needs to be solved. 
 
A factory class for generating StatusTestOutput objects. 
 
Iterated Modified Gram-Schmidt (IMGS) implementation of the Belos::OrthoManager class. 
 
ReturnType solve() override
This method performs possibly repeated calls to the underlying linear solver's iterate() routine unti...
 
Traits class which defines basic operations on multivectors. 
 
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > ¶ms) override
Set the parameters the solver manager should use to solve the linear problem. 
 
Belos::StatusTest for logically combining several status tests. 
 
bool isParameter(const std::string &name) const 
 
Classical Gram-Schmidt (with DGKS correction) implementation of the Belos::OrthoManager class...
 
A Belos::StatusTest class for specifying a maximum number of iterations. 
 
ResetType
How to reset the solver. 
 
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
 
Pure virtual base class which describes the basic interface for a solver manager. ...
 
int curDim
The current dimension of the reduction. 
 
static void summarize(Ptr< const Comm< int > > comm, std::ostream &out=std::cout, const bool alwaysWriteLocal=false, const bool writeGlobalStats=true, const bool writeZeroTimers=true, const ECounterSetOp setOp=Intersection, const std::string &filter="", const bool ignoreZeroTimers=false)
 
BlockGmresSolMgrLinearProblemFailure(const std::string &what_arg)
 
void reset(const ResetType type) override
Performs a reset of the solver manager specified by the ResetType. This informs the solver manager th...
 
A linear system to solve, and its associated information. 
 
Class which describes the linear problem to be solved by the iterative solver. 
 
This class implements the block flexible GMRES iteration, where a block Krylov subspace is constructe...
 
int getNumIters() const override
Get the iteration count for the most recent call to solve(). 
 
Teuchos::RCP< SolverManager< ScalarType, MV, OP > > clone() const override
clone for Inverted Injection (DII) 
 
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 ...
 
Belos concrete class for performing the block GMRES iteration. 
 
Iterated Classical Gram-Schmidt (ICGS) implementation of the Belos::OrthoManager class. 
 
void validateParameters(ParameterList const &validParamList, int const depth=1000, EValidateUsed const validateUsed=VALIDATE_USED_ENABLED, EValidateDefaults const validateDefaults=VALIDATE_DEFAULTS_ENABLED) const 
 
An implementation of the Belos::MatOrthoManager that performs orthogonalization using multiple steps ...
 
static const EVerbosityLevel verbLevel_default
 
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const override
Print the object with the given verbosity level to a FancyOStream. 
 
An implementation of the Belos::MatOrthoManager that performs orthogonalization using multiple steps ...
 
bool nonnull(const boost::shared_ptr< T > &p)
 
bool isLOADetected() const override
Return whether a loss of accuracy was detected by this solver during the most current solve...
 
Belos::StatusTestResNorm for specifying general residual norm stopping criteria. 
 
Belos::StatusTest for specifying an implicit residual norm stopping criteria that checks for loss of ...
 
BlockGmresSolMgr()
Empty constructor for BlockGmresSolMgr. This constructor takes no arguments and sets the default valu...
 
bool isType(const std::string &name) const 
 
An implementation of the Belos::MatOrthoManager that performs orthogonalization using (potentially) m...
 
BlockGmresSolMgrOrthoFailure(const std::string &what_arg)
 
A class for extending the status testing capabilities of Belos via logical combinations. 
 
std::string description() const override
Return a one-line description of this object. 
 
Class which defines basic traits for the operator type. 
 
Parent class to all Belos exceptions. 
 
BlockGmresSolMgrOrthoFailure is thrown when the orthogonalization manager is unable to generate ortho...
 
Default parameters common to most Belos solvers. 
 
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 ...
 
Belos header file which uses auto-configuration information to include necessary C++ headers...
 
Teuchos::RCP< const Teuchos::ParameterList > getCurrentParameters() const override
Get a parameter list containing the current parameters for this object. 
 
const LinearProblem< ScalarType, MV, OP > & getProblem() const override
Get current linear problem being solved for in this object. 
 
Belos concrete class for performing the block, flexible GMRES iteration.