Belos Package Browser (Single Doxygen Collection) Development
BelosBlockFGmresIter.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_BLOCK_FGMRES_ITER_HPP
00043 #define BELOS_BLOCK_FGMRES_ITER_HPP
00044 
00049 #include "BelosConfigDefs.hpp"
00050 #include "BelosTypes.hpp"
00051 #include "BelosGmresIteration.hpp"
00052 
00053 #include "BelosLinearProblem.hpp"
00054 #include "BelosMatOrthoManager.hpp"
00055 #include "BelosOutputManager.hpp"
00056 #include "BelosStatusTest.hpp"
00057 #include "BelosOperatorTraits.hpp"
00058 #include "BelosMultiVecTraits.hpp"
00059 
00060 #include "Teuchos_BLAS.hpp"
00061 #include "Teuchos_SerialDenseMatrix.hpp"
00062 #include "Teuchos_SerialDenseVector.hpp"
00063 #include "Teuchos_ScalarTraits.hpp"
00064 #include "Teuchos_ParameterList.hpp"
00065 #include "Teuchos_TimeMonitor.hpp"
00066 
00080 namespace Belos {
00081   
00082 template<class ScalarType, class MV, class OP>
00083 class BlockFGmresIter : virtual public GmresIteration<ScalarType,MV,OP> {
00084 
00085   public:
00086     
00087   //
00088   // Convenience typedefs
00089   //
00090   typedef MultiVecTraits<ScalarType,MV> MVT;
00091   typedef OperatorTraits<ScalarType,MV,OP> OPT;
00092   typedef Teuchos::ScalarTraits<ScalarType> SCT;
00093   typedef typename SCT::magnitudeType MagnitudeType;
00094 
00096 
00097 
00107   BlockFGmresIter( const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem, 
00108        const Teuchos::RCP<OutputManager<ScalarType> > &printer,
00109        const Teuchos::RCP<StatusTest<ScalarType,MV,OP> > &tester,
00110        const Teuchos::RCP<MatOrthoManager<ScalarType,MV,OP> > &ortho,
00111        Teuchos::ParameterList &params );
00112   
00114   virtual ~BlockFGmresIter() {};
00116   
00117   
00119 
00120   
00142   void iterate();
00143 
00165   void initializeGmres(GmresIterationState<ScalarType,MV> newstate);
00166 
00170   void initialize()
00171   {
00172     GmresIterationState<ScalarType,MV> empty;
00173     initializeGmres(empty);
00174   }
00175   
00183   GmresIterationState<ScalarType,MV> getState() const {
00184     GmresIterationState<ScalarType,MV> state;
00185     state.curDim = curDim_;
00186     state.V = V_;
00187     state.Z = Z_;
00188     state.H = H_;
00189     state.R = R_;
00190     state.z = z_;
00191     return state;
00192   }
00193 
00195 
00196   
00198 
00199 
00201   int getNumIters() const { return iter_; }
00202   
00204   void resetNumIters( int iter = 0 ) { iter_ = iter; }
00205 
00208   Teuchos::RCP<const MV> getNativeResiduals( std::vector<MagnitudeType> *norms ) const;
00209 
00211 
00216   Teuchos::RCP<MV> getCurrentUpdate() const;
00217 
00219 
00222   void updateLSQR( int dim = -1 );
00223 
00225   int getCurSubspaceDim() const { 
00226     if (!initialized_) return 0;
00227     return curDim_;
00228   };
00229 
00231   int getMaxSubspaceDim() const { return blockSize_*numBlocks_; }
00232 
00234 
00235   
00237 
00238 
00240   const LinearProblem<ScalarType,MV,OP>& getProblem() const { return *lp_; }
00241 
00243   int getBlockSize() const { return blockSize_; }
00244   
00246   void setBlockSize(int blockSize) { setSize( blockSize, numBlocks_ ); }
00247   
00249   int getNumBlocks() const { return numBlocks_; }
00250   
00252   void setNumBlocks(int numBlocks) { setSize( blockSize_, numBlocks ); }
00253   
00260   void setSize(int blockSize, int numBlocks);
00261 
00263   bool isInitialized() { return initialized_; }
00264 
00266 
00267   private:
00268 
00269   //
00270   // Internal methods
00271   //
00273   void setStateSize();
00274   
00275   //
00276   // Classes inputed through constructor that define the linear problem to be solved.
00277   //
00278   const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> >    lp_;
00279   const Teuchos::RCP<OutputManager<ScalarType> >          om_;
00280   const Teuchos::RCP<StatusTest<ScalarType,MV,OP> >       stest_;
00281   const Teuchos::RCP<OrthoManager<ScalarType,MV> >        ortho_;
00282 
00283   //
00284   // Algorithmic parameters
00285   //  
00286   // blockSize_ is the solver block size.
00287   // It controls the number of vectors added to the basis on each iteration.
00288   int blockSize_;
00289   // numBlocks_ is the size of the allocated space for the Krylov basis, in blocks.
00290   int numBlocks_; 
00291   
00292   // Storage for QR factorization of the least squares system.
00293   Teuchos::SerialDenseVector<int,ScalarType> beta, sn;
00294   Teuchos::SerialDenseVector<int,MagnitudeType> cs;
00295   
00296   // 
00297   // Current solver state
00298   //
00299   // initialized_ specifies that the basis vectors have been initialized and the iterate() routine
00300   // is capable of running; _initialize is controlled  by the initialize() member method
00301   // For the implications of the state of initialized_, please see documentation for initialize()
00302   bool initialized_;
00303 
00304   // stateStorageInitialized_ specified that the state storage has be initialized to the current
00305   // blockSize_ and numBlocks_.  This initialization may be postponed if the linear problem was
00306   // generated without the right-hand side or solution vectors.
00307   bool stateStorageInitialized_;
00308 
00309   // Current subspace dimension, and number of iterations performed.
00310   int curDim_, iter_;
00311   
00312   // 
00313   // State Storage
00314   //
00315   Teuchos::RCP<MV> V_;
00316   Teuchos::RCP<MV> Z_;
00317   //
00318   // Projected matrices
00319   // H_ : Projected matrix from the Krylov factorization AV = VH + FE^T
00320   //
00321   Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > H_;
00322   // 
00323   // QR decomposition of Projected matrices for solving the least squares system HY = B.
00324   // R_: Upper triangular reduction of H
00325   // z_: Q applied to right-hand side of the least squares system
00326   Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > R_;
00327   Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > z_;  
00328 };
00329 
00331   // Constructor.
00332   template<class ScalarType, class MV, class OP>
00333   BlockFGmresIter<ScalarType,MV,OP>::BlockFGmresIter(const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem, 
00334                const Teuchos::RCP<OutputManager<ScalarType> > &printer,
00335                const Teuchos::RCP<StatusTest<ScalarType,MV,OP> > &tester,
00336                const Teuchos::RCP<MatOrthoManager<ScalarType,MV,OP> > &ortho,
00337                Teuchos::ParameterList &params ):
00338     lp_(problem),
00339     om_(printer),
00340     stest_(tester),
00341     ortho_(ortho),
00342     blockSize_(0),
00343     numBlocks_(0),
00344     initialized_(false),
00345     stateStorageInitialized_(false),
00346     curDim_(0),
00347     iter_(0)
00348   {
00349     // Get the maximum number of blocks allowed for this Krylov subspace
00350     TEST_FOR_EXCEPTION(!params.isParameter("Num Blocks"), std::invalid_argument,
00351                        "Belos::BlockFGmresIter::constructor: mandatory parameter 'Num Blocks' is not specified.");
00352     int nb = Teuchos::getParameter<int>(params, "Num Blocks");
00353 
00354     // Set the block size and allocate data
00355     int bs = params.get("Block Size", 1);
00356     setSize( bs, nb );
00357   }
00358 
00360   // Set the block size and make necessary adjustments.
00361   template <class ScalarType, class MV, class OP>
00362   void BlockFGmresIter<ScalarType,MV,OP>::setSize (int blockSize, int numBlocks)
00363   {
00364     // This routine only allocates space; it doesn't not perform any computation
00365     // any change in size will invalidate the state of the solver.
00366 
00367     TEST_FOR_EXCEPTION(numBlocks <= 0 || blockSize <= 0, std::invalid_argument, "Belos::BlockFGmresIter::setSize was passed a non-positive argument.");
00368     if (blockSize == blockSize_ && numBlocks == numBlocks_) {
00369       // do nothing
00370       return;
00371     }
00372 
00373     if (blockSize!=blockSize_ || numBlocks!=numBlocks_)
00374       stateStorageInitialized_ = false;
00375 
00376     blockSize_ = blockSize;
00377     numBlocks_ = numBlocks;
00378 
00379     initialized_ = false;
00380     curDim_ = 0;
00381 
00382     // Use the current blockSize_ and numBlocks_ to initialize the state storage.    
00383     setStateSize();
00384     
00385   }
00386 
00388   // Setup the state storage.
00389   template <class ScalarType, class MV, class OP>
00390   void BlockFGmresIter<ScalarType,MV,OP>::setStateSize ()
00391   {
00392     if (!stateStorageInitialized_) {
00393 
00394       // Check if there is any multivector to clone from.
00395       Teuchos::RCP<const MV> lhsMV = lp_->getLHS();
00396       Teuchos::RCP<const MV> rhsMV = lp_->getRHS();
00397       if (lhsMV == Teuchos::null && rhsMV == Teuchos::null) {
00398   stateStorageInitialized_ = false;
00399   return;
00400       }
00401       else {
00402   
00404   // blockSize*numBlocks dependent
00405   //
00406   int newsd = blockSize_*(numBlocks_+1);
00407   
00408   if (blockSize_==1) {
00409     cs.resize( newsd );
00410     sn.resize( newsd );
00411   }
00412   else {
00413     beta.resize( newsd );
00414   }
00415   
00416   // Initialize the state storage
00417         TEST_FOR_EXCEPTION(blockSize_*numBlocks_ > MVT::GetVecLength(*rhsMV),std::invalid_argument,
00418                            "Belos::BlockFGmresIter::setStateSize(): Cannot generate a Krylov basis with dimension larger the operator!");
00419 
00420   // If the subspace has not be initialized before, generate it using the LHS or RHS from lp_.
00421   if (V_ == Teuchos::null) {
00422     // Get the multivector that is not null.
00423     Teuchos::RCP<const MV> tmp = ( (rhsMV!=Teuchos::null)? rhsMV: lhsMV );
00424     TEST_FOR_EXCEPTION(tmp == Teuchos::null,std::invalid_argument,
00425            "Belos::BlockFGmresIter::setStateSize(): linear problem does not specify multivectors to clone from.");
00426     V_ = MVT::Clone( *tmp, newsd );
00427   }
00428   else {
00429     // Generate V_ by cloning itself ONLY if more space is needed.
00430     if (MVT::GetNumberVecs(*V_) < newsd) {
00431       Teuchos::RCP<const MV> tmp = V_;
00432       V_ = MVT::Clone( *tmp, newsd );
00433     }
00434   }
00435 
00436   if (Z_ == Teuchos::null) {
00437     // Get the multivector that is not null.
00438     Teuchos::RCP<const MV> tmp = ( (rhsMV!=Teuchos::null)? rhsMV: lhsMV );
00439     TEST_FOR_EXCEPTION(tmp == Teuchos::null,std::invalid_argument,
00440            "Belos::BlockFGmresIter::setStateSize(): linear problem does not specify multivectors to clone from.");
00441     Z_ = MVT::Clone( *tmp, newsd );
00442   }
00443   else {
00444     // Generate Z_ by cloning itself ONLY if more space is needed.
00445     if (MVT::GetNumberVecs(*Z_) < newsd) {
00446       Teuchos::RCP<const MV> tmp = Z_;
00447       Z_ = MVT::Clone( *tmp, newsd );
00448     }
00449   } 
00450 
00451   // Generate H_ only if it doesn't exist, otherwise resize it.
00452   if (H_ == Teuchos::null)
00453     H_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( newsd, newsd-blockSize_ ) ); 
00454   else
00455     H_->shapeUninitialized( newsd, newsd-blockSize_ );
00456   
00457   // TODO:  Insert logic so that Hessenberg matrix can be saved and reduced matrix is stored in R_
00458   //R_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( newsd, newsd-blockSize_ ) );
00459   // Generate z_ only if it doesn't exist, otherwise resize it.
00460   if (z_ == Teuchos::null)
00461     z_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(newsd,blockSize_) );
00462   else
00463     z_->shapeUninitialized( newsd, blockSize_ );
00464   
00465   // State storage has now been initialized.
00466   stateStorageInitialized_ = true;
00467       }
00468     }
00469   }
00470 
00472   // Get the current update from this subspace.
00473   template <class ScalarType, class MV, class OP>
00474   Teuchos::RCP<MV> BlockFGmresIter<ScalarType,MV,OP>::getCurrentUpdate() const
00475   {
00476     //
00477     // If this is the first iteration of the Arnoldi factorization, 
00478     // there is no update, so return Teuchos::null. 
00479     //
00480     RCP<MV> currentUpdate = Teuchos::null;
00481     if (curDim_==0) { 
00482       return currentUpdate; 
00483     } else {
00484       const ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
00485       const ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
00486       Teuchos::BLAS<int,ScalarType> blas;
00487       currentUpdate = MVT::Clone( *Z_, blockSize_ );
00488       //
00489       //  Make a view and then copy the RHS of the least squares problem.  DON'T OVERWRITE IT!
00490       //
00491       Teuchos::SerialDenseMatrix<int,ScalarType> y( Teuchos::Copy, *z_, curDim_, blockSize_ );
00492       //
00493       //  Solve the least squares problem.
00494       //
00495       blas.TRSM( Teuchos::LEFT_SIDE, Teuchos::UPPER_TRI, Teuchos::NO_TRANS,
00496      Teuchos::NON_UNIT_DIAG, curDim_, blockSize_, one,  
00497      H_->values(), H_->stride(), y.values(), y.stride() );
00498       //
00499       //  Compute the current update.
00500       //
00501       std::vector<int> index(curDim_);
00502       for ( int i=0; i<curDim_; i++ ) {   
00503         index[i] = i;
00504       }
00505       RCP<const MV> Zjp1 = MVT::CloneView( *Z_, index );
00506       MVT::MvTimesMatAddMv( one, *Zjp1, y, zero, *currentUpdate );
00507     }
00508     return currentUpdate;
00509   }
00510 
00511 
00513   // Get the native residuals stored in this iteration.  
00514   // Note:  No residual std::vector will be returned by FGmres.
00515   template <class ScalarType, class MV, class OP>
00516   Teuchos::RCP<const MV> BlockFGmresIter<ScalarType,MV,OP>::getNativeResiduals( std::vector<MagnitudeType> *norms ) const 
00517   {
00518     //
00519     // NOTE: Make sure the incoming std::vector is the correct size!
00520     //
00521     if ( norms && (int)norms->size() < blockSize_ )                         
00522       norms->resize( blockSize_ );                                          
00523     
00524     if (norms) {
00525       Teuchos::BLAS<int,ScalarType> blas;
00526       for (int j=0; j<blockSize_; j++) {
00527         (*norms)[j] = blas.NRM2( blockSize_, &(*z_)(curDim_, j), 1);
00528       }
00529     }
00530     return Teuchos::null;
00531   }
00532 
00533   
00534 
00536   // Initialize this iteration object
00537   template <class ScalarType, class MV, class OP>
00538   void BlockFGmresIter<ScalarType,MV,OP>::initializeGmres(GmresIterationState<ScalarType,MV> newstate)
00539   {
00540     // Initialize the state storage if it isn't already.
00541     if (!stateStorageInitialized_) 
00542       setStateSize();
00543 
00544     TEST_FOR_EXCEPTION(!stateStorageInitialized_,std::invalid_argument,
00545            "Belos::BlockFGmresIter::initialize(): Cannot initialize state storage!");
00546     
00547     // NOTE:  In BlockFGmresIter, V and Z are required!!!  
00548     // inconsitent multivectors widths and lengths will not be tolerated, and
00549     // will be treated with exceptions.
00550     //
00551     std::string errstr("Belos::BlockFGmresIter::initialize(): Specified multivectors must have a consistent length and width.");
00552 
00553     if (newstate.V != Teuchos::null && newstate.z != Teuchos::null) {
00554 
00555       // initialize V_,z_, and curDim_
00556 
00557       TEST_FOR_EXCEPTION( MVT::GetVecLength(*newstate.V) != MVT::GetVecLength(*V_),
00558                           std::invalid_argument, errstr );
00559       TEST_FOR_EXCEPTION( MVT::GetNumberVecs(*newstate.V) < blockSize_,
00560                           std::invalid_argument, errstr );
00561       TEST_FOR_EXCEPTION( newstate.curDim > blockSize_*(numBlocks_+1),
00562                           std::invalid_argument, errstr );
00563 
00564       curDim_ = newstate.curDim;
00565       int lclDim = MVT::GetNumberVecs(*newstate.V);
00566 
00567       // check size of Z
00568       TEST_FOR_EXCEPTION(newstate.z->numRows() < curDim_ || newstate.z->numCols() < blockSize_, std::invalid_argument, errstr);
00569       
00570 
00571       // copy basis vectors from newstate into V
00572       if (newstate.V != V_) {
00573         // only copy over the first block and print a warning.
00574         if (curDim_ == 0 && lclDim > blockSize_) {
00575     om_->stream(Warnings) << "Belos::BlockFGmresIter::initialize(): the solver was initialized with a kernel of " << lclDim << std::endl
00576         << "The block size however is only " << blockSize_ << std::endl
00577         << "The last " << lclDim - blockSize_ << " vectors will be discarded." << std::endl;
00578   }
00579         std::vector<int> nevind(curDim_+blockSize_);
00580         for (int i=0; i<curDim_+blockSize_; i++) nevind[i] = i;
00581   Teuchos::RCP<const MV> newV = MVT::CloneView( *newstate.V, nevind );
00582   Teuchos::RCP<MV> lclV = MVT::CloneViewNonConst( *V_, nevind );
00583         MVT::MvAddMv( 1.0, *newV, 0.0, *newV, *lclV );
00584 
00585         // done with local pointers
00586         lclV = Teuchos::null;
00587       }
00588 
00589       // put data into z_, make sure old information is not still hanging around.
00590       if (newstate.z != z_) {
00591         z_->putScalar();
00592         Teuchos::SerialDenseMatrix<int,ScalarType> newZ(Teuchos::View,*newstate.z,curDim_+blockSize_,blockSize_);
00593         Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > lclz;
00594         lclz = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*z_,curDim_+blockSize_,blockSize_) );
00595         lclz->assign(newZ);
00596 
00597         // done with local pointers
00598         lclz = Teuchos::null;
00599       }
00600 
00601     }
00602     else {
00603 
00604       TEST_FOR_EXCEPTION(newstate.V == Teuchos::null,std::invalid_argument,
00605                          "Belos::BlockFGmresIter::initialize(): BlockFGmresStateIterState does not have initial kernel V_0.");
00606 
00607       TEST_FOR_EXCEPTION(newstate.z == Teuchos::null,std::invalid_argument,
00608                          "Belos::BlockFGmresIter::initialize(): BlockFGmresStateIterState does not have initial norms z_0.");
00609     }
00610 
00611     // the solver is initialized
00612     initialized_ = true;
00613 
00614     /*
00615     if (om_->isVerbosity( Debug ) ) {
00616       // Check almost everything here
00617       CheckList chk;
00618       chk.checkV = true;
00619       chk.checkArn = true;
00620       chk.checkAux = true;
00621       om_->print( Debug, accuracyCheck(chk, ": after initialize()") );
00622     }
00623     */
00624 
00625   }
00626 
00627 
00629   // Iterate until the status test informs us we should stop.
00630   template <class ScalarType, class MV, class OP>
00631   void BlockFGmresIter<ScalarType,MV,OP>::iterate()
00632   {
00633     //
00634     // Allocate/initialize data structures
00635     //
00636     if (initialized_ == false) {
00637       initialize();
00638     }
00639     
00640     // Compute the current search dimension. 
00641     int searchDim = blockSize_*numBlocks_;
00642 
00644     // iterate until the status test tells us to stop.
00645     //
00646     // also break if our basis is full
00647     //
00648     while (stest_->checkStatus(this) != Passed && curDim_+blockSize_ <= searchDim) {
00649 
00650       iter_++;
00651 
00652       // F can be found at the curDim_ block, but the next block is at curDim_ + blockSize_.
00653       int lclDim = curDim_ + blockSize_; 
00654 
00655       // Get the current part of the basis.
00656       std::vector<int> curind(blockSize_);
00657       for (int i=0; i<blockSize_; i++) { curind[i] = lclDim + i; }
00658       Teuchos::RCP<MV> Vnext = MVT::CloneViewNonConst(*V_,curind);
00659 
00660       // Get a view of the previous vectors.
00661       // This is used for orthogonalization and for computing V^H K H
00662       for (int i=0; i<blockSize_; i++) { curind[i] = curDim_ + i; }
00663       Teuchos::RCP<const MV> Vprev = MVT::CloneView(*V_,curind);
00664       Teuchos::RCP<MV> Znext = MVT::CloneViewNonConst(*Z_,curind);
00665 
00666       // Compute the next std::vector in the Krylov basis:  Znext = M*Vprev
00667       lp_->applyRightPrec(*Vprev,*Znext);
00668       Vprev = Teuchos::null;
00669 
00670       // Compute the next std::vector in the Krylov basis:  Vnext = A*Znext
00671       lp_->applyOp(*Znext,*Vnext);
00672       Znext = Teuchos::null;
00673       
00674       // Remove all previous Krylov basis vectors from Vnext      
00675       // Get a view of all the previous vectors
00676       std::vector<int> prevind(lclDim);
00677       for (int i=0; i<lclDim; i++) { prevind[i] = i; }
00678       Vprev = MVT::CloneView(*V_,prevind);
00679       Teuchos::Array<Teuchos::RCP<const MV> > AVprev(1, Vprev);
00680       
00681       // Get a view of the part of the Hessenberg matrix needed to hold the ortho coeffs.
00682       Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> >
00683   subH = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>
00684            ( Teuchos::View,*H_,lclDim,blockSize_,0,curDim_ ) );
00685       Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > AsubH;
00686       AsubH.append( subH );
00687       
00688       // Get a view of the part of the Hessenberg matrix needed to hold the norm coeffs.
00689       Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> >
00690   subR = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>
00691            ( Teuchos::View,*H_,blockSize_,blockSize_,lclDim,curDim_ ) );
00692       int rank = ortho_->projectAndNormalize(*Vnext,AsubH,subR,AVprev);
00693       TEST_FOR_EXCEPTION(rank != blockSize_,GmresIterationOrthoFailure,
00694        "Belos::BlockFGmresIter::iterate(): couldn't generate basis of full rank.");
00695       //
00696       // V has been extended, and H has been extended. 
00697       //
00698       // Update the QR factorization of the upper Hessenberg matrix
00699       //
00700       updateLSQR();
00701       //
00702       // Update basis dim and release all pointers.
00703       //
00704       Vnext = Teuchos::null;
00705       curDim_ += blockSize_;
00706       //        
00707       /*      
00708       // When required, monitor some orthogonalities
00709       if (om_->isVerbosity( Debug ) ) {
00710       // Check almost everything here
00711       CheckList chk;
00712       chk.checkV = true;
00713       chk.checkArn = true;
00714       om_->print( Debug, accuracyCheck(chk, ": after local update") );
00715       }
00716       else if (om_->isVerbosity( OrthoDetails ) ) {
00717         CheckList chk;
00718         chk.checkV = true;
00719         om_->print( OrthoDetails, accuracyCheck(chk, ": after local update") );
00720       }
00721       */ 
00722       
00723     } // end while (statusTest == false)
00724    
00725   }
00726 
00727   
00728   template<class ScalarType, class MV, class OP>
00729   void BlockFGmresIter<ScalarType,MV,OP>::updateLSQR( int dim )
00730   {
00731     int i, j, maxidx;
00732     ScalarType sigma, mu, vscale, maxelem;
00733     const ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
00734     
00735     // Get correct dimension based on input "dim"
00736     // Remember that ortho failures result in an exit before updateLSQR() is called.
00737     // Therefore, it is possible that dim == curDim_.
00738     int curDim = curDim_;
00739     if (dim >= curDim_ && dim < getMaxSubspaceDim()) {
00740       curDim = dim;
00741     }
00742     
00743     Teuchos::BLAS<int, ScalarType> blas;
00744     //
00745     // Apply previous transformations and compute new transformation to reduce upper-Hessenberg
00746     // system to upper-triangular form.
00747     //
00748     if (blockSize_ == 1) {
00749       //
00750       // QR factorization of Least-Squares system with Givens rotations
00751       //
00752       for (i=0; i<curDim; i++) {
00753   //
00754   // Apply previous Givens rotations to new column of Hessenberg matrix
00755   //
00756   blas.ROT( 1, &(*H_)(i,curDim), 1, &(*H_)(i+1, curDim), 1, &cs[i], &sn[i] );
00757       }
00758       //
00759       // Calculate new Givens rotation
00760       //
00761       blas.ROTG( &(*H_)(curDim,curDim), &(*H_)(curDim+1,curDim), &cs[curDim], &sn[curDim] );
00762       (*H_)(curDim+1,curDim) = zero;
00763       //
00764       // Update RHS w/ new transformation
00765       //
00766       blas.ROT( 1, &(*z_)(curDim,0), 1, &(*z_)(curDim+1,0), 1, &cs[curDim], &sn[curDim] );
00767     }
00768     else {
00769       //
00770       // QR factorization of Least-Squares system with Householder reflectors
00771       //
00772       for (j=0; j<blockSize_; j++) {
00773   //
00774   // Apply previous Householder reflectors to new block of Hessenberg matrix
00775   //
00776   for (i=0; i<curDim+j; i++) {
00777     sigma = blas.DOT( blockSize_, &(*H_)(i+1,i), 1, &(*H_)(i+1,curDim+j), 1);
00778     sigma += (*H_)(i,curDim+j);
00779     sigma *= beta[i];
00780     blas.AXPY(blockSize_, -sigma, &(*H_)(i+1,i), 1, &(*H_)(i+1,curDim+j), 1);
00781     (*H_)(i,curDim+j) -= sigma;
00782   }
00783   //
00784   // Compute new Householder reflector
00785   //
00786   maxidx = blas.IAMAX( blockSize_+1, &(*H_)(curDim+j,curDim+j), 1 );
00787   maxelem = (*H_)(curDim+j+maxidx-1,curDim+j);
00788   for (i=0; i<blockSize_+1; i++)
00789     (*H_)(curDim+j+i,curDim+j) /= maxelem;
00790   sigma = blas.DOT( blockSize_, &(*H_)(curDim+j+1,curDim+j), 1,
00791         &(*H_)(curDim+j+1,curDim+j), 1 );
00792   if (sigma == zero) {
00793     beta[curDim + j] = zero;
00794   } else {
00795     mu = Teuchos::ScalarTraits<ScalarType>::squareroot((*H_)(curDim+j,curDim+j)*(*H_)(curDim+j,curDim+j)+sigma);
00796     if ( Teuchos::ScalarTraits<ScalarType>::real((*H_)(curDim+j,curDim+j)) 
00797          < Teuchos::ScalarTraits<MagnitudeType>::zero() ) {
00798       vscale = (*H_)(curDim+j,curDim+j) - mu;
00799     } else {
00800       vscale = -sigma / ((*H_)(curDim+j,curDim+j) + mu);
00801     }
00802     beta[curDim+j] = 2.0*vscale*vscale/(sigma + vscale*vscale);
00803     (*H_)(curDim+j,curDim+j) = maxelem*mu;
00804     for (i=0; i<blockSize_; i++)
00805       (*H_)(curDim+j+1+i,curDim+j) /= vscale;
00806   }
00807   //
00808   // Apply new Householder reflector to rhs
00809   //
00810   for (i=0; i<blockSize_; i++) {
00811     sigma = blas.DOT( blockSize_, &(*H_)(curDim+j+1,curDim+j),
00812           1, &(*z_)(curDim+j+1,i), 1);
00813     sigma += (*z_)(curDim+j,i);
00814     sigma *= beta[curDim+j];
00815     blas.AXPY(blockSize_, -sigma, &(*H_)(curDim+j+1,curDim+j),
00816         1, &(*z_)(curDim+j+1,i), 1);
00817     (*z_)(curDim+j,i) -= sigma;
00818   }
00819       }
00820     } // end if (blockSize_ == 1)
00821 
00822     // If the least-squares problem is updated wrt "dim" then update the curDim_.
00823     if (dim >= curDim_ && dim < getMaxSubspaceDim()) {
00824       curDim_ = dim + blockSize_;
00825     }
00826 
00827   } // end updateLSQR()
00828 
00829 } // end Belos namespace
00830 
00831 #endif /* BELOS_BLOCK_FGMRES_ITER_HPP */
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines