Belos Version of the Day
BelosBlockGmresSolMgr.hpp
Go to the documentation of this file.
00001 //@HEADER
00002 // ************************************************************************
00003 //
00004 //                 Belos: Block Linear Solvers Package
00005 //                  Copyright 2004 Sandia Corporation
00006 //
00007 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
00008 // the U.S. Government retains certain rights in this software.
00009 //
00010 // Redistribution and use in source and binary forms, with or without
00011 // modification, are permitted provided that the following conditions are
00012 // met:
00013 //
00014 // 1. Redistributions of source code must retain the above copyright
00015 // notice, this list of conditions and the following disclaimer.
00016 //
00017 // 2. Redistributions in binary form must reproduce the above copyright
00018 // notice, this list of conditions and the following disclaimer in the
00019 // documentation and/or other materials provided with the distribution.
00020 //
00021 // 3. Neither the name of the Corporation nor the names of the
00022 // contributors may be used to endorse or promote products derived from
00023 // this software without specific prior written permission.
00024 //
00025 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
00026 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00027 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
00028 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
00029 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
00030 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
00031 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
00032 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
00033 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
00034 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
00035 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00036 //
00037 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
00038 //
00039 // ************************************************************************
00040 //@HEADER
00041 
00042 #ifndef BELOS_BLOCK_GMRES_SOLMGR_HPP
00043 #define BELOS_BLOCK_GMRES_SOLMGR_HPP
00044 
00049 #include "BelosConfigDefs.hpp"
00050 #include "BelosTypes.hpp"
00051 
00052 #include "BelosLinearProblem.hpp"
00053 #include "BelosSolverManager.hpp"
00054 
00055 #include "BelosGmresIteration.hpp"
00056 #include "BelosBlockGmresIter.hpp"
00057 #include "BelosBlockFGmresIter.hpp"
00058 #include "BelosDGKSOrthoManager.hpp"
00059 #include "BelosICGSOrthoManager.hpp"
00060 #include "BelosIMGSOrthoManager.hpp"
00061 #include "BelosStatusTestMaxIters.hpp"
00062 #include "BelosStatusTestGenResNorm.hpp"
00063 #include "BelosStatusTestImpResNorm.hpp"
00064 #include "BelosStatusTestCombo.hpp"
00065 #include "BelosStatusTestOutput.hpp"
00066 #include "BelosStatusTestOutputFactory.hpp"
00067 #include "BelosOutputManager.hpp"
00068 #include "Teuchos_BLAS.hpp"
00069 #include "Teuchos_TimeMonitor.hpp"
00070 
00090 namespace Belos {
00091   
00093 
00094   
00101 class BlockGmresSolMgrLinearProblemFailure : public BelosError {public:
00102   BlockGmresSolMgrLinearProblemFailure(const std::string& what_arg) : BelosError(what_arg)
00103     {}};
00104   
00111 class BlockGmresSolMgrOrthoFailure : public BelosError {public:
00112   BlockGmresSolMgrOrthoFailure(const std::string& what_arg) : BelosError(what_arg)
00113     {}};
00114   
00115 template<class ScalarType, class MV, class OP>
00116 class BlockGmresSolMgr : public SolverManager<ScalarType,MV,OP> {
00117     
00118 private:
00119   typedef MultiVecTraits<ScalarType,MV> MVT;
00120   typedef OperatorTraits<ScalarType,MV,OP> OPT;
00121   typedef Teuchos::ScalarTraits<ScalarType> SCT;
00122   typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType MagnitudeType;
00123   typedef Teuchos::ScalarTraits<MagnitudeType> MT;
00124     
00125 public:
00126     
00128 
00129    
00135   BlockGmresSolMgr();
00136  
00151   BlockGmresSolMgr( const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem,
00152     const Teuchos::RCP<Teuchos::ParameterList> &pl );
00153     
00155   virtual ~BlockGmresSolMgr() {};
00157     
00159 
00160     
00163   const LinearProblem<ScalarType,MV,OP>& getProblem() const {
00164     return *problem_;
00165   }
00166 
00169   Teuchos::RCP<const Teuchos::ParameterList> getValidParameters() const;
00170 
00173   Teuchos::RCP<const Teuchos::ParameterList> getCurrentParameters() const { return params_; }
00174  
00180   Teuchos::Array<Teuchos::RCP<Teuchos::Time> > getTimers() const {
00181     return tuple(timerSolve_);
00182   }
00183   
00185   int getNumIters() const {
00186     return numIters_;
00187   } 
00188  
00192   bool isLOADetected() const { return loaDetected_; }
00193  
00195     
00197 
00198   
00200   void setProblem( const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem ) { problem_ = problem; isSTSet_ = false; }
00201   
00203   void setParameters( const Teuchos::RCP<Teuchos::ParameterList> &params );
00204     
00206 
00208 
00209 
00213   void reset( const ResetType type ) { if ((type & Belos::Problem) && !Teuchos::is_null(problem_)) problem_->setProblem(); }
00215   
00217 
00218     
00236   ReturnType solve();
00237     
00239     
00242     
00244   std::string description() const;
00245     
00247     
00248 private:
00249 
00250   // Method to convert std::string to enumerated type for residual.
00251   Belos::ScaleType convertStringToScaleType( std::string& scaleType ) {
00252     if (scaleType == "Norm of Initial Residual") {
00253       return Belos::NormOfInitRes;
00254     } else if (scaleType == "Norm of Preconditioned Initial Residual") {
00255       return Belos::NormOfPrecInitRes;
00256     } else if (scaleType == "Norm of RHS") {
00257       return Belos::NormOfRHS;
00258     } else if (scaleType == "None") {
00259       return Belos::None;
00260     } else 
00261       TEST_FOR_EXCEPTION( true ,std::logic_error,
00262         "Belos::BlockGmresSolMgr(): Invalid residual scaling type.");
00263   }
00264   
00265   // Method for checking current status test against defined linear problem.
00266   bool checkStatusTest();
00267 
00268   // Linear problem.
00269   Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > problem_;
00270     
00271   // Output manager.
00272   Teuchos::RCP<OutputManager<ScalarType> > printer_;
00273   Teuchos::RCP<std::ostream> outputStream_;
00274 
00275   // Status test.
00276   Teuchos::RCP<StatusTest<ScalarType,MV,OP> > sTest_;
00277   Teuchos::RCP<StatusTestMaxIters<ScalarType,MV,OP> > maxIterTest_;
00278   Teuchos::RCP<StatusTest<ScalarType,MV,OP> > convTest_;
00279   Teuchos::RCP<StatusTestResNorm<ScalarType,MV,OP> > expConvTest_, impConvTest_;
00280   Teuchos::RCP<StatusTestOutput<ScalarType,MV,OP> > outputTest_;
00281 
00282   // Orthogonalization manager.
00283   Teuchos::RCP<MatOrthoManager<ScalarType,MV,OP> > ortho_; 
00284     
00285   // Current parameter list.
00286   Teuchos::RCP<ParameterList> params_;
00287 
00288   // Default solver values.
00289   static const MagnitudeType convtol_default_;
00290   static const MagnitudeType orthoKappa_default_;
00291   static const int maxRestarts_default_;
00292   static const int maxIters_default_;
00293   static const bool adaptiveBlockSize_default_;
00294   static const bool showMaxResNormOnly_default_;
00295   static const bool flexibleGmres_default_;
00296   static const bool expResTest_default_;
00297   static const int blockSize_default_;
00298   static const int numBlocks_default_;
00299   static const int verbosity_default_;
00300   static const int outputStyle_default_;
00301   static const int outputFreq_default_;
00302   static const std::string impResScale_default_; 
00303   static const std::string expResScale_default_; 
00304   static const std::string label_default_;
00305   static const std::string orthoType_default_;
00306   static const Teuchos::RCP<std::ostream> outputStream_default_;
00307 
00308   // Current solver values.
00309   MagnitudeType convtol_, orthoKappa_;
00310   int maxRestarts_, maxIters_, numIters_;
00311   int blockSize_, numBlocks_, verbosity_, outputStyle_, outputFreq_;
00312   bool adaptiveBlockSize_, showMaxResNormOnly_, isFlexible_, expResTest_;
00313   std::string orthoType_; 
00314   std::string impResScale_, expResScale_;
00315     
00316   // Timers.
00317   std::string label_;
00318   Teuchos::RCP<Teuchos::Time> timerSolve_;
00319 
00320   // Internal state variables.
00321   bool isSet_, isSTSet_;
00322   bool loaDetected_;
00323 };
00324 
00325 
00326 // Default solver values.
00327 template<class ScalarType, class MV, class OP>
00328 const typename BlockGmresSolMgr<ScalarType,MV,OP>::MagnitudeType BlockGmresSolMgr<ScalarType,MV,OP>::convtol_default_ = 1e-8;
00329 
00330 template<class ScalarType, class MV, class OP>
00331 const typename BlockGmresSolMgr<ScalarType,MV,OP>::MagnitudeType BlockGmresSolMgr<ScalarType,MV,OP>::orthoKappa_default_ = -1.0;
00332 
00333 template<class ScalarType, class MV, class OP>
00334 const int BlockGmresSolMgr<ScalarType,MV,OP>::maxRestarts_default_ = 20;
00335 
00336 template<class ScalarType, class MV, class OP>
00337 const int BlockGmresSolMgr<ScalarType,MV,OP>::maxIters_default_ = 1000;
00338 
00339 template<class ScalarType, class MV, class OP>
00340 const bool BlockGmresSolMgr<ScalarType,MV,OP>::adaptiveBlockSize_default_ = true;
00341 
00342 template<class ScalarType, class MV, class OP>
00343 const bool BlockGmresSolMgr<ScalarType,MV,OP>::showMaxResNormOnly_default_ = false;
00344 
00345 template<class ScalarType, class MV, class OP>
00346 const bool BlockGmresSolMgr<ScalarType,MV,OP>::flexibleGmres_default_ = false;
00347 
00348 template<class ScalarType, class MV, class OP>
00349 const bool BlockGmresSolMgr<ScalarType,MV,OP>::expResTest_default_ = false;
00350 
00351 template<class ScalarType, class MV, class OP>
00352 const int BlockGmresSolMgr<ScalarType,MV,OP>::blockSize_default_ = 1;
00353 
00354 template<class ScalarType, class MV, class OP>
00355 const int BlockGmresSolMgr<ScalarType,MV,OP>::numBlocks_default_ = 300;
00356 
00357 template<class ScalarType, class MV, class OP>
00358 const int BlockGmresSolMgr<ScalarType,MV,OP>::verbosity_default_ = Belos::Errors;
00359 
00360 template<class ScalarType, class MV, class OP>
00361 const int BlockGmresSolMgr<ScalarType,MV,OP>::outputStyle_default_ = Belos::General;
00362 
00363 template<class ScalarType, class MV, class OP>
00364 const int BlockGmresSolMgr<ScalarType,MV,OP>::outputFreq_default_ = -1;
00365 
00366 template<class ScalarType, class MV, class OP>
00367 const std::string BlockGmresSolMgr<ScalarType,MV,OP>::impResScale_default_ = "Norm of Preconditioned Initial Residual";
00368 
00369 template<class ScalarType, class MV, class OP>
00370 const std::string BlockGmresSolMgr<ScalarType,MV,OP>::expResScale_default_ = "Norm of Initial Residual";
00371 
00372 template<class ScalarType, class MV, class OP>
00373 const std::string BlockGmresSolMgr<ScalarType,MV,OP>::label_default_ = "Belos";
00374 
00375 template<class ScalarType, class MV, class OP>
00376 const std::string BlockGmresSolMgr<ScalarType,MV,OP>::orthoType_default_ = "DGKS";
00377 
00378 template<class ScalarType, class MV, class OP>
00379 const Teuchos::RCP<std::ostream> BlockGmresSolMgr<ScalarType,MV,OP>::outputStream_default_ = Teuchos::rcp(&std::cout,false);
00380 
00381 
00382 // Empty Constructor
00383 template<class ScalarType, class MV, class OP>
00384 BlockGmresSolMgr<ScalarType,MV,OP>::BlockGmresSolMgr() :
00385   outputStream_(outputStream_default_),
00386   convtol_(convtol_default_),
00387   orthoKappa_(orthoKappa_default_),
00388   maxRestarts_(maxRestarts_default_),
00389   maxIters_(maxIters_default_),
00390   blockSize_(blockSize_default_),
00391   numBlocks_(numBlocks_default_),
00392   verbosity_(verbosity_default_),
00393   outputStyle_(outputStyle_default_),
00394   outputFreq_(outputFreq_default_),
00395   adaptiveBlockSize_(adaptiveBlockSize_default_),
00396   showMaxResNormOnly_(showMaxResNormOnly_default_),
00397   isFlexible_(flexibleGmres_default_),
00398   expResTest_(expResTest_default_),
00399   orthoType_(orthoType_default_),
00400   impResScale_(impResScale_default_),
00401   expResScale_(expResScale_default_),
00402   label_(label_default_),
00403   isSet_(false),
00404   isSTSet_(false),
00405   loaDetected_(false)
00406 {}
00407 
00408 
00409 // Basic Constructor
00410 template<class ScalarType, class MV, class OP>
00411 BlockGmresSolMgr<ScalarType,MV,OP>::BlockGmresSolMgr( 
00412   const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem,
00413   const Teuchos::RCP<Teuchos::ParameterList> &pl ) : 
00414   problem_(problem),
00415   outputStream_(outputStream_default_),
00416   convtol_(convtol_default_),
00417   orthoKappa_(orthoKappa_default_),
00418   maxRestarts_(maxRestarts_default_),
00419   maxIters_(maxIters_default_),
00420   blockSize_(blockSize_default_),
00421   numBlocks_(numBlocks_default_),
00422   verbosity_(verbosity_default_),
00423   outputStyle_(outputStyle_default_),
00424   outputFreq_(outputFreq_default_),
00425   adaptiveBlockSize_(adaptiveBlockSize_default_),
00426   showMaxResNormOnly_(showMaxResNormOnly_default_),
00427   isFlexible_(flexibleGmres_default_),
00428   expResTest_(expResTest_default_),
00429   orthoType_(orthoType_default_),
00430   impResScale_(impResScale_default_),
00431   expResScale_(expResScale_default_),
00432   label_(label_default_),
00433   isSet_(false),
00434   isSTSet_(false),
00435   loaDetected_(false)
00436 {
00437 
00438   TEST_FOR_EXCEPTION(problem_ == Teuchos::null, std::invalid_argument, "Problem not given to solver manager.");
00439 
00440   // If the parameter list pointer is null, then set the current parameters to the default parameter list.
00441   if ( !is_null(pl) ) {
00442     setParameters( pl );  
00443   }
00444 
00445 }
00446 
00447 
00448 template<class ScalarType, class MV, class OP>
00449 Teuchos::RCP<const Teuchos::ParameterList>
00450 BlockGmresSolMgr<ScalarType,MV,OP>::getValidParameters() const
00451 {
00452   static Teuchos::RCP<const Teuchos::ParameterList> validPL;
00453   if (is_null(validPL)) {
00454     Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
00455     pl->set("Convergence Tolerance", convtol_default_,
00456       "The relative residual tolerance that needs to be achieved by the\n"
00457       "iterative solver in order for the linear system to be declared converged." );
00458     pl->set("Maximum Restarts", maxRestarts_default_,
00459       "The maximum number of restarts allowed for each\n"
00460       "set of RHS solved.");
00461     pl->set(
00462       "Maximum Iterations", maxIters_default_,
00463       "The maximum number of block iterations allowed for each\n"
00464       "set of RHS solved.");
00465     pl->set("Num Blocks", numBlocks_default_,
00466       "The maximum number of blocks allowed in the Krylov subspace\n"
00467       "for each set of RHS solved.");
00468     pl->set("Block Size", blockSize_default_,
00469       "The number of vectors in each block.  This number times the\n"
00470       "number of blocks is the total Krylov subspace dimension.");
00471     pl->set("Adaptive Block Size", adaptiveBlockSize_default_,
00472       "Whether the solver manager should adapt the block size\n"
00473       "based on the number of RHS to solve.");
00474     pl->set("Verbosity", verbosity_default_,
00475       "What type(s) of solver information should be outputted\n"
00476       "to the output stream.");
00477     pl->set("Output Style", outputStyle_default_,
00478       "What style is used for the solver information outputted\n"
00479       "to the output stream.");
00480     pl->set("Output Frequency", outputFreq_default_,
00481       "How often convergence information should be outputted\n"
00482       "to the output stream.");  
00483     pl->set("Output Stream", outputStream_default_,
00484       "A reference-counted pointer to the output stream where all\n"
00485       "solver output is sent.");
00486     pl->set("Show Maximum Residual Norm Only", showMaxResNormOnly_default_,
00487       "When convergence information is printed, only show the maximum\n"
00488       "relative residual norm when the block size is greater than one.");
00489     pl->set("Flexible Gmres", flexibleGmres_default_,
00490       "Whether the solver manager should use the flexible variant\n"
00491       "of GMRES.");
00492     pl->set("Explicit Residual Test", expResTest_default_,
00493       "Whether the explicitly computed residual should be used in the convergence test.");
00494     pl->set("Implicit Residual Scaling", impResScale_default_,
00495       "The type of scaling used in the implicit residual convergence test.");
00496     pl->set("Explicit Residual Scaling", expResScale_default_,
00497       "The type of scaling used in the explicit residual convergence test.");
00498     pl->set("Timer Label", label_default_,
00499       "The string to use as a prefix for the timer labels.");
00500     //  pl->set("Restart Timers", restartTimers_);
00501     pl->set("Orthogonalization", orthoType_default_,
00502       "The type of orthogonalization to use: DGKS, ICGS, or IMGS.");
00503     pl->set("Orthogonalization Constant",orthoKappa_default_,
00504       "The constant used by DGKS orthogonalization to determine\n"
00505       "whether another step of classical Gram-Schmidt is necessary.");
00506     validPL = pl;
00507   }
00508   return validPL;
00509 }
00510 
00511 
00512 template<class ScalarType, class MV, class OP>
00513 void BlockGmresSolMgr<ScalarType,MV,OP>::setParameters( const Teuchos::RCP<Teuchos::ParameterList> &params )
00514 {
00515 
00516   // Create the internal parameter list if ones doesn't already exist.
00517   if (params_ == Teuchos::null) {
00518     params_ = Teuchos::rcp( new Teuchos::ParameterList(*getValidParameters()) );
00519   }
00520   else {
00521     params->validateParameters(*getValidParameters());
00522   }
00523 
00524   // Check for maximum number of restarts
00525   if (params->isParameter("Maximum Restarts")) {
00526     maxRestarts_ = params->get("Maximum Restarts",maxRestarts_default_);
00527 
00528     // Update parameter in our list.
00529     params_->set("Maximum Restarts", maxRestarts_);
00530   }
00531 
00532   // Check for maximum number of iterations
00533   if (params->isParameter("Maximum Iterations")) {
00534     maxIters_ = params->get("Maximum Iterations",maxIters_default_);
00535 
00536     // Update parameter in our list and in status test.
00537     params_->set("Maximum Iterations", maxIters_);
00538     if (maxIterTest_!=Teuchos::null)
00539       maxIterTest_->setMaxIters( maxIters_ );
00540   }
00541 
00542   // Check for blocksize
00543   if (params->isParameter("Block Size")) {
00544     blockSize_ = params->get("Block Size",blockSize_default_);    
00545     TEST_FOR_EXCEPTION(blockSize_ <= 0, std::invalid_argument,
00546       "Belos::BlockGmresSolMgr: \"Block Size\" must be strictly positive.");
00547 
00548     // Update parameter in our list.
00549     params_->set("Block Size", blockSize_);
00550   }
00551 
00552   // Check if the blocksize should be adaptive
00553   if (params->isParameter("Adaptive Block Size")) {
00554     adaptiveBlockSize_ = params->get("Adaptive Block Size",adaptiveBlockSize_default_);
00555     
00556     // Update parameter in our list.
00557     params_->set("Adaptive Block Size", adaptiveBlockSize_);
00558   }
00559 
00560   // Check for the maximum number of blocks.
00561   if (params->isParameter("Num Blocks")) {
00562     numBlocks_ = params->get("Num Blocks",numBlocks_default_);
00563     TEST_FOR_EXCEPTION(numBlocks_ <= 0, std::invalid_argument,
00564       "Belos::BlockGmresSolMgr: \"Num Blocks\" must be strictly positive.");
00565 
00566     // Update parameter in our list.
00567     params_->set("Num Blocks", numBlocks_);
00568   }
00569 
00570   // Check to see if the timer label changed.
00571   if (params->isParameter("Timer Label")) {
00572     std::string tempLabel = params->get("Timer Label", label_default_);
00573 
00574     // Update parameter in our list and solver timer
00575     if (tempLabel != label_) {
00576       label_ = tempLabel;
00577       params_->set("Timer Label", label_);
00578       std::string solveLabel = label_ + ": BlockGmresSolMgr total solve time";
00579       timerSolve_ = Teuchos::TimeMonitor::getNewTimer(solveLabel);
00580     }
00581   }
00582 
00583   // Determine whether this solver should be "flexible".
00584   if (params->isParameter("Flexible Gmres")) {
00585     isFlexible_ = Teuchos::getParameter<bool>(*params,"Flexible Gmres");
00586     params_->set("Flexible Gmres", isFlexible_);
00587     if (isFlexible_ && expResTest_) {
00588       // Use an implicit convergence test if the Gmres solver is flexible
00589       isSTSet_ = false;
00590     }
00591   }
00592 
00593 
00594   // Check if the orthogonalization changed.
00595   if (params->isParameter("Orthogonalization")) {
00596     std::string tempOrthoType = params->get("Orthogonalization",orthoType_default_);
00597     TEST_FOR_EXCEPTION( tempOrthoType != "DGKS" && tempOrthoType != "ICGS" && tempOrthoType != "IMGS", 
00598                         std::invalid_argument,
00599       "Belos::BlockGmresSolMgr: \"Orthogonalization\" must be either \"DGKS\", \"ICGS\", or \"IMGS\".");
00600     if (tempOrthoType != orthoType_) {
00601       orthoType_ = tempOrthoType;
00602       // Create orthogonalization manager
00603       if (orthoType_=="DGKS") {
00604         if (orthoKappa_ <= 0) {
00605           ortho_ = Teuchos::rcp( new DGKSOrthoManager<ScalarType,MV,OP>( label_ ) );
00606         }
00607         else {
00608           ortho_ = Teuchos::rcp( new DGKSOrthoManager<ScalarType,MV,OP>( label_ ) );
00609           Teuchos::rcp_dynamic_cast<DGKSOrthoManager<ScalarType,MV,OP> >(ortho_)->setDepTol( orthoKappa_ );
00610         }
00611       }
00612       else if (orthoType_=="ICGS") {
00613         ortho_ = Teuchos::rcp( new ICGSOrthoManager<ScalarType,MV,OP>( label_ ) );
00614       } 
00615       else if (orthoType_=="IMGS") {
00616         ortho_ = Teuchos::rcp( new IMGSOrthoManager<ScalarType,MV,OP>( label_ ) );
00617       } 
00618     }  
00619   }
00620 
00621   // Check which orthogonalization constant to use.
00622   if (params->isParameter("Orthogonalization Constant")) {
00623     orthoKappa_ = params->get("Orthogonalization Constant",orthoKappa_default_);
00624 
00625     // Update parameter in our list.
00626     params_->set("Orthogonalization Constant",orthoKappa_);
00627     if (orthoType_=="DGKS") {
00628       if (orthoKappa_ > 0 && ortho_ != Teuchos::null) {
00629         Teuchos::rcp_dynamic_cast<DGKSOrthoManager<ScalarType,MV,OP> >(ortho_)->setDepTol( orthoKappa_ );
00630       }
00631     } 
00632   }
00633 
00634   // Check for a change in verbosity level
00635   if (params->isParameter("Verbosity")) {
00636     if (Teuchos::isParameterType<int>(*params,"Verbosity")) {
00637       verbosity_ = params->get("Verbosity", verbosity_default_);
00638     } else {
00639       verbosity_ = (int)Teuchos::getParameter<Belos::MsgType>(*params,"Verbosity");
00640     }
00641 
00642     // Update parameter in our list.
00643     params_->set("Verbosity", verbosity_);
00644     if (printer_ != Teuchos::null)
00645       printer_->setVerbosity(verbosity_);
00646   }
00647 
00648   // Check for a change in output style
00649   if (params->isParameter("Output Style")) {
00650     if (Teuchos::isParameterType<int>(*params,"Output Style")) {
00651       outputStyle_ = params->get("Output Style", outputStyle_default_);
00652     } else {
00653       outputStyle_ = (int)Teuchos::getParameter<Belos::OutputType>(*params,"Output Style");
00654     }
00655 
00656     // Reconstruct the convergence test if the explicit residual test is not being used.
00657     params_->set("Output Style", outputStyle_);
00658     if (outputTest_ != Teuchos::null) {
00659       isSTSet_ = false;
00660     }
00661   }
00662 
00663   // output stream
00664   if (params->isParameter("Output Stream")) {
00665     outputStream_ = Teuchos::getParameter<Teuchos::RCP<std::ostream> >(*params,"Output Stream");
00666 
00667     // Update parameter in our list.
00668     params_->set("Output Stream", outputStream_);
00669     if (printer_ != Teuchos::null)
00670       printer_->setOStream( outputStream_ );
00671   }
00672 
00673   // frequency level
00674   if (verbosity_ & Belos::StatusTestDetails) {
00675     if (params->isParameter("Output Frequency")) {
00676       outputFreq_ = params->get("Output Frequency", outputFreq_default_);
00677     }
00678 
00679     // Update parameter in out list and output status test.
00680     params_->set("Output Frequency", outputFreq_);
00681     if (outputTest_ != Teuchos::null)
00682       outputTest_->setOutputFrequency( outputFreq_ );
00683   }
00684 
00685   // Create output manager if we need to.
00686   if (printer_ == Teuchos::null) {
00687     printer_ = Teuchos::rcp( new OutputManager<ScalarType>(verbosity_, outputStream_) );
00688   }  
00689   
00690   // Check for convergence tolerance
00691   if (params->isParameter("Convergence Tolerance")) {
00692     convtol_ = params->get("Convergence Tolerance",convtol_default_);
00693 
00694     // Update parameter in our list and residual tests.
00695     params_->set("Convergence Tolerance", convtol_);
00696     if (impConvTest_ != Teuchos::null)
00697       impConvTest_->setTolerance( convtol_ );
00698     if (expConvTest_ != Teuchos::null)
00699       expConvTest_->setTolerance( convtol_ );
00700   }
00701  
00702   // Check for a change in scaling, if so we need to build new residual tests.
00703   if (params->isParameter("Implicit Residual Scaling")) {
00704     std::string tempImpResScale = Teuchos::getParameter<std::string>( *params, "Implicit Residual Scaling" );
00705 
00706     // Only update the scaling if it's different.
00707     if (impResScale_ != tempImpResScale) {
00708       Belos::ScaleType impResScaleType = convertStringToScaleType( tempImpResScale );
00709       impResScale_ = tempImpResScale;
00710 
00711       // Update parameter in our list and residual tests
00712       params_->set("Implicit Residual Scaling", impResScale_);
00713       if (impConvTest_ != Teuchos::null) {
00714         try { 
00715           impConvTest_->defineScaleForm( impResScaleType, Belos::TwoNorm );
00716         }
00717         catch (std::exception& e) { 
00718           // Make sure the convergence test gets constructed again.
00719           isSTSet_ = false;
00720         }
00721       }
00722     }      
00723   }
00724   
00725   if (params->isParameter("Explicit Residual Scaling")) {
00726     std::string tempExpResScale = Teuchos::getParameter<std::string>( *params, "Explicit Residual Scaling" );
00727 
00728     // Only update the scaling if it's different.
00729     if (expResScale_ != tempExpResScale) {
00730       Belos::ScaleType expResScaleType = convertStringToScaleType( tempExpResScale );
00731       expResScale_ = tempExpResScale;
00732 
00733       // Update parameter in our list and residual tests
00734       params_->set("Explicit Residual Scaling", expResScale_);
00735       if (expConvTest_ != Teuchos::null) {
00736         try { 
00737           expConvTest_->defineScaleForm( expResScaleType, Belos::TwoNorm );
00738         }
00739         catch (std::exception& e) {
00740           // Make sure the convergence test gets constructed again.
00741           isSTSet_ = false;
00742         }
00743       }
00744     }      
00745   }
00746 
00747   if (params->isParameter("Explicit Residual Test")) {
00748     expResTest_ = Teuchos::getParameter<bool>( *params,"Explicit Residual Test" );
00749 
00750     // Reconstruct the convergence test if the explicit residual test is not being used.
00751     params_->set("Explicit Residual Test", expResTest_);
00752     if (expConvTest_ == Teuchos::null) {
00753       isSTSet_ = false;
00754     }
00755   }
00756 
00757 
00758   if (params->isParameter("Show Maximum Residual Norm Only")) {
00759     showMaxResNormOnly_ = Teuchos::getParameter<bool>(*params,"Show Maximum Residual Norm Only");
00760 
00761     // Update parameter in our list and residual tests
00762     params_->set("Show Maximum Residual Norm Only", showMaxResNormOnly_);
00763     if (impConvTest_ != Teuchos::null)
00764       impConvTest_->setShowMaxResNormOnly( showMaxResNormOnly_ );
00765     if (expConvTest_ != Teuchos::null)
00766       expConvTest_->setShowMaxResNormOnly( showMaxResNormOnly_ );
00767   }
00768 
00769   // Create orthogonalization manager if we need to.
00770   if (ortho_ == Teuchos::null) {
00771     if (orthoType_=="DGKS") {
00772       if (orthoKappa_ <= 0) {
00773         ortho_ = Teuchos::rcp( new DGKSOrthoManager<ScalarType,MV,OP>( label_ ) );
00774       }
00775       else {
00776         ortho_ = Teuchos::rcp( new DGKSOrthoManager<ScalarType,MV,OP>( label_ ) );
00777         Teuchos::rcp_dynamic_cast<DGKSOrthoManager<ScalarType,MV,OP> >(ortho_)->setDepTol( orthoKappa_ );
00778       }
00779     }
00780     else if (orthoType_=="ICGS") {
00781       ortho_ = Teuchos::rcp( new ICGSOrthoManager<ScalarType,MV,OP>( label_ ) );
00782     } 
00783     else if (orthoType_=="IMGS") {
00784       ortho_ = Teuchos::rcp( new IMGSOrthoManager<ScalarType,MV,OP>( label_ ) );
00785     } 
00786     else {
00787       TEST_FOR_EXCEPTION(orthoType_!="ICGS"&&orthoType_!="DGKS"&&orthoType_!="IMGS",std::logic_error,
00788         "Belos::BlockGmresSolMgr(): Invalid orthogonalization type.");
00789     }  
00790   }
00791 
00792   // Create the timer if we need to.
00793   if (timerSolve_ == Teuchos::null) {
00794     std::string solveLabel = label_ + ": BlockGmresSolMgr total solve time";
00795     timerSolve_ = Teuchos::TimeMonitor::getNewTimer(solveLabel);
00796   }
00797 
00798   // Inform the solver manager that the current parameters were set.
00799   isSet_ = true;
00800 }
00801 
00802 // Check the status test versus the defined linear problem
00803 template<class ScalarType, class MV, class OP>
00804 bool BlockGmresSolMgr<ScalarType,MV,OP>::checkStatusTest() {
00805 
00806   typedef Belos::StatusTestCombo<ScalarType,MV,OP>  StatusTestCombo_t;
00807   typedef Belos::StatusTestGenResNorm<ScalarType,MV,OP>  StatusTestGenResNorm_t;
00808   typedef Belos::StatusTestImpResNorm<ScalarType,MV,OP>  StatusTestImpResNorm_t;
00809 
00810   // Basic test checks maximum iterations and native residual.
00811   maxIterTest_ = Teuchos::rcp( new StatusTestMaxIters<ScalarType,MV,OP>( maxIters_ ) );
00812 
00813   // If there is a left preconditioner, we create a combined status test that checks the implicit
00814   // and then explicit residual norm to see if we have convergence.
00815   if (!Teuchos::is_null(problem_->getLeftPrec()) && !isFlexible_) {
00816     expResTest_ = true;
00817   }
00818 
00819   if (expResTest_) {
00820    
00821     // Implicit residual test, using the native residual to determine if convergence was achieved.
00822     Teuchos::RCP<StatusTestGenResNorm_t> tmpImpConvTest =
00823       Teuchos::rcp( new StatusTestGenResNorm_t( convtol_ ) );
00824     tmpImpConvTest->defineScaleForm( convertStringToScaleType(impResScale_), Belos::TwoNorm );
00825     tmpImpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
00826     impConvTest_ = tmpImpConvTest;
00827 
00828     // Explicit residual test once the native residual is below the tolerance
00829     Teuchos::RCP<StatusTestGenResNorm_t> tmpExpConvTest =
00830       Teuchos::rcp( new StatusTestGenResNorm_t( convtol_ ) );
00831     tmpExpConvTest->defineResForm( StatusTestGenResNorm_t::Explicit, Belos::TwoNorm );
00832     tmpExpConvTest->defineScaleForm( convertStringToScaleType(expResScale_), Belos::TwoNorm );
00833     tmpExpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
00834     expConvTest_ = tmpExpConvTest;
00835 
00836     // The convergence test is a combination of the "cheap" implicit test and explicit test.
00837     convTest_ = Teuchos::rcp( new StatusTestCombo_t( StatusTestCombo_t::SEQ, impConvTest_, expConvTest_ ) );
00838   }
00839   else {
00840 
00841     if (isFlexible_) {
00842       // Implicit residual test, using the native residual to determine if convergence was achieved.
00843       Teuchos::RCP<StatusTestGenResNorm_t> tmpImpConvTest =
00844         Teuchos::rcp( new StatusTestGenResNorm_t( convtol_ ) );
00845       tmpImpConvTest->defineScaleForm( convertStringToScaleType(impResScale_), Belos::TwoNorm );
00846       tmpImpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
00847       impConvTest_ = tmpImpConvTest;
00848     }
00849     else {
00850       // Implicit residual test, using the native residual to determine if convergence was achieved.
00851       // Use test that checks for loss of accuracy.
00852       Teuchos::RCP<StatusTestImpResNorm_t> tmpImpConvTest =
00853         Teuchos::rcp( new StatusTestImpResNorm_t( convtol_ ) );
00854       tmpImpConvTest->defineScaleForm( convertStringToScaleType(impResScale_), Belos::TwoNorm );
00855       tmpImpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
00856       impConvTest_ = tmpImpConvTest;
00857     }
00858 
00859     // Set the explicit and total convergence test to this implicit test that checks for accuracy loss.
00860     expConvTest_ = impConvTest_;
00861     convTest_ = impConvTest_;
00862   }
00863 
00864   // Create the status test.
00865   sTest_ = Teuchos::rcp( new StatusTestCombo_t( StatusTestCombo_t::OR, maxIterTest_, convTest_ ) );
00866 
00867   // Create the status test output class.
00868   // This class manages and formats the output from the status test.
00869   StatusTestOutputFactory<ScalarType,MV,OP> stoFactory( outputStyle_ );
00870   outputTest_ = stoFactory.create( printer_, sTest_, outputFreq_, Passed+Failed+Undefined );
00871 
00872   // Set the solver string for the output test
00873   std::string solverDesc = " Block Gmres ";
00874   if (isFlexible_)
00875     solverDesc = "Flexible" + solverDesc;
00876   outputTest_->setSolverDesc( solverDesc );
00877 
00878   // The status test is now set.
00879   isSTSet_ = true;
00880 
00881   return false;
00882 }
00883   
00884 // solve()
00885 template<class ScalarType, class MV, class OP>
00886 ReturnType BlockGmresSolMgr<ScalarType,MV,OP>::solve() {
00887 
00888   // Set the current parameters if they were not set before.
00889   // NOTE:  This may occur if the user generated the solver manager with the default constructor and 
00890   // then didn't set any parameters using setParameters().
00891   if (!isSet_) {
00892     setParameters(Teuchos::parameterList(*getValidParameters()));
00893   }
00894 
00895   Teuchos::BLAS<int,ScalarType> blas;
00896   
00897   TEST_FOR_EXCEPTION(problem_ == Teuchos::null,BlockGmresSolMgrLinearProblemFailure,
00898     "Belos::BlockGmresSolMgr::solve(): Linear problem is not a valid object.");
00899 
00900   TEST_FOR_EXCEPTION(!problem_->isProblemSet(),BlockGmresSolMgrLinearProblemFailure,
00901     "Belos::BlockGmresSolMgr::solve(): Linear problem is not ready, setProblem() has not been called.");
00902 
00903   if (isFlexible_) {
00904     TEST_FOR_EXCEPTION(problem_->getRightPrec()==Teuchos::null,BlockGmresSolMgrLinearProblemFailure,
00905       "Belos::BlockGmresSolMgr::solve(): Linear problem does not have a preconditioner required for flexible GMRES, call setRightPrec().");
00906   }
00907  
00908   if (!isSTSet_ || (!expResTest_ && !Teuchos::is_null(problem_->getLeftPrec())) ) {
00909     TEST_FOR_EXCEPTION( checkStatusTest(),BlockGmresSolMgrLinearProblemFailure,
00910       "Belos::BlockGmresSolMgr::solve(): Linear problem and requested status tests are incompatible.");
00911   }
00912 
00913   // Create indices for the linear systems to be solved.
00914   int startPtr = 0;
00915   int numRHS2Solve = MVT::GetNumberVecs( *(problem_->getRHS()) );
00916   int numCurrRHS = ( numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;
00917 
00918   std::vector<int> currIdx;
00919   //  If an adaptive block size is allowed then only the linear systems that need to be solved are solved.
00920   //  Otherwise, the index set is generated that informs the linear problem that some linear systems are augmented.
00921   if ( adaptiveBlockSize_ ) {
00922     blockSize_ = numCurrRHS;
00923     currIdx.resize( numCurrRHS  );
00924     for (int i=0; i<numCurrRHS; ++i) 
00925     { currIdx[i] = startPtr+i; }
00926     
00927   }
00928   else {
00929     currIdx.resize( blockSize_ );
00930     for (int i=0; i<numCurrRHS; ++i) 
00931     { currIdx[i] = startPtr+i; }
00932     for (int i=numCurrRHS; i<blockSize_; ++i) 
00933     { currIdx[i] = -1; }
00934   }
00935 
00936   // Inform the linear problem of the current linear system to solve.
00937   problem_->setLSIndex( currIdx );
00938 
00940   // Parameter list
00941   Teuchos::ParameterList plist;
00942   plist.set("Block Size",blockSize_);
00943 
00944   int dim = MVT::GetVecLength( *(problem_->getRHS()) );  
00945   if (blockSize_*numBlocks_ > dim) {
00946     int tmpNumBlocks = 0;
00947     if (blockSize_ == 1)
00948       tmpNumBlocks = dim / blockSize_;  // Allow for a good breakdown.
00949     else
00950       tmpNumBlocks = ( dim - blockSize_) / blockSize_;  // Allow for restarting.
00951     printer_->stream(Warnings) << 
00952       "Warning! Requested Krylov subspace dimension is larger than operator dimension!" << std::endl <<
00953       " The maximum number of blocks allowed for the Krylov subspace will be adjusted to " << tmpNumBlocks << std::endl;
00954     plist.set("Num Blocks",tmpNumBlocks);
00955   } 
00956   else 
00957     plist.set("Num Blocks",numBlocks_);
00958   
00959   // Reset the status test.  
00960   outputTest_->reset();
00961   loaDetected_ = false;
00962 
00963   // Assume convergence is achieved, then let any failed convergence set this to false.
00964   bool isConverged = true;  
00965 
00967   // BlockGmres solver
00968 
00969   Teuchos::RCP<GmresIteration<ScalarType,MV,OP> > block_gmres_iter;
00970 
00971   if (isFlexible_)
00972     block_gmres_iter = Teuchos::rcp( new BlockFGmresIter<ScalarType,MV,OP>(problem_,printer_,outputTest_,ortho_,plist) );
00973   else
00974     block_gmres_iter = Teuchos::rcp( new BlockGmresIter<ScalarType,MV,OP>(problem_,printer_,outputTest_,ortho_,plist) );
00975   
00976   // Enter solve() iterations
00977   {
00978     Teuchos::TimeMonitor slvtimer(*timerSolve_);
00979 
00980     while ( numRHS2Solve > 0 ) {
00981 
00982       // Set the current number of blocks and blocksize with the Gmres iteration.
00983       if (blockSize_*numBlocks_ > dim) {
00984         int tmpNumBlocks = 0;
00985         if (blockSize_ == 1)
00986           tmpNumBlocks = dim / blockSize_;  // Allow for a good breakdown.
00987         else
00988           tmpNumBlocks = ( dim - blockSize_) / blockSize_;  // Allow for restarting.
00989         block_gmres_iter->setSize( blockSize_, tmpNumBlocks );
00990       }
00991       else
00992         block_gmres_iter->setSize( blockSize_, numBlocks_ );
00993 
00994       // Reset the number of iterations.
00995       block_gmres_iter->resetNumIters();
00996 
00997       // Reset the number of calls that the status test output knows about.
00998       outputTest_->resetNumCalls();
00999 
01000       // Create the first block in the current Krylov basis.
01001       Teuchos::RCP<MV> V_0;
01002       if (isFlexible_) {
01003         // Load the correct residual if the system is augmented
01004         if (currIdx[blockSize_-1] == -1) {
01005           V_0 = MVT::Clone( *(problem_->getInitResVec()), blockSize_ );
01006           problem_->computeCurrResVec( &*V_0 );
01007         }
01008         else {
01009           V_0 = MVT::CloneCopy( *(problem_->getInitResVec()), currIdx );
01010         }
01011       }
01012       else { 
01013         // Load the correct residual if the system is augmented
01014         if (currIdx[blockSize_-1] == -1) {
01015           V_0 = MVT::Clone( *(problem_->getInitPrecResVec()), blockSize_ );
01016           problem_->computeCurrPrecResVec( &*V_0 );
01017         }
01018         else {
01019           V_0 = MVT::CloneCopy( *(problem_->getInitPrecResVec()), currIdx );
01020         }
01021       }
01022 
01023       // Get a matrix to hold the orthonormalization coefficients.
01024       Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > z_0 = 
01025         rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( blockSize_, blockSize_ ) );
01026  
01027       // Orthonormalize the new V_0
01028       int rank = ortho_->normalize( *V_0, z_0 );
01029       TEST_FOR_EXCEPTION(rank != blockSize_,BlockGmresSolMgrOrthoFailure,
01030         "Belos::BlockGmresSolMgr::solve(): Failed to compute initial block of orthonormal vectors.");
01031       
01032       // Set the new state and initialize the solver.
01033       GmresIterationState<ScalarType,MV> newstate;
01034       newstate.V = V_0;
01035       newstate.z = z_0;
01036       newstate.curDim = 0;
01037       block_gmres_iter->initializeGmres(newstate);
01038       int numRestarts = 0;
01039 
01040       while(1) {
01041         // tell block_gmres_iter to iterate
01042         try {
01043           block_gmres_iter->iterate();
01044     
01046           //
01047           // check convergence first
01048           //
01050           if ( convTest_->getStatus() == Passed ) {
01051             if ( expConvTest_->getLOADetected() ) {
01052               // we don't have convergence
01053               loaDetected_ = true;
01054               isConverged = false;
01055             }
01056             break;  // break from while(1){block_gmres_iter->iterate()}
01057           }
01059           //
01060           // check for maximum iterations
01061           //
01063           else if ( maxIterTest_->getStatus() == Passed ) {
01064             // we don't have convergence
01065             isConverged = false;
01066             break;  // break from while(1){block_gmres_iter->iterate()}
01067           }
01069           //
01070           // check for restarting, i.e. the subspace is full
01071           //
01073           else if ( block_gmres_iter->getCurSubspaceDim() == block_gmres_iter->getMaxSubspaceDim() ) {
01074       
01075             if ( numRestarts >= maxRestarts_ ) {
01076               isConverged = false;
01077               break; // break from while(1){block_gmres_iter->iterate()}
01078             }
01079             numRestarts++;
01080       
01081             printer_->stream(Debug) << " Performing restart number " << numRestarts << " of " << maxRestarts_ << std::endl << std::endl;
01082       
01083             // Update the linear problem.
01084             Teuchos::RCP<MV> update = block_gmres_iter->getCurrentUpdate();
01085             if (isFlexible_) {
01086               // Update the solution manually, since the preconditioning doesn't need to be undone.
01087               Teuchos::RCP<MV> curX = problem_->getCurrLHSVec();
01088               MVT::MvAddMv( 1.0, *curX, 1.0, *update, *curX );
01089             }
01090             else 
01091               problem_->updateSolution( update, true );
01092 
01093             // Get the state.
01094             GmresIterationState<ScalarType,MV> oldState = block_gmres_iter->getState();
01095       
01096             // Compute the restart std::vector.
01097             // Get a view of the current Krylov basis.
01098             V_0  = MVT::Clone( *(oldState.V), blockSize_ );
01099             if (isFlexible_)
01100               problem_->computeCurrResVec( &*V_0 );
01101             else
01102               problem_->computeCurrPrecResVec( &*V_0 );
01103 
01104             // Get a view of the first block of the Krylov basis.
01105             z_0 = rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( blockSize_, blockSize_ ) );
01106       
01107             // Orthonormalize the new V_0
01108             rank = ortho_->normalize( *V_0, z_0 );
01109             TEST_FOR_EXCEPTION(rank != blockSize_,BlockGmresSolMgrOrthoFailure,
01110               "Belos::BlockGmresSolMgr::solve(): Failed to compute initial block of orthonormal vectors after restart.");
01111 
01112             // Set the new state and initialize the solver.
01113             newstate.V = V_0;
01114             newstate.z = z_0;
01115             newstate.curDim = 0;
01116             block_gmres_iter->initializeGmres(newstate);
01117 
01118           } // end of restarting
01119 
01121           //
01122           // we returned from iterate(), but none of our status tests Passed.
01123           // something is wrong, and it is probably our fault.
01124           //
01126 
01127           else {
01128             TEST_FOR_EXCEPTION(true,std::logic_error,
01129               "Belos::BlockGmresSolMgr::solve(): Invalid return from BlockGmresIter::iterate().");
01130           }
01131         }
01132         catch (const GmresIterationOrthoFailure &e) {
01133           // If the block size is not one, it's not considered a lucky breakdown.
01134           if (blockSize_ != 1) {
01135             printer_->stream(Errors) << "Error! Caught std::exception in BlockGmresIter::iterate() at iteration " 
01136                                      << block_gmres_iter->getNumIters() << std::endl 
01137                                      << e.what() << std::endl;
01138             if (convTest_->getStatus() != Passed)
01139               isConverged = false;
01140             break;
01141           } 
01142           else {
01143             // If the block size is one, try to recover the most recent least-squares solution
01144             block_gmres_iter->updateLSQR( block_gmres_iter->getCurSubspaceDim() );
01145 
01146             // Check to see if the most recent least-squares solution yielded convergence.
01147             sTest_->checkStatus( &*block_gmres_iter );
01148             if (convTest_->getStatus() != Passed)
01149               isConverged = false;
01150             break;
01151           }
01152         }
01153         catch (const std::exception &e) {
01154           printer_->stream(Errors) << "Error! Caught std::exception in BlockGmresIter::iterate() at iteration " 
01155                                    << block_gmres_iter->getNumIters() << std::endl 
01156                                    << e.what() << std::endl;
01157           throw;
01158         }
01159       }
01160       
01161       // Compute the current solution.
01162       // Update the linear problem.
01163       if (isFlexible_) {
01164         // Update the solution manually, since the preconditioning doesn't need to be undone.
01165         Teuchos::RCP<MV> update = block_gmres_iter->getCurrentUpdate();
01166         Teuchos::RCP<MV> curX = problem_->getCurrLHSVec();
01167         MVT::MvAddMv( 1.0, *curX, 1.0, *update, *curX );
01168       }
01169       else {
01170         // Attempt to get the current solution from the residual status test, if it has one.
01171         if ( !Teuchos::is_null(expConvTest_->getSolution()) ) {
01172           Teuchos::RCP<MV> newX = expConvTest_->getSolution();
01173           Teuchos::RCP<MV> curX = problem_->getCurrLHSVec();
01174           MVT::MvAddMv( 0.0, *newX, 1.0, *newX, *curX );
01175         }
01176         else {
01177           Teuchos::RCP<MV> update = block_gmres_iter->getCurrentUpdate();
01178           problem_->updateSolution( update, true );
01179         }  
01180       }
01181 
01182       // Inform the linear problem that we are finished with this block linear system.
01183       problem_->setCurrLS();
01184       
01185       // Update indices for the linear systems to be solved.
01186       startPtr += numCurrRHS;
01187       numRHS2Solve -= numCurrRHS;
01188       if ( numRHS2Solve > 0 ) {
01189         numCurrRHS = ( numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;
01190 
01191         if ( adaptiveBlockSize_ ) {
01192           blockSize_ = numCurrRHS;
01193           currIdx.resize( numCurrRHS  );
01194           for (int i=0; i<numCurrRHS; ++i) 
01195           { currIdx[i] = startPtr+i; }    
01196         }
01197         else {
01198           currIdx.resize( blockSize_ );
01199           for (int i=0; i<numCurrRHS; ++i) 
01200           { currIdx[i] = startPtr+i; }
01201           for (int i=numCurrRHS; i<blockSize_; ++i) 
01202           { currIdx[i] = -1; }
01203         }
01204         // Set the next indices.
01205         problem_->setLSIndex( currIdx );
01206       }
01207       else {
01208         currIdx.resize( numRHS2Solve );
01209       }
01210       
01211     }// while ( numRHS2Solve > 0 )
01212     
01213   }
01214 
01215   // print final summary
01216   sTest_->print( printer_->stream(FinalSummary) );
01217  
01218   // print timing information
01219   Teuchos::TimeMonitor::summarize( printer_->stream(TimingDetails) );
01220 
01221   // get iteration information for this solve
01222   numIters_ = maxIterTest_->getNumIters();
01223   
01224   if (!isConverged || loaDetected_) {
01225     return Unconverged; // return from BlockGmresSolMgr::solve() 
01226   }
01227   return Converged; // return from BlockGmresSolMgr::solve() 
01228 }
01229 
01230 //  This method requires the solver manager to return a std::string that describes itself.
01231 template<class ScalarType, class MV, class OP>
01232 std::string BlockGmresSolMgr<ScalarType,MV,OP>::description() const
01233 {
01234   std::ostringstream oss;
01235   oss << "Belos::BlockGmresSolMgr<...,"<<Teuchos::ScalarTraits<ScalarType>::name()<<">";
01236   oss << "{";
01237   if (isFlexible_) {
01238     oss << "Variant=\'Flexible\', ";
01239   }
01240   oss << "Ortho Type='"<<orthoType_<<"\', Block Size=" <<blockSize_;
01241   oss << ", Num Blocks=" <<numBlocks_<< ", Max Restarts=" << maxRestarts_;
01242   oss << "}";
01243   return oss.str();
01244 }
01245   
01246 } // end Belos namespace
01247 
01248 #endif /* BELOS_BLOCK_GMRES_SOLMGR_HPP */
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines