BelosPseudoBlockGmresSolMgr.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 terms of Contract DE-AC04-94AL85000, there is a non-exclusive
00008 // license for use of this work by or on behalf of the U.S. Government.
00009 //
00010 // This library is free software; you can redistribute it and/or modify
00011 // it under the terms of the GNU Lesser General Public License as
00012 // published by the Free Software Foundation; either version 2.1 of the
00013 // License, or (at your option) any later version.
00014 //
00015 // This library is distributed in the hope that it will be useful, but
00016 // WITHOUT ANY WARRANTY; without even the implied warranty of
00017 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00018 // Lesser General Public License for more details.
00019 //
00020 // You should have received a copy of the GNU Lesser General Public
00021 // License along with this library; if not, write to the Free Software
00022 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307
00023 // USA
00024 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
00025 //
00026 // ***********************************************************************
00027 // @HEADER
00028 
00029 #ifndef BELOS_PSEUDO_BLOCK_GMRES_SOLMGR_HPP
00030 #define BELOS_PSEUDO_BLOCK_GMRES_SOLMGR_HPP
00031 
00036 #include "BelosConfigDefs.hpp"
00037 #include "BelosTypes.hpp"
00038 
00039 #include "BelosLinearProblem.hpp"
00040 #include "BelosSolverManager.hpp"
00041 
00042 #include "BelosPseudoBlockGmresIter.hpp"
00043 #include "BelosDGKSOrthoManager.hpp"
00044 #include "BelosICGSOrthoManager.hpp"
00045 #include "BelosIMGSOrthoManager.hpp"
00046 #include "BelosStatusTestMaxIters.hpp"
00047 #include "BelosStatusTestGenResNorm.hpp"
00048 #include "BelosStatusTestImpResNorm.hpp"
00049 #include "BelosStatusTestCombo.hpp"
00050 #include "BelosStatusTestOutput.hpp"
00051 #include "BelosOutputManager.hpp"
00052 #include "Teuchos_BLAS.hpp"
00053 #include "Teuchos_TimeMonitor.hpp"
00054 
00070 namespace Belos {
00071   
00073 
00074   
00081   class PseudoBlockGmresSolMgrLinearProblemFailure : public BelosError {public:
00082     PseudoBlockGmresSolMgrLinearProblemFailure(const std::string& what_arg) : BelosError(what_arg)
00083     {}};
00084   
00091   class PseudoBlockGmresSolMgrOrthoFailure : public BelosError {public:
00092     PseudoBlockGmresSolMgrOrthoFailure(const std::string& what_arg) : BelosError(what_arg)
00093     {}};
00094   
00095   template<class ScalarType, class MV, class OP>
00096   class PseudoBlockGmresSolMgr : public SolverManager<ScalarType,MV,OP> {
00097     
00098   private:
00099     typedef MultiVecTraits<ScalarType,MV> MVT;
00100     typedef OperatorTraits<ScalarType,MV,OP> OPT;
00101     typedef Teuchos::ScalarTraits<ScalarType> SCT;
00102     typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType MagnitudeType;
00103     typedef Teuchos::ScalarTraits<MagnitudeType> MT;
00104     
00105   public:
00106     
00108 
00109 
00115     PseudoBlockGmresSolMgr();
00116     
00128     PseudoBlockGmresSolMgr( const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem,
00129                 const Teuchos::RCP<Teuchos::ParameterList> &pl );
00130     
00132     virtual ~PseudoBlockGmresSolMgr() {};
00134     
00136 
00137     
00138     const LinearProblem<ScalarType,MV,OP>& getProblem() const {
00139       return *problem_;
00140     }
00141     
00144     Teuchos::RCP<const Teuchos::ParameterList> getValidParameters() const;
00145    
00148     Teuchos::RCP<const Teuchos::ParameterList> getCurrentParameters() const { return params_; }
00149  
00155     Teuchos::Array<Teuchos::RCP<Teuchos::Time> > getTimers() const {
00156       return tuple(timerSolve_);
00157     }
00158 
00160     int getNumIters() const {
00161       return numIters_;
00162     }
00163   
00167     bool isLOADetected() const { return loaDetected_; }
00168   
00170     
00172 
00173    
00175     void setProblem( const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem ) { problem_ = problem; }
00176    
00178     void setParameters( const Teuchos::RCP<Teuchos::ParameterList> &params );
00179     
00181 
00183 
00184 
00188     void reset( const ResetType type ) { if ((type & Belos::Problem) && !Teuchos::is_null(problem_)) problem_->setProblem(); }
00190 
00192 
00193     
00211     ReturnType solve();
00212     
00214     
00217     
00219     std::string description() const;
00220     
00222     
00223   private:
00224 
00225     // Method to convert std::string to enumerated type for residual.
00226     Belos::ScaleType convertStringToScaleType( const std::string& scaleType ) {
00227       if (scaleType == "Norm of Initial Residual") {
00228         return Belos::NormOfInitRes;
00229       } else if (scaleType == "Norm of Preconditioned Initial Residual") {
00230         return Belos::NormOfPrecInitRes;
00231       } else if (scaleType == "Norm of RHS") {
00232         return Belos::NormOfRHS;
00233       } else if (scaleType == "None") {
00234         return Belos::None;
00235       } else 
00236         TEST_FOR_EXCEPTION( true ,std::logic_error,
00237           "Belos::PseudoBlockGmresSolMgr(): Invalid residual scaling type.");
00238     }
00239 
00240     // Method for checking current status test against defined linear problem.
00241     bool checkStatusTest();
00242     
00243     // Linear problem.
00244     Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > problem_;
00245     
00246     // Output manager.
00247     Teuchos::RCP<OutputManager<ScalarType> > printer_;
00248     Teuchos::RCP<std::ostream> outputStream_;
00249 
00250     // Status test.
00251     Teuchos::RCP<StatusTest<ScalarType,MV,OP> > sTest_;
00252     Teuchos::RCP<StatusTestMaxIters<ScalarType,MV,OP> > maxIterTest_;
00253     Teuchos::RCP<StatusTest<ScalarType,MV,OP> > convTest_;
00254     Teuchos::RCP<StatusTestResNorm<ScalarType,MV,OP> > impConvTest_, expConvTest_;
00255     Teuchos::RCP<StatusTestOutput<ScalarType,MV,OP> > outputTest_;
00256 
00257     // Orthogonalization manager.
00258     Teuchos::RCP<MatOrthoManager<ScalarType,MV,OP> > ortho_; 
00259 
00260      // Current parameter list.
00261     Teuchos::RCP<ParameterList> params_;
00262    
00263     // Default solver values.
00264     static const MagnitudeType convtol_default_;
00265     static const MagnitudeType orthoKappa_default_;
00266     static const int maxRestarts_default_;
00267     static const int maxIters_default_;
00268     static const bool adaptiveBlockSize_default_;
00269     static const bool showMaxResNormOnly_default_;
00270     static const int blockSize_default_;
00271     static const int numBlocks_default_;
00272     static const int verbosity_default_;
00273     static const int outputFreq_default_;
00274     static const int defQuorum_default_;
00275     static const std::string impResScale_default_; 
00276     static const std::string expResScale_default_; 
00277     static const std::string label_default_;
00278     static const std::string orthoType_default_;
00279     static const Teuchos::RCP<std::ostream> outputStream_default_;
00280 
00281     // Current solver values.
00282     MagnitudeType convtol_, orthoKappa_;
00283     int maxRestarts_, maxIters_, numIters_;
00284     int blockSize_, numBlocks_, verbosity_, outputFreq_, defQuorum_;
00285     bool adaptiveBlockSize_, showMaxResNormOnly_;
00286     std::string orthoType_; 
00287     std::string impResScale_, expResScale_;       
00288  
00289     // Timers.
00290     std::string label_;
00291     Teuchos::RCP<Teuchos::Time> timerSolve_;
00292 
00293     // Internal state variables.
00294     bool isSet_, isSTSet_, expResTest_;
00295     bool loaDetected_;
00296   };
00297   
00298   
00299 // Default solver values.
00300 template<class ScalarType, class MV, class OP>
00301 const typename Teuchos::ScalarTraits<ScalarType>::magnitudeType PseudoBlockGmresSolMgr<ScalarType,MV,OP>::convtol_default_ = 1e-8;
00302 
00303 template<class ScalarType, class MV, class OP>
00304 const typename Teuchos::ScalarTraits<ScalarType>::magnitudeType PseudoBlockGmresSolMgr<ScalarType,MV,OP>::orthoKappa_default_ = -1.0;
00305 
00306 template<class ScalarType, class MV, class OP>
00307 const int PseudoBlockGmresSolMgr<ScalarType,MV,OP>::maxRestarts_default_ = 20;
00308 
00309 template<class ScalarType, class MV, class OP>
00310 const int PseudoBlockGmresSolMgr<ScalarType,MV,OP>::maxIters_default_ = 1000;
00311 
00312 template<class ScalarType, class MV, class OP>
00313 const bool PseudoBlockGmresSolMgr<ScalarType,MV,OP>::adaptiveBlockSize_default_ = true;
00314 
00315 template<class ScalarType, class MV, class OP>
00316 const bool PseudoBlockGmresSolMgr<ScalarType,MV,OP>::showMaxResNormOnly_default_ = false;
00317 
00318 template<class ScalarType, class MV, class OP>
00319 const int PseudoBlockGmresSolMgr<ScalarType,MV,OP>::blockSize_default_ = 1;
00320 
00321 template<class ScalarType, class MV, class OP>
00322 const int PseudoBlockGmresSolMgr<ScalarType,MV,OP>::numBlocks_default_ = 300;
00323 
00324 template<class ScalarType, class MV, class OP>
00325 const int PseudoBlockGmresSolMgr<ScalarType,MV,OP>::verbosity_default_ = Belos::Errors;
00326 
00327 template<class ScalarType, class MV, class OP>
00328 const int PseudoBlockGmresSolMgr<ScalarType,MV,OP>::outputFreq_default_ = -1;
00329 
00330 template<class ScalarType, class MV, class OP>
00331 const int PseudoBlockGmresSolMgr<ScalarType,MV,OP>::defQuorum_default_ = 1;
00332 
00333 template<class ScalarType, class MV, class OP>
00334 const std::string PseudoBlockGmresSolMgr<ScalarType,MV,OP>::impResScale_default_ = "Norm of Preconditioned Initial Residual";
00335 
00336 template<class ScalarType, class MV, class OP>
00337 const std::string PseudoBlockGmresSolMgr<ScalarType,MV,OP>::expResScale_default_ = "Norm of Initial Residual";
00338 
00339 template<class ScalarType, class MV, class OP>
00340 const std::string PseudoBlockGmresSolMgr<ScalarType,MV,OP>::label_default_ = "Belos";
00341 
00342 template<class ScalarType, class MV, class OP>
00343 const std::string PseudoBlockGmresSolMgr<ScalarType,MV,OP>::orthoType_default_ = "DGKS";
00344 
00345 template<class ScalarType, class MV, class OP>
00346 const Teuchos::RCP<std::ostream> PseudoBlockGmresSolMgr<ScalarType,MV,OP>::outputStream_default_ = Teuchos::rcp(&std::cout,false);
00347 
00348 
00349 // Empty Constructor
00350 template<class ScalarType, class MV, class OP>
00351 PseudoBlockGmresSolMgr<ScalarType,MV,OP>::PseudoBlockGmresSolMgr() :
00352   outputStream_(outputStream_default_),
00353   convtol_(convtol_default_),
00354   orthoKappa_(orthoKappa_default_),
00355   maxRestarts_(maxRestarts_default_),
00356   maxIters_(maxIters_default_),
00357   blockSize_(blockSize_default_),
00358   numBlocks_(numBlocks_default_),
00359   verbosity_(verbosity_default_),
00360   outputFreq_(outputFreq_default_),
00361   defQuorum_(defQuorum_default_),
00362   adaptiveBlockSize_(adaptiveBlockSize_default_),
00363   showMaxResNormOnly_(showMaxResNormOnly_default_),
00364   orthoType_(orthoType_default_),
00365   impResScale_(impResScale_default_),
00366   expResScale_(expResScale_default_),
00367   label_(label_default_),
00368   isSet_(false),
00369   isSTSet_(false),
00370   expResTest_(false),
00371   loaDetected_(false)
00372 {}
00373 
00374 // Basic Constructor
00375 template<class ScalarType, class MV, class OP>
00376 PseudoBlockGmresSolMgr<ScalarType,MV,OP>::PseudoBlockGmresSolMgr( 
00377                  const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem,
00378                  const Teuchos::RCP<Teuchos::ParameterList> &pl ) : 
00379   problem_(problem),
00380   outputStream_(outputStream_default_),
00381   convtol_(convtol_default_),
00382   orthoKappa_(orthoKappa_default_),
00383   maxRestarts_(maxRestarts_default_),
00384   maxIters_(maxIters_default_),
00385   blockSize_(blockSize_default_),
00386   numBlocks_(numBlocks_default_),
00387   verbosity_(verbosity_default_),
00388   outputFreq_(outputFreq_default_),
00389   defQuorum_(defQuorum_default_),
00390   adaptiveBlockSize_(adaptiveBlockSize_default_),
00391   showMaxResNormOnly_(showMaxResNormOnly_default_),
00392   orthoType_(orthoType_default_),
00393   impResScale_(impResScale_default_),
00394   expResScale_(expResScale_default_),
00395   label_(label_default_),
00396   isSet_(false),
00397   isSTSet_(false),
00398   expResTest_(false),
00399   loaDetected_(false)
00400 {
00401   TEST_FOR_EXCEPTION(problem_ == Teuchos::null, std::invalid_argument, "Problem not given to solver manager.");
00402 
00403   // If the parameter list pointer is null, then set the current parameters to the default parameter list.
00404   if (!is_null(pl)) {
00405     // Set the parameters using the list that was passed in.
00406     setParameters( pl );  
00407   }
00408 }
00409 
00410 template<class ScalarType, class MV, class OP>
00411 void PseudoBlockGmresSolMgr<ScalarType,MV,OP>::setParameters( const Teuchos::RCP<Teuchos::ParameterList> &params )
00412 {
00413   // Create the internal parameter list if ones doesn't already exist.
00414   if (params_ == Teuchos::null) {
00415     params_ = Teuchos::rcp( new Teuchos::ParameterList(*getValidParameters()) );
00416   }
00417   else {
00418     params->validateParameters(*getValidParameters());
00419   }
00420 
00421   // Check for maximum number of restarts
00422   if (params->isParameter("Maximum Restarts")) {
00423     maxRestarts_ = params->get("Maximum Restarts",maxRestarts_default_);
00424 
00425     // Update parameter in our list.
00426     params_->set("Maximum Restarts", maxRestarts_);
00427   }
00428 
00429   // Check for maximum number of iterations
00430   if (params->isParameter("Maximum Iterations")) {
00431     maxIters_ = params->get("Maximum Iterations",maxIters_default_);
00432 
00433     // Update parameter in our list and in status test.
00434     params_->set("Maximum Iterations", maxIters_);
00435     if (maxIterTest_!=Teuchos::null)
00436       maxIterTest_->setMaxIters( maxIters_ );
00437   }
00438 
00439   // Check for blocksize
00440   if (params->isParameter("Block Size")) {
00441     blockSize_ = params->get("Block Size",blockSize_default_);    
00442     TEST_FOR_EXCEPTION(blockSize_ <= 0, std::invalid_argument,
00443            "Belos::PseudoBlockGmresSolMgr: \"Block Size\" must be strictly positive.");
00444 
00445     // Update parameter in our list.
00446     params_->set("Block Size", blockSize_);
00447   }
00448 
00449   // Check if the blocksize should be adaptive
00450   if (params->isParameter("Adapative Block Size")) {
00451     adaptiveBlockSize_ = params->get("Adaptive Block Size",adaptiveBlockSize_default_);
00452     
00453     // Update parameter in our list.
00454     params_->set("Adaptive Block Size", adaptiveBlockSize_);
00455   }
00456 
00457   // Check for the maximum number of blocks.
00458   if (params->isParameter("Num Blocks")) {
00459     numBlocks_ = params->get("Num Blocks",numBlocks_default_);
00460     TEST_FOR_EXCEPTION(numBlocks_ <= 0, std::invalid_argument,
00461            "Belos::PseudoBlockGmresSolMgr: \"Num Blocks\" must be strictly positive.");
00462 
00463     // Update parameter in our list.
00464     params_->set("Num Blocks", numBlocks_);
00465   }
00466 
00467   // Check to see if the timer label changed.
00468   if (params->isParameter("Timer Label")) {
00469     std::string tempLabel = params->get("Timer Label", label_default_);
00470 
00471     // Update parameter in our list and solver timer
00472     if (tempLabel != label_) {
00473       label_ = tempLabel;
00474       params_->set("Timer Label", label_);
00475       std::string solveLabel = label_ + ": PseudoBlockGmresSolMgr total solve time";
00476       timerSolve_ = Teuchos::TimeMonitor::getNewTimer(solveLabel);
00477     }
00478   }
00479 
00480   // Check if the orthogonalization changed.
00481   if (params->isParameter("Orthogonalization")) {
00482     std::string tempOrthoType = params->get("Orthogonalization",orthoType_default_);
00483     TEST_FOR_EXCEPTION( tempOrthoType != "DGKS" && tempOrthoType != "ICGS" && tempOrthoType != "IMGS", 
00484                         std::invalid_argument,
00485       "Belos::PseudoBlockGmresSolMgr: \"Orthogonalization\" must be either \"DGKS\",\"ICGS\", or \"IMGS\".");
00486     if (tempOrthoType != orthoType_) {
00487       orthoType_ = tempOrthoType;
00488       // Create orthogonalization manager
00489       if (orthoType_=="DGKS") {
00490   if (orthoKappa_ <= 0) {
00491     ortho_ = Teuchos::rcp( new DGKSOrthoManager<ScalarType,MV,OP>( label_ ) );
00492   }
00493   else {
00494     ortho_ = Teuchos::rcp( new DGKSOrthoManager<ScalarType,MV,OP>( label_ ) );
00495     Teuchos::rcp_dynamic_cast<DGKSOrthoManager<ScalarType,MV,OP> >(ortho_)->setDepTol( orthoKappa_ );
00496   }
00497       }
00498       else if (orthoType_=="ICGS") {
00499   ortho_ = Teuchos::rcp( new ICGSOrthoManager<ScalarType,MV,OP>( label_ ) );
00500       } 
00501       else if (orthoType_=="IMGS") {
00502   ortho_ = Teuchos::rcp( new IMGSOrthoManager<ScalarType,MV,OP>( label_ ) );
00503       } 
00504     }  
00505   }
00506 
00507   // Check which orthogonalization constant to use.
00508   if (params->isParameter("Orthogonalization Constant")) {
00509     orthoKappa_ = params->get("Orthogonalization Constant",orthoKappa_default_);
00510 
00511     // Update parameter in our list.
00512     params_->set("Orthogonalization Constant",orthoKappa_);
00513     if (orthoType_=="DGKS") {
00514       if (orthoKappa_ > 0 && ortho_ != Teuchos::null) {
00515   Teuchos::rcp_dynamic_cast<DGKSOrthoManager<ScalarType,MV,OP> >(ortho_)->setDepTol( orthoKappa_ );
00516       }
00517     } 
00518   }
00519 
00520   // Check for a change in verbosity level
00521   if (params->isParameter("Verbosity")) {
00522     if (Teuchos::isParameterType<int>(*params,"Verbosity")) {
00523       verbosity_ = params->get("Verbosity", verbosity_default_);
00524     } else {
00525       verbosity_ = (int)Teuchos::getParameter<Belos::MsgType>(*params,"Verbosity");
00526     }
00527 
00528     // Update parameter in our list.
00529     params_->set("Verbosity", verbosity_);
00530     if (printer_ != Teuchos::null)
00531       printer_->setVerbosity(verbosity_);
00532   }
00533 
00534   // output stream
00535   if (params->isParameter("Output Stream")) {
00536     outputStream_ = Teuchos::getParameter<Teuchos::RCP<std::ostream> >(*params,"Output Stream");
00537 
00538     // Update parameter in our list.
00539     params_->set("Output Stream", outputStream_);
00540     if (printer_ != Teuchos::null)
00541       printer_->setOStream( outputStream_ );
00542   }
00543 
00544   // frequency level
00545   if (verbosity_ & Belos::StatusTestDetails) {
00546     if (params->isParameter("Output Frequency")) {
00547       outputFreq_ = params->get("Output Frequency", outputFreq_default_);
00548     }
00549 
00550     // Update parameter in out list and output status test.
00551     params_->set("Output Frequency", outputFreq_);
00552     if (outputTest_ != Teuchos::null)
00553       outputTest_->setOutputFrequency( outputFreq_ );
00554   }
00555 
00556   // Create output manager if we need to.
00557   if (printer_ == Teuchos::null) {
00558     printer_ = Teuchos::rcp( new OutputManager<ScalarType>(verbosity_, outputStream_) );
00559   }
00560   
00561   // Convergence
00562   typedef Belos::StatusTestCombo<ScalarType,MV,OP>  StatusTestCombo_t;
00563   typedef Belos::StatusTestGenResNorm<ScalarType,MV,OP>  StatusTestResNorm_t;
00564 
00565   // Check for convergence tolerance
00566   if (params->isParameter("Convergence Tolerance")) {
00567     convtol_ = params->get("Convergence Tolerance",convtol_default_);
00568 
00569     // Update parameter in our list and residual tests.
00570     params_->set("Convergence Tolerance", convtol_);
00571     if (impConvTest_ != Teuchos::null)
00572       impConvTest_->setTolerance( convtol_ );
00573     if (expConvTest_ != Teuchos::null)
00574       expConvTest_->setTolerance( convtol_ );
00575   }
00576 
00577   // Check for a change in scaling, if so we need to build new residual tests.
00578   bool newImpResTest = false, newExpResTest = false; 
00579   if (params->isParameter("Implicit Residual Scaling")) {
00580     std::string tempImpResScale = Teuchos::getParameter<std::string>( *params, "Implicit Residual Scaling" );
00581 
00582     // Only update the scaling if it's different.
00583     if (impResScale_ != tempImpResScale) {
00584       Belos::ScaleType impResScaleType = convertStringToScaleType( tempImpResScale );
00585       impResScale_ = tempImpResScale;
00586 
00587       // Update parameter in our list and residual tests
00588       params_->set("Implicit Residual Scaling", impResScale_);
00589       if (impConvTest_ != Teuchos::null) {
00590         try {
00591           impConvTest_->defineScaleForm( impResScaleType, Belos::TwoNorm );
00592         }
00593         catch (std::exception& e) {
00594     // Make sure the convergence test gets constructed again.
00595           isSTSet_ = false;
00596     newImpResTest = true;
00597         }
00598       }
00599     }      
00600   }
00601   
00602   if (params->isParameter("Explicit Residual Scaling")) {
00603     std::string tempExpResScale = Teuchos::getParameter<std::string>( *params, "Explicit Residual Scaling" );
00604 
00605     // Only update the scaling if it's different.
00606     if (expResScale_ != tempExpResScale) {
00607       Belos::ScaleType expResScaleType = convertStringToScaleType( tempExpResScale );
00608       expResScale_ = tempExpResScale;
00609 
00610       // Update parameter in our list and residual tests
00611       params_->set("Explicit Residual Scaling", expResScale_);
00612       if (expConvTest_ != Teuchos::null) {
00613         try {
00614           expConvTest_->defineScaleForm( expResScaleType, Belos::TwoNorm );
00615         }
00616         catch (std::exception& e) {
00617     // Make sure the convergence test gets constructed again.
00618           isSTSet_ = false;
00619     newExpResTest = true;
00620         }
00621       }
00622     }      
00623   }
00624 
00625 
00626   if (params->isParameter("Show Maximum Residual Norm Only")) {
00627     showMaxResNormOnly_ = Teuchos::getParameter<bool>(*params,"Show Maximum Residual Norm Only");
00628 
00629     // Update parameter in our list and residual tests
00630     params_->set("Show Maximum Residual Norm Only", showMaxResNormOnly_);
00631     if (impConvTest_ != Teuchos::null)
00632       impConvTest_->setShowMaxResNormOnly( showMaxResNormOnly_ );
00633     if (expConvTest_ != Teuchos::null)
00634       expConvTest_->setShowMaxResNormOnly( showMaxResNormOnly_ );
00635   }
00636 
00637   // Create status tests if we need to.
00638 
00639   // Get the deflation quorum, or number of converged systems before deflation is allowed
00640   if (params->isParameter("Deflation Quorum")) {
00641     defQuorum_ = params->get("Deflation Quorum", defQuorum_);
00642     TEST_FOR_EXCEPTION(defQuorum_ > blockSize_, std::invalid_argument,
00643            "Belos::PseudoBlockGmresSolMgr: \"Deflation Quorum\" cannot be larger than \"Block Size\".");
00644     params_->set("Deflation Quorum", defQuorum_);
00645     if (impConvTest_ != Teuchos::null)
00646       impConvTest_->setQuorum( defQuorum_ );
00647     if (expConvTest_ != Teuchos::null)
00648       expConvTest_->setQuorum( defQuorum_ );
00649   }
00650 
00651   // Create orthogonalization manager if we need to.
00652   if (ortho_ == Teuchos::null) {
00653     if (orthoType_=="DGKS") {
00654       if (orthoKappa_ <= 0) {
00655   ortho_ = Teuchos::rcp( new DGKSOrthoManager<ScalarType,MV,OP>( label_ ) );
00656       }
00657       else {
00658   ortho_ = Teuchos::rcp( new DGKSOrthoManager<ScalarType,MV,OP>( label_ ) );
00659   Teuchos::rcp_dynamic_cast<DGKSOrthoManager<ScalarType,MV,OP> >(ortho_)->setDepTol( orthoKappa_ );
00660       }
00661     }
00662     else if (orthoType_=="ICGS") {
00663       ortho_ = Teuchos::rcp( new ICGSOrthoManager<ScalarType,MV,OP>( label_ ) );
00664     } 
00665     else if (orthoType_=="IMGS") {
00666       ortho_ = Teuchos::rcp( new IMGSOrthoManager<ScalarType,MV,OP>( label_ ) );
00667     } 
00668     else {
00669       TEST_FOR_EXCEPTION(orthoType_!="ICGS"&&orthoType_!="DGKS"&&orthoType_!="IMGS",std::logic_error,
00670        "Belos::PseudoBlockGmresSolMgr(): Invalid orthogonalization type.");
00671     }  
00672   }
00673 
00674   // Create the timer if we need to.
00675   if (timerSolve_ == Teuchos::null) {
00676     std::string solveLabel = label_ + ": PseudoBlockGmresSolMgr total solve time";
00677     timerSolve_ = Teuchos::TimeMonitor::getNewTimer(solveLabel);
00678   }
00679 
00680   // Inform the solver manager that the current parameters were set.
00681   isSet_ = true;
00682 }
00683 
00684     
00685 template<class ScalarType, class MV, class OP>
00686 Teuchos::RCP<const Teuchos::ParameterList>
00687 PseudoBlockGmresSolMgr<ScalarType,MV,OP>::getValidParameters() const
00688 {
00689   static Teuchos::RCP<const Teuchos::ParameterList> validPL;
00690   if (is_null(validPL)) {
00691     Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
00692   // Set all the valid parameters and their default values.
00693     pl= Teuchos::rcp( new Teuchos::ParameterList() );
00694     pl->set("Convergence Tolerance", convtol_default_,
00695       "The relative residual tolerance that needs to be achieved by the\n"
00696       "iterative solver in order for the linera system to be declared converged.");
00697     pl->set("Maximum Restarts", maxRestarts_default_,
00698       "The maximum number of restarts allowed for each\n"
00699       "set of RHS solved.");
00700     pl->set("Maximum Iterations", maxIters_default_,
00701       "The maximum number of block iterations allowed for each\n"
00702       "set of RHS solved.");
00703     pl->set("Num Blocks", numBlocks_default_,
00704       "The maximum number of vectors allowed in the Krylov subspace\n"
00705       "for each set of RHS solved.");
00706     pl->set("Block Size", blockSize_default_,
00707       "The number of RHS solved simultaneously.");
00708     pl->set("Adaptive Block Size", adaptiveBlockSize_default_,
00709       "Whether the solver manager should adapt the block size\n"
00710       "based on the number of RHS to solve.");
00711     pl->set("Verbosity", verbosity_default_,
00712       "What type(s) of solver information should be outputted\n"
00713       "to the output stream.");
00714     pl->set("Output Frequency", outputFreq_default_,
00715       "How often convergence information should be outputted\n"
00716       "to the output stream.");  
00717     pl->set("Deflation Quorum", defQuorum_default_,
00718       "The number of linear systems that need to converge before\n"
00719       "they are deflated.  This number should be <= block size.");
00720     pl->set("Output Stream", outputStream_default_,
00721       "A reference-counted pointer to the output stream where all\n"
00722       "solver output is sent.");
00723     pl->set("Show Maximum Residual Norm Only", showMaxResNormOnly_default_,
00724       "When convergence information is printed, only show the maximum\n"
00725       "relative residual norm when the block size is greater than one.");
00726     pl->set("Implicit Residual Scaling", impResScale_default_,
00727       "The type of scaling used in the implicit residual convergence test.");
00728     pl->set("Explicit Residual Scaling", expResScale_default_,
00729       "The type of scaling used in the explicit residual convergence test.");
00730     pl->set("Timer Label", label_default_,
00731       "The string to use as a prefix for the timer labels.");
00732     //  defaultParams_->set("Restart Timers", restartTimers_);
00733     pl->set("Orthogonalization", orthoType_default_,
00734       "The type of orthogonalization to use: DGKS, ICGS, IMGS.");
00735     pl->set("Orthogonalization Constant",orthoKappa_default_,
00736       "The constant used by DGKS orthogonalization to determine\n"
00737       "whether another step of classical Gram-Schmidt is necessary.");
00738     validPL = pl;
00739   }
00740   return validPL;
00741 }
00742 
00743 // Check the status test versus the defined linear problem
00744 template<class ScalarType, class MV, class OP>
00745 bool PseudoBlockGmresSolMgr<ScalarType,MV,OP>::checkStatusTest() {
00746 
00747   typedef Belos::StatusTestCombo<ScalarType,MV,OP>  StatusTestCombo_t;
00748   typedef Belos::StatusTestGenResNorm<ScalarType,MV,OP>  StatusTestGenResNorm_t;
00749   typedef Belos::StatusTestImpResNorm<ScalarType,MV,OP>  StatusTestImpResNorm_t;
00750 
00751   // Basic test checks maximum iterations and native residual.
00752   maxIterTest_ = Teuchos::rcp( new StatusTestMaxIters<ScalarType,MV,OP>( maxIters_ ) );
00753 
00754   // If there is a left preconditioner, we create a combined status test that checks the implicit
00755   // and then explicit residual norm to see if we have convergence.
00756   if ( !Teuchos::is_null(problem_->getLeftPrec()) ) {
00757     expResTest_ = true;
00758   }
00759 
00760   if (expResTest_) {
00761    
00762     // Implicit residual test, using the native residual to determine if convergence was achieved.
00763     Teuchos::RCP<StatusTestGenResNorm_t> tmpImpConvTest =
00764       Teuchos::rcp( new StatusTestGenResNorm_t( convtol_, defQuorum_ ) );
00765     tmpImpConvTest->defineScaleForm( convertStringToScaleType(impResScale_), Belos::TwoNorm );
00766     tmpImpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
00767     impConvTest_ = tmpImpConvTest;
00768 
00769     // Explicit residual test once the native residual is below the tolerance
00770     Teuchos::RCP<StatusTestGenResNorm_t> tmpExpConvTest =
00771       Teuchos::rcp( new StatusTestGenResNorm_t( convtol_, defQuorum_ ) );
00772     tmpExpConvTest->defineResForm( StatusTestGenResNorm_t::Explicit, Belos::TwoNorm );
00773     tmpExpConvTest->defineScaleForm( convertStringToScaleType(expResScale_), Belos::TwoNorm );
00774     tmpExpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
00775     expConvTest_ = tmpExpConvTest;
00776 
00777     // The convergence test is a combination of the "cheap" implicit test and explicit test.
00778     convTest_ = Teuchos::rcp( new StatusTestCombo_t( StatusTestCombo_t::SEQ, impConvTest_, expConvTest_ ) );
00779   }
00780   else {
00781 
00782     // Implicit residual test, using the native residual to determine if convergence was achieved.
00783     // Use test that checks for loss of accuracy.
00784     Teuchos::RCP<StatusTestImpResNorm_t> tmpImpConvTest =
00785       Teuchos::rcp( new StatusTestImpResNorm_t( convtol_, defQuorum_ ) );
00786     tmpImpConvTest->defineScaleForm( convertStringToScaleType(impResScale_), Belos::TwoNorm );
00787     tmpImpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
00788     impConvTest_ = tmpImpConvTest;
00789     
00790     // Set the explicit and total convergence test to this implicit test that checks for accuracy loss.
00791     expConvTest_ = impConvTest_;
00792     convTest_ = impConvTest_;
00793   }
00794 
00795   sTest_ = Teuchos::rcp( new StatusTestCombo_t( StatusTestCombo_t::OR, maxIterTest_, convTest_ ) );
00796 
00797   if (outputFreq_ > 0) {
00798     outputTest_ = Teuchos::rcp( new StatusTestOutput<ScalarType,MV,OP>( printer_,
00799         sTest_,
00800         outputFreq_,
00801         Passed+Failed+Undefined ) );
00802   }
00803   else {
00804     outputTest_ = Teuchos::rcp( new StatusTestOutput<ScalarType,MV,OP>( printer_,
00805         sTest_, 1 ) );
00806   }
00807 
00808   // The status test is now set.
00809   isSTSet_ = true;
00810 
00811   return false;
00812 }
00813 
00814 
00815 // solve()
00816 template<class ScalarType, class MV, class OP>
00817 ReturnType PseudoBlockGmresSolMgr<ScalarType,MV,OP>::solve() {
00818 
00819   // Set the current parameters if they were not set before.
00820   // NOTE:  This may occur if the user generated the solver manager with the default constructor and 
00821   // then didn't set any parameters using setParameters().
00822   if (!isSet_) { setParameters( params_ ); }
00823   
00824   Teuchos::BLAS<int,ScalarType> blas;
00825   
00826   TEST_FOR_EXCEPTION(!problem_->isProblemSet(),PseudoBlockGmresSolMgrLinearProblemFailure,
00827                      "Belos::PseudoBlockGmresSolMgr::solve(): Linear problem is not ready, setProblem() has not been called.");
00828 
00829   // Check if we have to create the status tests.
00830   if (!isSTSet_ || (!expResTest_ && !Teuchos::is_null(problem_->getLeftPrec())) ) {
00831     TEST_FOR_EXCEPTION( checkStatusTest(),PseudoBlockGmresSolMgrLinearProblemFailure,
00832       "Belos::BlockGmresSolMgr::solve(): Linear problem and requested status tests are incompatible.");
00833   }
00834 
00835   // Create indices for the linear systems to be solved.
00836   int startPtr = 0;
00837   int numRHS2Solve = MVT::GetNumberVecs( *(problem_->getRHS()) );
00838   int numCurrRHS = ( numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;
00839 
00840   std::vector<int> currIdx;
00841   //  If an adaptive block size is allowed then only the linear systems that need to be solved are solved.
00842   //  Otherwise, the index set is generated that informs the linear problem that some linear systems are augmented.
00843   if ( adaptiveBlockSize_ ) {
00844     blockSize_ = numCurrRHS;
00845     currIdx.resize( numCurrRHS  );
00846     for (int i=0; i<numCurrRHS; ++i) 
00847       { currIdx[i] = startPtr+i; }
00848     
00849   }
00850   else {
00851     currIdx.resize( blockSize_ );
00852     for (int i=0; i<numCurrRHS; ++i) 
00853       { currIdx[i] = startPtr+i; }
00854     for (int i=numCurrRHS; i<blockSize_; ++i) 
00855       { currIdx[i] = -1; }
00856   }
00857 
00858   // Inform the linear problem of the current linear system to solve.
00859   problem_->setLSIndex( currIdx );
00860 
00862   // Parameter list
00863   Teuchos::ParameterList plist;
00864   plist.set("Num Blocks",numBlocks_);
00865   
00866   // Reset the status test.  
00867   outputTest_->reset();
00868   loaDetected_ = false;
00869 
00870   // Assume convergence is achieved, then let any failed convergence set this to false.
00871   bool isConverged = true;  
00872 
00874   // BlockGmres solver
00875 
00876   Teuchos::RCP<PseudoBlockGmresIter<ScalarType,MV,OP> > block_gmres_iter
00877     = Teuchos::rcp( new PseudoBlockGmresIter<ScalarType,MV,OP>(problem_,printer_,outputTest_,ortho_,plist) );  
00878 
00879   // Enter solve() iterations
00880   {
00881     Teuchos::TimeMonitor slvtimer(*timerSolve_);
00882 
00883     while ( numRHS2Solve > 0 ) {
00884 
00885       // Reset the active / converged vectors from this block
00886       std::vector<int> convRHSIdx;
00887       std::vector<int> currRHSIdx( currIdx );
00888       currRHSIdx.resize(numCurrRHS);
00889 
00890       // Set the current number of blocks with the pseudo Gmres iteration.
00891       block_gmres_iter->setNumBlocks( numBlocks_ );
00892 
00893       // Reset the number of iterations.
00894       block_gmres_iter->resetNumIters();
00895 
00896       // Reset the number of calls that the status test output knows about.
00897       outputTest_->resetNumCalls();
00898 
00899       // Get a new state struct and initialize the solver.
00900       PseudoBlockGmresIterState<ScalarType,MV> newState;
00901 
00902       // Create the first block in the current Krylov basis for each right-hand side.
00903       std::vector<int> index(1);
00904       Teuchos::RCP<MV> tmpV, R_0 = MVT::CloneCopy( *(problem_->getInitPrecResVec()), currIdx );
00905       newState.V.resize( blockSize_ );
00906       newState.Z.resize( blockSize_ );
00907       for (int i=0; i<blockSize_; ++i) {
00908   index[0]=i;
00909   tmpV = MVT::CloneView( *R_0, index );
00910   
00911   // Get a matrix to hold the orthonormalization coefficients.
00912   Teuchos::RCP<Teuchos::SerialDenseVector<int,ScalarType> > tmpZ
00913     = Teuchos::rcp( new Teuchos::SerialDenseVector<int,ScalarType>( 1 ));
00914       
00915   // Orthonormalize the new V_0
00916   int rank = ortho_->normalize( *tmpV, tmpZ );
00917   TEST_FOR_EXCEPTION(rank != 1, PseudoBlockGmresSolMgrOrthoFailure,
00918          "Belos::PseudoBlockGmresSolMgr::solve(): Failed to compute initial block of orthonormal vectors.");
00919 
00920   newState.V[i] = tmpV;
00921   newState.Z[i] = tmpZ;
00922       }
00923 
00924       newState.curDim = 0;
00925       block_gmres_iter->initialize(newState);
00926       int numRestarts = 0;
00927 
00928       while(1) {
00929   
00930   // tell block_gmres_iter to iterate
00931   try {
00932     block_gmres_iter->iterate();
00933     
00935     //
00936     // check convergence first
00937     //
00939     if ( convTest_->getStatus() == Passed ) {
00940 
00941       if ( expConvTest_->getLOADetected() ) {
00942               // we don't have convergence but we will deflate out the linear systems and move on.
00943               loaDetected_ = true;
00944               isConverged = false;
00945             }
00946 
00947       // Figure out which linear systems converged.
00948       std::vector<int> convIdx = expConvTest_->convIndices();
00949 
00950       // If the number of converged linear systems is equal to the
00951             // number of current linear systems, then we are done with this block.
00952       if (convIdx.size() == currRHSIdx.size())
00953         break;  // break from while(1){block_gmres_iter->iterate()}
00954 
00955       // Get a new state struct and initialize the solver.
00956       PseudoBlockGmresIterState<ScalarType,MV> defState;
00957 
00958       // Inform the linear problem that we are finished with this current linear system.
00959       problem_->setCurrLS();
00960 
00961       // Get the state.
00962       PseudoBlockGmresIterState<ScalarType,MV> oldState = block_gmres_iter->getState();
00963       int curDim = oldState.curDim;
00964 
00965       // Get a new state struct and reset currRHSIdx to have the right-hand sides that 
00966       // are left to converge for this block.
00967       int have = 0;
00968       std::vector<int> oldRHSIdx( currRHSIdx );
00969       std::vector<int> defRHSIdx;
00970       for (unsigned int i=0; i<currRHSIdx.size(); ++i) {
00971         bool found = false;
00972         for (unsigned int j=0; j<convIdx.size(); ++j) {
00973     if (currRHSIdx[i] == convIdx[j]) {
00974       found = true;
00975       break;
00976     }
00977         }
00978         if (found) {
00979     defRHSIdx.push_back( i );
00980         }
00981         else {
00982     defState.V.push_back( Teuchos::rcp_const_cast<MV>( oldState.V[i] ) );
00983     defState.Z.push_back( Teuchos::rcp_const_cast<Teuchos::SerialDenseVector<int,ScalarType> >( oldState.Z[i] ) );
00984     defState.H.push_back( Teuchos::rcp_const_cast<Teuchos::SerialDenseMatrix<int,ScalarType> >( oldState.H[i] ) );
00985     defState.sn.push_back( Teuchos::rcp_const_cast<Teuchos::SerialDenseVector<int,ScalarType> >( oldState.sn[i] ) );
00986     defState.cs.push_back( Teuchos::rcp_const_cast<Teuchos::SerialDenseVector<int,MagnitudeType> >(oldState.cs[i] ) );
00987     currRHSIdx[have] = currRHSIdx[i];
00988     have++;
00989         }
00990       }
00991       defRHSIdx.resize(currRHSIdx.size()-have);
00992       currRHSIdx.resize(have);
00993 
00994       // Compute the current solution that needs to be deflated if this solver has taken any steps.
00995       if (curDim) {
00996         Teuchos::RCP<MV> update = block_gmres_iter->getCurrentUpdate();
00997         Teuchos::RCP<MV> defUpdate = MVT::CloneView( *update, defRHSIdx );
00998         
00999         // Set the deflated indices so we can update the solution.
01000         problem_->setLSIndex( convIdx );
01001         
01002         // Update the linear problem.
01003         problem_->updateSolution( defUpdate, true );
01004       }
01005       
01006       // Set the remaining indices after deflation.
01007       problem_->setLSIndex( currRHSIdx );
01008       
01009       // Set the dimension of the subspace, which is the same as the old subspace size.
01010       defState.curDim = curDim;
01011       
01012       // Initialize the solver with the deflated system.
01013       block_gmres_iter->initialize(defState);
01014     }
01016     //
01017     // check for maximum iterations
01018     //
01020     else if ( maxIterTest_->getStatus() == Passed ) {
01021       // we don't have convergence
01022       isConverged = false;
01023       break;  // break from while(1){block_gmres_iter->iterate()}
01024     }
01026     //
01027     // check for restarting, i.e. the subspace is full
01028     //
01030     else if ( block_gmres_iter->getCurSubspaceDim() == block_gmres_iter->getMaxSubspaceDim() ) {
01031   
01032       if ( numRestarts >= maxRestarts_ ) {
01033         isConverged = false;
01034         break; // break from while(1){block_gmres_iter->iterate()}
01035       }
01036       numRestarts++;
01037 
01038       printer_->stream(Debug) << " Performing restart number " << numRestarts << " of " << maxRestarts_ << std::endl << std::endl;
01039       
01040       // Update the linear problem.
01041       Teuchos::RCP<MV> update = block_gmres_iter->getCurrentUpdate();
01042       problem_->updateSolution( update, true );
01043       
01044       // Get the state.
01045       PseudoBlockGmresIterState<ScalarType,MV> oldState = block_gmres_iter->getState();
01046       
01047       // Set the new state.
01048       PseudoBlockGmresIterState<ScalarType,MV> newstate;
01049       newstate.V.resize(currRHSIdx.size());
01050       newstate.Z.resize(currRHSIdx.size());
01051 
01052       // Compute the restart vectors
01053       // NOTE: Force the linear problem to update the current residual since the solution was updated.
01054       R_0 = MVT::Clone( *(problem_->getInitPrecResVec()), currRHSIdx.size() );
01055       problem_->computeCurrPrecResVec( &*R_0 );
01056       for (unsigned int i=0; i<currRHSIdx.size(); ++i) {
01057         index[0] = i;  // index(1) vector declared on line 891
01058   
01059         tmpV = MVT::CloneView( *R_0, index );
01060   
01061         // Get a matrix to hold the orthonormalization coefficients.
01062         Teuchos::RCP<Teuchos::SerialDenseVector<int,ScalarType> > tmpZ
01063     = Teuchos::rcp( new Teuchos::SerialDenseVector<int,ScalarType>( 1 ));
01064         
01065         // Orthonormalize the new V_0
01066         int rank = ortho_->normalize( *tmpV, tmpZ );
01067         TEST_FOR_EXCEPTION(rank != 1 ,PseudoBlockGmresSolMgrOrthoFailure,
01068          "Belos::PseudoBlockGmresSolMgr::solve(): Failed to compute initial block of orthonormal vectors after the restart.");
01069         
01070         newstate.V[i] = tmpV;
01071         newstate.Z[i] = tmpZ;
01072       }
01073 
01074       // Initialize the solver.
01075       newstate.curDim = 0;
01076       block_gmres_iter->initialize(newstate);
01077 
01078     } // end of restarting
01079 
01081     //
01082     // we returned from iterate(), but none of our status tests Passed.
01083     // something is wrong, and it is probably our fault.
01084     //
01086 
01087     else {
01088       TEST_FOR_EXCEPTION(true,std::logic_error,
01089              "Belos::PseudoBlockGmresSolMgr::solve(): Invalid return from PseudoBlockGmresIter::iterate().");
01090     }
01091   }
01092         catch (const PseudoBlockGmresIterOrthoFailure &e) {
01093      
01094     // Try to recover the most recent least-squares solution
01095     block_gmres_iter->updateLSQR( block_gmres_iter->getCurSubspaceDim() );
01096 
01097     // Check to see if the most recent least-squares solution yielded convergence.
01098     sTest_->checkStatus( &*block_gmres_iter );
01099     if (convTest_->getStatus() != Passed)
01100       isConverged = false;
01101     break;
01102         }
01103   catch (const std::exception &e) {
01104     printer_->stream(Errors) << "Error! Caught std::exception in PseudoBlockGmresIter::iterate() at iteration " 
01105           << block_gmres_iter->getNumIters() << std::endl 
01106           << e.what() << std::endl;
01107     throw;
01108   }
01109       }
01110       
01111       // Compute the current solution.
01112       // Update the linear problem.
01113       if ( !Teuchos::is_null(expConvTest_->getSolution()) ) {
01114         Teuchos::RCP<MV> newX = expConvTest_->getSolution();
01115         Teuchos::RCP<MV> curX = problem_->getCurrLHSVec();
01116         MVT::MvAddMv( 0.0, *newX, 1.0, *newX, *curX );
01117       }
01118       else {
01119         Teuchos::RCP<MV> update = block_gmres_iter->getCurrentUpdate();
01120         problem_->updateSolution( update, true );
01121       }
01122 
01123       // Inform the linear problem that we are finished with this block linear system.
01124       problem_->setCurrLS();
01125       
01126       // Update indices for the linear systems to be solved.
01127       startPtr += numCurrRHS;
01128       numRHS2Solve -= numCurrRHS;
01129       if ( numRHS2Solve > 0 ) {
01130   numCurrRHS = ( numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;
01131 
01132   if ( adaptiveBlockSize_ ) {
01133           blockSize_ = numCurrRHS;
01134     currIdx.resize( numCurrRHS  );
01135     for (int i=0; i<numCurrRHS; ++i) 
01136       { currIdx[i] = startPtr+i; }    
01137     
01138     // Adapt the status test quorum if we need to.
01139     if (defQuorum_ > blockSize_) {
01140       if (impConvTest_ != Teuchos::null)
01141         impConvTest_->setQuorum( blockSize_ );
01142       if (expConvTest_ != Teuchos::null)
01143         expConvTest_->setQuorum( blockSize_ );
01144     }
01145   }
01146   else {
01147     currIdx.resize( blockSize_ );
01148     for (int i=0; i<numCurrRHS; ++i) 
01149       { currIdx[i] = startPtr+i; }
01150     for (int i=numCurrRHS; i<blockSize_; ++i) 
01151       { currIdx[i] = -1; }
01152   }
01153   // Set the next indices.
01154   problem_->setLSIndex( currIdx );
01155       }
01156       else {
01157         currIdx.resize( numRHS2Solve );
01158       }
01159       
01160     }// while ( numRHS2Solve > 0 )
01161     
01162   }
01163 
01164   // print final summary
01165   sTest_->print( printer_->stream(FinalSummary) );
01166  
01167   // print timing information
01168   Teuchos::TimeMonitor::summarize( printer_->stream(TimingDetails) );
01169  
01170   // get iteration information for this solve
01171   numIters_ = maxIterTest_->getNumIters();
01172  
01173   if (!isConverged || loaDetected_) {
01174     return Unconverged; // return from PseudoBlockGmresSolMgr::solve() 
01175   }
01176   return Converged; // return from PseudoBlockGmresSolMgr::solve() 
01177 }
01178 
01179 //  This method requires the solver manager to return a std::string that describes itself.
01180 template<class ScalarType, class MV, class OP>
01181 std::string PseudoBlockGmresSolMgr<ScalarType,MV,OP>::description() const
01182 {
01183   std::ostringstream oss;
01184   oss << "Belos::PseudoBlockGmresSolMgr<...,"<<Teuchos::ScalarTraits<ScalarType>::name()<<">";
01185   oss << "{";
01186   oss << "Ortho Type='"<<orthoType_<<"\', Block Size=" <<blockSize_;
01187   oss << ", Num Blocks="<<numBlocks_<< ", Max Restarts=" << maxRestarts_;
01188   oss << "}";
01189   return oss.str();
01190 }
01191   
01192 } // end Belos namespace
01193 
01194 #endif /* BELOS_PSEUDO_BLOCK_GMRES_SOLMGR_HPP */

Generated on Wed May 12 21:30:08 2010 for Belos by  doxygen 1.4.7