Belos Package Browser (Single Doxygen Collection) Development
BelosRCGSolMgr.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_RCG_SOLMGR_HPP
00043 #define BELOS_RCG_SOLMGR_HPP
00044 
00049 #include "BelosConfigDefs.hpp"
00050 #include "BelosTypes.hpp"
00051 
00052 #include "BelosLinearProblem.hpp"
00053 #include "BelosSolverManager.hpp"
00054 
00055 #include "BelosRCGIter.hpp"
00056 #include "BelosDGKSOrthoManager.hpp"
00057 #include "BelosICGSOrthoManager.hpp"
00058 #include "BelosIMGSOrthoManager.hpp"
00059 #include "BelosStatusTestMaxIters.hpp"
00060 #include "BelosStatusTestGenResNorm.hpp"
00061 #include "BelosStatusTestCombo.hpp"
00062 #include "BelosStatusTestOutputFactory.hpp"
00063 #include "BelosOutputManager.hpp"
00064 #include "Teuchos_BLAS.hpp"
00065 #include "Teuchos_LAPACK.hpp"
00066 #include "Teuchos_as.hpp"
00067 #ifdef BELOS_TEUCHOS_TIME_MONITOR
00068 #include "Teuchos_TimeMonitor.hpp"
00069 #endif
00070 
00110 namespace Belos {
00111 
00113 
00114 
00121   class RCGSolMgrLinearProblemFailure : public BelosError {public:
00122     RCGSolMgrLinearProblemFailure(const std::string& what_arg) : BelosError(what_arg)
00123     {}};
00124 
00131   class RCGSolMgrLAPACKFailure : public BelosError {public:
00132     RCGSolMgrLAPACKFailure(const std::string& what_arg) : BelosError(what_arg)
00133     {}};
00134 
00141   class RCGSolMgrRecyclingFailure : public BelosError {public:
00142     RCGSolMgrRecyclingFailure(const std::string& what_arg) : BelosError(what_arg)
00143     {}};
00144 
00146 
00147 
00148   // Partial specialization for complex ScalarType.
00149   // This contains a trivial implementation.
00150   // See discussion in the class documentation above.
00151   template<class ScalarType, class MV, class OP,
00152            const bool scalarTypeIsComplex = Teuchos::ScalarTraits<ScalarType>::isComplex>
00153   class RCGSolMgr :
00154     public Details::RealSolverManager<ScalarType, MV, OP,
00155                                       Teuchos::ScalarTraits<ScalarType>::isComplex>
00156   {
00157     static const bool isComplex = Teuchos::ScalarTraits<ScalarType>::isComplex;
00158     typedef Details::RealSolverManager<ScalarType, MV, OP, isComplex> base_type;
00159 
00160   public:
00161     RCGSolMgr () :
00162       base_type ()
00163     {}
00164     RCGSolMgr (const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem,
00165                const Teuchos::RCP<Teuchos::ParameterList> &pl) :
00166       base_type ()
00167     {}
00168     virtual ~RCGSolMgr () {}
00169   };
00170 
00171 
00172   // Partial specialization for real ScalarType.
00173   // This contains the actual working implementation of RCG.
00174   // See discussion in the class documentation above.
00175   template<class ScalarType, class MV, class OP>
00176   class RCGSolMgr<ScalarType, MV, OP, false> :
00177     public Details::RealSolverManager<ScalarType, MV, OP, false> {
00178   private:
00179     typedef MultiVecTraits<ScalarType,MV> MVT;
00180     typedef MultiVecTraitsExt<ScalarType,MV> MVText;
00181     typedef OperatorTraits<ScalarType,MV,OP> OPT;
00182     typedef Teuchos::ScalarTraits<ScalarType> SCT;
00183     typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType MagnitudeType;
00184     typedef Teuchos::ScalarTraits<MagnitudeType> MT;
00185 
00186   public:
00187 
00189 
00190 
00196      RCGSolMgr();
00197 
00219     RCGSolMgr( const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem,
00220                    const Teuchos::RCP<Teuchos::ParameterList> &pl );
00221 
00223     virtual ~RCGSolMgr() {};
00225 
00227 
00228 
00229     const LinearProblem<ScalarType,MV,OP>& getProblem() const {
00230       return *problem_;
00231     }
00232 
00234     Teuchos::RCP<const Teuchos::ParameterList> getValidParameters() const;
00235 
00237     Teuchos::RCP<const Teuchos::ParameterList> getCurrentParameters() const { return params_; }
00238 
00244     Teuchos::Array<Teuchos::RCP<Teuchos::Time> > getTimers() const {
00245       return Teuchos::tuple(timerSolve_);
00246     }
00247 
00252     MagnitudeType achievedTol() const {
00253       return achievedTol_;
00254     }
00255 
00257     int getNumIters() const {
00258       return numIters_;
00259     }
00260 
00262     bool isLOADetected() const { return false; }
00263 
00265 
00267 
00268 
00270     void setProblem( const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem ) { problem_ = problem; }
00271 
00273     void setParameters( const Teuchos::RCP<Teuchos::ParameterList> &params );
00274 
00276 
00278 
00279 
00285     void reset( const ResetType type ) {
00286       if ((type & Belos::Problem) && !Teuchos::is_null(problem_)) problem_->setProblem();
00287       else if (type & Belos::RecycleSubspace) existU_ = false;
00288     }
00290 
00292 
00293 
00311     ReturnType solve();
00312 
00314 
00317 
00319     std::string description() const;
00320 
00322 
00323   private:
00324 
00325     // Called by all constructors; Contains init instructions common to all constructors
00326     void init();
00327 
00328     //  Computes harmonic eigenpairs of projected matrix created during one cycle.
00329     //  Y contains the harmonic Ritz vectors corresponding to the recycleBlocks eigenvalues of smallest magnitude.
00330     void getHarmonicVecs(const Teuchos::SerialDenseMatrix<int,ScalarType> &F,
00331                          const Teuchos::SerialDenseMatrix<int,ScalarType> &G,
00332                          Teuchos::SerialDenseMatrix<int,ScalarType>& Y);
00333 
00334     // Sort list of n floating-point numbers and return permutation vector
00335     void sort(std::vector<ScalarType>& dlist, int n, std::vector<int>& iperm);
00336 
00337     // Initialize solver state storage
00338     void initializeStateStorage();
00339 
00340     // Linear problem.
00341     Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > problem_;
00342 
00343     // Output manager.
00344     Teuchos::RCP<OutputManager<ScalarType> > printer_;
00345     Teuchos::RCP<std::ostream> outputStream_;
00346 
00347     // Status test.
00348     Teuchos::RCP<StatusTest<ScalarType,MV,OP> > sTest_;
00349     Teuchos::RCP<StatusTestMaxIters<ScalarType,MV,OP> > maxIterTest_;
00350     Teuchos::RCP<StatusTestGenResNorm<ScalarType,MV,OP> > convTest_;
00351     Teuchos::RCP<StatusTestOutput<ScalarType,MV,OP> > outputTest_;
00352 
00353     // Current parameter list.
00354     Teuchos::RCP<Teuchos::ParameterList> params_;
00355 
00356     // Default solver values.
00357     static const MagnitudeType convtol_default_;
00358     static const int maxIters_default_;
00359     static const int blockSize_default_;
00360     static const int numBlocks_default_;
00361     static const int recycleBlocks_default_;
00362     static const bool showMaxResNormOnly_default_;
00363     static const int verbosity_default_;
00364     static const int outputStyle_default_;
00365     static const int outputFreq_default_;
00366     static const std::string label_default_;
00367     static const Teuchos::RCP<std::ostream> outputStream_default_;
00368 
00369     //
00370     // Current solver values.
00371     //
00372 
00374     MagnitudeType convtol_;
00375 
00380     MagnitudeType achievedTol_;
00381 
00383     int maxIters_;
00384 
00386     int numIters_;
00387 
00388     int numBlocks_, recycleBlocks_;
00389     bool showMaxResNormOnly_;
00390     int verbosity_, outputStyle_, outputFreq_;
00391 
00393     // Solver State Storage
00395     // Search vectors
00396     Teuchos::RCP<MV> P_;
00397     //
00398     // A times current search direction
00399     Teuchos::RCP<MV> Ap_;
00400     //
00401     // Residual vector
00402     Teuchos::RCP<MV> r_;
00403     //
00404     // Preconditioned residual
00405     Teuchos::RCP<MV> z_;
00406     //
00407     // Flag indicating that the recycle space should be used
00408     bool existU_;
00409     //
00410     // Flag indicating that the updated recycle space has been created
00411     bool existU1_;
00412     //
00413     // Recycled subspace and its image
00414     Teuchos::RCP<MV> U_, AU_;
00415     //
00416     // Recycled subspace for next system and its image
00417     Teuchos::RCP<MV> U1_;
00418     //
00419     // Coefficients arising in RCG iteration
00420     Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > Alpha_;
00421     Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > Beta_;
00422     Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > D_;
00423     //
00424     // Solutions to local least-squares problems
00425     Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > Delta_;
00426     //
00427     // The matrix U^T A U
00428     Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > UTAU_;
00429     //
00430     // LU factorization of U^T A U
00431     Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > LUUTAU_;
00432     //
00433     // Data from LU factorization of UTAU
00434     Teuchos::RCP<std::vector<int> > ipiv_;
00435     //
00436     // The matrix (AU)^T AU
00437     Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > AUTAU_;
00438     //
00439     // The scalar r'*z
00440     Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > rTz_old_;
00441     //
00442     // Matrices needed for calculation of harmonic Ritz eigenproblem
00443     Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > F_,G_,Y_;
00444     //
00445     // Matrices needed for updating recycle space
00446     Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > L2_,DeltaL2_,AU1TUDeltaL2_;
00447     Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > AU1TAU1_, AU1TU1_, AU1TAP_;
00448     Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > FY_,GY_;
00449     Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > APTAP_;
00450     Teuchos::RCP<MV> U1Y1_, PY2_;
00451     Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > AUTAP_, AU1TU_;
00452     ScalarType dold;
00454 
00455     // Timers.
00456     std::string label_;
00457     Teuchos::RCP<Teuchos::Time> timerSolve_;
00458 
00459     // Internal state variables.
00460     bool params_Set_;
00461   };
00462 
00463 
00464 // Default solver values.
00465 template<class ScalarType, class MV, class OP>
00466 const typename RCGSolMgr<ScalarType,MV,OP,false>::MagnitudeType
00467 RCGSolMgr<ScalarType,MV,OP,false>::convtol_default_ = 1e-8;
00468 
00469 template<class ScalarType, class MV, class OP>
00470 const int RCGSolMgr<ScalarType,MV,OP,false>::maxIters_default_ = 1000;
00471 
00472 template<class ScalarType, class MV, class OP>
00473 const int RCGSolMgr<ScalarType,MV,OP,false>::numBlocks_default_ = 25;
00474 
00475 template<class ScalarType, class MV, class OP>
00476 const int RCGSolMgr<ScalarType,MV,OP,false>::blockSize_default_ = 1;
00477 
00478 template<class ScalarType, class MV, class OP>
00479 const int RCGSolMgr<ScalarType,MV,OP,false>::recycleBlocks_default_ = 3;
00480 
00481 template<class ScalarType, class MV, class OP>
00482 const bool RCGSolMgr<ScalarType,MV,OP,false>::showMaxResNormOnly_default_ = false;
00483 
00484 template<class ScalarType, class MV, class OP>
00485 const int RCGSolMgr<ScalarType,MV,OP,false>::verbosity_default_ = Belos::Errors;
00486 
00487 template<class ScalarType, class MV, class OP>
00488 const int RCGSolMgr<ScalarType,MV,OP,false>::outputStyle_default_ = Belos::General;
00489 
00490 template<class ScalarType, class MV, class OP>
00491 const int RCGSolMgr<ScalarType,MV,OP,false>::outputFreq_default_ = -1;
00492 
00493 template<class ScalarType, class MV, class OP>
00494 const std::string RCGSolMgr<ScalarType,MV,OP,false>::label_default_ = "Belos";
00495 
00496 template<class ScalarType, class MV, class OP>
00497 const Teuchos::RCP<std::ostream> RCGSolMgr<ScalarType,MV,OP,false>::outputStream_default_ = Teuchos::rcp(&std::cout,false);
00498 
00499 // Empty Constructor
00500 template<class ScalarType, class MV, class OP>
00501 RCGSolMgr<ScalarType,MV,OP,false>::RCGSolMgr() {
00502   init();
00503 }
00504 
00505 // Basic Constructor
00506 template<class ScalarType, class MV, class OP>
00507 RCGSolMgr<ScalarType,MV,OP,false>::RCGSolMgr(
00508                                                      const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem,
00509                                                      const Teuchos::RCP<Teuchos::ParameterList> &pl ) :
00510   problem_(problem)
00511 {
00512   init();
00513   TEUCHOS_TEST_FOR_EXCEPTION(problem_ == Teuchos::null, std::invalid_argument, "Problem not given to solver manager.");
00514 
00515   // If the parameter list pointer is null, then set the current parameters to the default parameter list.
00516   if ( !is_null(pl) ) {
00517     setParameters( pl );
00518   }
00519 }
00520 
00521 // Common instructions executed in all constructors
00522 template<class ScalarType, class MV, class OP>
00523 void RCGSolMgr<ScalarType,MV,OP,false>::init()
00524 {
00525   outputStream_ = outputStream_default_;
00526   convtol_ = convtol_default_;
00527   maxIters_ = maxIters_default_;
00528   numBlocks_ = numBlocks_default_;
00529   recycleBlocks_ = recycleBlocks_default_;
00530   verbosity_ = verbosity_default_;
00531   outputStyle_= outputStyle_default_;
00532   outputFreq_= outputFreq_default_;
00533   showMaxResNormOnly_ = showMaxResNormOnly_default_;
00534   label_ = label_default_;
00535   params_Set_ = false;
00536   P_ = Teuchos::null;
00537   Ap_ = Teuchos::null;
00538   r_ = Teuchos::null;
00539   z_ = Teuchos::null;
00540   existU_ = false;
00541   existU1_ = false;
00542   U_ = Teuchos::null;
00543   AU_ = Teuchos::null;
00544   U1_ = Teuchos::null;
00545   Alpha_ = Teuchos::null;
00546   Beta_ = Teuchos::null;
00547   D_ = Teuchos::null;
00548   Delta_ = Teuchos::null;
00549   UTAU_ = Teuchos::null;
00550   LUUTAU_ = Teuchos::null;
00551   ipiv_ = Teuchos::null;
00552   AUTAU_ = Teuchos::null;
00553   rTz_old_ = Teuchos::null;
00554   F_ = Teuchos::null;
00555   G_ = Teuchos::null;
00556   Y_ = Teuchos::null;
00557   L2_ = Teuchos::null;
00558   DeltaL2_ = Teuchos::null;
00559   AU1TUDeltaL2_ = Teuchos::null;
00560   AU1TAU1_ = Teuchos::null;
00561   AU1TU1_ = Teuchos::null;
00562   AU1TAP_ = Teuchos::null;
00563   FY_ = Teuchos::null;
00564   GY_ = Teuchos::null;
00565   APTAP_ = Teuchos::null;
00566   U1Y1_ = Teuchos::null;
00567   PY2_ = Teuchos::null;
00568   AUTAP_ = Teuchos::null;
00569   AU1TU_ = Teuchos::null;
00570   dold = 0.;
00571 }
00572 
00573 template<class ScalarType, class MV, class OP>
00574 void RCGSolMgr<ScalarType,MV,OP,false>::setParameters( const Teuchos::RCP<Teuchos::ParameterList> &params )
00575 {
00576   // Create the internal parameter list if ones doesn't already exist.
00577   if (params_ == Teuchos::null) {
00578     params_ = Teuchos::rcp( new Teuchos::ParameterList(*getValidParameters()) );
00579   }
00580   else {
00581     params->validateParameters(*getValidParameters());
00582   }
00583 
00584   // Check for maximum number of iterations
00585   if (params->isParameter("Maximum Iterations")) {
00586     maxIters_ = params->get("Maximum Iterations",maxIters_default_);
00587 
00588     // Update parameter in our list and in status test.
00589     params_->set("Maximum Iterations", maxIters_);
00590     if (maxIterTest_!=Teuchos::null)
00591       maxIterTest_->setMaxIters( maxIters_ );
00592   }
00593 
00594   // Check for the maximum number of blocks.
00595   if (params->isParameter("Num Blocks")) {
00596     numBlocks_ = params->get("Num Blocks",numBlocks_default_);
00597     TEUCHOS_TEST_FOR_EXCEPTION(numBlocks_ <= 0, std::invalid_argument,
00598                        "Belos::RCGSolMgr: \"Num Blocks\" must be strictly positive.");
00599 
00600     // Update parameter in our list.
00601     params_->set("Num Blocks", numBlocks_);
00602   }
00603 
00604   // Check for the maximum number of blocks.
00605   if (params->isParameter("Num Recycled Blocks")) {
00606     recycleBlocks_ = params->get("Num Recycled Blocks",recycleBlocks_default_);
00607     TEUCHOS_TEST_FOR_EXCEPTION(recycleBlocks_ <= 0, std::invalid_argument,
00608                        "Belos::RCGSolMgr: \"Num Recycled Blocks\" must be strictly positive.");
00609 
00610     TEUCHOS_TEST_FOR_EXCEPTION(recycleBlocks_ >= numBlocks_, std::invalid_argument,
00611                        "Belos::RCGSolMgr: \"Num Recycled Blocks\" must be less than \"Num Blocks\".");
00612 
00613     // Update parameter in our list.
00614     params_->set("Num Recycled Blocks", recycleBlocks_);
00615   }
00616 
00617   // Check to see if the timer label changed.
00618   if (params->isParameter("Timer Label")) {
00619     std::string tempLabel = params->get("Timer Label", label_default_);
00620 
00621     // Update parameter in our list and solver timer
00622     if (tempLabel != label_) {
00623       label_ = tempLabel;
00624       params_->set("Timer Label", label_);
00625       std::string solveLabel = label_ + ": RCGSolMgr total solve time";
00626 #ifdef BELOS_TEUCHOS_TIME_MONITOR
00627       timerSolve_ = Teuchos::TimeMonitor::getNewCounter(solveLabel);
00628 #endif
00629     }
00630   }
00631 
00632   // Check for a change in verbosity level
00633   if (params->isParameter("Verbosity")) {
00634     if (Teuchos::isParameterType<int>(*params,"Verbosity")) {
00635       verbosity_ = params->get("Verbosity", verbosity_default_);
00636     } else {
00637       verbosity_ = (int)Teuchos::getParameter<Belos::MsgType>(*params,"Verbosity");
00638     }
00639 
00640     // Update parameter in our list.
00641     params_->set("Verbosity", verbosity_);
00642     if (printer_ != Teuchos::null)
00643       printer_->setVerbosity(verbosity_);
00644   }
00645 
00646   // Check for a change in output style
00647   if (params->isParameter("Output Style")) {
00648     if (Teuchos::isParameterType<int>(*params,"Output Style")) {
00649       outputStyle_ = params->get("Output Style", outputStyle_default_);
00650     } else {
00651       outputStyle_ = (int)Teuchos::getParameter<Belos::OutputType>(*params,"Output Style");
00652     }
00653 
00654     // Reconstruct the convergence test if the explicit residual test is not being used.
00655     params_->set("Output Style", outputStyle_);
00656     outputTest_ = Teuchos::null;
00657   }
00658 
00659   // output stream
00660   if (params->isParameter("Output Stream")) {
00661     outputStream_ = Teuchos::getParameter<Teuchos::RCP<std::ostream> >(*params,"Output Stream");
00662 
00663     // Update parameter in our list.
00664     params_->set("Output Stream", outputStream_);
00665     if (printer_ != Teuchos::null)
00666       printer_->setOStream( outputStream_ );
00667   }
00668 
00669   // frequency level
00670   if (verbosity_ & Belos::StatusTestDetails) {
00671     if (params->isParameter("Output Frequency")) {
00672       outputFreq_ = params->get("Output Frequency", outputFreq_default_);
00673     }
00674 
00675     // Update parameter in out list and output status test.
00676     params_->set("Output Frequency", outputFreq_);
00677     if (outputTest_ != Teuchos::null)
00678       outputTest_->setOutputFrequency( outputFreq_ );
00679   }
00680 
00681   // Create output manager if we need to.
00682   if (printer_ == Teuchos::null) {
00683     printer_ = Teuchos::rcp( new OutputManager<ScalarType>(verbosity_, outputStream_) );
00684   }
00685 
00686   // Convergence
00687   typedef Belos::StatusTestCombo<ScalarType,MV,OP>  StatusTestCombo_t;
00688   typedef Belos::StatusTestGenResNorm<ScalarType,MV,OP>  StatusTestResNorm_t;
00689 
00690   // Check for convergence tolerance
00691   if (params->isParameter("Convergence Tolerance")) {
00692     convtol_ = params->get("Convergence Tolerance",convtol_default_);
00693 
00694     // Update parameter in our list and residual tests.
00695     params_->set("Convergence Tolerance", convtol_);
00696     if (convTest_ != Teuchos::null)
00697       convTest_->setTolerance( convtol_ );
00698   }
00699 
00700   if (params->isParameter("Show Maximum Residual Norm Only")) {
00701     showMaxResNormOnly_ = Teuchos::getParameter<bool>(*params,"Show Maximum Residual Norm Only");
00702 
00703     // Update parameter in our list and residual tests
00704     params_->set("Show Maximum Residual Norm Only", showMaxResNormOnly_);
00705     if (convTest_ != Teuchos::null)
00706       convTest_->setShowMaxResNormOnly( showMaxResNormOnly_ );
00707   }
00708 
00709   // Create status tests if we need to.
00710 
00711   // Basic test checks maximum iterations and native residual.
00712   if (maxIterTest_ == Teuchos::null)
00713     maxIterTest_ = Teuchos::rcp( new StatusTestMaxIters<ScalarType,MV,OP>( maxIters_ ) );
00714 
00715   // Implicit residual test, using the native residual to determine if convergence was achieved.
00716   if (convTest_ == Teuchos::null)
00717     convTest_ = Teuchos::rcp( new StatusTestResNorm_t( convtol_, 1 ) );
00718 
00719   if (sTest_ == Teuchos::null)
00720     sTest_ = Teuchos::rcp( new StatusTestCombo_t( StatusTestCombo_t::OR, maxIterTest_, convTest_ ) );
00721 
00722   if (outputTest_ == Teuchos::null) {
00723 
00724     // Create the status test output class.
00725     // This class manages and formats the output from the status test.
00726     StatusTestOutputFactory<ScalarType,MV,OP> stoFactory( outputStyle_ );
00727     outputTest_ = stoFactory.create( printer_, sTest_, outputFreq_, Passed+Failed+Undefined );
00728 
00729     // Set the solver string for the output test
00730     std::string solverDesc = " Recycling CG ";
00731     outputTest_->setSolverDesc( solverDesc );
00732   }
00733 
00734   // Create the timer if we need to.
00735   if (timerSolve_ == Teuchos::null) {
00736     std::string solveLabel = label_ + ": RCGSolMgr total solve time";
00737 #ifdef BELOS_TEUCHOS_TIME_MONITOR
00738     timerSolve_ = Teuchos::TimeMonitor::getNewCounter(solveLabel);
00739 #endif
00740   }
00741 
00742   // Inform the solver manager that the current parameters were set.
00743   params_Set_ = true;
00744 }
00745 
00746 
00747 template<class ScalarType, class MV, class OP>
00748 Teuchos::RCP<const Teuchos::ParameterList>
00749 RCGSolMgr<ScalarType,MV,OP,false>::getValidParameters() const
00750 {
00751   static Teuchos::RCP<const Teuchos::ParameterList> validPL;
00752 
00753   // Set all the valid parameters and their default values.
00754   if(is_null(validPL)) {
00755     Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
00756     pl->set("Convergence Tolerance", convtol_default_,
00757       "The relative residual tolerance that needs to be achieved by the\n"
00758       "iterative solver in order for the linear system to be declared converged.");
00759     pl->set("Maximum Iterations", maxIters_default_,
00760       "The maximum number of block iterations allowed for each\n"
00761       "set of RHS solved.");
00762     pl->set("Block Size", blockSize_default_,
00763       "Block Size Parameter -- currently must be 1 for RCG");
00764     pl->set("Num Blocks", numBlocks_default_,
00765       "The length of a cycle (and this max number of search vectors kept)\n");
00766     pl->set("Num Recycled Blocks", recycleBlocks_default_,
00767       "The number of vectors in the recycle subspace.");
00768     pl->set("Verbosity", verbosity_default_,
00769       "What type(s) of solver information should be outputted\n"
00770       "to the output stream.");
00771     pl->set("Output Style", outputStyle_default_,
00772       "What style is used for the solver information outputted\n"
00773       "to the output stream.");
00774     pl->set("Output Frequency", outputFreq_default_,
00775       "How often convergence information should be outputted\n"
00776       "to the output stream.");
00777     pl->set("Output Stream", outputStream_default_,
00778       "A reference-counted pointer to the output stream where all\n"
00779       "solver output is sent.");
00780     pl->set("Show Maximum Residual Norm Only", showMaxResNormOnly_default_,
00781       "When convergence information is printed, only show the maximum\n"
00782       "relative residual norm when the block size is greater than one.");
00783     pl->set("Timer Label", label_default_,
00784       "The string to use as a prefix for the timer labels.");
00785     //  pl->set("Restart Timers", restartTimers_);
00786     validPL = pl;
00787   }
00788   return validPL;
00789 }
00790 
00791 // initializeStateStorage
00792 template<class ScalarType, class MV, class OP>
00793 void RCGSolMgr<ScalarType,MV,OP,false>::initializeStateStorage() {
00794 
00795     // Check if there is any multivector to clone from.
00796     Teuchos::RCP<const MV> rhsMV = problem_->getRHS();
00797     if (rhsMV == Teuchos::null) {
00798       // Nothing to do
00799       return;
00800     }
00801     else {
00802 
00803       // Initialize the state storage
00804       TEUCHOS_TEST_FOR_EXCEPTION(static_cast<ptrdiff_t>(numBlocks_) > MVText::GetGlobalLength(*rhsMV),std::invalid_argument,
00805                          "Belos::RCGSolMgr::initializeStateStorage(): Cannot generate a Krylov basis with dimension larger the operator!");
00806 
00807       // If the subspace has not been initialized before, generate it using the RHS from lp_.
00808       if (P_ == Teuchos::null) {
00809         P_ = MVT::Clone( *rhsMV, numBlocks_+2 );
00810       }
00811       else {
00812         // Generate P_ by cloning itself ONLY if more space is needed.
00813         if (MVT::GetNumberVecs(*P_) < numBlocks_+2) {
00814           Teuchos::RCP<const MV> tmp = P_;
00815           P_ = MVT::Clone( *tmp, numBlocks_+2 );
00816         }
00817       }
00818 
00819       // Generate Ap_ only if it doesn't exist
00820       if (Ap_ == Teuchos::null)
00821         Ap_ = MVT::Clone( *rhsMV, 1 );
00822 
00823       // Generate r_ only if it doesn't exist
00824       if (r_ == Teuchos::null)
00825         r_ = MVT::Clone( *rhsMV, 1 );
00826 
00827       // Generate z_ only if it doesn't exist
00828       if (z_ == Teuchos::null)
00829         z_ = MVT::Clone( *rhsMV, 1 );
00830 
00831       // If the recycle space has not been initialized before, generate it using the RHS from lp_.
00832       if (U_ == Teuchos::null) {
00833         U_ = MVT::Clone( *rhsMV, recycleBlocks_ );
00834       }
00835       else {
00836         // Generate U_ by cloning itself ONLY if more space is needed.
00837         if (MVT::GetNumberVecs(*U_) < recycleBlocks_) {
00838           Teuchos::RCP<const MV> tmp = U_;
00839           U_ = MVT::Clone( *tmp, recycleBlocks_ );
00840         }
00841       }
00842 
00843       // If the recycle space has not be initialized before, generate it using the RHS from lp_.
00844       if (AU_ == Teuchos::null) {
00845         AU_ = MVT::Clone( *rhsMV, recycleBlocks_ );
00846       }
00847       else {
00848         // Generate AU_ by cloning itself ONLY if more space is needed.
00849         if (MVT::GetNumberVecs(*AU_) < recycleBlocks_) {
00850           Teuchos::RCP<const MV> tmp = AU_;
00851           AU_ = MVT::Clone( *tmp, recycleBlocks_ );
00852         }
00853       }
00854 
00855       // If the recycle space has not been initialized before, generate it using the RHS from lp_.
00856       if (U1_ == Teuchos::null) {
00857         U1_ = MVT::Clone( *rhsMV, recycleBlocks_ );
00858       }
00859       else {
00860         // Generate U1_ by cloning itself ONLY if more space is needed.
00861         if (MVT::GetNumberVecs(*U1_) < recycleBlocks_) {
00862           Teuchos::RCP<const MV> tmp = U1_;
00863           U1_ = MVT::Clone( *tmp, recycleBlocks_ );
00864         }
00865       }
00866 
00867       // Generate Alpha_ only if it doesn't exist, otherwise resize it.
00868       if (Alpha_ == Teuchos::null)
00869         Alpha_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( numBlocks_, 1 ) );
00870       else {
00871         if ( (Alpha_->numRows() != numBlocks_) || (Alpha_->numCols() != 1) )
00872           Alpha_->reshape( numBlocks_, 1 );
00873       }
00874 
00875       // Generate Beta_ only if it doesn't exist, otherwise resize it.
00876       if (Beta_ == Teuchos::null)
00877         Beta_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( numBlocks_ + 1, 1 ) );
00878       else {
00879         if ( (Beta_->numRows() != (numBlocks_+1)) || (Beta_->numCols() != 1) )
00880           Beta_->reshape( numBlocks_ + 1, 1 );
00881       }
00882 
00883       // Generate D_ only if it doesn't exist, otherwise resize it.
00884       if (D_ == Teuchos::null)
00885         D_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( numBlocks_ , 1 ) );
00886       else {
00887         if ( (D_->numRows() != numBlocks_) || (D_->numCols() != 1) )
00888           D_->reshape( numBlocks_, 1 );
00889       }
00890 
00891       // Generate Delta_ only if it doesn't exist, otherwise resize it.
00892       if (Delta_ == Teuchos::null)
00893         Delta_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( recycleBlocks_, numBlocks_ + 1 ) );
00894       else {
00895         if ( (Delta_->numRows() != recycleBlocks_) || (Delta_->numCols() != (numBlocks_ + 1)) )
00896           Delta_->reshape( recycleBlocks_, numBlocks_ + 1 );
00897       }
00898 
00899       // Generate UTAU_ only if it doesn't exist, otherwise resize it.
00900       if (UTAU_ == Teuchos::null)
00901         UTAU_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( recycleBlocks_, recycleBlocks_ ) );
00902       else {
00903         if ( (UTAU_->numRows() != recycleBlocks_) || (UTAU_->numCols() != recycleBlocks_) )
00904           UTAU_->reshape( recycleBlocks_, recycleBlocks_ );
00905       }
00906 
00907       // Generate LUUTAU_ only if it doesn't exist, otherwise resize it.
00908       if (LUUTAU_ == Teuchos::null)
00909         LUUTAU_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( recycleBlocks_, recycleBlocks_ ) );
00910       else {
00911         if ( (LUUTAU_->numRows() != recycleBlocks_) || (LUUTAU_->numCols() != recycleBlocks_) )
00912           LUUTAU_->reshape( recycleBlocks_, recycleBlocks_ );
00913       }
00914 
00915       // Generate ipiv_ only if it doesn't exist, otherwise resize it.
00916       if (ipiv_ == Teuchos::null)
00917         ipiv_ = Teuchos::rcp( new std::vector<int>(recycleBlocks_) );
00918       else {
00919         if ( (int)ipiv_->size() != recycleBlocks_ ) // if ipiv not correct size, always resize it
00920           ipiv_->resize(recycleBlocks_);
00921       }
00922 
00923       // Generate AUTAU_ only if it doesn't exist, otherwise resize it.
00924       if (AUTAU_ == Teuchos::null)
00925         AUTAU_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( recycleBlocks_, recycleBlocks_ ) );
00926       else {
00927         if ( (AUTAU_->numRows() != recycleBlocks_) || (AUTAU_->numCols() != recycleBlocks_) )
00928           AUTAU_->reshape( recycleBlocks_, recycleBlocks_ );
00929       }
00930 
00931       // Generate rTz_old_ only if it doesn't exist
00932       if (rTz_old_ == Teuchos::null)
00933         rTz_old_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( 1, 1 ) );
00934       else {
00935         if ( (rTz_old_->numRows() != 1) || (rTz_old_->numCols() != 1) )
00936           rTz_old_->reshape( 1, 1 );
00937       }
00938 
00939       // Generate F_ only if it doesn't exist
00940       if (F_ == Teuchos::null)
00941         F_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( numBlocks_+recycleBlocks_, numBlocks_+recycleBlocks_ ) );
00942       else {
00943         if ( (F_->numRows() != (numBlocks_+recycleBlocks_)) || (F_->numCols() != numBlocks_+recycleBlocks_) )
00944           F_->reshape( numBlocks_+recycleBlocks_, numBlocks_+recycleBlocks_ );
00945       }
00946 
00947       // Generate G_ only if it doesn't exist
00948       if (G_ == Teuchos::null)
00949         G_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( numBlocks_+recycleBlocks_, numBlocks_+recycleBlocks_ ) );
00950       else {
00951         if ( (G_->numRows() != (numBlocks_+recycleBlocks_)) || (G_->numCols() != numBlocks_+recycleBlocks_) )
00952           G_->reshape( numBlocks_+recycleBlocks_, numBlocks_+recycleBlocks_ );
00953       }
00954 
00955       // Generate Y_ only if it doesn't exist
00956       if (Y_ == Teuchos::null)
00957         Y_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( numBlocks_+recycleBlocks_, recycleBlocks_ ) );
00958       else {
00959         if ( (Y_->numRows() != (numBlocks_+recycleBlocks_)) || (Y_->numCols() != numBlocks_+recycleBlocks_) )
00960           Y_->reshape( numBlocks_+recycleBlocks_, numBlocks_+recycleBlocks_ );
00961       }
00962 
00963       // Generate L2_ only if it doesn't exist
00964       if (L2_ == Teuchos::null)
00965         L2_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( numBlocks_+1, numBlocks_ ) );
00966       else {
00967         if ( (L2_->numRows() != (numBlocks_+1)) || (L2_->numCols() != numBlocks_) )
00968           L2_->reshape( numBlocks_+1, numBlocks_ );
00969       }
00970 
00971       // Generate DeltaL2_ only if it doesn't exist
00972       if (DeltaL2_ == Teuchos::null)
00973         DeltaL2_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( recycleBlocks_, numBlocks_ ) );
00974       else {
00975         if ( (DeltaL2_->numRows() != (recycleBlocks_)) || (DeltaL2_->numCols() != (numBlocks_) ) )
00976           DeltaL2_->reshape( recycleBlocks_, numBlocks_ );
00977       }
00978 
00979       // Generate AU1TUDeltaL2_ only if it doesn't exist
00980       if (AU1TUDeltaL2_ == Teuchos::null)
00981         AU1TUDeltaL2_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( recycleBlocks_, numBlocks_ ) );
00982       else {
00983         if ( (AU1TUDeltaL2_->numRows() != (recycleBlocks_)) || (AU1TUDeltaL2_->numCols() != (numBlocks_) ) )
00984           AU1TUDeltaL2_->reshape( recycleBlocks_, numBlocks_ );
00985       }
00986 
00987       // Generate AU1TAU1_ only if it doesn't exist
00988       if (AU1TAU1_ == Teuchos::null)
00989         AU1TAU1_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( recycleBlocks_, recycleBlocks_ ) );
00990       else {
00991         if ( (AU1TAU1_->numRows() != (recycleBlocks_)) || (AU1TAU1_->numCols() != (recycleBlocks_) ) )
00992           AU1TAU1_->reshape( recycleBlocks_, recycleBlocks_ );
00993       }
00994 
00995       // Generate GY_ only if it doesn't exist
00996       if (GY_ == Teuchos::null)
00997         GY_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( numBlocks_ + recycleBlocks_, recycleBlocks_ ) );
00998       else {
00999         if ( (GY_->numRows() != (numBlocks_ + recycleBlocks_)) || (GY_->numCols() != (recycleBlocks_) ) )
01000           GY_->reshape( numBlocks_+recycleBlocks_, recycleBlocks_ );
01001       }
01002 
01003       // Generate AU1TU1_ only if it doesn't exist
01004       if (AU1TU1_ == Teuchos::null)
01005         AU1TU1_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( recycleBlocks_, recycleBlocks_ ) );
01006       else {
01007         if ( (AU1TU1_->numRows() != (recycleBlocks_)) || (AU1TU1_->numCols() != (recycleBlocks_) ) )
01008           AU1TU1_->reshape( recycleBlocks_, recycleBlocks_ );
01009       }
01010 
01011       // Generate FY_ only if it doesn't exist
01012       if (FY_ == Teuchos::null)
01013         FY_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( numBlocks_ + recycleBlocks_, recycleBlocks_ ) );
01014       else {
01015         if ( (FY_->numRows() != (numBlocks_ + recycleBlocks_)) || (FY_->numCols() != (recycleBlocks_) ) )
01016           FY_->reshape( numBlocks_+recycleBlocks_, recycleBlocks_ );
01017       }
01018 
01019       // Generate AU1TAP_ only if it doesn't exist
01020       if (AU1TAP_ == Teuchos::null)
01021         AU1TAP_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( recycleBlocks_, numBlocks_ ) );
01022       else {
01023         if ( (AU1TAP_->numRows() != (recycleBlocks_)) || (AU1TAP_->numCols() != (numBlocks_) ) )
01024           AU1TAP_->reshape( recycleBlocks_, numBlocks_ );
01025       }
01026 
01027       // Generate APTAP_ only if it doesn't exist
01028       if (APTAP_ == Teuchos::null)
01029         APTAP_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( numBlocks_, numBlocks_ ) );
01030       else {
01031         if ( (APTAP_->numRows() != (numBlocks_)) || (APTAP_->numCols() != (numBlocks_) ) )
01032           APTAP_->reshape( numBlocks_, numBlocks_ );
01033       }
01034 
01035       // If the subspace has not been initialized before, generate it using the RHS from lp_.
01036       if (U1Y1_ == Teuchos::null) {
01037         U1Y1_ = MVT::Clone( *rhsMV, recycleBlocks_ );
01038       }
01039       else {
01040         // Generate U1Y1_ by cloning itself ONLY if more space is needed.
01041         if (MVT::GetNumberVecs(*U1Y1_) < recycleBlocks_) {
01042           Teuchos::RCP<const MV> tmp = U1Y1_;
01043           U1Y1_ = MVT::Clone( *tmp, recycleBlocks_ );
01044         }
01045       }
01046 
01047       // If the subspace has not been initialized before, generate it using the RHS from lp_.
01048       if (PY2_ == Teuchos::null) {
01049         PY2_ = MVT::Clone( *rhsMV, recycleBlocks_ );
01050       }
01051       else {
01052         // Generate PY2_ by cloning itself ONLY if more space is needed.
01053         if (MVT::GetNumberVecs(*PY2_) < recycleBlocks_) {
01054           Teuchos::RCP<const MV> tmp = PY2_;
01055           PY2_ = MVT::Clone( *tmp, recycleBlocks_ );
01056         }
01057       }
01058 
01059       // Generate AUTAP_ only if it doesn't exist
01060       if (AUTAP_ == Teuchos::null)
01061         AUTAP_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( recycleBlocks_, numBlocks_ ) );
01062       else {
01063         if ( (AUTAP_->numRows() != (recycleBlocks_)) || (AUTAP_->numCols() != (numBlocks_) ) )
01064           AUTAP_->reshape( recycleBlocks_, numBlocks_ );
01065       }
01066 
01067       // Generate AU1TU_ only if it doesn't exist
01068       if (AU1TU_ == Teuchos::null)
01069         AU1TU_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( recycleBlocks_, recycleBlocks_ ) );
01070       else {
01071         if ( (AU1TU_->numRows() != (recycleBlocks_)) || (AU1TU_->numCols() != (recycleBlocks_) ) )
01072           AU1TU_->reshape( recycleBlocks_, recycleBlocks_ );
01073       }
01074 
01075 
01076     }
01077 }
01078 
01079 template<class ScalarType, class MV, class OP>
01080 ReturnType RCGSolMgr<ScalarType,MV,OP,false>::solve() {
01081 
01082   Teuchos::LAPACK<int,ScalarType> lapack;
01083   std::vector<int> index(recycleBlocks_);
01084   ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
01085   ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
01086 
01087   // Count of number of cycles performed on current rhs
01088   int cycle = 0;
01089 
01090   // Set the current parameters if they were not set before.
01091   // NOTE:  This may occur if the user generated the solver manager with the default constructor and
01092   // then didn't set any parameters using setParameters().
01093   if (!params_Set_) {
01094     setParameters(Teuchos::parameterList(*getValidParameters()));
01095   }
01096 
01097   TEUCHOS_TEST_FOR_EXCEPTION(problem_ == Teuchos::null,RCGSolMgrLinearProblemFailure,
01098                      "Belos::RCGSolMgr::solve(): Linear problem is not a valid object.");
01099   TEUCHOS_TEST_FOR_EXCEPTION(!problem_->isProblemSet(),RCGSolMgrLinearProblemFailure,
01100                      "Belos::RCGSolMgr::solve(): Linear problem is not ready, setProblem() has not been called.");
01101   TEUCHOS_TEST_FOR_EXCEPTION((problem_->getLeftPrec() != Teuchos::null)&&(problem_->getRightPrec() != Teuchos::null),
01102                      RCGSolMgrLinearProblemFailure,
01103                      "Belos::RCGSolMgr::solve(): RCG does not support split preconditioning, only set left or right preconditioner.");
01104 
01105   // Grab the preconditioning object
01106   Teuchos::RCP<OP> precObj;
01107   if (problem_->getLeftPrec() != Teuchos::null) {
01108     precObj = Teuchos::rcp_const_cast<OP>(problem_->getLeftPrec());
01109   }
01110   else if (problem_->getRightPrec() != Teuchos::null) {
01111     precObj = Teuchos::rcp_const_cast<OP>(problem_->getRightPrec());
01112   }
01113 
01114   // Create indices for the linear systems to be solved.
01115   int numRHS2Solve = MVT::GetNumberVecs( *(problem_->getRHS()) );
01116   std::vector<int> currIdx(1);
01117   currIdx[0] = 0;
01118 
01119   // Inform the linear problem of the current linear system to solve.
01120   problem_->setLSIndex( currIdx );
01121 
01122   // Check the number of blocks and change them if necessary.
01123   ptrdiff_t dim = MVText::GetGlobalLength( *(problem_->getRHS()) );
01124   if (numBlocks_ > dim) {
01125     numBlocks_ = Teuchos::asSafe<int>(dim);
01126     params_->set("Num Blocks", numBlocks_);
01127     printer_->stream(Warnings) <<
01128       "Warning! Requested Krylov subspace dimension is larger than operator dimension!" << std::endl <<
01129       " The maximum number of blocks allowed for the Krylov subspace will be adjusted to " << numBlocks_ << std::endl;
01130   }
01131 
01132   // Initialize storage for all state variables
01133   initializeStateStorage();
01134 
01135   // Parameter list
01136   Teuchos::ParameterList plist;
01137   plist.set("Num Blocks",numBlocks_);
01138   plist.set("Recycled Blocks",recycleBlocks_);
01139 
01140   // Reset the status test.
01141   outputTest_->reset();
01142 
01143   // Assume convergence is achieved, then let any failed convergence set this to false.
01144   bool isConverged = true;
01145 
01146   // Compute AU = A*U, UTAU = U'*AU, AUTAU = (AU)'*(AU)
01147   if (existU_) {
01148     index.resize(recycleBlocks_);
01149     for (int i=0; i<recycleBlocks_; ++i) { index[i] = i; }
01150     Teuchos::RCP<const MV> Utmp  = MVT::CloneView( *U_,  index );
01151     index.resize(recycleBlocks_);
01152     for (int i=0; i<recycleBlocks_; ++i) { index[i] = i; }
01153     Teuchos::RCP<MV> AUtmp = MVT::CloneViewNonConst( *AU_, index );
01154     // Initialize AU
01155     problem_->applyOp( *Utmp, *AUtmp );
01156     // Initialize UTAU
01157     Teuchos::SerialDenseMatrix<int,ScalarType> UTAUtmp( Teuchos::View, *UTAU_, recycleBlocks_, recycleBlocks_ );
01158     MVT::MvTransMv( one, *Utmp, *AUtmp, UTAUtmp );
01159     // Initialize AUTAU  ( AUTAU = AU'*(M\AU) )
01160     Teuchos::SerialDenseMatrix<int,ScalarType> AUTAUtmp( Teuchos::View, *AUTAU_, recycleBlocks_, recycleBlocks_ );
01161     if ( precObj != Teuchos::null ) {
01162       index.resize(recycleBlocks_);
01163       for (int i=0; i<recycleBlocks_; ++i) { index[i] = i; }
01164       index.resize(recycleBlocks_);
01165       for (int ii=0; ii<recycleBlocks_; ++ii) { index[ii] = ii; }
01166       Teuchos::RCP<MV> PCAU = MVT::CloneViewNonConst( *U1_, index ); // use U1 as temp storage
01167       OPT::Apply( *precObj, *AUtmp, *PCAU );
01168       MVT::MvTransMv( one, *AUtmp, *PCAU, AUTAUtmp );
01169     } else {
01170       MVT::MvTransMv( one, *AUtmp, *AUtmp, AUTAUtmp );
01171     }
01172   }
01173 
01175   // RCG solver
01176 
01177   Teuchos::RCP<RCGIter<ScalarType,MV,OP> > rcg_iter;
01178   rcg_iter = Teuchos::rcp( new RCGIter<ScalarType,MV,OP>(problem_,printer_,outputTest_,plist) );
01179 
01180   // Enter solve() iterations
01181   {
01182 #ifdef BELOS_TEUCHOS_TIME_MONITOR
01183     Teuchos::TimeMonitor slvtimer(*timerSolve_);
01184 #endif
01185 
01186     while ( numRHS2Solve > 0 ) {
01187 
01188       // Debugging output to tell use if recycle space exists and will be used
01189       if (printer_->isVerbosity( Debug ) ) {
01190         if (existU_) printer_->print( Debug, "Using recycle space generated from previous call to solve()." );
01191         else printer_->print( Debug, "No recycle space exists." );
01192       }
01193 
01194       // Reset the number of iterations.
01195       rcg_iter->resetNumIters();
01196 
01197       // Set the current number of recycle blocks and subspace dimension with the RCG iteration.
01198       rcg_iter->setSize( recycleBlocks_, numBlocks_ );
01199 
01200       // Reset the number of calls that the status test output knows about.
01201       outputTest_->resetNumCalls();
01202 
01203       // indicate that updated recycle space has not yet been generated for this linear system
01204       existU1_ = false;
01205 
01206       // reset cycle count
01207       cycle = 0;
01208 
01209       // Get the current residual
01210       problem_->computeCurrResVec( &*r_ );
01211 
01212       // If U exists, find best soln over this space first
01213       if (existU_) {
01214         // Solve linear system UTAU * y = (U'*r)
01215         Teuchos::SerialDenseMatrix<int,ScalarType> Utr(recycleBlocks_,1);
01216         index.resize(recycleBlocks_);
01217         for (int ii=0; ii<recycleBlocks_; ++ii) { index[ii] = ii; }
01218         Teuchos::RCP<const MV> Utmp  = MVT::CloneView( *U_,  index );
01219         MVT::MvTransMv( one, *Utmp, *r_, Utr );
01220         Teuchos::SerialDenseMatrix<int,ScalarType> UTAUtmp( Teuchos::View, *UTAU_, recycleBlocks_, recycleBlocks_ );
01221         Teuchos::SerialDenseMatrix<int,ScalarType> LUUTAUtmp( Teuchos::View, *LUUTAU_, recycleBlocks_, recycleBlocks_ );
01222         LUUTAUtmp.assign(UTAUtmp);
01223         int info = 0;
01224         lapack.GESV(recycleBlocks_, 1, LUUTAUtmp.values(), LUUTAUtmp.stride(), &(*ipiv_)[0], Utr.values(), Utr.stride(), &info);
01225         TEUCHOS_TEST_FOR_EXCEPTION(info != 0, RCGSolMgrLAPACKFailure,
01226                            "Belos::RCGSolMgr::solve(): LAPACK GESV failed to compute a solution.");
01227 
01228         // Update solution (x = x + U*y)
01229         MVT::MvTimesMatAddMv( one, *Utmp, Utr, one, *problem_->getCurrLHSVec() );
01230 
01231         // Update residual ( r = r - AU*y )
01232         index.resize(recycleBlocks_);
01233         for (int ii=0; ii<recycleBlocks_; ++ii) { index[ii] = ii; }
01234         Teuchos::RCP<const MV> AUtmp  = MVT::CloneView( *AU_, index );
01235         MVT::MvTimesMatAddMv( -one, *AUtmp, Utr, one, *r_ );
01236       }
01237 
01238       if ( precObj != Teuchos::null ) {
01239         OPT::Apply( *precObj, *r_, *z_ );
01240       } else {
01241         z_ = r_;
01242       }
01243 
01244       // rTz_old = r'*z
01245       MVT::MvTransMv( one, *r_, *z_, *rTz_old_ );
01246 
01247       if ( existU_ ) {
01248         // mu = UTAU\(AU'*z);
01249         Teuchos::SerialDenseMatrix<int,ScalarType> mu( Teuchos::View, *Delta_, recycleBlocks_, 1);
01250         index.resize(recycleBlocks_);
01251         for (int ii=0; ii<recycleBlocks_; ++ii) { index[ii] = ii; }
01252         Teuchos::RCP<const MV> AUtmp  = MVT::CloneView( *AU_, index );
01253         MVT::MvTransMv( one, *AUtmp, *z_, mu );
01254         Teuchos::SerialDenseMatrix<int,ScalarType> LUUTAUtmp( Teuchos::View, *LUUTAU_, recycleBlocks_, recycleBlocks_ );
01255         char TRANS = 'N';
01256         int info;
01257         lapack.GETRS( TRANS, recycleBlocks_, 1, LUUTAUtmp.values(), LUUTAUtmp.stride(), &(*ipiv_)[0], mu.values(), mu.stride(), &info );
01258         TEUCHOS_TEST_FOR_EXCEPTION(info != 0, RCGSolMgrLAPACKFailure,
01259                            "Belos::RCGSolMgr::solve(): LAPACK GETRS failed to compute a solution.");
01260         // p  = z - U*mu;
01261         index.resize( 1 );
01262         index[0] = 0;
01263         Teuchos::RCP<MV> Ptmp  = MVT::CloneViewNonConst( *P_, index );
01264         MVT::MvAddMv(one,*z_,zero,*z_,*Ptmp);
01265         MVT::MvTimesMatAddMv( -one, *U_, mu, one, *Ptmp );
01266       } else {
01267         // p = z;
01268         index.resize( 1 );
01269         index[0] = 0;
01270         Teuchos::RCP<MV> Ptmp  = MVT::CloneViewNonConst( *P_, index );
01271         MVT::MvAddMv(one,*z_,zero,*z_,*Ptmp);
01272       }
01273 
01274       // Set the new state and initialize the solver.
01275       RCGIterState<ScalarType,MV> newstate;
01276 
01277       // Create RCP views here
01278       index.resize( numBlocks_+1 );
01279       for (int ii=0; ii<(numBlocks_+1); ++ii) { index[ii] = ii; }
01280       newstate.P  = MVT::CloneViewNonConst( *P_,  index );
01281       index.resize( recycleBlocks_ );
01282       for (int ii=0; ii<recycleBlocks_; ++ii) { index[ii] = ii; }
01283       newstate.U  = MVT::CloneViewNonConst( *U_,  index );
01284       index.resize( recycleBlocks_ );
01285       for (int ii=0; ii<recycleBlocks_; ++ii) { index[ii] = ii; }
01286       newstate.AU  = MVT::CloneViewNonConst( *AU_,  index );
01287       newstate.Alpha = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( Teuchos::View, *Alpha_, numBlocks_, 1 ) );
01288       newstate.Beta = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( Teuchos::View, *Beta_, numBlocks_, 1 ) );
01289       newstate.D = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( Teuchos::View, *D_, numBlocks_, 1 ) );
01290       newstate.Delta = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( Teuchos::View, *Delta_, recycleBlocks_, numBlocks_, 0, 1 ) );
01291       newstate.LUUTAU = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( Teuchos::View, *LUUTAU_, recycleBlocks_, recycleBlocks_ ) );
01292       // assign the rest of the values in the struct
01293       newstate.curDim = 1; // We have initialized the first search vector
01294       newstate.Ap = Ap_;
01295       newstate.r = r_;
01296       newstate.z = z_;
01297       newstate.existU = existU_;
01298       newstate.ipiv = ipiv_;
01299       newstate.rTz_old = rTz_old_;
01300 
01301       rcg_iter->initialize(newstate);
01302 
01303       while(1) {
01304 
01305         // tell rcg_iter to iterate
01306         try {
01307           rcg_iter->iterate();
01308 
01310           //
01311           // check convergence first
01312           //
01314           if ( convTest_->getStatus() == Passed ) {
01315             // We have convergence
01316             break; // break from while(1){rcg_iter->iterate()}
01317           }
01319           //
01320           // check for maximum iterations
01321           //
01323           else if ( maxIterTest_->getStatus() == Passed ) {
01324             // we don't have convergence
01325             isConverged = false;
01326             break; // break from while(1){rcg_iter->iterate()}
01327           }
01329           //
01330           // check if cycle complete; update for next cycle
01331           //
01333           else if ( rcg_iter->getCurSubspaceDim() == rcg_iter->getMaxSubspaceDim() ) {
01334             // index into P_ of last search vector generated this cycle
01335             int lastp = -1;
01336             // index into Beta_ of last entry generated this cycle
01337             int lastBeta = -1;
01338             if (recycleBlocks_ > 0) {
01339               if (!existU_) {
01340                 if (cycle == 0) { // No U, no U1
01341 
01342                    Teuchos::SerialDenseMatrix<int,ScalarType> Ftmp( Teuchos::View, *F_, numBlocks_, numBlocks_ );
01343                    Teuchos::SerialDenseMatrix<int,ScalarType> Gtmp( Teuchos::View, *G_, numBlocks_, numBlocks_ );
01344                    Teuchos::SerialDenseMatrix<int,ScalarType> Dtmp( Teuchos::View, *D_, numBlocks_, 1 );
01345                    Teuchos::SerialDenseMatrix<int,ScalarType> Alphatmp( Teuchos::View, *Alpha_, numBlocks_, 1 );
01346                    Teuchos::SerialDenseMatrix<int,ScalarType> Betatmp( Teuchos::View, *Beta_, numBlocks_, 1 );
01347                    Ftmp.putScalar(zero);
01348                    Gtmp.putScalar(zero);
01349                    for (int ii=0;ii<numBlocks_;ii++) {
01350                      Gtmp(ii,ii) = (Dtmp(ii,0) / Alphatmp(ii,0))*(1 + Betatmp(ii,0));
01351                      if (ii > 0) {
01352                        Gtmp(ii-1,ii) = -Dtmp(ii,0)/Alphatmp(ii-1,0);
01353                        Gtmp(ii,ii-1) = -Dtmp(ii,0)/Alphatmp(ii-1,0);
01354                      }
01355                      Ftmp(ii,ii) = Dtmp(ii,0);
01356                    }
01357 
01358                    // compute harmonic Ritz vectors
01359                    Teuchos::SerialDenseMatrix<int,ScalarType> Ytmp( Teuchos::View, *Y_, numBlocks_, recycleBlocks_ );
01360                    getHarmonicVecs(Ftmp,Gtmp,Ytmp);
01361 
01362                    // U1 = [P(:,1:end-1)*Y];
01363                    index.resize( numBlocks_ );
01364                    for (int ii=0; ii<numBlocks_; ++ii) { index[ii] = ii; }
01365                    Teuchos::RCP<const MV> Ptmp = MVT::CloneView( *P_,  index );
01366                    index.resize( recycleBlocks_ );
01367                    for (int ii=0; ii<recycleBlocks_; ++ii) { index[ii] = ii; }
01368                    Teuchos::RCP<MV> U1tmp  = MVT::CloneViewNonConst( *U1_,  index );
01369                    MVT::MvTimesMatAddMv( one, *Ptmp, Ytmp, zero, *U1tmp );
01370 
01371                    // Precompute some variables for next cycle
01372 
01373                    // AU1TAU1     = Y'*G*Y;
01374                    Teuchos::SerialDenseMatrix<int,ScalarType> GYtmp( Teuchos::View, *GY_, numBlocks_, recycleBlocks_ );
01375                    GYtmp.multiply(Teuchos::NO_TRANS,Teuchos::NO_TRANS,one,Gtmp,Ytmp,zero);
01376                    Teuchos::SerialDenseMatrix<int,ScalarType> AU1TAU1tmp( Teuchos::View, *AU1TAU1_, recycleBlocks_, recycleBlocks_ );
01377                    AU1TAU1tmp.multiply(Teuchos::TRANS,Teuchos::NO_TRANS,one,Ytmp,GYtmp,zero);
01378 
01379                    // AU1TU1      = Y'*F*Y;
01380                    Teuchos::SerialDenseMatrix<int,ScalarType> FYtmp( Teuchos::View, *FY_, numBlocks_, recycleBlocks_ );
01381                    FYtmp.multiply(Teuchos::NO_TRANS,Teuchos::NO_TRANS,one,Ftmp,Ytmp,zero);
01382                    Teuchos::SerialDenseMatrix<int,ScalarType> AU1TU1tmp( Teuchos::View, *AU1TU1_, recycleBlocks_, recycleBlocks_ );
01383                    AU1TU1tmp.multiply(Teuchos::TRANS,Teuchos::NO_TRANS,one,Ytmp,FYtmp,zero);
01384 
01385                    Teuchos::SerialDenseMatrix<int,ScalarType> AU1TAPtmp( Teuchos::View, *AU1TAP_, recycleBlocks_, 1 );
01386                    // Must reinitialize AU1TAP; can become dense later
01387                    AU1TAPtmp.putScalar(zero);
01388                    // AU1TAP(:,1) = Y(end,:)' * (-1/Alpha(end));
01389                    ScalarType alphatmp = -1.0 / Alphatmp(numBlocks_-1,0);
01390                    for (int ii=0; ii<recycleBlocks_; ++ii) {
01391                       AU1TAPtmp(ii,0) = Ytmp(numBlocks_-1,ii) * alphatmp;
01392                    }
01393 
01394                    // indicate that updated recycle space now defined
01395                    existU1_ = true;
01396 
01397                    // Indicate the size of the P, Beta structures generated this cycle
01398                    lastp = numBlocks_;
01399                    lastBeta = numBlocks_-1;
01400 
01401                 } // if (cycle == 0)
01402                 else { // No U, but U1 guaranteed to exist now
01403 
01404                    // Finish computation of subblocks
01405                    // AU1TAP = AU1TAP * D(1);
01406                    Teuchos::SerialDenseMatrix<int,ScalarType> Dtmp( Teuchos::View, *D_, numBlocks_, 1 );
01407                    Teuchos::SerialDenseMatrix<int,ScalarType> AU1TAPtmp( Teuchos::View, *AU1TAP_, recycleBlocks_, numBlocks_ );
01408                    AU1TAPtmp.scale(Dtmp(0,0));
01409 
01410                    Teuchos::SerialDenseMatrix<int,ScalarType> Alphatmp( Teuchos::View, *Alpha_, numBlocks_, 1 );
01411                    Teuchos::SerialDenseMatrix<int,ScalarType> Betatmp( Teuchos::View, *Beta_, numBlocks_+1, 1 );
01412                    Teuchos::SerialDenseMatrix<int,ScalarType> APTAPtmp( Teuchos::View, *APTAP_, numBlocks_, numBlocks_ );
01413                    APTAPtmp.putScalar(zero);
01414                    for (int ii=0; ii<numBlocks_; ii++) {
01415                      APTAPtmp(ii,ii) = (Dtmp(ii,0) / Alphatmp(ii,0))*(1 + Betatmp(ii+1,0));
01416                      if (ii > 0) {
01417                        APTAPtmp(ii-1,ii) = -Dtmp(ii,0)/Alphatmp(ii-1,0);
01418                        APTAPtmp(ii,ii-1) = -Dtmp(ii,0)/Alphatmp(ii-1,0);
01419                      }
01420                    }
01421 
01422                    // F = [AU1TU1 zeros(k,m); zeros(m,k) diag(D)];
01423                    Teuchos::SerialDenseMatrix<int,ScalarType> Ftmp( Teuchos::View, *F_, (numBlocks_+recycleBlocks_), (numBlocks_+recycleBlocks_) );
01424                    Teuchos::SerialDenseMatrix<int,ScalarType> F11( Teuchos::View, *F_, recycleBlocks_, recycleBlocks_ );
01425                    Teuchos::SerialDenseMatrix<int,ScalarType> F12( Teuchos::View, *F_, recycleBlocks_, numBlocks_, 0, recycleBlocks_ );
01426                    Teuchos::SerialDenseMatrix<int,ScalarType> F21( Teuchos::View, *F_, numBlocks_, recycleBlocks_, recycleBlocks_, 0 );
01427                    Teuchos::SerialDenseMatrix<int,ScalarType> F22( Teuchos::View, *F_, numBlocks_, numBlocks_, recycleBlocks_, recycleBlocks_ );
01428                    Teuchos::SerialDenseMatrix<int,ScalarType> AU1TU1tmp( Teuchos::View, *AU1TU1_, recycleBlocks_, recycleBlocks_ );
01429                    F11.assign(AU1TU1tmp);
01430                    F12.putScalar(zero);
01431                    F21.putScalar(zero);
01432                    F22.putScalar(zero);
01433                    for(int ii=0;ii<numBlocks_;ii++) {
01434                      F22(ii,ii) = Dtmp(ii,0);
01435                    }
01436 
01437                    // G = [AU1TAU1 AU1TAP; AU1TAP' APTAP];
01438                    Teuchos::SerialDenseMatrix<int,ScalarType> Gtmp( Teuchos::View, *G_, (numBlocks_+recycleBlocks_), (numBlocks_+recycleBlocks_) );
01439                    Teuchos::SerialDenseMatrix<int,ScalarType> AU1TAU1tmp( Teuchos::View, *AU1TAU1_, recycleBlocks_, recycleBlocks_ );
01440                    Teuchos::SerialDenseMatrix<int,ScalarType> G11( Teuchos::View, *G_, recycleBlocks_, recycleBlocks_ );
01441                    Teuchos::SerialDenseMatrix<int,ScalarType> G12( Teuchos::View, *G_, recycleBlocks_, numBlocks_, 0, recycleBlocks_ );
01442                    Teuchos::SerialDenseMatrix<int,ScalarType> G21( Teuchos::View, *G_, numBlocks_, recycleBlocks_, recycleBlocks_, 0 );
01443                    Teuchos::SerialDenseMatrix<int,ScalarType> G22( Teuchos::View, *G_, numBlocks_, numBlocks_, recycleBlocks_, recycleBlocks_ );
01444                    G11.assign(AU1TAU1tmp);
01445                    G12.assign(AU1TAPtmp);
01446                    // G21 = G12'; (no transpose operator exists for SerialDenseMatrix; Do copy manually)
01447                    for (int ii=0;ii<recycleBlocks_;++ii)
01448                      for (int jj=0;jj<numBlocks_;++jj)
01449                        G21(jj,ii) = G12(ii,jj);
01450                    G22.assign(APTAPtmp);
01451 
01452                    // compute harmonic Ritz vectors
01453                    Teuchos::SerialDenseMatrix<int,ScalarType> Ytmp( Teuchos::View, *Y_, (recycleBlocks_+numBlocks_), recycleBlocks_ );
01454                    getHarmonicVecs(Ftmp,Gtmp,Ytmp);
01455 
01456                    // U1 = [U1 P(:,2:end-1)]*Y;
01457                    index.resize( numBlocks_ );
01458                    for (int ii=0; ii<numBlocks_; ++ii) { index[ii] = ii+1; }
01459                    Teuchos::RCP<const MV> Ptmp = MVT::CloneView( *P_,  index );
01460                    index.resize( recycleBlocks_ );
01461                    for (int ii=0; ii<recycleBlocks_; ++ii) { index[ii] = ii; }
01462                    Teuchos::RCP<MV> PY2tmp = MVT::CloneViewNonConst( *PY2_,  index );
01463                    Teuchos::SerialDenseMatrix<int,ScalarType> Y2( Teuchos::View, *Y_, numBlocks_, recycleBlocks_, recycleBlocks_, 0 );
01464                    index.resize( recycleBlocks_ );
01465                    for (int ii=0; ii<recycleBlocks_; ++ii) { index[ii] = ii; }
01466                    Teuchos::RCP<MV> U1tmp  = MVT::CloneViewNonConst( *U1_,  index );
01467                    index.resize( recycleBlocks_ );
01468                    for (int ii=0; ii<recycleBlocks_; ++ii) { index[ii] = ii; }
01469                    Teuchos::RCP<MV> U1Y1tmp  = MVT::CloneViewNonConst( *U1Y1_,  index );
01470                    Teuchos::SerialDenseMatrix<int,ScalarType> Y1( Teuchos::View, *Y_, recycleBlocks_, recycleBlocks_ );
01471                    MVT::MvTimesMatAddMv( one, *Ptmp, Y2, zero, *PY2tmp );
01472                    MVT::MvTimesMatAddMv( one, *U1tmp, Y1, zero, *U1Y1tmp );
01473                    MVT::MvAddMv(one,*U1Y1tmp, one, *PY2tmp, *U1tmp);
01474 
01475                    // Precompute some variables for next cycle
01476 
01477                    // AU1TAU1     = Y'*G*Y;
01478                    Teuchos::SerialDenseMatrix<int,ScalarType> GYtmp( Teuchos::View, *GY_, (numBlocks_+recycleBlocks_), recycleBlocks_ );
01479                    GYtmp.multiply(Teuchos::NO_TRANS,Teuchos::NO_TRANS,one,Gtmp,Ytmp,zero);
01480                    AU1TAU1tmp.multiply(Teuchos::TRANS,Teuchos::NO_TRANS,one,Ytmp,GYtmp,zero);
01481 
01482                    // AU1TAP      = zeros(k,m);
01483                    // AU1TAP(:,1) = Y(end,:)' * (-1/Alpha(end));
01484                    AU1TAPtmp.putScalar(zero);
01485                    ScalarType alphatmp = -1.0 / Alphatmp(numBlocks_-1,0);
01486                    for (int ii=0; ii<recycleBlocks_; ++ii) {
01487                       AU1TAPtmp(ii,0) = Ytmp(numBlocks_+recycleBlocks_-1,ii) * alphatmp;
01488                    }
01489 
01490                    // AU1TU1      = Y'*F*Y;
01491                    Teuchos::SerialDenseMatrix<int,ScalarType> FYtmp( Teuchos::View, *FY_, (numBlocks_+recycleBlocks_), recycleBlocks_ );
01492                    FYtmp.multiply(Teuchos::NO_TRANS,Teuchos::NO_TRANS,one,Ftmp,Ytmp,zero);
01493                    //Teuchos::SerialDenseMatrix<int,ScalarType> AU1TU1tmp( Teuchos::View, *AU1TU1_, recycleBlocks_, recycleBlocks_ );
01494                    AU1TU1tmp.multiply(Teuchos::TRANS,Teuchos::NO_TRANS,one,Ytmp,FYtmp,zero);
01495 
01496                    // Indicate the size of the P, Beta structures generated this cycle
01497                    lastp = numBlocks_+1;
01498                    lastBeta = numBlocks_;
01499 
01500                 } // if (cycle != 1)
01501               } // if (!existU_)
01502               else { // U exists
01503                 if (cycle == 0) { // No U1, but U exists
01504                    Teuchos::SerialDenseMatrix<int,ScalarType> Alphatmp( Teuchos::View, *Alpha_, numBlocks_, 1 );
01505                    Teuchos::SerialDenseMatrix<int,ScalarType> Betatmp( Teuchos::View, *Beta_, numBlocks_, 1 );
01506                    Teuchos::SerialDenseMatrix<int,ScalarType> Dtmp( Teuchos::View, *D_, numBlocks_, 1 );
01507                    Teuchos::SerialDenseMatrix<int,ScalarType> APTAPtmp( Teuchos::View, *APTAP_, numBlocks_, numBlocks_ );
01508                    APTAPtmp.putScalar(zero);
01509                    for (int ii=0; ii<numBlocks_; ii++) {
01510                      APTAPtmp(ii,ii) = (Dtmp(ii,0) / Alphatmp(ii,0))*(1 + Betatmp(ii,0));
01511                      if (ii > 0) {
01512                        APTAPtmp(ii-1,ii) = -Dtmp(ii,0)/Alphatmp(ii-1,0);
01513                        APTAPtmp(ii,ii-1) = -Dtmp(ii,0)/Alphatmp(ii-1,0);
01514                      }
01515                    }
01516 
01517                    Teuchos::SerialDenseMatrix<int,ScalarType> L2tmp( Teuchos::View, *L2_, numBlocks_+1, numBlocks_ );
01518                    L2tmp.putScalar(zero);
01519                    for(int ii=0;ii<numBlocks_;ii++) {
01520                      L2tmp(ii,ii)   = 1./Alphatmp(ii,0);
01521                      L2tmp(ii+1,ii) = -1./Alphatmp(ii,0);
01522                    }
01523 
01524                    // AUTAP = UTAU*Delta*L2;
01525                    Teuchos::SerialDenseMatrix<int,ScalarType> AUTAPtmp( Teuchos::View, *AUTAP_, recycleBlocks_, numBlocks_ );
01526                    Teuchos::SerialDenseMatrix<int,ScalarType> UTAUtmp( Teuchos::View, *UTAU_, recycleBlocks_, recycleBlocks_ );
01527                    Teuchos::SerialDenseMatrix<int,ScalarType> Deltatmp( Teuchos::View, *Delta_, recycleBlocks_, numBlocks_+1 );
01528                    Teuchos::SerialDenseMatrix<int,ScalarType> DeltaL2tmp( Teuchos::View, *DeltaL2_, recycleBlocks_, numBlocks_ );
01529                    DeltaL2tmp.multiply(Teuchos::NO_TRANS,Teuchos::NO_TRANS,one,Deltatmp,L2tmp,zero);
01530                    AUTAPtmp.multiply(Teuchos::NO_TRANS,Teuchos::NO_TRANS,one,UTAUtmp,DeltaL2tmp,zero);
01531 
01532                    // F = [UTAU zeros(k,m); zeros(m,k) diag(D)];
01533                    Teuchos::SerialDenseMatrix<int,ScalarType> Ftmp( Teuchos::View, *F_, (numBlocks_+recycleBlocks_), (numBlocks_+recycleBlocks_) );
01534                    Teuchos::SerialDenseMatrix<int,ScalarType> F11( Teuchos::View, *F_, recycleBlocks_, recycleBlocks_ );
01535                    Teuchos::SerialDenseMatrix<int,ScalarType> F12( Teuchos::View, *F_, recycleBlocks_, numBlocks_, 0, recycleBlocks_ );
01536                    Teuchos::SerialDenseMatrix<int,ScalarType> F21( Teuchos::View, *F_, numBlocks_, recycleBlocks_, recycleBlocks_, 0 );
01537                    Teuchos::SerialDenseMatrix<int,ScalarType> F22( Teuchos::View, *F_, numBlocks_, numBlocks_, recycleBlocks_, recycleBlocks_ );
01538                    F11.assign(UTAUtmp);
01539                    F12.putScalar(zero);
01540                    F21.putScalar(zero);
01541                    F22.putScalar(zero);
01542                    for(int ii=0;ii<numBlocks_;ii++) {
01543                      F22(ii,ii) = Dtmp(ii,0);
01544                    }
01545 
01546                    // G = [AUTAU AUTAP; AUTAP' APTAP];
01547                    Teuchos::SerialDenseMatrix<int,ScalarType> Gtmp( Teuchos::View, *G_, (numBlocks_+recycleBlocks_), (numBlocks_+recycleBlocks_) );
01548                    Teuchos::SerialDenseMatrix<int,ScalarType> G11( Teuchos::View, *G_, recycleBlocks_, recycleBlocks_ );
01549                    Teuchos::SerialDenseMatrix<int,ScalarType> G12( Teuchos::View, *G_, recycleBlocks_, numBlocks_, 0, recycleBlocks_ );
01550                    Teuchos::SerialDenseMatrix<int,ScalarType> G21( Teuchos::View, *G_, numBlocks_, recycleBlocks_, recycleBlocks_, 0 );
01551                    Teuchos::SerialDenseMatrix<int,ScalarType> G22( Teuchos::View, *G_, numBlocks_, numBlocks_, recycleBlocks_, recycleBlocks_ );
01552                    Teuchos::SerialDenseMatrix<int,ScalarType> AUTAUtmp( Teuchos::View, *AUTAU_, recycleBlocks_, recycleBlocks_ );
01553                    G11.assign(AUTAUtmp);
01554                    G12.assign(AUTAPtmp);
01555                    // G21 = G12'; (no transpose operator exists for SerialDenseMatrix; Do copy manually)
01556                    for (int ii=0;ii<recycleBlocks_;++ii)
01557                      for (int jj=0;jj<numBlocks_;++jj)
01558                        G21(jj,ii) = G12(ii,jj);
01559                    G22.assign(APTAPtmp);
01560 
01561                    // compute harmonic Ritz vectors
01562                    Teuchos::SerialDenseMatrix<int,ScalarType> Ytmp( Teuchos::View, *Y_, (recycleBlocks_+numBlocks_), recycleBlocks_ );
01563                    getHarmonicVecs(Ftmp,Gtmp,Ytmp);
01564 
01565                    // U1 = [U P(:,1:end-1)]*Y;
01566                    index.resize( recycleBlocks_ );
01567                    for (int ii=0; ii<(recycleBlocks_); ++ii) { index[ii] = ii; }
01568                    Teuchos::RCP<const MV> Utmp = MVT::CloneView( *U_,  index );
01569                    index.resize( numBlocks_ );
01570                    for (int ii=0; ii<numBlocks_; ++ii) { index[ii] = ii; }
01571                    Teuchos::RCP<const MV> Ptmp = MVT::CloneView( *P_,  index );
01572                    index.resize( recycleBlocks_ );
01573                    for (int ii=0; ii<recycleBlocks_; ++ii) { index[ii] = ii; }
01574                    Teuchos::RCP<MV> PY2tmp = MVT::CloneViewNonConst( *PY2_,  index );
01575                    Teuchos::SerialDenseMatrix<int,ScalarType> Y2( Teuchos::View, *Y_, numBlocks_, recycleBlocks_, recycleBlocks_, 0 );
01576                    index.resize( recycleBlocks_ );
01577                    for (int ii=0; ii<recycleBlocks_; ++ii) { index[ii] = ii; }
01578                    Teuchos::RCP<MV> UY1tmp  = MVT::CloneViewNonConst( *U1Y1_,  index );
01579                    Teuchos::SerialDenseMatrix<int,ScalarType> Y1( Teuchos::View, *Y_, recycleBlocks_, recycleBlocks_ );
01580                    index.resize( recycleBlocks_ );
01581                    for (int ii=0; ii<recycleBlocks_; ++ii) { index[ii] = ii; }
01582                    Teuchos::RCP<MV> U1tmp  = MVT::CloneViewNonConst( *U1_,  index );
01583                    MVT::MvTimesMatAddMv( one, *Ptmp, Y2, zero, *PY2tmp );
01584                    MVT::MvTimesMatAddMv( one, *Utmp, Y1, zero, *UY1tmp );
01585                    MVT::MvAddMv(one,*UY1tmp, one, *PY2tmp, *U1tmp);
01586 
01587                    // Precompute some variables for next cycle
01588 
01589                    // AU1TAU1     = Y'*G*Y;
01590                    Teuchos::SerialDenseMatrix<int,ScalarType> GYtmp( Teuchos::View, *GY_, (numBlocks_+recycleBlocks_), recycleBlocks_ );
01591                    GYtmp.multiply(Teuchos::NO_TRANS,Teuchos::NO_TRANS,one,Gtmp,Ytmp,zero);
01592                    Teuchos::SerialDenseMatrix<int,ScalarType> AU1TAU1tmp( Teuchos::View, *AU1TAU1_, recycleBlocks_, recycleBlocks_ );
01593                    AU1TAU1tmp.multiply(Teuchos::TRANS,Teuchos::NO_TRANS,one,Ytmp,GYtmp,zero);
01594 
01595                    // AU1TU1      = Y'*F*Y;
01596                    Teuchos::SerialDenseMatrix<int,ScalarType> FYtmp( Teuchos::View, *FY_, (numBlocks_+recycleBlocks_), recycleBlocks_ );
01597                    FYtmp.multiply(Teuchos::NO_TRANS,Teuchos::NO_TRANS,one,Ftmp,Ytmp,zero);
01598                    Teuchos::SerialDenseMatrix<int,ScalarType> AU1TU1tmp( Teuchos::View, *AU1TU1_, recycleBlocks_, recycleBlocks_ );
01599                    AU1TU1tmp.multiply(Teuchos::TRANS,Teuchos::NO_TRANS,one,Ytmp,FYtmp,zero);
01600 
01601                    // AU1TU   = UTAU;
01602                    Teuchos::SerialDenseMatrix<int,ScalarType> AU1TUtmp( Teuchos::View, *AU1TU_, recycleBlocks_, recycleBlocks_ );
01603                    AU1TUtmp.assign(UTAUtmp);
01604 
01605                    // dold    = D(end);
01606                    dold = Dtmp(numBlocks_-1,0);
01607 
01608                    // indicate that updated recycle space now defined
01609                    existU1_ = true;
01610 
01611                    // Indicate the size of the P, Beta structures generated this cycle
01612                    lastp = numBlocks_;
01613                    lastBeta = numBlocks_-1;
01614                 }
01615                 else { // Have U and U1
01616                    Teuchos::SerialDenseMatrix<int,ScalarType> Alphatmp( Teuchos::View, *Alpha_, numBlocks_, 1 );
01617                    Teuchos::SerialDenseMatrix<int,ScalarType> Betatmp( Teuchos::View, *Beta_, numBlocks_+1, 1 );
01618                    Teuchos::SerialDenseMatrix<int,ScalarType> Dtmp( Teuchos::View, *D_, numBlocks_, 1 );
01619                    Teuchos::SerialDenseMatrix<int,ScalarType> APTAPtmp( Teuchos::View, *APTAP_, numBlocks_, numBlocks_ );
01620                    for (int ii=0; ii<numBlocks_; ii++) {
01621                      APTAPtmp(ii,ii) = (Dtmp(ii,0) / Alphatmp(ii,0))*(1 + Betatmp(ii+1,0));
01622                      if (ii > 0) {
01623                        APTAPtmp(ii-1,ii) = -Dtmp(ii,0)/Alphatmp(ii-1,0);
01624                        APTAPtmp(ii,ii-1) = -Dtmp(ii,0)/Alphatmp(ii-1,0);
01625                      }
01626                    }
01627 
01628                    Teuchos::SerialDenseMatrix<int,ScalarType> L2tmp( Teuchos::View, *L2_, numBlocks_+1, numBlocks_ );
01629                    for(int ii=0;ii<numBlocks_;ii++) {
01630                      L2tmp(ii,ii)   = 1./Alphatmp(ii,0);
01631                      L2tmp(ii+1,ii) = -1./Alphatmp(ii,0);
01632                    }
01633 
01634                    // M(end,1) = dold*(-Beta(1)/Alpha(1));
01635                    // AU1TAP = Y'*[AU1TU*Delta*L2; M];
01636                    Teuchos::SerialDenseMatrix<int,ScalarType> DeltaL2( Teuchos::View, *DeltaL2_, recycleBlocks_, numBlocks_ );
01637                    Teuchos::SerialDenseMatrix<int,ScalarType> Deltatmp( Teuchos::View, *Delta_, recycleBlocks_, numBlocks_+1 );
01638                    DeltaL2.multiply(Teuchos::NO_TRANS,Teuchos::NO_TRANS,one,Deltatmp,L2tmp,zero);
01639                    Teuchos::SerialDenseMatrix<int,ScalarType> AU1TUDeltaL2( Teuchos::View, *AU1TUDeltaL2_, recycleBlocks_, numBlocks_ );
01640                    Teuchos::SerialDenseMatrix<int,ScalarType> AU1TUtmp( Teuchos::View, *AU1TU_, recycleBlocks_, recycleBlocks_ );
01641                    AU1TUDeltaL2.multiply(Teuchos::NO_TRANS,Teuchos::NO_TRANS,one,AU1TUtmp,DeltaL2,zero);
01642                    Teuchos::SerialDenseMatrix<int,ScalarType> Y1( Teuchos::View, *Y_, recycleBlocks_, recycleBlocks_ );
01643                    Teuchos::SerialDenseMatrix<int,ScalarType> AU1TAPtmp( Teuchos::View, *AU1TAP_, recycleBlocks_, numBlocks_ );
01644                    AU1TAPtmp.putScalar(zero);
01645                    AU1TAPtmp.multiply(Teuchos::TRANS,Teuchos::NO_TRANS,one,Y1,AU1TUDeltaL2,zero);
01646                    Teuchos::SerialDenseMatrix<int,ScalarType> Y2( Teuchos::View, *Y_, numBlocks_, recycleBlocks_, recycleBlocks_, 0 );
01647                    ScalarType val = dold * (-Betatmp(0,0)/Alphatmp(0,0));
01648                    for(int ii=0;ii<recycleBlocks_;ii++) {
01649                      AU1TAPtmp(ii,0) += Y2(numBlocks_-1,ii)*val;
01650                    }
01651 
01652                    // AU1TU = Y1'*AU1TU
01653                    Teuchos::SerialDenseMatrix<int,ScalarType> Y1TAU1TU( Teuchos::View, *GY_, recycleBlocks_, recycleBlocks_ );
01654                    Y1TAU1TU.multiply(Teuchos::TRANS,Teuchos::NO_TRANS,one,Y1,AU1TUtmp,zero);
01655                    AU1TUtmp.assign(Y1TAU1TU);
01656 
01657                    // F = [AU1TU1 zeros(k,m); zeros(m,k) diag(D)];
01658                    Teuchos::SerialDenseMatrix<int,ScalarType> Ftmp( Teuchos::View, *F_, (numBlocks_+recycleBlocks_), (numBlocks_+recycleBlocks_) );
01659                    Teuchos::SerialDenseMatrix<int,ScalarType> F11( Teuchos::View, *F_, recycleBlocks_, recycleBlocks_ );
01660                    Teuchos::SerialDenseMatrix<int,ScalarType> F12( Teuchos::View, *F_, recycleBlocks_, numBlocks_, 0, recycleBlocks_ );
01661                    Teuchos::SerialDenseMatrix<int,ScalarType> F21( Teuchos::View, *F_, numBlocks_, recycleBlocks_, recycleBlocks_, 0 );
01662                    Teuchos::SerialDenseMatrix<int,ScalarType> F22( Teuchos::View, *F_, numBlocks_, numBlocks_, recycleBlocks_, recycleBlocks_ );
01663                    Teuchos::SerialDenseMatrix<int,ScalarType> AU1TU1tmp( Teuchos::View, *AU1TU1_, recycleBlocks_, recycleBlocks_ );
01664                    F11.assign(AU1TU1tmp);
01665                    for(int ii=0;ii<numBlocks_;ii++) {
01666                      F22(ii,ii) = Dtmp(ii,0);
01667                    }
01668 
01669                    // G = [AU1TAU1 AU1TAP; AU1TAP' APTAP];
01670                    Teuchos::SerialDenseMatrix<int,ScalarType> Gtmp( Teuchos::View, *G_, (numBlocks_+recycleBlocks_), (numBlocks_+recycleBlocks_) );
01671                    Teuchos::SerialDenseMatrix<int,ScalarType> G11( Teuchos::View, *G_, recycleBlocks_, recycleBlocks_ );
01672                    Teuchos::SerialDenseMatrix<int,ScalarType> G12( Teuchos::View, *G_, recycleBlocks_, numBlocks_, 0, recycleBlocks_ );
01673                    Teuchos::SerialDenseMatrix<int,ScalarType> G21( Teuchos::View, *G_, numBlocks_, recycleBlocks_, recycleBlocks_, 0 );
01674                    Teuchos::SerialDenseMatrix<int,ScalarType> G22( Teuchos::View, *G_, numBlocks_, numBlocks_, recycleBlocks_, recycleBlocks_ );
01675                    Teuchos::SerialDenseMatrix<int,ScalarType> AU1TAU1tmp( Teuchos::View, *AU1TAU1_, recycleBlocks_, recycleBlocks_ );
01676                    G11.assign(AU1TAU1tmp);
01677                    G12.assign(AU1TAPtmp);
01678                    // G21 = G12'; (no transpose operator exists for SerialDenseMatrix; Do copy manually)
01679                    for (int ii=0;ii<recycleBlocks_;++ii)
01680                      for (int jj=0;jj<numBlocks_;++jj)
01681                        G21(jj,ii) = G12(ii,jj);
01682                    G22.assign(APTAPtmp);
01683 
01684                    // compute harmonic Ritz vectors
01685                    Teuchos::SerialDenseMatrix<int,ScalarType> Ytmp( Teuchos::View, *Y_, (recycleBlocks_+numBlocks_), recycleBlocks_ );
01686                    getHarmonicVecs(Ftmp,Gtmp,Ytmp);
01687 
01688                    // U1 = [U1 P(:,2:end-1)]*Y;
01689                    index.resize( numBlocks_ );
01690                    for (int ii=0; ii<numBlocks_; ++ii) { index[ii] = ii+1; }
01691                    Teuchos::RCP<const MV> Ptmp = MVT::CloneView( *P_,  index );
01692                    index.resize( recycleBlocks_ );
01693                    for (int ii=0; ii<recycleBlocks_; ++ii) { index[ii] = ii; }
01694                    Teuchos::RCP<MV> PY2tmp = MVT::CloneViewNonConst( *PY2_,  index );
01695                    index.resize( recycleBlocks_ );
01696                    for (int ii=0; ii<recycleBlocks_; ++ii) { index[ii] = ii; }
01697                    Teuchos::RCP<MV> U1tmp  = MVT::CloneViewNonConst( *U1_,  index );
01698                    index.resize( recycleBlocks_ );
01699                    for (int ii=0; ii<recycleBlocks_; ++ii) { index[ii] = ii; }
01700                    Teuchos::RCP<MV> U1Y1tmp  = MVT::CloneViewNonConst( *U1Y1_,  index );
01701                    MVT::MvTimesMatAddMv( one, *Ptmp, Y2, zero, *PY2tmp );
01702                    MVT::MvTimesMatAddMv( one, *U1tmp, Y1, zero, *U1Y1tmp );
01703                    MVT::MvAddMv(one,*U1Y1tmp, one, *PY2tmp, *U1tmp);
01704 
01705                    // Precompute some variables for next cycle
01706 
01707                    // AU1TAU1     = Y'*G*Y;
01708                    Teuchos::SerialDenseMatrix<int,ScalarType> GYtmp( Teuchos::View, *GY_, (numBlocks_+recycleBlocks_), recycleBlocks_ );
01709                    GYtmp.multiply(Teuchos::NO_TRANS,Teuchos::NO_TRANS,one,Gtmp,Ytmp,zero);
01710                    AU1TAU1tmp.multiply(Teuchos::TRANS,Teuchos::NO_TRANS,one,Ytmp,GYtmp,zero);
01711 
01712                    // AU1TU1      = Y'*F*Y;
01713                    Teuchos::SerialDenseMatrix<int,ScalarType> FYtmp( Teuchos::View, *FY_, (numBlocks_+recycleBlocks_), recycleBlocks_ );
01714                    FYtmp.multiply(Teuchos::NO_TRANS,Teuchos::NO_TRANS,one,Ftmp,Ytmp,zero);
01715                    AU1TU1tmp.multiply(Teuchos::TRANS,Teuchos::NO_TRANS,one,Ytmp,FYtmp,zero);
01716 
01717                    // dold    = D(end);
01718                    dold = Dtmp(numBlocks_-1,0);
01719 
01720                    // Indicate the size of the P, Beta structures generated this cycle
01721                    lastp = numBlocks_+1;
01722                    lastBeta = numBlocks_;
01723 
01724                 }
01725               }
01726             } // if (recycleBlocks_ > 0)
01727 
01728             // Cleanup after end of cycle
01729 
01730             // P = P(:,end-1:end);
01731             index.resize( 2 );
01732             index[0] = lastp-1; index[1] = lastp;
01733             Teuchos::RCP<const MV> Ptmp2 = MVT::CloneView( *P_,  index );
01734             index[0] = 0;       index[1] = 1;
01735             MVT::SetBlock(*Ptmp2,index,*P_);
01736 
01737             // Beta = Beta(end);
01738             (*Beta_)(0,0) = (*Beta_)(lastBeta,0);
01739 
01740             // Delta = Delta(:,end);
01741             if (existU_) { // Delta only initialized if U exists
01742               Teuchos::SerialDenseMatrix<int,ScalarType> mu1( Teuchos::View, *Delta_, recycleBlocks_, 1, 0, 0 );
01743               Teuchos::SerialDenseMatrix<int,ScalarType> mu2( Teuchos::View, *Delta_, recycleBlocks_, 1, 0, numBlocks_ );
01744               mu1.assign(mu2);
01745             }
01746 
01747             // Now reinitialize state variables for next cycle
01748             newstate.P = Teuchos::null;
01749             index.resize( numBlocks_+1 );
01750             for (int ii=0; ii<(numBlocks_+1); ++ii) { index[ii] = ii+1; }
01751             newstate.P  = MVT::CloneViewNonConst( *P_,  index );
01752 
01753             newstate.Beta = Teuchos::null;
01754             newstate.Beta = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( Teuchos::View, *Beta_, numBlocks_, 1, 1, 0 ) );
01755 
01756             newstate.Delta = Teuchos::null;
01757             newstate.Delta = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( Teuchos::View, *Delta_, recycleBlocks_, numBlocks_, 0, 1 ) );
01758 
01759             newstate.curDim = 1; // We have initialized the first search vector
01760 
01761             // Pass to iteration object
01762             rcg_iter->initialize(newstate);
01763 
01764             // increment cycle count
01765             cycle = cycle + 1;
01766 
01767           }
01769           //
01770           // we returned from iterate(), but none of our status tests Passed.
01771           // something is wrong, and it is probably our fault.
01772           //
01774           else {
01775             TEUCHOS_TEST_FOR_EXCEPTION(true,std::logic_error,
01776                                "Belos::RCGSolMgr::solve(): Invalid return from RCGIter::iterate().");
01777           }
01778         }
01779         catch (const std::exception &e) {
01780           printer_->stream(Errors) << "Error! Caught std::exception in RCGIter::iterate() at iteration "
01781                                    << rcg_iter->getNumIters() << std::endl
01782                                    << e.what() << std::endl;
01783           throw;
01784         }
01785       }
01786 
01787       // Inform the linear problem that we are finished with this block linear system.
01788       problem_->setCurrLS();
01789 
01790       // Update indices for the linear systems to be solved.
01791       numRHS2Solve -= 1;
01792       if ( numRHS2Solve > 0 ) {
01793         currIdx[0]++;
01794         // Set the next indices.
01795         problem_->setLSIndex( currIdx );
01796       }
01797       else {
01798         currIdx.resize( numRHS2Solve );
01799       }
01800 
01801       // Update the recycle space for the next linear system
01802       if (existU1_) { // be sure updated recycle space was created
01803         // U = U1
01804         index.resize(recycleBlocks_);
01805         for (int ii=0; ii<recycleBlocks_; ++ii) { index[ii] = ii; }
01806         MVT::SetBlock(*U1_,index,*U_);
01807         // Set flag indicating recycle space is now defined
01808         existU_ = true;
01809         if (numRHS2Solve > 0) { // also update AU, UTAU, and AUTAU
01810           // Free pointers in newstate
01811           newstate.P = Teuchos::null;
01812           newstate.Ap = Teuchos::null;
01813           newstate.r = Teuchos::null;
01814           newstate.z = Teuchos::null;
01815           newstate.U = Teuchos::null;
01816           newstate.AU = Teuchos::null;
01817           newstate.Alpha = Teuchos::null;
01818           newstate.Beta = Teuchos::null;
01819           newstate.D = Teuchos::null;
01820           newstate.Delta = Teuchos::null;
01821           newstate.LUUTAU = Teuchos::null;
01822           newstate.ipiv = Teuchos::null;
01823           newstate.rTz_old = Teuchos::null;
01824 
01825           // Reinitialize AU, UTAU, AUTAU
01826           index.resize(recycleBlocks_);
01827           for (int i=0; i<recycleBlocks_; ++i) { index[i] = i; }
01828           Teuchos::RCP<const MV> Utmp  = MVT::CloneView( *U_,  index );
01829           index.resize(recycleBlocks_);
01830           for (int i=0; i<recycleBlocks_; ++i) { index[i] = i; }
01831           Teuchos::RCP<MV> AUtmp = MVT::CloneViewNonConst( *AU_, index );
01832           // Initialize AU
01833           problem_->applyOp( *Utmp, *AUtmp );
01834           // Initialize UTAU
01835           Teuchos::SerialDenseMatrix<int,ScalarType> UTAUtmp( Teuchos::View, *UTAU_, recycleBlocks_, recycleBlocks_ );
01836           MVT::MvTransMv( one, *Utmp, *AUtmp, UTAUtmp );
01837           // Initialize AUTAU  ( AUTAU = AU'*(M\AU) )
01838           Teuchos::SerialDenseMatrix<int,ScalarType> AUTAUtmp( Teuchos::View, *AUTAU_, recycleBlocks_, recycleBlocks_ );
01839           if ( precObj != Teuchos::null ) {
01840             index.resize(recycleBlocks_);
01841             for (int i=0; i<recycleBlocks_; ++i) { index[i] = i; }
01842             index.resize(recycleBlocks_);
01843             for (int ii=0; ii<recycleBlocks_; ++ii) { index[ii] = ii; }
01844             Teuchos::RCP<MV> LeftPCAU = MVT::CloneViewNonConst( *U1_, index ); // use U1 as temp storage
01845             OPT::Apply( *precObj, *AUtmp, *LeftPCAU );
01846             MVT::MvTransMv( one, *AUtmp, *LeftPCAU, AUTAUtmp );
01847           } else {
01848             MVT::MvTransMv( one, *AUtmp, *AUtmp, AUTAUtmp );
01849           }
01850         } // if (numRHS2Solve > 0)
01851 
01852       } // if (existU1)
01853     }// while ( numRHS2Solve > 0 )
01854 
01855   }
01856 
01857   // print final summary
01858   sTest_->print( printer_->stream(FinalSummary) );
01859 
01860   // print timing information
01861 #ifdef BELOS_TEUCHOS_TIME_MONITOR
01862   // Calling summarize() can be expensive, so don't call unless the
01863   // user wants to print out timing details.  summarize() will do all
01864   // the work even if it's passed a "black hole" output stream.
01865   if (verbosity_ & TimingDetails)
01866     Teuchos::TimeMonitor::summarize( printer_->stream(TimingDetails) );
01867 #endif
01868 
01869   // get iteration information for this solve
01870   numIters_ = maxIterTest_->getNumIters();
01871 
01872   // Save the convergence test value ("achieved tolerance") for this solve.
01873   {
01874     using Teuchos::rcp_dynamic_cast;
01875     typedef StatusTestGenResNorm<ScalarType,MV,OP> conv_test_type;
01876     // testValues is nonnull and not persistent.
01877     const std::vector<MagnitudeType>* pTestValues =
01878       rcp_dynamic_cast<conv_test_type>(convTest_)->getTestValue();
01879 
01880     TEUCHOS_TEST_FOR_EXCEPTION(pTestValues == NULL, std::logic_error,
01881       "Belos::RCGSolMgr::solve(): The convergence test's getTestValue() "
01882       "method returned NULL.  Please report this bug to the Belos developers.");
01883 
01884     TEUCHOS_TEST_FOR_EXCEPTION(pTestValues->size() < 1, std::logic_error,
01885       "Belos::RCGSolMgr::solve(): The convergence test's getTestValue() "
01886       "method returned a vector of length zero.  Please report this bug to the "
01887       "Belos developers.");
01888 
01889     // FIXME (mfh 12 Dec 2011) Does pTestValues really contain the
01890     // achieved tolerances for all vectors in the current solve(), or
01891     // just for the vectors from the last deflation?
01892     achievedTol_ = *std::max_element (pTestValues->begin(), pTestValues->end());
01893   }
01894 
01895   if (!isConverged) {
01896     return Unconverged; // return from RCGSolMgr::solve()
01897   }
01898   return Converged; // return from RCGSolMgr::solve()
01899 }
01900 
01901 //  Compute the harmonic eigenpairs of the projected, dense system.
01902 template<class ScalarType, class MV, class OP>
01903 void RCGSolMgr<ScalarType,MV,OP,false>::getHarmonicVecs(const Teuchos::SerialDenseMatrix<int,ScalarType>& F,
01904                                                   const Teuchos::SerialDenseMatrix<int,ScalarType>& G,
01905                                                         Teuchos::SerialDenseMatrix<int,ScalarType>& Y) {
01906   // order of F,G
01907   int n = F.numCols();
01908 
01909   // The LAPACK interface
01910   Teuchos::LAPACK<int,ScalarType> lapack;
01911 
01912   // Magnitude of harmonic Ritz values
01913   std::vector<MagnitudeType> w(n);
01914 
01915   // Sorted order of harmonic Ritz values
01916   std::vector<int> iperm(n);
01917 
01918   // Compute k smallest harmonic Ritz pairs
01919   // SUBROUTINE DSYGV( ITYPE, JOBZ, UPLO, N, A, LDA, B, LDB, W, WORK, LWORK, INFO )
01920   int itype = 1; // solve A*x = (lambda)*B*x
01921   char jobz='V'; // compute eigenvalues and eigenvectors
01922   char uplo='U'; // since F,G symmetric, reference only their upper triangular data
01923   std::vector<ScalarType> work(1);
01924   int lwork = -1;
01925   int info = 0;
01926   // since SYGV destroys workspace, create copies of F,G
01927   Teuchos::SerialDenseMatrix<int,ScalarType> F2( Teuchos::Copy, *F_, F_->numRows(), F_->numCols() );
01928   Teuchos::SerialDenseMatrix<int,ScalarType> G2( Teuchos::Copy, *G_, G_->numRows(), G_->numCols() );
01929 
01930   // query for optimal workspace size
01931   lapack.SYGV(itype, jobz, uplo, n, G2.values(), G2.stride(), F2.values(), F2.stride(), &w[0], &work[0], lwork, &info);
01932   TEUCHOS_TEST_FOR_EXCEPTION(info != 0, RCGSolMgrLAPACKFailure,
01933                      "Belos::RCGSolMgr::solve(): LAPACK SYGV failed to query optimal work size.");
01934   lwork = (int)work[0];
01935   work.resize(lwork);
01936   lapack.SYGV(itype, jobz, uplo, n, G2.values(), G2.stride(), F2.values(), F2.stride(), &w[0], &work[0], lwork, &info);
01937   TEUCHOS_TEST_FOR_EXCEPTION(info != 0, RCGSolMgrLAPACKFailure,
01938                      "Belos::RCGSolMgr::solve(): LAPACK SYGV failed to compute eigensolutions.");
01939 
01940 
01941   // Construct magnitude of each harmonic Ritz value
01942   this->sort(w,n,iperm);
01943 
01944   // Select recycledBlocks_ smallest eigenvectors
01945   for( int i=0; i<recycleBlocks_; i++ ) {
01946     for( int j=0; j<n; j++ ) {
01947       Y(j,i) = G2(j,iperm[i]);
01948     }
01949   }
01950 
01951 }
01952 
01953 // This method sorts list of n floating-point numbers and return permutation vector
01954 template<class ScalarType, class MV, class OP>
01955 void RCGSolMgr<ScalarType,MV,OP,false>::sort(std::vector<ScalarType>& dlist, int n, std::vector<int>& iperm)
01956 {
01957   int l, r, j, i, flag;
01958   int    RR2;
01959   double dRR, dK;
01960 
01961   // Initialize the permutation vector.
01962   for(j=0;j<n;j++)
01963     iperm[j] = j;
01964 
01965   if (n <= 1) return;
01966 
01967   l    = n / 2 + 1;
01968   r    = n - 1;
01969   l    = l - 1;
01970   dRR  = dlist[l - 1];
01971   dK   = dlist[l - 1];
01972 
01973   RR2 = iperm[l - 1];
01974   while (r != 0) {
01975     j = l;
01976     flag = 1;
01977 
01978     while (flag == 1) {
01979       i = j;
01980       j = j + j;
01981 
01982       if (j > r + 1)
01983         flag = 0;
01984       else {
01985         if (j < r + 1)
01986           if (dlist[j] > dlist[j - 1]) j = j + 1;
01987 
01988         if (dlist[j - 1] > dK) {
01989           dlist[i - 1] = dlist[j - 1];
01990           iperm[i - 1] = iperm[j - 1];
01991         }
01992         else {
01993           flag = 0;
01994         }
01995       }
01996     }
01997     dlist[i - 1] = dRR;
01998     iperm[i - 1] = RR2;
01999     if (l == 1) {
02000       dRR  = dlist [r];
02001       RR2 = iperm[r];
02002       dK = dlist[r];
02003       dlist[r] = dlist[0];
02004       iperm[r] = iperm[0];
02005       r = r - 1;
02006     }
02007     else {
02008       l   = l - 1;
02009       dRR  = dlist[l - 1];
02010       RR2  = iperm[l - 1];
02011       dK   = dlist[l - 1];
02012     }
02013   }
02014   dlist[0] = dRR;
02015   iperm[0] = RR2;
02016 }
02017 
02018 //  This method requires the solver manager to return a std::string that describes itself.
02019 template<class ScalarType, class MV, class OP>
02020 std::string RCGSolMgr<ScalarType,MV,OP,false>::description() const
02021 {
02022   std::ostringstream oss;
02023   oss << "Belos::RCGSolMgr<...,"<<Teuchos::ScalarTraits<ScalarType>::name()<<">";
02024   return oss.str();
02025 }
02026 
02027 } // end Belos namespace
02028 
02029 #endif /* BELOS_RCG_SOLMGR_HPP */
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines