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

Generated on Tue Oct 20 12:48:34 2009 for Belos by doxygen 1.4.7