BelosBlockGmresIter.hpp

Go to the documentation of this file.
00001 // @HEADER
00002 // ***********************************************************************
00003 //
00004 //                 Belos: Block Linear Solvers Package
00005 //                 Copyright (2004) Sandia Corporation
00006 //
00007 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
00008 // license for use of this work by or on behalf of the U.S. Government.
00009 //
00010 // This library is free software; you can redistribute it and/or modify
00011 // it under the terms of the GNU Lesser General Public License as
00012 // published by the Free Software Foundation; either version 2.1 of the
00013 // License, or (at your option) any later version.
00014 //
00015 // This library is distributed in the hope that it will be useful, but
00016 // WITHOUT ANY WARRANTY; without even the implied warranty of
00017 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00018 // Lesser General Public License for more details.
00019 //
00020 // You should have received a copy of the GNU Lesser General Public
00021 // License along with this library; if not, write to the Free Software
00022 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307
00023 // USA
00024 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
00025 //
00026 // ***********************************************************************
00027 // @HEADER
00028 
00029 #ifndef BELOS_BLOCK_GMRES_ITER_HPP
00030 #define BELOS_BLOCK_GMRES_ITER_HPP
00031 
00036 #include "BelosConfigDefs.hpp"
00037 #include "BelosTypes.hpp"
00038 #include "BelosGmresIteration.hpp"
00039 
00040 #include "BelosLinearProblem.hpp"
00041 #include "BelosMatOrthoManager.hpp"
00042 #include "BelosOutputManager.hpp"
00043 #include "BelosStatusTest.hpp"
00044 #include "BelosOperatorTraits.hpp"
00045 #include "BelosMultiVecTraits.hpp"
00046 
00047 #include "Teuchos_BLAS.hpp"
00048 #include "Teuchos_SerialDenseMatrix.hpp"
00049 #include "Teuchos_SerialDenseVector.hpp"
00050 #include "Teuchos_ScalarTraits.hpp"
00051 #include "Teuchos_ParameterList.hpp"
00052 #include "Teuchos_TimeMonitor.hpp"
00053 
00067 namespace Belos {
00068   
00069 template<class ScalarType, class MV, class OP>
00070 class BlockGmresIter : virtual public GmresIteration<ScalarType,MV,OP> {
00071 
00072   public:
00073     
00074   //
00075   // Convenience typedefs
00076   //
00077   typedef MultiVecTraits<ScalarType,MV> MVT;
00078   typedef OperatorTraits<ScalarType,MV,OP> OPT;
00079   typedef Teuchos::ScalarTraits<ScalarType> SCT;
00080   typedef typename SCT::magnitudeType MagnitudeType;
00081 
00083 
00084 
00094   BlockGmresIter( const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem, 
00095       const Teuchos::RCP<OutputManager<ScalarType> > &printer,
00096       const Teuchos::RCP<StatusTest<ScalarType,MV,OP> > &tester,
00097       const Teuchos::RCP<MatOrthoManager<ScalarType,MV,OP> > &ortho,
00098       Teuchos::ParameterList &params );
00099 
00101   virtual ~BlockGmresIter() {};
00103 
00104 
00106 
00107   
00129   void iterate();
00130 
00152   void initializeGmres(GmresIterationState<ScalarType,MV> newstate);
00153 
00157   void initialize()
00158   {
00159     GmresIterationState<ScalarType,MV> empty;
00160     initializeGmres(empty);
00161   }
00162   
00170   GmresIterationState<ScalarType,MV> getState() const {
00171     GmresIterationState<ScalarType,MV> state;
00172     state.curDim = curDim_;
00173     state.V = V_;
00174     state.H = H_;
00175     state.R = R_;
00176     state.z = z_;
00177     return state;
00178   }
00179 
00181 
00182   
00184 
00185 
00187   int getNumIters() const { return iter_; }
00188   
00190   void resetNumIters( int iter = 0 ) { iter_ = iter; }
00191 
00194   Teuchos::RCP<const MV> getNativeResiduals( std::vector<MagnitudeType> *norms ) const;
00195 
00197 
00202   Teuchos::RCP<MV> getCurrentUpdate() const;
00203 
00205 
00208   void updateLSQR( int dim = -1 );
00209 
00211   int getCurSubspaceDim() const { 
00212     if (!initialized_) return 0;
00213     return curDim_;
00214   };
00215 
00217   int getMaxSubspaceDim() const { return blockSize_*numBlocks_; }
00218 
00220 
00221   
00223 
00224 
00226   const LinearProblem<ScalarType,MV,OP>& getProblem() const { return *lp_; }
00227 
00229   int getBlockSize() const { return blockSize_; }
00230   
00232   void setBlockSize(int blockSize) { setSize( blockSize, numBlocks_ ); }
00233   
00235   int getNumBlocks() const { return numBlocks_; }
00236   
00238   void setNumBlocks(int numBlocks) { setSize( blockSize_, numBlocks ); }
00239   
00246   void setSize(int blockSize, int numBlocks);
00247 
00249   bool isInitialized() { return initialized_; }
00250 
00252 
00253   private:
00254 
00255   //
00256   // Internal methods
00257   //
00259   void setStateSize();
00260   
00261   //
00262   // Classes inputed through constructor that define the linear problem to be solved.
00263   //
00264   const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> >    lp_;
00265   const Teuchos::RCP<OutputManager<ScalarType> >          om_;
00266   const Teuchos::RCP<StatusTest<ScalarType,MV,OP> >       stest_;
00267   const Teuchos::RCP<OrthoManager<ScalarType,MV> >        ortho_;
00268 
00269   //
00270   // Algorithmic parameters
00271   //  
00272   // blockSize_ is the solver block size.
00273   // It controls the number of vectors added to the basis on each iteration.
00274   int blockSize_;
00275   // numBlocks_ is the size of the allocated space for the Krylov basis, in blocks.
00276   int numBlocks_; 
00277   
00278   // Storage for QR factorization of the least squares system.
00279   Teuchos::SerialDenseVector<int,ScalarType> beta, sn;
00280   Teuchos::SerialDenseVector<int,MagnitudeType> cs;
00281   
00282   // 
00283   // Current solver state
00284   //
00285   // initialized_ specifies that the basis vectors have been initialized and the iterate() routine
00286   // is capable of running; _initialize is controlled  by the initialize() member method
00287   // For the implications of the state of initialized_, please see documentation for initialize()
00288   bool initialized_;
00289 
00290   // stateStorageInitialized_ specifies that the state storage has be initialized to the current
00291   // blockSize_ and numBlocks_.  This initialization may be postponed if the linear problem was
00292   // generated without the right-hand side or solution vectors.
00293   bool stateStorageInitialized_;
00294 
00295   // keepHessenberg_ specifies that the iteration must keep the Hessenberg matrix formed via the 
00296   // Arnoldi factorization and the upper triangular matrix that is the Hessenberg matrix reduced via
00297   // QR factorization separate.
00298   bool keepHessenberg_;
00299 
00300   // initHessenberg_ specifies that the iteration should reinitialize the Hessenberg matrix by zeroing
00301   // out all entries before an iteration is started.
00302   bool initHessenberg_;
00303  
00304   // Current subspace dimension, and number of iterations performed.
00305   int curDim_, iter_;
00306   
00307   // 
00308   // State Storage
00309   //
00310   Teuchos::RCP<MV> V_;
00311   //
00312   // Projected matrices
00313   // H_ : Projected matrix from the Krylov factorization AV = VH + FE^T
00314   //
00315   Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > H_;
00316   // 
00317   // QR decomposition of Projected matrices for solving the least squares system HY = B.
00318   // R_: Upper triangular reduction of H
00319   // z_: Q applied to right-hand side of the least squares system
00320   Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > R_;
00321   Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > z_;  
00322 };
00323 
00325   // Constructor.
00326   template<class ScalarType, class MV, class OP>
00327   BlockGmresIter<ScalarType,MV,OP>::BlockGmresIter(const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem, 
00328                const Teuchos::RCP<OutputManager<ScalarType> > &printer,
00329                const Teuchos::RCP<StatusTest<ScalarType,MV,OP> > &tester,
00330                const Teuchos::RCP<MatOrthoManager<ScalarType,MV,OP> > &ortho,
00331                Teuchos::ParameterList &params ):
00332     lp_(problem),
00333     om_(printer),
00334     stest_(tester),
00335     ortho_(ortho),
00336     blockSize_(0),
00337     numBlocks_(0),
00338     initialized_(false),
00339     stateStorageInitialized_(false),
00340     keepHessenberg_(false),
00341     initHessenberg_(false),
00342     curDim_(0),
00343     iter_(0)
00344   {
00345     // Find out whether we are saving the Hessenberg matrix.
00346     keepHessenberg_ = params.get("Keep Hessenberg", false);
00347 
00348     // Find out whether we are initializing the Hessenberg matrix.
00349     initHessenberg_ = params.get("Initialize Hessenberg", false);
00350 
00351     // Get the maximum number of blocks allowed for this Krylov subspace
00352     TEST_FOR_EXCEPTION(!params.isParameter("Num Blocks"), std::invalid_argument,
00353                        "Belos::BlockGmresIter::constructor: mandatory parameter 'Num Blocks' is not specified.");
00354     int nb = Teuchos::getParameter<int>(params, "Num Blocks");
00355 
00356     // Set the block size and allocate data
00357     int bs = params.get("Block Size", 1);
00358     setSize( bs, nb );
00359   }
00360 
00362   // Set the block size and make necessary adjustments.
00363   template <class ScalarType, class MV, class OP>
00364   void BlockGmresIter<ScalarType,MV,OP>::setSize (int blockSize, int numBlocks)
00365   {
00366     // This routine only allocates space; it doesn't not perform any computation
00367     // any change in size will invalidate the state of the solver.
00368 
00369     TEST_FOR_EXCEPTION(numBlocks <= 0 || blockSize <= 0, std::invalid_argument, "Belos::BlockGmresIter::setSize was passed a non-positive argument.");
00370     if (blockSize == blockSize_ && numBlocks == numBlocks_) {
00371       // do nothing
00372       return;
00373     }
00374 
00375     if (blockSize!=blockSize_ || numBlocks!=numBlocks_)
00376       stateStorageInitialized_ = false;
00377 
00378     blockSize_ = blockSize;
00379     numBlocks_ = numBlocks;
00380 
00381     initialized_ = false;
00382     curDim_ = 0;
00383 
00384     // Use the current blockSize_ and numBlocks_ to initialize the state storage.    
00385     setStateSize();
00386     
00387   }
00388 
00390   // Setup the state storage.
00391   template <class ScalarType, class MV, class OP>
00392   void BlockGmresIter<ScalarType,MV,OP>::setStateSize ()
00393   {
00394     if (!stateStorageInitialized_) {
00395 
00396       // Check if there is any multivector to clone from.
00397       Teuchos::RCP<const MV> lhsMV = lp_->getLHS();
00398       Teuchos::RCP<const MV> rhsMV = lp_->getRHS();
00399       if (lhsMV == Teuchos::null && rhsMV == Teuchos::null) {
00400   stateStorageInitialized_ = false;
00401   return;
00402       }
00403       else {
00404   
00406   // blockSize*numBlocks dependent
00407   //
00408   int newsd = blockSize_*(numBlocks_+1);
00409   
00410   if (blockSize_==1) {
00411     cs.resize( newsd );
00412     sn.resize( newsd );
00413   }
00414   else {
00415     beta.resize( newsd );
00416   }
00417   
00418   // Initialize the state storage
00419         TEST_FOR_EXCEPTION(blockSize_*numBlocks_ > MVT::GetVecLength(*rhsMV),std::invalid_argument,
00420                            "Belos::BlockGmresIter::setStateSize(): Cannot generate a Krylov basis with dimension larger the operator!");
00421 
00422   // If the subspace has not be initialized before, generate it using the LHS or RHS from lp_.
00423   if (V_ == Teuchos::null) {
00424     // Get the multivector that is not null.
00425     Teuchos::RCP<const MV> tmp = ( (rhsMV!=Teuchos::null)? rhsMV: lhsMV );
00426     TEST_FOR_EXCEPTION(tmp == Teuchos::null,std::invalid_argument,
00427            "Belos::BlockGmresIter::setStateSize(): linear problem does not specify multivectors to clone from.");
00428     V_ = MVT::Clone( *tmp, newsd );
00429   }
00430   else {
00431     // Generate V_ by cloning itself ONLY if more space is needed.
00432     if (MVT::GetNumberVecs(*V_) < newsd) {
00433       Teuchos::RCP<const MV> tmp = V_;
00434       V_ = MVT::Clone( *tmp, newsd );
00435     }
00436   }
00437   
00438   // Generate R_ only if it doesn't exist, otherwise resize it.
00439   if (R_ == Teuchos::null) {
00440     R_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>() );
00441         }
00442         if (initHessenberg_) {
00443           R_->shape( newsd, newsd-blockSize_ );
00444         }
00445         else {
00446           if (R_->numRows() < newsd || R_->numCols() < newsd-blockSize_) {
00447       R_->shapeUninitialized( newsd, newsd-blockSize_ );
00448     }
00449         }
00450   
00451   // Generate H_ only if it doesn't exist, and we are keeping the upper Hessenberg matrix.
00452   if (keepHessenberg_) {
00453     if (H_ == Teuchos::null) {
00454       H_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>() );
00455           }
00456           if (initHessenberg_) {
00457             H_->shape( newsd, newsd-blockSize_ );
00458           }
00459           else {
00460             if (H_->numRows() < newsd || H_->numCols() < newsd-blockSize_) {
00461         H_->shapeUninitialized( newsd, newsd-blockSize_ );
00462       }
00463           }
00464   }
00465   else {
00466     // Point H_ and R_ at the same object.
00467     H_ = R_;
00468   }
00469   
00470   // Generate z_ only if it doesn't exist, otherwise resize it.
00471   if (z_ == Teuchos::null) {
00472     z_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>() );
00473         }
00474         if (z_-> numRows() < newsd || z_->numCols() < blockSize_) {
00475     z_->shapeUninitialized( newsd, blockSize_ );
00476   }
00477   
00478   // State storage has now been initialized.
00479   stateStorageInitialized_ = true;
00480       }
00481     }
00482   }
00483 
00485   // Get the current update from this subspace.
00486   template <class ScalarType, class MV, class OP>
00487   Teuchos::RCP<MV> BlockGmresIter<ScalarType,MV,OP>::getCurrentUpdate() const
00488   {
00489     //
00490     // If this is the first iteration of the Arnoldi factorization, 
00491     // there is no update, so return Teuchos::null. 
00492     //
00493     RCP<MV> currentUpdate = Teuchos::null;
00494     if (curDim_==0) { 
00495       return currentUpdate; 
00496     } else {
00497       const ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
00498       const ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
00499       Teuchos::BLAS<int,ScalarType> blas;
00500       currentUpdate = MVT::Clone( *V_, blockSize_ );
00501       //
00502       //  Make a view and then copy the RHS of the least squares problem.  DON'T OVERWRITE IT!
00503       //
00504       Teuchos::SerialDenseMatrix<int,ScalarType> y( Teuchos::Copy, *z_, curDim_, blockSize_ );
00505       //
00506       //  Solve the least squares problem.
00507       //
00508       blas.TRSM( Teuchos::LEFT_SIDE, Teuchos::UPPER_TRI, Teuchos::NO_TRANS,
00509      Teuchos::NON_UNIT_DIAG, curDim_, blockSize_, one,  
00510      R_->values(), R_->stride(), y.values(), y.stride() );
00511       //
00512       //  Compute the current update.
00513       //
00514       std::vector<int> index(curDim_);
00515       for ( int i=0; i<curDim_; i++ ) {
00516         index[i] = i;
00517       }
00518       RCP<const MV> Vjp1 = MVT::CloneView( *V_, index );
00519       MVT::MvTimesMatAddMv( one, *Vjp1, y, zero, *currentUpdate );
00520     }
00521     return currentUpdate;
00522   }
00523 
00524 
00526   // Get the native residuals stored in this iteration.  
00527   // Note:  No residual std::vector will be returned by Gmres.
00528   template <class ScalarType, class MV, class OP>
00529   Teuchos::RCP<const MV> BlockGmresIter<ScalarType,MV,OP>::getNativeResiduals( std::vector<MagnitudeType> *norms ) const 
00530   {
00531     //
00532     // NOTE: Make sure the incoming std::vector is the correct size!
00533     //
00534     if ( norms && (int)norms->size() < blockSize_ )                         
00535       norms->resize( blockSize_ );                                          
00536     
00537     if (norms) {
00538       Teuchos::BLAS<int,ScalarType> blas;
00539       for (int j=0; j<blockSize_; j++) {
00540         (*norms)[j] = blas.NRM2( blockSize_, &(*z_)(curDim_, j), 1);
00541       }
00542     }
00543     return Teuchos::null;
00544   }
00545 
00546   
00547 
00549   // Initialize this iteration object
00550   template <class ScalarType, class MV, class OP>
00551   void BlockGmresIter<ScalarType,MV,OP>::initializeGmres(GmresIterationState<ScalarType,MV> newstate)
00552   {
00553     // Initialize the state storage if it isn't already.
00554     if (!stateStorageInitialized_) 
00555       setStateSize();
00556 
00557     TEST_FOR_EXCEPTION(!stateStorageInitialized_,std::invalid_argument,
00558            "Belos::BlockGmresIter::initialize(): Cannot initialize state storage!");
00559     
00560     // NOTE:  In BlockGmresIter, V and Z are required!!!  
00561     // inconsitent multivectors widths and lengths will not be tolerated, and
00562     // will be treated with exceptions.
00563     //
00564     std::string errstr("Belos::BlockGmresIter::initialize(): Specified multivectors must have a consistent length and width.");
00565 
00566     if (newstate.V != Teuchos::null && newstate.z != Teuchos::null) {
00567 
00568       // initialize V_,z_, and curDim_
00569 
00570       TEST_FOR_EXCEPTION( MVT::GetVecLength(*newstate.V) != MVT::GetVecLength(*V_),
00571                           std::invalid_argument, errstr );
00572       TEST_FOR_EXCEPTION( MVT::GetNumberVecs(*newstate.V) < blockSize_,
00573                           std::invalid_argument, errstr );
00574       TEST_FOR_EXCEPTION( newstate.curDim > blockSize_*(numBlocks_+1),
00575                           std::invalid_argument, errstr );
00576 
00577       curDim_ = newstate.curDim;
00578       int lclDim = MVT::GetNumberVecs(*newstate.V);
00579 
00580       // check size of Z
00581       TEST_FOR_EXCEPTION(newstate.z->numRows() < curDim_ || newstate.z->numCols() < blockSize_, std::invalid_argument, errstr);
00582       
00583 
00584       // copy basis vectors from newstate into V
00585       if (newstate.V != V_) {
00586         // only copy over the first block and print a warning.
00587         if (curDim_ == 0 && lclDim > blockSize_) {
00588     om_->stream(Warnings) << "Belos::BlockGmresIter::initialize(): the solver was initialized with a kernel of " << lclDim << std::endl
00589         << "The block size however is only " << blockSize_ << std::endl
00590         << "The last " << lclDim - blockSize_ << " vectors will be discarded." << std::endl;
00591   }
00592         std::vector<int> nevind(curDim_+blockSize_);
00593         for (int i=0; i<curDim_+blockSize_; i++) nevind[i] = i;
00594   Teuchos::RCP<const MV> newV = MVT::CloneView( *newstate.V, nevind );
00595   Teuchos::RCP<MV> lclV = MVT::CloneViewNonConst( *V_, nevind );
00596         MVT::MvAddMv( 1.0, *newV, 0.0, *newV, *lclV );
00597 
00598         // done with local pointers
00599         lclV = Teuchos::null;
00600       }
00601 
00602       // put data into z_, make sure old information is not still hanging around.
00603       if (newstate.z != z_) {
00604         z_->putScalar();
00605         Teuchos::SerialDenseMatrix<int,ScalarType> newZ(Teuchos::View,*newstate.z,curDim_+blockSize_,blockSize_);
00606         Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > lclZ;
00607         lclZ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*z_,curDim_+blockSize_,blockSize_) );
00608         lclZ->assign(newZ);
00609 
00610         // done with local pointers
00611         lclZ = Teuchos::null;
00612       }
00613 
00614     }
00615     else {
00616 
00617       TEST_FOR_EXCEPTION(newstate.V == Teuchos::null,std::invalid_argument,
00618                          "Belos::BlockGmresIter::initialize(): BlockGmresStateIterState does not have initial kernel V_0.");
00619 
00620       TEST_FOR_EXCEPTION(newstate.z == Teuchos::null,std::invalid_argument,
00621                          "Belos::BlockGmresIter::initialize(): BlockGmresStateIterState does not have initial norms z_0.");
00622     }
00623 
00624     // the solver is initialized
00625     initialized_ = true;
00626 
00627     /*
00628     if (om_->isVerbosity( Debug ) ) {
00629       // Check almost everything here
00630       CheckList chk;
00631       chk.checkV = true;
00632       chk.checkArn = true;
00633       chk.checkAux = true;
00634       om_->print( Debug, accuracyCheck(chk, ": after initialize()") );
00635     }
00636     */
00637 
00638   }
00639 
00640 
00642   // Iterate until the status test informs us we should stop.
00643   template <class ScalarType, class MV, class OP>
00644   void BlockGmresIter<ScalarType,MV,OP>::iterate()
00645   {
00646     //
00647     // Allocate/initialize data structures
00648     //
00649     if (initialized_ == false) {
00650       initialize();
00651     }
00652     
00653     // Compute the current search dimension. 
00654     int searchDim = blockSize_*numBlocks_;
00655 
00657     // iterate until the status test tells us to stop.
00658     //
00659     // also break if our basis is full
00660     //
00661     while (stest_->checkStatus(this) != Passed && curDim_+blockSize_ <= searchDim) {
00662 
00663       iter_++;
00664 
00665       // F can be found at the curDim_ block, but the next block is at curDim_ + blockSize_.
00666       int lclDim = curDim_ + blockSize_; 
00667 
00668       // Get the current part of the basis.
00669       std::vector<int> curind(blockSize_);
00670       for (int i=0; i<blockSize_; i++) { curind[i] = lclDim + i; }
00671       Teuchos::RCP<MV> Vnext = MVT::CloneViewNonConst(*V_,curind);
00672 
00673       // Get a view of the previous vectors.
00674       // This is used for orthogonalization and for computing V^H K H
00675       for (int i=0; i<blockSize_; i++) { curind[i] = curDim_ + i; }
00676       Teuchos::RCP<const MV> Vprev = MVT::CloneView(*V_,curind);
00677 
00678       // Compute the next std::vector in the Krylov basis:  Vnext = Op*Vprev
00679       lp_->apply(*Vprev,*Vnext);
00680       Vprev = Teuchos::null;
00681       
00682       // Remove all previous Krylov basis vectors from Vnext      
00683       // Get a view of all the previous vectors
00684       std::vector<int> prevind(lclDim);
00685       for (int i=0; i<lclDim; i++) { prevind[i] = i; }
00686       Vprev = MVT::CloneView(*V_,prevind);
00687       Teuchos::Array<Teuchos::RCP<const MV> > AVprev(1, Vprev);
00688       
00689       // Get a view of the part of the Hessenberg matrix needed to hold the ortho coeffs.
00690       Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> >
00691   subH = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>
00692            ( Teuchos::View,*H_,lclDim,blockSize_,0,curDim_ ) );
00693       Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > AsubH;
00694       AsubH.append( subH );
00695       
00696       // Get a view of the part of the Hessenberg matrix needed to hold the norm coeffs.
00697       Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> >
00698   subH2 = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>
00699             ( Teuchos::View,*H_,blockSize_,blockSize_,lclDim,curDim_ ) );
00700       subH2->putScalar();  // Initialize subdiagonal to zero
00701       int rank = ortho_->projectAndNormalize(*Vnext,AsubH,subH2,AVprev);
00702 
00703       // Copy over the coefficients if we are saving the upper Hessenberg matrix,
00704       // just in case we run into an error.
00705       if (keepHessenberg_) {
00706   // Copy over the orthogonalization coefficients.
00707   Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> >
00708     subR = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>
00709              ( Teuchos::View,*R_,lclDim,blockSize_,0,curDim_ ) );
00710   subR->assign(*subH);
00711   
00712   // Copy over the lower diagonal block of the Hessenberg matrix.
00713   Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> >
00714     subR2 = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>
00715         ( Teuchos::View,*R_,blockSize_,blockSize_,lclDim,curDim_ ) );
00716   subR2->assign(*subH2);
00717       }
00718 
00719       TEST_FOR_EXCEPTION(rank != blockSize_,GmresIterationOrthoFailure,
00720        "Belos::BlockGmresIter::iterate(): couldn't generate basis of full rank.");
00721       //
00722       // V has been extended, and H has been extended. 
00723       //
00724       // Update the QR factorization of the upper Hessenberg matrix
00725       //
00726       updateLSQR();
00727       //
00728       // Update basis dim and release all pointers.
00729       //
00730       Vnext = Teuchos::null;
00731       curDim_ += blockSize_;
00732       //        
00733       /*      
00734       // When required, monitor some orthogonalities
00735       if (om_->isVerbosity( Debug ) ) {
00736       // Check almost everything here
00737       CheckList chk;
00738       chk.checkV = true;
00739       chk.checkArn = true;
00740       om_->print( Debug, accuracyCheck(chk, ": after local update") );
00741       }
00742       else if (om_->isVerbosity( OrthoDetails ) ) {
00743         CheckList chk;
00744         chk.checkV = true;
00745         om_->print( OrthoDetails, accuracyCheck(chk, ": after local update") );
00746       }
00747       */ 
00748       
00749     } // end while (statusTest == false)
00750    
00751   }
00752 
00753   
00754   template<class ScalarType, class MV, class OP>
00755   void BlockGmresIter<ScalarType,MV,OP>::updateLSQR( int dim )
00756   {
00757     int i, j, maxidx;
00758     ScalarType sigma, mu, vscale, maxelem;
00759     const ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
00760     
00761     // Get correct dimension based on input "dim"
00762     // Remember that ortho failures result in an exit before updateLSQR() is called.
00763     // Therefore, it is possible that dim == curDim_.
00764     int curDim = curDim_;
00765     if (dim >= curDim_ && dim < getMaxSubspaceDim()) {
00766       curDim = dim;
00767     }
00768     
00769     Teuchos::BLAS<int, ScalarType> blas;
00770     //
00771     // Apply previous transformations and compute new transformation to reduce upper-Hessenberg
00772     // system to upper-triangular form.
00773     //
00774     if (blockSize_ == 1) {
00775       //
00776       // QR factorization of Least-Squares system with Givens rotations
00777       //
00778       for (i=0; i<curDim; i++) {
00779   //
00780   // Apply previous Givens rotations to new column of Hessenberg matrix
00781   //
00782   blas.ROT( 1, &(*R_)(i,curDim), 1, &(*R_)(i+1, curDim), 1, &cs[i], &sn[i] );
00783       }
00784       //
00785       // Calculate new Givens rotation
00786       //
00787       blas.ROTG( &(*R_)(curDim,curDim), &(*R_)(curDim+1,curDim), &cs[curDim], &sn[curDim] );
00788       (*R_)(curDim+1,curDim) = zero;
00789       //
00790       // Update RHS w/ new transformation
00791       //
00792       blas.ROT( 1, &(*z_)(curDim,0), 1, &(*z_)(curDim+1,0), 1, &cs[curDim], &sn[curDim] );
00793     }
00794     else {
00795       //
00796       // QR factorization of Least-Squares system with Householder reflectors
00797       //
00798       for (j=0; j<blockSize_; j++) {
00799   //
00800   // Apply previous Householder reflectors to new block of Hessenberg matrix
00801   //
00802   for (i=0; i<curDim+j; i++) {
00803     sigma = blas.DOT( blockSize_, &(*R_)(i+1,i), 1, &(*R_)(i+1,curDim+j), 1);
00804     sigma += (*R_)(i,curDim+j);
00805     sigma *= beta[i];
00806     blas.AXPY(blockSize_, -sigma, &(*R_)(i+1,i), 1, &(*R_)(i+1,curDim+j), 1);
00807     (*R_)(i,curDim+j) -= sigma;
00808   }
00809   //
00810   // Compute new Householder reflector
00811   //
00812   maxidx = blas.IAMAX( blockSize_+1, &(*R_)(curDim+j,curDim+j), 1 );
00813   maxelem = (*R_)(curDim+j+maxidx-1,curDim+j);
00814   for (i=0; i<blockSize_+1; i++)
00815     (*R_)(curDim+j+i,curDim+j) /= maxelem;
00816   sigma = blas.DOT( blockSize_, &(*R_)(curDim+j+1,curDim+j), 1,
00817         &(*R_)(curDim+j+1,curDim+j), 1 );
00818   if (sigma == zero) {
00819     beta[curDim + j] = zero;
00820   } else {
00821     mu = Teuchos::ScalarTraits<ScalarType>::squareroot((*R_)(curDim+j,curDim+j)*(*R_)(curDim+j,curDim+j)+sigma);
00822     if ( Teuchos::ScalarTraits<ScalarType>::real((*R_)(curDim+j,curDim+j)) 
00823          < Teuchos::ScalarTraits<MagnitudeType>::zero() ) {
00824       vscale = (*R_)(curDim+j,curDim+j) - mu;
00825     } else {
00826       vscale = -sigma / ((*R_)(curDim+j,curDim+j) + mu);
00827     }
00828     beta[curDim+j] = 2.0*vscale*vscale/(sigma + vscale*vscale);
00829     (*R_)(curDim+j,curDim+j) = maxelem*mu;
00830     for (i=0; i<blockSize_; i++)
00831       (*R_)(curDim+j+1+i,curDim+j) /= vscale;
00832   }
00833   //
00834   // Apply new Householder reflector to rhs
00835   //
00836   for (i=0; i<blockSize_; i++) {
00837     sigma = blas.DOT( blockSize_, &(*R_)(curDim+j+1,curDim+j),
00838           1, &(*z_)(curDim+j+1,i), 1);
00839     sigma += (*z_)(curDim+j,i);
00840     sigma *= beta[curDim+j];
00841     blas.AXPY(blockSize_, -sigma, &(*R_)(curDim+j+1,curDim+j),
00842         1, &(*z_)(curDim+j+1,i), 1);
00843     (*z_)(curDim+j,i) -= sigma;
00844   }
00845       }
00846     } // end if (blockSize_ == 1)
00847 
00848     // If the least-squares problem is updated wrt "dim" then update the curDim_.
00849     if (dim >= curDim_ && dim < getMaxSubspaceDim()) {
00850       curDim_ = dim + blockSize_;
00851     }
00852 
00853   } // end updateLSQR()
00854 
00855 } // end Belos namespace
00856 
00857 #endif /* BELOS_BLOCK_GMRES_ITER_HPP */
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines
Generated on Wed Apr 13 10:05:05 2011 for Belos by  doxygen 1.6.3