10#ifndef BELOS_BLOCK_GMRES_ITER_HPP
11#define BELOS_BLOCK_GMRES_ITER_HPP
28#include "Teuchos_BLAS.hpp"
29#include "Teuchos_LAPACK.hpp"
30#include "Teuchos_SerialDenseMatrix.hpp"
31#include "Teuchos_SerialDenseVector.hpp"
32#include "Teuchos_ScalarTraits.hpp"
33#include "Teuchos_ParameterList.hpp"
34#include "Teuchos_TimeMonitor.hpp"
51template<
class ScalarType,
class MV,
class OP>
61 typedef Teuchos::ScalarTraits<ScalarType>
SCT;
80 Teuchos::ParameterList ¶ms );
194 if (!initialized_)
return 0;
228 void setSize(
int blockSize,
int numBlocks);
243 CheckList() : checkV(false), checkArn(false) {};
249 std::string accuracyCheck(
const CheckList &chk,
const std::string &where)
const;
257 const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > lp_;
258 const Teuchos::RCP<OutputManager<ScalarType> > om_;
259 const Teuchos::RCP<StatusTest<ScalarType,MV,OP> > stest_;
260 const Teuchos::RCP<OrthoManager<ScalarType,MV> > ortho_;
272 Teuchos::SerialDenseVector<int,ScalarType> beta, sn;
273 Teuchos::SerialDenseVector<int,MagnitudeType> cs;
286 bool stateStorageInitialized_;
291 bool keepHessenberg_;
295 bool initHessenberg_;
308 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > H_;
313 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > R_;
314 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > z_;
319 template<
class ScalarType,
class MV,
class OP>
324 Teuchos::ParameterList ¶ms ):
332 stateStorageInitialized_(false),
333 keepHessenberg_(false),
334 initHessenberg_(false),
339 if ( om_->isVerbosity(
Debug ) )
340 keepHessenberg_ =
true;
342 keepHessenberg_ = params.get(
"Keep Hessenberg",
false);
345 initHessenberg_ = params.get(
"Initialize Hessenberg",
false);
348 TEUCHOS_TEST_FOR_EXCEPTION(!params.isParameter(
"Num Blocks"), std::invalid_argument,
349 "Belos::BlockGmresIter::constructor: mandatory parameter 'Num Blocks' is not specified.");
350 int nb = Teuchos::getParameter<int>(params,
"Num Blocks");
353 int bs = params.get(
"Block Size", 1);
359 template <
class ScalarType,
class MV,
class OP>
365 TEUCHOS_TEST_FOR_EXCEPTION(numBlocks <= 0 || blockSize <= 0, std::invalid_argument,
"Belos::BlockGmresIter::setSize was passed a non-positive argument.");
366 if (blockSize == blockSize_ && numBlocks == numBlocks_) {
371 if (blockSize!=blockSize_ || numBlocks!=numBlocks_)
372 stateStorageInitialized_ =
false;
374 blockSize_ = blockSize;
375 numBlocks_ = numBlocks;
377 initialized_ =
false;
387 template <
class ScalarType,
class MV,
class OP>
388 void BlockGmresIter<ScalarType,MV,OP>::setStateSize ()
390 if (!stateStorageInitialized_) {
393 Teuchos::RCP<const MV> lhsMV = lp_->getLHS();
394 Teuchos::RCP<const MV> rhsMV = lp_->getRHS();
395 if (lhsMV == Teuchos::null && rhsMV == Teuchos::null) {
396 stateStorageInitialized_ =
false;
404 int newsd = blockSize_*(numBlocks_+1);
411 beta.resize( newsd );
415 TEUCHOS_TEST_FOR_EXCEPTION(blockSize_*
static_cast<ptrdiff_t
>(numBlocks_) > MVT::GetGlobalLength(*rhsMV),std::invalid_argument,
416 "Belos::BlockGmresIter::setStateSize(): Cannot generate a Krylov basis with dimension larger the operator!");
419 if (V_ == Teuchos::null) {
421 Teuchos::RCP<const MV> tmp = ( (rhsMV!=Teuchos::null)? rhsMV: lhsMV );
422 TEUCHOS_TEST_FOR_EXCEPTION(tmp == Teuchos::null,std::invalid_argument,
423 "Belos::BlockGmresIter::setStateSize(): linear problem does not specify multivectors to clone from.");
424 V_ = MVT::Clone( *tmp, newsd );
428 if (MVT::GetNumberVecs(*V_) < newsd) {
429 Teuchos::RCP<const MV> tmp = V_;
430 V_ = MVT::Clone( *tmp, newsd );
435 if (R_ == Teuchos::null) {
436 R_ = Teuchos::rcp(
new Teuchos::SerialDenseMatrix<int,ScalarType>() );
438 if (initHessenberg_) {
439 R_->shape( newsd, newsd-blockSize_ );
442 if (R_->numRows() < newsd || R_->numCols() < newsd-blockSize_) {
443 R_->shapeUninitialized( newsd, newsd-blockSize_ );
448 if (keepHessenberg_) {
449 if (H_ == Teuchos::null) {
450 H_ = Teuchos::rcp(
new Teuchos::SerialDenseMatrix<int,ScalarType>() );
452 if (initHessenberg_) {
453 H_->shape( newsd, newsd-blockSize_ );
456 if (H_->numRows() < newsd || H_->numCols() < newsd-blockSize_) {
457 H_->shapeUninitialized( newsd, newsd-blockSize_ );
467 if (z_ == Teuchos::null) {
468 z_ = Teuchos::rcp(
new Teuchos::SerialDenseMatrix<int,ScalarType>() );
470 if (z_-> numRows() < newsd || z_->numCols() < blockSize_) {
471 z_->shapeUninitialized( newsd, blockSize_ );
475 stateStorageInitialized_ =
true;
482 template <
class ScalarType,
class MV,
class OP>
489 Teuchos::RCP<MV> currentUpdate = Teuchos::null;
491 return currentUpdate;
493 const ScalarType one = SCT::one();
494 const ScalarType zero = SCT::zero();
495 Teuchos::BLAS<int,ScalarType> blas;
496 currentUpdate =
MVT::Clone( *V_, blockSize_ );
500 Teuchos::SerialDenseMatrix<int,ScalarType> y( Teuchos::Copy, *z_, curDim_, blockSize_ );
504 blas.TRSM( Teuchos::LEFT_SIDE, Teuchos::UPPER_TRI, Teuchos::NO_TRANS,
505 Teuchos::NON_UNIT_DIAG, curDim_, blockSize_, one,
506 R_->values(), R_->stride(), y.values(), y.stride() );
510 std::vector<int> index(curDim_);
511 for (
int i=0; i<curDim_; i++ ) {
517 return currentUpdate;
524 template <
class ScalarType,
class MV,
class OP>
530 if ( norms && (
int)norms->size() < blockSize_ )
531 norms->resize( blockSize_ );
534 Teuchos::BLAS<int,ScalarType> blas;
535 for (
int j=0; j<blockSize_; j++) {
536 (*norms)[j] = blas.NRM2( blockSize_, &(*z_)(curDim_, j), 1);
539 return Teuchos::null;
546 template <
class ScalarType,
class MV,
class OP>
550 if (!stateStorageInitialized_)
553 TEUCHOS_TEST_FOR_EXCEPTION(!stateStorageInitialized_,std::invalid_argument,
554 "Belos::BlockGmresIter::initialize(): Cannot initialize state storage!");
556 const ScalarType zero = SCT::zero(), one = SCT::one();
562 std::string errstr(
"Belos::BlockGmresIter::initialize(): Specified multivectors must have a consistent length and width.");
564 if (newstate.
V != Teuchos::null && newstate.
z != Teuchos::null) {
569 std::invalid_argument, errstr );
571 std::invalid_argument, errstr );
572 TEUCHOS_TEST_FOR_EXCEPTION( newstate.
curDim > blockSize_*(numBlocks_+1),
573 std::invalid_argument, errstr );
575 curDim_ = newstate.
curDim;
579 TEUCHOS_TEST_FOR_EXCEPTION(newstate.
z->numRows() < curDim_ || newstate.
z->numCols() < blockSize_, std::invalid_argument, errstr);
583 if (newstate.
V != V_) {
585 if (curDim_ == 0 && lclDim > blockSize_) {
586 om_->stream(
Warnings) <<
"Belos::BlockGmresIter::initialize(): the solver was initialized with a kernel of " << lclDim << std::endl
587 <<
"The block size however is only " << blockSize_ << std::endl
588 <<
"The last " << lclDim - blockSize_ <<
" vectors will be discarded." << std::endl;
590 std::vector<int> nevind(curDim_+blockSize_);
591 for (
int i=0; i<curDim_+blockSize_; i++) nevind[i] = i;
597 lclV = Teuchos::null;
601 if (newstate.
z != z_) {
603 Teuchos::SerialDenseMatrix<int,ScalarType> newZ(Teuchos::View,*newstate.
z,curDim_+blockSize_,blockSize_);
604 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > lclZ;
605 lclZ = Teuchos::rcp(
new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*z_,curDim_+blockSize_,blockSize_) );
609 lclZ = Teuchos::null;
615 TEUCHOS_TEST_FOR_EXCEPTION(newstate.
V == Teuchos::null,std::invalid_argument,
616 "Belos::BlockGmresIter::initialize(): BlockGmresStateIterState does not have initial kernel V_0.");
618 TEUCHOS_TEST_FOR_EXCEPTION(newstate.
z == Teuchos::null,std::invalid_argument,
619 "Belos::BlockGmresIter::initialize(): BlockGmresStateIterState does not have initial norms z_0.");
625 if (om_->isVerbosity(
Debug ) ) {
630 om_->print(
Debug, accuracyCheck(chk,
": after initialize()") );
638 template <
class ScalarType,
class MV,
class OP>
644 if (initialized_ ==
false) {
649 int searchDim = blockSize_*numBlocks_;
656 while (stest_->checkStatus(
this) !=
Passed && curDim_+blockSize_ <= searchDim) {
661 int lclDim = curDim_ + blockSize_;
664 std::vector<int> curind(blockSize_);
665 for (
int i=0; i<blockSize_; i++) { curind[i] = lclDim + i; }
670 for (
int i=0; i<blockSize_; i++) { curind[i] = curDim_ + i; }
674 lp_->apply(*Vprev,*Vnext);
675 Vprev = Teuchos::null;
679 std::vector<int> prevind(lclDim);
680 for (
int i=0; i<lclDim; i++) { prevind[i] = i; }
682 Teuchos::Array<Teuchos::RCP<const MV> > AVprev(1, Vprev);
685 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> >
686 subH = Teuchos::rcp(
new Teuchos::SerialDenseMatrix<int,ScalarType>
687 ( Teuchos::View,*H_,lclDim,blockSize_,0,curDim_ ) );
688 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > AsubH;
689 AsubH.append( subH );
692 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> >
693 subH2 = Teuchos::rcp(
new Teuchos::SerialDenseMatrix<int,ScalarType>
694 ( Teuchos::View,*H_,blockSize_,blockSize_,lclDim,curDim_ ) );
696 int rank = ortho_->projectAndNormalize(*Vnext,AsubH,subH2,AVprev);
700 if (keepHessenberg_) {
702 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> >
703 subR = Teuchos::rcp(
new Teuchos::SerialDenseMatrix<int,ScalarType>
704 ( Teuchos::View,*R_,lclDim,blockSize_,0,curDim_ ) );
708 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> >
709 subR2 = Teuchos::rcp(
new Teuchos::SerialDenseMatrix<int,ScalarType>
710 ( Teuchos::View,*R_,blockSize_,blockSize_,lclDim,curDim_ ) );
711 subR2->assign(*subH2);
715 "Belos::BlockGmresIter::iterate(): couldn't generate basis of full rank.");
725 Vnext = Teuchos::null;
726 curDim_ += blockSize_;
729 if (om_->isVerbosity(
Debug ) ) {
734 om_->print(
Debug, accuracyCheck(chk,
": after local update") );
739 om_->print(
OrthoDetails, accuracyCheck(chk,
": after local update") );
747 template<
class ScalarType,
class MV,
class OP>
751 ScalarType sigma, mu, vscale, maxelem;
752 const ScalarType zero = SCT::zero();
757 int curDim = curDim_;
762 Teuchos::BLAS<int, ScalarType> blas;
767 if (blockSize_ == 1) {
771 for (i=0; i<curDim; i++) {
775 blas.ROT( 1, &(*R_)(i,curDim), 1, &(*R_)(i+1, curDim), 1, &cs[i], &sn[i] );
780 blas.ROTG( &(*R_)(curDim,curDim), &(*R_)(curDim+1,curDim), &cs[curDim], &sn[curDim] );
781 (*R_)(curDim+1,curDim) = zero;
785 blas.ROT( 1, &(*z_)(curDim,0), 1, &(*z_)(curDim+1,0), 1, &cs[curDim], &sn[curDim] );
791 for (j=0; j<blockSize_; j++) {
795 for (i=0; i<curDim+j; i++) {
796 sigma = blas.DOT( blockSize_, &(*R_)(i+1,i), 1, &(*R_)(i+1,curDim+j), 1);
797 sigma += (*R_)(i,curDim+j);
798 sigma *= SCT::conjugate(beta[i]);
799 blas.AXPY(blockSize_, ScalarType(-sigma), &(*R_)(i+1,i), 1, &(*R_)(i+1,curDim+j), 1);
800 (*R_)(i,curDim+j) -= sigma;
805 maxidx = blas.IAMAX( blockSize_+1, &(*R_)(curDim+j,curDim+j), 1 );
806 maxelem = SCT::magnitude((*R_)(curDim+j+maxidx-1,curDim+j));
807 for (i=0; i<blockSize_+1; i++)
808 (*R_)(curDim+j+i,curDim+j) /= maxelem;
809 sigma = blas.DOT( blockSize_, &(*R_)(curDim+j+1,curDim+j), 1,
810 &(*R_)(curDim+j+1,curDim+j), 1 );
811 MagnitudeType sign_Rjj = -SCT::real((*R_)(curDim+j,curDim+j)) /
812 SCT::magnitude(SCT::real(((*R_)(curDim+j,curDim+j))));
814 beta[curDim + j] = zero;
816 mu = SCT::squareroot(SCT::conjugate((*R_)(curDim+j,curDim+j))*(*R_)(curDim+j,curDim+j)+sigma);
817 vscale = (*R_)(curDim+j,curDim+j) - Teuchos::as<ScalarType>(sign_Rjj)*mu;
818 beta[curDim+j] = -Teuchos::as<ScalarType>(sign_Rjj) * vscale / mu;
819 (*R_)(curDim+j,curDim+j) = Teuchos::as<ScalarType>(sign_Rjj)*maxelem*mu;
820 for (i=0; i<blockSize_; i++)
821 (*R_)(curDim+j+1+i,curDim+j) /= vscale;
826 for (i=0; i<blockSize_; i++) {
827 sigma = blas.DOT( blockSize_, &(*R_)(curDim+j+1,curDim+j),
828 1, &(*z_)(curDim+j+1,i), 1);
829 sigma += (*z_)(curDim+j,i);
830 sigma *= SCT::conjugate(beta[curDim+j]);
831 blas.AXPY(blockSize_, ScalarType(-sigma), &(*R_)(curDim+j+1,curDim+j),
832 1, &(*z_)(curDim+j+1,i), 1);
833 (*z_)(curDim+j,i) -= sigma;
840 curDim_ = dim + blockSize_;
856 template <
class ScalarType,
class MV,
class OP>
857 std::string BlockGmresIter<ScalarType,MV,OP>::accuracyCheck(
const CheckList &chk,
const std::string &where )
const
859 std::stringstream os;
861 os.setf(std::ios::scientific, std::ios::floatfield);
864 os <<
" Debugging checks: iteration " << iter_ << where << std::endl;
867 std::vector<int> lclind(curDim_);
868 for (
int i=0; i<curDim_; i++) lclind[i] = i;
869 std::vector<int> bsind(blockSize_);
870 for (
int i=0; i<blockSize_; i++) { bsind[i] = curDim_ + i; }
872 Teuchos::RCP<const MV> lclV,lclF;
873 Teuchos::RCP<MV> lclAV;
875 lclV = MVT::CloneView(*V_,lclind);
876 lclF = MVT::CloneView(*V_,bsind);
880 tmp = ortho_->orthonormError(*lclV);
881 os <<
" >> Error in V^H M V == I : " << tmp << std::endl;
883 tmp = ortho_->orthonormError(*lclF);
884 os <<
" >> Error in F^H M F == I : " << tmp << std::endl;
886 tmp = ortho_->orthogError(*lclV,*lclF);
887 os <<
" >> Error in V^H M F == 0 : " << tmp << std::endl;
895 lclAV = MVT::Clone(*V_,curDim_);
896 lp_->apply(*lclV,*lclAV);
899 const ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
900 Teuchos::SerialDenseMatrix<int,ScalarType> subH(Teuchos::View,*H_,curDim_,curDim_);
901 MVT::MvTimesMatAddMv( -one, *lclV, subH, one, *lclAV );
904 Teuchos::SerialDenseMatrix<int,ScalarType> curB(Teuchos::View,*H_,
905 blockSize_,curDim_, curDim_ );
906 MVT::MvTimesMatAddMv( -one, *lclF, curB, one, *lclAV );
909 std::vector<MagnitudeType> arnNorms( curDim_ );
910 ortho_->norm( *lclAV, arnNorms );
912 for (
int i=0; i<curDim_; i++) {
913 os <<
" >> Error in Krylov factorization (R = AV-VH-FB^H), ||R[" << i <<
"]|| : " << arnNorms[i] << std::endl;
Belos header file which uses auto-configuration information to include necessary C++ headers.
Pure virtual base class which augments the basic interface for a Gmres linear solver iteration.
Class which describes the linear problem to be solved by the iterative solver.
Templated virtual class for providing orthogonalization/orthonormalization methods with matrix-based ...
Declaration of basic traits for the multivector type.
Class which defines basic traits for the operator type.
Class which manages the output and verbosity of the Belos solvers.
Pure virtual base class for defining the status testing capabilities of Belos.
Collection of types and exceptions used within the Belos solvers.
void setNumBlocks(int numBlocks)
Set the maximum number of blocks used by the iterative solver.
GmresIterationState< ScalarType, MV > getState() const
Get the current state of the linear solver.
int getMaxSubspaceDim() const
Get the maximum dimension allocated for the search subspace.
Teuchos::ScalarTraits< ScalarType > SCT
int getBlockSize() const
Get the blocksize to be used by the iterative solver in solving this linear problem.
const LinearProblem< ScalarType, MV, OP > & getProblem() const
Get a constant reference to the linear problem.
SCT::magnitudeType MagnitudeType
virtual ~BlockGmresIter()
Destructor.
Teuchos::RCP< const MV > getNativeResiduals(std::vector< MagnitudeType > *norms) const
void setSize(int blockSize, int numBlocks)
Set the blocksize and number of blocks to be used by the iterative solver in solving this linear prob...
BlockGmresIter(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem, const Teuchos::RCP< OutputManager< ScalarType > > &printer, const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &tester, const Teuchos::RCP< MatOrthoManager< ScalarType, MV, OP > > &ortho, Teuchos::ParameterList ¶ms)
BlockGmresIter constructor with linear problem, solver utilities, and parameter list of solver option...
void iterate()
This method performs block Gmres iterations until the status test indicates the need to stop or an er...
MultiVecTraits< ScalarType, MV > MVT
void updateLSQR(int dim=-1)
Method for updating QR factorization of upper Hessenberg matrix.
bool isInitialized()
States whether the solver has been initialized or not.
int getCurSubspaceDim() const
Get the dimension of the search subspace used to generate the current solution to the linear problem.
void resetNumIters(int iter=0)
Reset the iteration count.
void setBlockSize(int blockSize)
Set the blocksize.
int getNumBlocks() const
Get the maximum number of blocks used by the iterative solver in solving this linear problem.
void initializeGmres(GmresIterationState< ScalarType, MV > &newstate)
Initialize the solver to an iterate, providing a complete state.
Teuchos::RCP< MV > getCurrentUpdate() const
Get the current update to the linear system.
void initialize()
Initialize the solver with the initial vectors from the linear problem or random data.
OperatorTraits< ScalarType, MV, OP > OPT
int getNumIters() const
Get the current iteration count.
GmresIterationOrthoFailure is thrown when the GmresIteration object is unable to compute independent ...
Belos's templated virtual class for providing routines for orthogonalization and orthonormzalition of...
Traits class which defines basic operations on multivectors.
static void MvTimesMatAddMv(const ScalarType alpha, const MV &A, const Teuchos::SerialDenseMatrix< int, ScalarType > &B, const ScalarType beta, MV &mv)
Update mv with .
static ptrdiff_t GetGlobalLength(const MV &mv)
Return the number of rows in the given multivector mv.
static Teuchos::RCP< MV > Clone(const MV &mv, const int numvecs)
Creates a new empty MV containing numvecs columns.
static Teuchos::RCP< const MV > CloneView(const MV &mv, const std::vector< int > &index)
Creates a new const MV that shares the selected contents of mv (shallow copy).
static Teuchos::RCP< MV > CloneViewNonConst(MV &mv, const std::vector< int > &index)
Creates a new MV that shares the selected contents of mv (shallow copy).
static void MvAddMv(const ScalarType alpha, const MV &A, const ScalarType beta, const MV &B, MV &mv)
Replace mv with .
static int GetNumberVecs(const MV &mv)
Obtain the number of vectors in mv.
Class which defines basic traits for the operator type.
Belos's basic output manager for sending information of select verbosity levels to the appropriate ou...
A pure virtual class for defining the status tests for the Belos iterative solvers.
Structure to contain pointers to GmresIteration state variables.
Teuchos::RCP< const MV > V
The current Krylov basis.
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > z
The current right-hand side of the least squares system RY = Z.
int curDim
The current dimension of the reduction.
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > H
The current Hessenberg matrix.
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > R
The current upper-triangular matrix from the QR reduction of H.