Belos Version of the Day
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 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_PSEUDO_BLOCK_GMRES_SOLMGR_HPP
00043 #define BELOS_PSEUDO_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 "BelosPseudoBlockGmresIter.hpp"
00056 #include "BelosDGKSOrthoManager.hpp"
00057 #include "BelosICGSOrthoManager.hpp"
00058 #include "BelosIMGSOrthoManager.hpp"
00059 #include "BelosStatusTestMaxIters.hpp"
00060 #include "BelosStatusTestGenResNorm.hpp"
00061 #include "BelosStatusTestImpResNorm.hpp"
00062 #include "BelosStatusTestCombo.hpp"
00063 #include "BelosStatusTestOutputFactory.hpp"
00064 #include "BelosOutputManager.hpp"
00065 #include "Teuchos_BLAS.hpp"
00066 #ifdef BELOS_TEUCHOS_TIME_MONITOR
00067 #include "Teuchos_TimeMonitor.hpp"
00068 #endif
00069 
00086 namespace Belos {
00087   
00089 
00090   
00097   class PseudoBlockGmresSolMgrLinearProblemFailure : public BelosError {public:
00098     PseudoBlockGmresSolMgrLinearProblemFailure(const std::string& what_arg) : BelosError(what_arg)
00099     {}};
00100   
00107   class PseudoBlockGmresSolMgrOrthoFailure : public BelosError {public:
00108     PseudoBlockGmresSolMgrOrthoFailure(const std::string& what_arg) : BelosError(what_arg)
00109     {}};
00110   
00111   template<class ScalarType, class MV, class OP>
00112   class PseudoBlockGmresSolMgr : public SolverManager<ScalarType,MV,OP> {
00113     
00114   private:
00115     typedef MultiVecTraits<ScalarType,MV> MVT;
00116     typedef OperatorTraits<ScalarType,MV,OP> OPT;
00117     typedef Teuchos::ScalarTraits<ScalarType> SCT;
00118     typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType MagnitudeType;
00119     typedef Teuchos::ScalarTraits<MagnitudeType> MT;
00120     
00121   public:
00122     
00124 
00125 
00131     PseudoBlockGmresSolMgr();
00132     
00145     PseudoBlockGmresSolMgr( const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem,
00146                             const Teuchos::RCP<Teuchos::ParameterList> &pl );
00147     
00149     virtual ~PseudoBlockGmresSolMgr() {};
00151     
00153 
00154     
00155     const LinearProblem<ScalarType,MV,OP>& getProblem() const {
00156       return *problem_;
00157     }
00158     
00161     Teuchos::RCP<const Teuchos::ParameterList> getValidParameters() const;
00162    
00165     Teuchos::RCP<const Teuchos::ParameterList> getCurrentParameters() const { return params_; }
00166  
00172     Teuchos::Array<Teuchos::RCP<Teuchos::Time> > getTimers() const {
00173       return Teuchos::tuple(timerSolve_);
00174     }
00175 
00177     int getNumIters() const {
00178       return numIters_;
00179     }
00180   
00237     bool isLOADetected() const { return loaDetected_; }
00238   
00240     
00242 
00243    
00245     void setProblem (const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem) {
00246       problem_ = problem; 
00247     }
00248    
00250     void setParameters (const Teuchos::RCP<Teuchos::ParameterList> &params);
00251 
00258     virtual void setUserConvStatusTest(
00259       const Teuchos::RCP<StatusTest<ScalarType,MV,OP> > &userConvStatusTest
00260       );
00261     
00263 
00265 
00266 
00270     void reset( const ResetType type ) { if ((type & Belos::Problem) && !Teuchos::is_null(problem_)) problem_->setProblem(); }
00272 
00274 
00275     
00293     ReturnType solve();
00294     
00296     
00299     
00301     std::string description() const;
00302     
00304     
00305   private:
00306 
00308     Belos::ScaleType convertStringToScaleType (const std::string& scaleType) {
00309       if (scaleType == "Norm of Initial Residual") {
00310         return Belos::NormOfInitRes;
00311       } else if (scaleType == "Norm of Preconditioned Initial Residual") {
00312         return Belos::NormOfPrecInitRes;
00313       } else if (scaleType == "Norm of RHS") {
00314         return Belos::NormOfRHS;
00315       } else if (scaleType == "None") {
00316         return Belos::None;
00317       } else 
00318         TEST_FOR_EXCEPTION( true ,std::logic_error,
00319             "Belos::PseudoBlockGmresSolMgr(): Invalid residual scaling type.");
00320     }
00321 
00336     bool checkStatusTest();
00337     
00339     Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > problem_;
00340     
00341     // Output manager.
00342     Teuchos::RCP<OutputManager<ScalarType> > printer_;
00343     Teuchos::RCP<std::ostream> outputStream_;
00344 
00345     // Status tests.
00346     Teuchos::RCP<StatusTest<ScalarType,MV,OP> > userConvStatusTest_;
00347     Teuchos::RCP<StatusTest<ScalarType,MV,OP> > sTest_;
00348     Teuchos::RCP<StatusTestMaxIters<ScalarType,MV,OP> > maxIterTest_;
00349     Teuchos::RCP<StatusTest<ScalarType,MV,OP> > convTest_;
00350     Teuchos::RCP<StatusTestResNorm<ScalarType,MV,OP> > impConvTest_, expConvTest_;
00351     Teuchos::RCP<StatusTestOutput<ScalarType,MV,OP> > outputTest_;
00352 
00353     // Orthogonalization manager.
00354     Teuchos::RCP<MatOrthoManager<ScalarType,MV,OP> > ortho_; 
00355 
00356      // Current parameter list.
00357     Teuchos::RCP<Teuchos::ParameterList> params_;
00358    
00359     // Default solver values.
00360     static const MagnitudeType convtol_default_;
00361     static const MagnitudeType orthoKappa_default_;
00362     static const int maxRestarts_default_;
00363     static const int maxIters_default_;
00364     static const bool showMaxResNormOnly_default_;
00365     static const int blockSize_default_;
00366     static const int numBlocks_default_;
00367     static const int verbosity_default_;
00368     static const int outputStyle_default_;
00369     static const int outputFreq_default_;
00370     static const int defQuorum_default_;
00371     static const std::string impResScale_default_; 
00372     static const std::string expResScale_default_; 
00373     static const std::string label_default_;
00374     static const std::string orthoType_default_;
00375     static const Teuchos::RCP<std::ostream> outputStream_default_;
00376 
00377     // Current solver values.
00378     MagnitudeType convtol_, orthoKappa_;
00379     int maxRestarts_, maxIters_, numIters_;
00380     int blockSize_, numBlocks_, verbosity_, outputStyle_, outputFreq_, defQuorum_;
00381     bool showMaxResNormOnly_;
00382     std::string orthoType_; 
00383     std::string impResScale_, expResScale_;       
00384  
00385     // Timers.
00386     std::string label_;
00387     Teuchos::RCP<Teuchos::Time> timerSolve_;
00388 
00389     // Internal state variables.
00390     bool isSet_, isSTSet_, expResTest_;
00391     bool loaDetected_;
00392   };
00393   
00394   
00395 // Default solver values.
00396 template<class ScalarType, class MV, class OP>
00397 const typename PseudoBlockGmresSolMgr<ScalarType,MV,OP>::MagnitudeType PseudoBlockGmresSolMgr<ScalarType,MV,OP>::convtol_default_ = 1e-8;
00398 
00399 template<class ScalarType, class MV, class OP>
00400 const typename PseudoBlockGmresSolMgr<ScalarType,MV,OP>::MagnitudeType PseudoBlockGmresSolMgr<ScalarType,MV,OP>::orthoKappa_default_ = -1.0;
00401 
00402 template<class ScalarType, class MV, class OP>
00403 const int PseudoBlockGmresSolMgr<ScalarType,MV,OP>::maxRestarts_default_ = 20;
00404 
00405 template<class ScalarType, class MV, class OP>
00406 const int PseudoBlockGmresSolMgr<ScalarType,MV,OP>::maxIters_default_ = 1000;
00407 
00408 template<class ScalarType, class MV, class OP>
00409 const bool PseudoBlockGmresSolMgr<ScalarType,MV,OP>::showMaxResNormOnly_default_ = false;
00410 
00411 template<class ScalarType, class MV, class OP>
00412 const int PseudoBlockGmresSolMgr<ScalarType,MV,OP>::blockSize_default_ = 1;
00413 
00414 template<class ScalarType, class MV, class OP>
00415 const int PseudoBlockGmresSolMgr<ScalarType,MV,OP>::numBlocks_default_ = 300;
00416 
00417 template<class ScalarType, class MV, class OP>
00418 const int PseudoBlockGmresSolMgr<ScalarType,MV,OP>::verbosity_default_ = Belos::Errors;
00419 
00420 template<class ScalarType, class MV, class OP>
00421 const int PseudoBlockGmresSolMgr<ScalarType,MV,OP>::outputStyle_default_ = Belos::General;
00422 
00423 template<class ScalarType, class MV, class OP>
00424 const int PseudoBlockGmresSolMgr<ScalarType,MV,OP>::outputFreq_default_ = -1;
00425 
00426 template<class ScalarType, class MV, class OP>
00427 const int PseudoBlockGmresSolMgr<ScalarType,MV,OP>::defQuorum_default_ = 1;
00428 
00429 template<class ScalarType, class MV, class OP>
00430 const std::string PseudoBlockGmresSolMgr<ScalarType,MV,OP>::impResScale_default_ = "Norm of Preconditioned Initial Residual";
00431 
00432 template<class ScalarType, class MV, class OP>
00433 const std::string PseudoBlockGmresSolMgr<ScalarType,MV,OP>::expResScale_default_ = "Norm of Initial Residual";
00434 
00435 template<class ScalarType, class MV, class OP>
00436 const std::string PseudoBlockGmresSolMgr<ScalarType,MV,OP>::label_default_ = "Belos";
00437 
00438 template<class ScalarType, class MV, class OP>
00439 const std::string PseudoBlockGmresSolMgr<ScalarType,MV,OP>::orthoType_default_ = "DGKS";
00440 
00441 template<class ScalarType, class MV, class OP>
00442 const Teuchos::RCP<std::ostream> PseudoBlockGmresSolMgr<ScalarType,MV,OP>::outputStream_default_ = Teuchos::rcp(&std::cout,false);
00443 
00444 
00445 // Empty Constructor
00446 template<class ScalarType, class MV, class OP>
00447 PseudoBlockGmresSolMgr<ScalarType,MV,OP>::PseudoBlockGmresSolMgr() :
00448   outputStream_(outputStream_default_),
00449   convtol_(convtol_default_),
00450   orthoKappa_(orthoKappa_default_),
00451   maxRestarts_(maxRestarts_default_),
00452   maxIters_(maxIters_default_),
00453   blockSize_(blockSize_default_),
00454   numBlocks_(numBlocks_default_),
00455   verbosity_(verbosity_default_),
00456   outputStyle_(outputStyle_default_),
00457   outputFreq_(outputFreq_default_),
00458   defQuorum_(defQuorum_default_),
00459   showMaxResNormOnly_(showMaxResNormOnly_default_),
00460   orthoType_(orthoType_default_),
00461   impResScale_(impResScale_default_),
00462   expResScale_(expResScale_default_),
00463   label_(label_default_),
00464   isSet_(false),
00465   isSTSet_(false),
00466   expResTest_(false),
00467   loaDetected_(false)
00468 {}
00469 
00470 // Basic Constructor
00471 template<class ScalarType, class MV, class OP>
00472 PseudoBlockGmresSolMgr<ScalarType,MV,OP>::PseudoBlockGmresSolMgr( 
00473                                         const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem,
00474                                         const Teuchos::RCP<Teuchos::ParameterList> &pl ) : 
00475   problem_(problem),
00476   outputStream_(outputStream_default_),
00477   convtol_(convtol_default_),
00478   orthoKappa_(orthoKappa_default_),
00479   maxRestarts_(maxRestarts_default_),
00480   maxIters_(maxIters_default_),
00481   blockSize_(blockSize_default_),
00482   numBlocks_(numBlocks_default_),
00483   verbosity_(verbosity_default_),
00484   outputStyle_(outputStyle_default_),
00485   outputFreq_(outputFreq_default_),
00486   defQuorum_(defQuorum_default_),
00487   showMaxResNormOnly_(showMaxResNormOnly_default_),
00488   orthoType_(orthoType_default_),
00489   impResScale_(impResScale_default_),
00490   expResScale_(expResScale_default_),
00491   label_(label_default_),
00492   isSet_(false),
00493   isSTSet_(false),
00494   expResTest_(false),
00495   loaDetected_(false)
00496 {
00497   TEST_FOR_EXCEPTION(problem_ == Teuchos::null, std::invalid_argument, "Problem not given to solver manager.");
00498 
00499   // If the parameter list pointer is null, then set the current parameters to the default parameter list.
00500   if (!is_null(pl)) {
00501     // Set the parameters using the list that was passed in.
00502     setParameters( pl );  
00503   }
00504 }
00505 
00506 template<class ScalarType, class MV, class OP>
00507 void PseudoBlockGmresSolMgr<ScalarType,MV,OP>::setParameters( const Teuchos::RCP<Teuchos::ParameterList> &params )
00508 {
00509   // Create the internal parameter list if ones doesn't already exist.
00510   if (params_ == Teuchos::null) {
00511     params_ = Teuchos::rcp( new Teuchos::ParameterList(*getValidParameters()) );
00512   }
00513   else {
00514     params->validateParameters(*getValidParameters());
00515   }
00516 
00517   // Check for maximum number of restarts
00518   if (params->isParameter("Maximum Restarts")) {
00519     maxRestarts_ = params->get("Maximum Restarts",maxRestarts_default_);
00520 
00521     // Update parameter in our list.
00522     params_->set("Maximum Restarts", maxRestarts_);
00523   }
00524 
00525   // Check for maximum number of iterations
00526   if (params->isParameter("Maximum Iterations")) {
00527     maxIters_ = params->get("Maximum Iterations",maxIters_default_);
00528 
00529     // Update parameter in our list and in status test.
00530     params_->set("Maximum Iterations", maxIters_);
00531     if (maxIterTest_!=Teuchos::null)
00532       maxIterTest_->setMaxIters( maxIters_ );
00533   }
00534 
00535   // Check for blocksize
00536   if (params->isParameter("Block Size")) {
00537     blockSize_ = params->get("Block Size",blockSize_default_);    
00538     TEST_FOR_EXCEPTION(blockSize_ <= 0, std::invalid_argument,
00539         "Belos::PseudoBlockGmresSolMgr: \"Block Size\" must be strictly positive.");
00540 
00541     // Update parameter in our list.
00542     params_->set("Block Size", blockSize_);
00543   }
00544 
00545   // Check for the maximum number of blocks.
00546   if (params->isParameter("Num Blocks")) {
00547     numBlocks_ = params->get("Num Blocks",numBlocks_default_);
00548     TEST_FOR_EXCEPTION(numBlocks_ <= 0, std::invalid_argument,
00549         "Belos::PseudoBlockGmresSolMgr: \"Num Blocks\" must be strictly positive.");
00550 
00551     // Update parameter in our list.
00552     params_->set("Num Blocks", numBlocks_);
00553   }
00554 
00555   // Check to see if the timer label changed.
00556   if (params->isParameter("Timer Label")) {
00557     std::string tempLabel = params->get("Timer Label", label_default_);
00558 
00559     // Update parameter in our list and solver timer
00560     if (tempLabel != label_) {
00561       label_ = tempLabel;
00562       params_->set("Timer Label", label_);
00563       std::string solveLabel = label_ + ": PseudoBlockGmresSolMgr total solve time";
00564 #ifdef BELOS_TEUCHOS_TIME_MONITOR
00565       timerSolve_ = Teuchos::TimeMonitor::getNewTimer(solveLabel);
00566 #endif
00567     }
00568   }
00569 
00570   // Check if the orthogonalization changed.
00571   if (params->isParameter("Orthogonalization")) {
00572     std::string tempOrthoType = params->get("Orthogonalization",orthoType_default_);
00573     TEST_FOR_EXCEPTION( tempOrthoType != "DGKS" && tempOrthoType != "ICGS" && tempOrthoType != "IMGS", 
00574                         std::invalid_argument,
00575                         "Belos::PseudoBlockGmresSolMgr: \"Orthogonalization\" must be either \"DGKS\",\"ICGS\", or \"IMGS\".");
00576     if (tempOrthoType != orthoType_) {
00577       orthoType_ = tempOrthoType;
00578       // Create orthogonalization manager
00579       if (orthoType_=="DGKS") {
00580         if (orthoKappa_ <= 0) {
00581           ortho_ = Teuchos::rcp( new DGKSOrthoManager<ScalarType,MV,OP>( label_ ) );
00582         }
00583         else {
00584           ortho_ = Teuchos::rcp( new DGKSOrthoManager<ScalarType,MV,OP>( label_ ) );
00585           Teuchos::rcp_dynamic_cast<DGKSOrthoManager<ScalarType,MV,OP> >(ortho_)->setDepTol( orthoKappa_ );
00586         }
00587       }
00588       else if (orthoType_=="ICGS") {
00589         ortho_ = Teuchos::rcp( new ICGSOrthoManager<ScalarType,MV,OP>( label_ ) );
00590       } 
00591       else if (orthoType_=="IMGS") {
00592         ortho_ = Teuchos::rcp( new IMGSOrthoManager<ScalarType,MV,OP>( label_ ) );
00593       } 
00594     }  
00595   }
00596 
00597   // Check which orthogonalization constant to use.
00598   if (params->isParameter("Orthogonalization Constant")) {
00599     orthoKappa_ = params->get("Orthogonalization Constant",orthoKappa_default_);
00600 
00601     // Update parameter in our list.
00602     params_->set("Orthogonalization Constant",orthoKappa_);
00603     if (orthoType_=="DGKS") {
00604       if (orthoKappa_ > 0 && ortho_ != Teuchos::null) {
00605         Teuchos::rcp_dynamic_cast<DGKSOrthoManager<ScalarType,MV,OP> >(ortho_)->setDepTol( orthoKappa_ );
00606       }
00607     } 
00608   }
00609 
00610   // Check for a change in verbosity level
00611   if (params->isParameter("Verbosity")) {
00612     if (Teuchos::isParameterType<int>(*params,"Verbosity")) {
00613       verbosity_ = params->get("Verbosity", verbosity_default_);
00614     } else {
00615       verbosity_ = (int)Teuchos::getParameter<Belos::MsgType>(*params,"Verbosity");
00616     }
00617 
00618     // Update parameter in our list.
00619     params_->set("Verbosity", verbosity_);
00620     if (printer_ != Teuchos::null)
00621       printer_->setVerbosity(verbosity_);
00622   }
00623 
00624   // Check for a change in output style.
00625   if (params->isParameter("Output Style")) {
00626     if (Teuchos::isParameterType<int>(*params,"Output Style")) {
00627       outputStyle_ = params->get("Output Style", outputStyle_default_);
00628     } else {
00629       outputStyle_ = (int)Teuchos::getParameter<Belos::OutputType>(*params,"Output Style");
00630     }
00631 
00632     // Update parameter in our list.
00633     params_->set("Output Style", verbosity_);
00634     if (outputTest_ != Teuchos::null)
00635       isSTSet_ = false;
00636   }
00637 
00638   // output stream
00639   if (params->isParameter("Output Stream")) {
00640     outputStream_ = Teuchos::getParameter<Teuchos::RCP<std::ostream> >(*params,"Output Stream");
00641 
00642     // Update parameter in our list.
00643     params_->set("Output Stream", outputStream_);
00644     if (printer_ != Teuchos::null)
00645       printer_->setOStream( outputStream_ );
00646   }
00647 
00648   // frequency level
00649   if (verbosity_ & Belos::StatusTestDetails) {
00650     if (params->isParameter("Output Frequency")) {
00651       outputFreq_ = params->get("Output Frequency", outputFreq_default_);
00652     }
00653 
00654     // Update parameter in out list and output status test.
00655     params_->set("Output Frequency", outputFreq_);
00656     if (outputTest_ != Teuchos::null)
00657       outputTest_->setOutputFrequency( outputFreq_ );
00658   }
00659 
00660   // Create output manager if we need to.
00661   if (printer_ == Teuchos::null) {
00662     printer_ = Teuchos::rcp( new OutputManager<ScalarType>(verbosity_, outputStream_) );
00663   }
00664   
00665   // Convergence
00666   typedef Belos::StatusTestCombo<ScalarType,MV,OP>  StatusTestCombo_t;
00667   typedef Belos::StatusTestGenResNorm<ScalarType,MV,OP>  StatusTestResNorm_t;
00668 
00669   // Check for convergence tolerance
00670   if (params->isParameter("Convergence Tolerance")) {
00671     convtol_ = params->get("Convergence Tolerance",convtol_default_);
00672 
00673     // Update parameter in our list and residual tests.
00674     params_->set("Convergence Tolerance", convtol_);
00675     if (impConvTest_ != Teuchos::null)
00676       impConvTest_->setTolerance( convtol_ );
00677     if (expConvTest_ != Teuchos::null)
00678       expConvTest_->setTolerance( convtol_ );
00679   }
00680 
00681   // Check for a change in scaling, if so we need to build new residual tests.
00682   bool newImpResTest = false, newExpResTest = false; 
00683   if (params->isParameter("Implicit Residual Scaling")) {
00684     std::string tempImpResScale = Teuchos::getParameter<std::string>( *params, "Implicit Residual Scaling" );
00685 
00686     // Only update the scaling if it's different.
00687     if (impResScale_ != tempImpResScale) {
00688       Belos::ScaleType impResScaleType = convertStringToScaleType( tempImpResScale );
00689       impResScale_ = tempImpResScale;
00690 
00691       // Update parameter in our list and residual tests
00692       params_->set("Implicit Residual Scaling", impResScale_);
00693       if (impConvTest_ != Teuchos::null) {
00694         try {
00695           impConvTest_->defineScaleForm( impResScaleType, Belos::TwoNorm );
00696         }
00697         catch (std::exception& e) {
00698           // Make sure the convergence test gets constructed again.
00699           isSTSet_ = false;
00700           newImpResTest = true;
00701         }
00702       }
00703     }      
00704   }
00705 
00706   if (params->isParameter("Explicit Residual Scaling")) {
00707     std::string tempExpResScale = Teuchos::getParameter<std::string>( *params, "Explicit Residual Scaling" );
00708 
00709     // Only update the scaling if it's different.
00710     if (expResScale_ != tempExpResScale) {
00711       Belos::ScaleType expResScaleType = convertStringToScaleType( tempExpResScale );
00712       expResScale_ = tempExpResScale;
00713 
00714       // Update parameter in our list and residual tests
00715       params_->set("Explicit Residual Scaling", expResScale_);
00716       if (expConvTest_ != Teuchos::null) {
00717         try {
00718           expConvTest_->defineScaleForm( expResScaleType, Belos::TwoNorm );
00719         }
00720         catch (std::exception& e) {
00721           // Make sure the convergence test gets constructed again.
00722           isSTSet_ = false;
00723           newExpResTest = true;
00724         }
00725       }
00726     }      
00727   }
00728 
00729 
00730   if (params->isParameter("Show Maximum Residual Norm Only")) {
00731     showMaxResNormOnly_ = Teuchos::getParameter<bool>(*params,"Show Maximum Residual Norm Only");
00732 
00733     // Update parameter in our list and residual tests
00734     params_->set("Show Maximum Residual Norm Only", showMaxResNormOnly_);
00735     if (impConvTest_ != Teuchos::null)
00736       impConvTest_->setShowMaxResNormOnly( showMaxResNormOnly_ );
00737     if (expConvTest_ != Teuchos::null)
00738       expConvTest_->setShowMaxResNormOnly( showMaxResNormOnly_ );
00739   }
00740 
00741   // Create status tests if we need to.
00742 
00743   // Get the deflation quorum, or number of converged systems before deflation is allowed
00744   if (params->isParameter("Deflation Quorum")) {
00745     defQuorum_ = params->get("Deflation Quorum", defQuorum_);
00746     TEST_FOR_EXCEPTION(defQuorum_ > blockSize_, std::invalid_argument,
00747         "Belos::PseudoBlockGmresSolMgr: \"Deflation Quorum\" cannot be larger than \"Block Size\".");
00748     params_->set("Deflation Quorum", defQuorum_);
00749     if (impConvTest_ != Teuchos::null)
00750       impConvTest_->setQuorum( defQuorum_ );
00751     if (expConvTest_ != Teuchos::null)
00752       expConvTest_->setQuorum( defQuorum_ );
00753   }
00754 
00755   // Create orthogonalization manager if we need to.
00756   if (ortho_ == Teuchos::null) {
00757     if (orthoType_=="DGKS") {
00758       if (orthoKappa_ <= 0) {
00759         ortho_ = Teuchos::rcp( new DGKSOrthoManager<ScalarType,MV,OP>( label_ ) );
00760       }
00761       else {
00762         ortho_ = Teuchos::rcp( new DGKSOrthoManager<ScalarType,MV,OP>( label_ ) );
00763         Teuchos::rcp_dynamic_cast<DGKSOrthoManager<ScalarType,MV,OP> >(ortho_)->setDepTol( orthoKappa_ );
00764       }
00765     }
00766     else if (orthoType_=="ICGS") {
00767       ortho_ = Teuchos::rcp( new ICGSOrthoManager<ScalarType,MV,OP>( label_ ) );
00768     } 
00769     else if (orthoType_=="IMGS") {
00770       ortho_ = Teuchos::rcp( new IMGSOrthoManager<ScalarType,MV,OP>( label_ ) );
00771     } 
00772     else {
00773       TEST_FOR_EXCEPTION(orthoType_!="ICGS"&&orthoType_!="DGKS"&&orthoType_!="IMGS",std::logic_error,
00774           "Belos::PseudoBlockGmresSolMgr(): Invalid orthogonalization type.");
00775     }  
00776   }
00777 
00778   // Create the timer if we need to.
00779   if (timerSolve_ == Teuchos::null) {
00780     std::string solveLabel = label_ + ": PseudoBlockGmresSolMgr total solve time";
00781 #ifdef BELOS_TEUCHOS_TIME_MONITOR
00782     timerSolve_ = Teuchos::TimeMonitor::getNewTimer(solveLabel);
00783 #endif
00784   }
00785 
00786   // Inform the solver manager that the current parameters were set.
00787   isSet_ = true;
00788 }
00789 
00790     
00791 template<class ScalarType, class MV, class OP>
00792 void 
00793 PseudoBlockGmresSolMgr<ScalarType,MV,OP>::setUserConvStatusTest(
00794   const Teuchos::RCP<StatusTest<ScalarType,MV,OP> > &userConvStatusTest
00795   )
00796 {
00797   userConvStatusTest_ = userConvStatusTest;
00798 }
00799 
00800     
00801 template<class ScalarType, class MV, class OP>
00802 Teuchos::RCP<const Teuchos::ParameterList>
00803 PseudoBlockGmresSolMgr<ScalarType,MV,OP>::getValidParameters() const
00804 {
00805   static Teuchos::RCP<const Teuchos::ParameterList> validPL;
00806   if (is_null(validPL)) {
00807     Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
00808   // Set all the valid parameters and their default values.
00809     pl= Teuchos::rcp( new Teuchos::ParameterList() );
00810     pl->set("Convergence Tolerance", convtol_default_,
00811       "The relative residual tolerance that needs to be achieved by the\n"
00812       "iterative solver in order for the linear system to be declared converged.");
00813     pl->set("Maximum Restarts", maxRestarts_default_,
00814       "The maximum number of restarts allowed for each\n"
00815       "set of RHS solved.");
00816     pl->set("Maximum Iterations", maxIters_default_,
00817       "The maximum number of block iterations allowed for each\n"
00818       "set of RHS solved.");
00819     pl->set("Num Blocks", numBlocks_default_,
00820       "The maximum number of vectors allowed in the Krylov subspace\n"
00821       "for each set of RHS solved.");
00822     pl->set("Block Size", blockSize_default_,
00823       "The number of RHS solved simultaneously.");
00824     pl->set("Verbosity", verbosity_default_,
00825       "What type(s) of solver information should be outputted\n"
00826       "to the output stream.");
00827     pl->set("Output Style", outputStyle_default_,
00828       "What style is used for the solver information outputted\n"
00829       "to the output stream.");
00830     pl->set("Output Frequency", outputFreq_default_,
00831       "How often convergence information should be outputted\n"
00832       "to the output stream.");  
00833     pl->set("Deflation Quorum", defQuorum_default_,
00834       "The number of linear systems that need to converge before\n"
00835       "they are deflated.  This number should be <= block size.");
00836     pl->set("Output Stream", outputStream_default_,
00837       "A reference-counted pointer to the output stream where all\n"
00838       "solver output is sent.");
00839     pl->set("Show Maximum Residual Norm Only", showMaxResNormOnly_default_,
00840       "When convergence information is printed, only show the maximum\n"
00841       "relative residual norm when the block size is greater than one.");
00842     pl->set("Implicit Residual Scaling", impResScale_default_,
00843       "The type of scaling used in the implicit residual convergence test.");
00844     pl->set("Explicit Residual Scaling", expResScale_default_,
00845       "The type of scaling used in the explicit residual convergence test.");
00846     pl->set("Timer Label", label_default_,
00847       "The string to use as a prefix for the timer labels.");
00848     //  defaultParams_->set("Restart Timers", restartTimers_);
00849     pl->set("Orthogonalization", orthoType_default_,
00850       "The type of orthogonalization to use: DGKS, ICGS, IMGS.");
00851     pl->set("Orthogonalization Constant",orthoKappa_default_,
00852       "The constant used by DGKS orthogonalization to determine\n"
00853       "whether another step of classical Gram-Schmidt is necessary.");
00854     validPL = pl;
00855   }
00856   return validPL;
00857 }
00858 
00859 // Check the status test versus the defined linear problem
00860 template<class ScalarType, class MV, class OP>
00861 bool PseudoBlockGmresSolMgr<ScalarType,MV,OP>::checkStatusTest() {
00862 
00863   typedef Belos::StatusTestCombo<ScalarType,MV,OP>  StatusTestCombo_t;
00864   typedef Belos::StatusTestGenResNorm<ScalarType,MV,OP>  StatusTestGenResNorm_t;
00865   typedef Belos::StatusTestImpResNorm<ScalarType,MV,OP>  StatusTestImpResNorm_t;
00866 
00867   // Basic test checks maximum iterations and native residual.
00868   maxIterTest_ = Teuchos::rcp( new StatusTestMaxIters<ScalarType,MV,OP>( maxIters_ ) );
00869 
00870   // If there is a left preconditioner, we create a combined status test that checks the implicit
00871   // and then explicit residual norm to see if we have convergence.
00872   if ( !Teuchos::is_null(problem_->getLeftPrec()) ) {
00873     expResTest_ = true;
00874   }
00875 
00876   if (expResTest_) {
00877    
00878     // Implicit residual test, using the native residual to determine if convergence was achieved.
00879     Teuchos::RCP<StatusTestGenResNorm_t> tmpImpConvTest =
00880       Teuchos::rcp( new StatusTestGenResNorm_t( convtol_, defQuorum_ ) );
00881     tmpImpConvTest->defineScaleForm( convertStringToScaleType(impResScale_), Belos::TwoNorm );
00882     tmpImpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
00883     impConvTest_ = tmpImpConvTest;
00884 
00885     // Explicit residual test once the native residual is below the tolerance
00886     Teuchos::RCP<StatusTestGenResNorm_t> tmpExpConvTest =
00887       Teuchos::rcp( new StatusTestGenResNorm_t( convtol_, defQuorum_ ) );
00888     tmpExpConvTest->defineResForm( StatusTestGenResNorm_t::Explicit, Belos::TwoNorm );
00889     tmpExpConvTest->defineScaleForm( convertStringToScaleType(expResScale_), Belos::TwoNorm );
00890     tmpExpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
00891     expConvTest_ = tmpExpConvTest;
00892 
00893     // The convergence test is a combination of the "cheap" implicit test and explicit test.
00894     convTest_ = Teuchos::rcp( new StatusTestCombo_t( StatusTestCombo_t::SEQ, impConvTest_, expConvTest_ ) );
00895   }
00896   else {
00897 
00898     // Implicit residual test, using the native residual to determine if convergence was achieved.
00899     // Use test that checks for loss of accuracy.
00900     Teuchos::RCP<StatusTestImpResNorm_t> tmpImpConvTest =
00901       Teuchos::rcp( new StatusTestImpResNorm_t( convtol_, defQuorum_ ) );
00902     tmpImpConvTest->defineScaleForm( convertStringToScaleType(impResScale_), Belos::TwoNorm );
00903     tmpImpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
00904     impConvTest_ = tmpImpConvTest;
00905     
00906     // Set the explicit and total convergence test to this implicit test that checks for accuracy loss.
00907     expConvTest_ = impConvTest_;
00908     convTest_ = impConvTest_;
00909   }
00910 
00911   if (nonnull(userConvStatusTest_) ) {
00912     // Override the overall convergence test with the users convergence test
00913     convTest_ = Teuchos::rcp(
00914       new StatusTestCombo_t( StatusTestCombo_t::SEQ, convTest_, userConvStatusTest_ ) );
00915     // NOTE: Above, you have to run the other convergence tests also because
00916     // the logic in this class depends on it.  This is very unfortunate.
00917   }
00918 
00919   sTest_ = Teuchos::rcp( new StatusTestCombo_t( StatusTestCombo_t::OR, maxIterTest_, convTest_ ) );
00920 
00921   // Create the status test output class.
00922   // This class manages and formats the output from the status test.
00923   StatusTestOutputFactory<ScalarType,MV,OP> stoFactory( outputStyle_ );
00924   outputTest_ = stoFactory.create( printer_, sTest_, outputFreq_, Passed+Failed+Undefined );
00925   
00926   // Set the solver string for the output test
00927   std::string solverDesc = " Pseudo Block Gmres ";
00928   outputTest_->setSolverDesc( solverDesc );
00929 
00930 
00931   // The status test is now set.
00932   isSTSet_ = true;
00933 
00934   return false;
00935 }
00936 
00937 
00938 // solve()
00939 template<class ScalarType, class MV, class OP>
00940 ReturnType PseudoBlockGmresSolMgr<ScalarType,MV,OP>::solve() {
00941 
00942   // Set the current parameters if they were not set before.
00943   // NOTE:  This may occur if the user generated the solver manager with the default constructor and 
00944   // then didn't set any parameters using setParameters().
00945   if (!isSet_) { setParameters( params_ ); }
00946   
00947   Teuchos::BLAS<int,ScalarType> blas;
00948   
00949   TEST_FOR_EXCEPTION(!problem_->isProblemSet(),PseudoBlockGmresSolMgrLinearProblemFailure,
00950                      "Belos::PseudoBlockGmresSolMgr::solve(): Linear problem is not ready, setProblem() has not been called.");
00951 
00952   // Check if we have to create the status tests.
00953   if (!isSTSet_ || (!expResTest_ && !Teuchos::is_null(problem_->getLeftPrec())) ) {
00954     TEST_FOR_EXCEPTION( checkStatusTest(),PseudoBlockGmresSolMgrLinearProblemFailure,
00955       "Belos::BlockGmresSolMgr::solve(): Linear problem and requested status tests are incompatible.");
00956   }
00957 
00958   // Create indices for the linear systems to be solved.
00959   int startPtr = 0;
00960   int numRHS2Solve = MVT::GetNumberVecs( *(problem_->getRHS()) );
00961   int numCurrRHS = ( numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;
00962 
00963   std::vector<int> currIdx( numCurrRHS );
00964   blockSize_ = numCurrRHS;
00965   for (int i=0; i<numCurrRHS; ++i) 
00966     { currIdx[i] = startPtr+i; }
00967 
00968   // Inform the linear problem of the current linear system to solve.
00969   problem_->setLSIndex( currIdx );
00970 
00972   // Parameter list
00973   Teuchos::ParameterList plist;
00974   plist.set("Num Blocks",numBlocks_);
00975   
00976   // Reset the status test.  
00977   outputTest_->reset();
00978   loaDetected_ = false;
00979 
00980   // Assume convergence is achieved, then let any failed convergence set this to false.
00981   bool isConverged = true;
00982 
00984   // BlockGmres solver
00985 
00986   Teuchos::RCP<PseudoBlockGmresIter<ScalarType,MV,OP> > block_gmres_iter
00987     = Teuchos::rcp( new PseudoBlockGmresIter<ScalarType,MV,OP>(problem_,printer_,outputTest_,ortho_,plist) );  
00988 
00989   // Enter solve() iterations
00990   {
00991 #ifdef BELOS_TEUCHOS_TIME_MONITOR
00992     Teuchos::TimeMonitor slvtimer(*timerSolve_);
00993 #endif
00994 
00995     while ( numRHS2Solve > 0 ) {
00996 
00997       // Reset the active / converged vectors from this block
00998       std::vector<int> convRHSIdx;
00999       std::vector<int> currRHSIdx( currIdx );
01000       currRHSIdx.resize(numCurrRHS);
01001 
01002       // Set the current number of blocks with the pseudo Gmres iteration.
01003       block_gmres_iter->setNumBlocks( numBlocks_ );
01004 
01005       // Reset the number of iterations.
01006       block_gmres_iter->resetNumIters();
01007 
01008       // Reset the number of calls that the status test output knows about.
01009       outputTest_->resetNumCalls();
01010 
01011       // Get a new state struct and initialize the solver.
01012       PseudoBlockGmresIterState<ScalarType,MV> newState;
01013 
01014       // Create the first block in the current Krylov basis for each right-hand side.
01015       std::vector<int> index(1);
01016       Teuchos::RCP<MV> tmpV, R_0 = MVT::CloneCopy( *(problem_->getInitPrecResVec()), currIdx );
01017       newState.V.resize( blockSize_ );
01018       newState.Z.resize( blockSize_ );
01019       for (int i=0; i<blockSize_; ++i) {
01020         index[0]=i;
01021         tmpV = MVT::CloneViewNonConst( *R_0, index );
01022 
01023         // Get a matrix to hold the orthonormalization coefficients.
01024         Teuchos::RCP<Teuchos::SerialDenseVector<int,ScalarType> > tmpZ
01025           = Teuchos::rcp( new Teuchos::SerialDenseVector<int,ScalarType>( 1 ));
01026 
01027         // Orthonormalize the new V_0
01028         int rank = ortho_->normalize( *tmpV, tmpZ );
01029         TEST_FOR_EXCEPTION(rank != 1, PseudoBlockGmresSolMgrOrthoFailure,
01030             "Belos::PseudoBlockGmresSolMgr::solve(): Failed to compute initial block of orthonormal vectors.");
01031 
01032         newState.V[i] = tmpV;
01033         newState.Z[i] = tmpZ;
01034       }
01035 
01036       newState.curDim = 0;
01037       block_gmres_iter->initialize(newState);
01038       int numRestarts = 0;
01039 
01040       while(1) {
01041         
01042         // tell block_gmres_iter to iterate
01043         try {
01044           block_gmres_iter->iterate();
01045           
01047           //
01048           // check convergence first
01049           //
01051           if ( convTest_->getStatus() == Passed ) {
01052         
01053             if ( expConvTest_->getLOADetected() ) {
01054         // Oops!  There was a loss of accuracy (LOA) for one or
01055         // more right-hand sides.  That means the implicit
01056         // (a.k.a. "native") residuals claim convergence,
01057         // whereas the explicit (hence expConvTest_, i.e.,
01058         // "explicit convergence test") residuals have not
01059         // converged.
01060         //
01061         // We choose to deal with this situation by deflating
01062         // out the affected right-hand sides and moving on.
01063               loaDetected_ = true;
01064               printer_->stream(Warnings) <<
01065                 "Belos::PseudoBlockGmresSolMgr::solve(): Warning! Solver has experienced a loss of accuracy!" << std::endl;
01066               isConverged = false;
01067             }
01068         
01069             // Figure out which linear systems converged.
01070             std::vector<int> convIdx = expConvTest_->convIndices();
01071         
01072             // If the number of converged linear systems is equal to the
01073             // number of current linear systems, then we are done with this block.
01074             if (convIdx.size() == currRHSIdx.size())
01075               break;  // break from while(1){block_gmres_iter->iterate()}
01076         
01077             // Get a new state struct and initialize the solver.
01078             PseudoBlockGmresIterState<ScalarType,MV> defState;
01079         
01080             // Inform the linear problem that we are finished with this current linear system.
01081             problem_->setCurrLS();
01082         
01083             // Get the state.
01084             PseudoBlockGmresIterState<ScalarType,MV> oldState = block_gmres_iter->getState();
01085             int curDim = oldState.curDim;
01086         
01087             // Get a new state struct and reset currRHSIdx to have the right-hand sides that 
01088             // are left to converge for this block.
01089             int have = 0;
01090             std::vector<int> oldRHSIdx( currRHSIdx );
01091             std::vector<int> defRHSIdx;
01092             for (unsigned int i=0; i<currRHSIdx.size(); ++i) {
01093               bool found = false;
01094               for (unsigned int j=0; j<convIdx.size(); ++j) {
01095                 if (currRHSIdx[i] == convIdx[j]) {
01096                   found = true;
01097                   break;
01098                 }
01099               }
01100               if (found) {
01101                 defRHSIdx.push_back( i );
01102               }
01103               else {
01104                 defState.V.push_back( Teuchos::rcp_const_cast<MV>( oldState.V[i] ) );
01105                 defState.Z.push_back( Teuchos::rcp_const_cast<Teuchos::SerialDenseVector<int,ScalarType> >( oldState.Z[i] ) );
01106                 defState.H.push_back( Teuchos::rcp_const_cast<Teuchos::SerialDenseMatrix<int,ScalarType> >( oldState.H[i] ) );
01107                 defState.sn.push_back( Teuchos::rcp_const_cast<Teuchos::SerialDenseVector<int,ScalarType> >( oldState.sn[i] ) );
01108                 defState.cs.push_back( Teuchos::rcp_const_cast<Teuchos::SerialDenseVector<int,MagnitudeType> >(oldState.cs[i] ) );
01109                 currRHSIdx[have] = currRHSIdx[i];
01110                 have++;
01111               }
01112             }
01113             defRHSIdx.resize(currRHSIdx.size()-have);
01114             currRHSIdx.resize(have);
01115         
01116             // Compute the current solution that needs to be deflated if this solver has taken any steps.
01117             if (curDim) {
01118               Teuchos::RCP<MV> update = block_gmres_iter->getCurrentUpdate();
01119               Teuchos::RCP<MV> defUpdate = MVT::CloneViewNonConst( *update, defRHSIdx );
01120               
01121               // Set the deflated indices so we can update the solution.
01122               problem_->setLSIndex( convIdx );
01123               
01124               // Update the linear problem.
01125               problem_->updateSolution( defUpdate, true );
01126             }
01127             
01128             // Set the remaining indices after deflation.
01129             problem_->setLSIndex( currRHSIdx );
01130             
01131             // Set the dimension of the subspace, which is the same as the old subspace size.
01132             defState.curDim = curDim;
01133             
01134             // Initialize the solver with the deflated system.
01135             block_gmres_iter->initialize(defState);
01136           }
01138           //
01139           // check for maximum iterations
01140           //
01142           else if ( maxIterTest_->getStatus() == Passed ) {
01143             // we don't have convergence
01144             isConverged = false;
01145             break;  // break from while(1){block_gmres_iter->iterate()}
01146           }
01148           //
01149           // check for restarting, i.e. the subspace is full
01150           //
01152           else if ( block_gmres_iter->getCurSubspaceDim() == block_gmres_iter->getMaxSubspaceDim() ) {
01153         
01154             if ( numRestarts >= maxRestarts_ ) {
01155               isConverged = false;
01156               break; // break from while(1){block_gmres_iter->iterate()}
01157             }
01158             numRestarts++;
01159         
01160             printer_->stream(Debug) << " Performing restart number " << numRestarts << " of " << maxRestarts_ << std::endl << std::endl;
01161             
01162             // Update the linear problem.
01163             Teuchos::RCP<MV> update = block_gmres_iter->getCurrentUpdate();
01164             problem_->updateSolution( update, true );
01165             
01166             // Get the state.
01167             PseudoBlockGmresIterState<ScalarType,MV> oldState = block_gmres_iter->getState();
01168             
01169             // Set the new state.
01170             PseudoBlockGmresIterState<ScalarType,MV> newstate;
01171             newstate.V.resize(currRHSIdx.size());
01172             newstate.Z.resize(currRHSIdx.size());
01173         
01174             // Compute the restart vectors
01175             // NOTE: Force the linear problem to update the current residual since the solution was updated.
01176             R_0 = MVT::Clone( *(problem_->getInitPrecResVec()), currRHSIdx.size() );
01177             problem_->computeCurrPrecResVec( &*R_0 );
01178             for (unsigned int i=0; i<currRHSIdx.size(); ++i) {
01179               index[0] = i;  // index(1) vector declared on line 891
01180         
01181               tmpV = MVT::CloneViewNonConst( *R_0, index );
01182         
01183               // Get a matrix to hold the orthonormalization coefficients.
01184               Teuchos::RCP<Teuchos::SerialDenseVector<int,ScalarType> > tmpZ
01185                 = Teuchos::rcp( new Teuchos::SerialDenseVector<int,ScalarType>( 1 ));
01186               
01187               // Orthonormalize the new V_0
01188               int rank = ortho_->normalize( *tmpV, tmpZ );
01189               TEST_FOR_EXCEPTION(rank != 1 ,PseudoBlockGmresSolMgrOrthoFailure,
01190                   "Belos::PseudoBlockGmresSolMgr::solve(): Failed to compute initial block of orthonormal vectors after the restart.");
01191               
01192               newstate.V[i] = tmpV;
01193               newstate.Z[i] = tmpZ;
01194             }
01195         
01196             // Initialize the solver.
01197             newstate.curDim = 0;
01198             block_gmres_iter->initialize(newstate);
01199         
01200           } // end of restarting
01201         
01203           //
01204           // we returned from iterate(), but none of our status tests Passed.
01205           // something is wrong, and it is probably our fault.
01206           //
01208         
01209           else {
01210             TEST_FOR_EXCEPTION(true,std::logic_error,
01211                 "Belos::PseudoBlockGmresSolMgr::solve(): Invalid return from PseudoBlockGmresIter::iterate().");
01212           }
01213         }
01214         catch (const PseudoBlockGmresIterOrthoFailure &e) {
01215         
01216           // Try to recover the most recent least-squares solution
01217           block_gmres_iter->updateLSQR( block_gmres_iter->getCurSubspaceDim() );
01218         
01219           // Check to see if the most recent least-squares solution yielded convergence.
01220           sTest_->checkStatus( &*block_gmres_iter );
01221           if (convTest_->getStatus() != Passed)
01222             isConverged = false;
01223           break;
01224         }
01225         catch (const std::exception &e) {
01226           printer_->stream(Errors) << "Error! Caught std::exception in PseudoBlockGmresIter::iterate() at iteration " 
01227                                    << block_gmres_iter->getNumIters() << std::endl 
01228                                    << e.what() << std::endl;
01229           throw;
01230         }
01231       }
01232       
01233       // Compute the current solution.
01234       // Update the linear problem.
01235       if (nonnull(userConvStatusTest_)) {
01236         //std::cout << "\nnonnull(userConvStatusTest_)\n";
01237         Teuchos::RCP<MV> update = block_gmres_iter->getCurrentUpdate();
01238         problem_->updateSolution( update, true );
01239       }
01240       else if (nonnull(expConvTest_->getSolution())) {
01241         //std::cout << "\nexpConvTest_->getSolution()\n";
01242         Teuchos::RCP<MV> newX = expConvTest_->getSolution();
01243         Teuchos::RCP<MV> curX = problem_->getCurrLHSVec();
01244         MVT::MvAddMv( 0.0, *newX, 1.0, *newX, *curX );
01245       }
01246       else {
01247         //std::cout << "\nblock_gmres_iter->getCurrentUpdate()\n";
01248         Teuchos::RCP<MV> update = block_gmres_iter->getCurrentUpdate();
01249         problem_->updateSolution( update, true );
01250       }
01251 
01252       // Inform the linear problem that we are finished with this block linear system.
01253       problem_->setCurrLS();
01254       
01255       // Update indices for the linear systems to be solved.
01256       startPtr += numCurrRHS;
01257       numRHS2Solve -= numCurrRHS;
01258       if ( numRHS2Solve > 0 ) {
01259         numCurrRHS = ( numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;
01260 
01261         blockSize_ = numCurrRHS;
01262         currIdx.resize( numCurrRHS  );
01263         for (int i=0; i<numCurrRHS; ++i) 
01264         { currIdx[i] = startPtr+i; }
01265 
01266         // Adapt the status test quorum if we need to.
01267         if (defQuorum_ > blockSize_) {
01268           if (impConvTest_ != Teuchos::null)
01269             impConvTest_->setQuorum( blockSize_ );
01270           if (expConvTest_ != Teuchos::null)
01271             expConvTest_->setQuorum( blockSize_ );
01272         }
01273 
01274         // Set the next indices.
01275         problem_->setLSIndex( currIdx );
01276       }
01277       else {
01278         currIdx.resize( numRHS2Solve );
01279       }
01280 
01281     }// while ( numRHS2Solve > 0 )
01282 
01283   }
01284 
01285   // print final summary
01286   sTest_->print( printer_->stream(FinalSummary) );
01287  
01288   // print timing information
01289 #ifdef BELOS_TEUCHOS_TIME_MONITOR
01290   // Calling summarize() can be expensive, so don't call unless the
01291   // user wants to print out timing details.  summarize() will do all
01292   // the work even if it's passed a "black hole" output stream.
01293   if (verbosity_ & TimingDetails)
01294     Teuchos::TimeMonitor::summarize( printer_->stream(TimingDetails) );
01295 #endif
01296  
01297   // get iteration information for this solve
01298   numIters_ = maxIterTest_->getNumIters();
01299  
01300   if (!isConverged || loaDetected_) {
01301     return Unconverged; // return from PseudoBlockGmresSolMgr::solve() 
01302   }
01303   return Converged; // return from PseudoBlockGmresSolMgr::solve() 
01304 }
01305 
01306 //  This method requires the solver manager to return a std::string that describes itself.
01307 template<class ScalarType, class MV, class OP>
01308 std::string PseudoBlockGmresSolMgr<ScalarType,MV,OP>::description() const
01309 {
01310   std::ostringstream oss;
01311   oss << "Belos::PseudoBlockGmresSolMgr<...,"<<Teuchos::ScalarTraits<ScalarType>::name()<<">";
01312   oss << "{";
01313   oss << "Ortho Type='"<<orthoType_<<"\', Block Size=" <<blockSize_;
01314   oss << ", Num Blocks="<<numBlocks_<< ", Max Restarts=" << maxRestarts_;
01315   oss << "}";
01316   return oss.str();
01317 }
01318   
01319 } // end Belos namespace
01320 
01321 #endif /* BELOS_PSEUDO_BLOCK_GMRES_SOLMGR_HPP */
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines