10#ifndef BELOS_BLOCK_FGMRES_ITER_HPP
11#define BELOS_BLOCK_FGMRES_ITER_HPP
28#include "Teuchos_BLAS.hpp"
29#include "Teuchos_SerialDenseMatrix.hpp"
30#include "Teuchos_SerialDenseVector.hpp"
31#include "Teuchos_ScalarTraits.hpp"
32#include "Teuchos_ParameterList.hpp"
33#include "Teuchos_TimeMonitor.hpp"
50template<
class ScalarType,
class MV,
class OP>
60 typedef Teuchos::ScalarTraits<ScalarType>
SCT;
79 Teuchos::ParameterList &
params );
153 state.curDim = curDim_;
194 if (!initialized_)
return 0;
246 const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > lp_;
247 const Teuchos::RCP<OutputManager<ScalarType> > om_;
248 const Teuchos::RCP<StatusTest<ScalarType,MV,OP> > stest_;
249 const Teuchos::RCP<OrthoManager<ScalarType,MV> > ortho_;
261 Teuchos::SerialDenseVector<int,ScalarType> beta, sn;
262 Teuchos::SerialDenseVector<int,MagnitudeType> cs;
275 bool stateStorageInitialized_;
280 bool keepHessenberg_;
294 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > H_;
299 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > R_;
300 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > z_;
305 template<
class ScalarType,
class MV,
class OP>
311 Teuchos::ParameterList &
params ):
319 stateStorageInitialized_(
false),
320 keepHessenberg_(
false),
325 if (om_->isVerbosity(
Debug))
326 keepHessenberg_ =
true;
328 keepHessenberg_ =
params.get(
"Keep Hessenberg",
false);
332 !
params.isParameter (
"Num Blocks"), std::invalid_argument,
333 "Belos::BlockFGmresIter::constructor: mandatory parameter 'Num Blocks' is not specified.");
334 const int nb =
params.get<
int> (
"Num Blocks");
337 const int bs =
params.get (
"Block Size", 1);
343 template <
class ScalarType,
class MV,
class OP>
356 stateStorageInitialized_ =
false;
361 initialized_ =
false;
371 template <
class ScalarType,
class MV,
class OP>
376 typedef Teuchos::SerialDenseMatrix<int, ScalarType> SDM;
378 if (! stateStorageInitialized_) {
382 if (
lhsMV == Teuchos::null &&
rhsMV == Teuchos::null) {
383 stateStorageInitialized_ =
false;
390 int newsd = blockSize_*(numBlocks_+1);
402 blockSize_ *
static_cast<ptrdiff_t> (numBlocks_) > MVT::GetGlobalLength (*
rhsMV),
403 std::invalid_argument,
"Belos::BlockFGmresIter::setStateSize(): "
404 "Cannot generate a Krylov basis with dimension larger the operator!");
407 if (V_ == Teuchos::null) {
411 tmp == Teuchos::null, std::invalid_argument,
412 "Belos::BlockFGmresIter::setStateSize(): "
413 "linear problem does not specify multivectors to clone from.");
418 if (MVT::GetNumberVecs (*V_) <
newsd) {
424 if (Z_ == Teuchos::null) {
428 tmp == Teuchos::null, std::invalid_argument,
429 "Belos::BlockFGmresIter::setStateSize(): "
430 "linear problem does not specify multivectors to clone from.");
435 if (MVT::GetNumberVecs (*Z_) <
newsd) {
442 if (R_ == Teuchos::null) {
446 R_->shapeUninitialized (
newsd,
newsd - blockSize_);
452 if (keepHessenberg_) {
453 if (H_ == Teuchos::null) {
457 H_->shapeUninitialized (
newsd,
newsd - blockSize_);
465 if (z_ == Teuchos::null) {
466 z_ =
rcp (
new SDM (
newsd, blockSize_));
469 z_->shapeUninitialized (
newsd, blockSize_);
473 stateStorageInitialized_ =
true;
479 template <
class ScalarType,
class MV,
class OP>
483 typedef Teuchos::SerialDenseMatrix<int, ScalarType> SDM;
492 const ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero ();
493 const ScalarType one = Teuchos::ScalarTraits<ScalarType>::one ();
494 Teuchos::BLAS<int,ScalarType>
blas;
499 SDM y (Teuchos::Copy, *z_, curDim_, blockSize_);
502 blas.TRSM (Teuchos::LEFT_SIDE, Teuchos::UPPER_TRI, Teuchos::NO_TRANS,
503 Teuchos::NON_UNIT_DIAG, curDim_, blockSize_,
one,
504 R_->values (), R_->stride (), y.values (), y.stride ());
507 std::vector<int> index (curDim_);
508 for (
int i = 0;
i < curDim_; ++
i) {
511 Teuchos::RCP<const MV>
Zjp1 = MVT::CloneView (*Z_, index);
518 template <
class ScalarType,
class MV,
class OP>
519 Teuchos::RCP<const MV>
525 norms->resize (blockSize_);
529 Teuchos::BLAS<int, ScalarType>
blas;
530 for (
int j = 0;
j < blockSize_; ++
j) {
531 (*norms)[
j] =
blas.NRM2 (blockSize_, &(*z_)(curDim_,
j), 1);
536 return Teuchos::null;
540 template <
class ScalarType,
class MV,
class OP>
547 typedef Teuchos::ScalarTraits<ScalarType> STS;
548 typedef Teuchos::SerialDenseMatrix<int, ScalarType> SDM;
553 if (! stateStorageInitialized_) {
558 ! stateStorageInitialized_, std::invalid_argument,
559 "Belos::BlockFGmresIter::initialize(): Cannot initialize state storage!");
564 const char errstr[] =
"Belos::BlockFGmresIter::initialize(): The given "
565 "multivectors must have a consistent length and width.";
572 MVT::GetGlobalLength(*
newstate.V) != MVT::GetGlobalLength(*V_),
573 std::invalid_argument,
errstr );
575 MVT::GetNumberVecs(*
newstate.V) < blockSize_,
576 std::invalid_argument,
errstr );
578 newstate.curDim > blockSize_*(numBlocks_+1),
579 std::invalid_argument,
errstr );
587 std::invalid_argument,
errstr);
592 if (curDim_ == 0 &&
lclDim > blockSize_) {
594 warn <<
"Belos::BlockFGmresIter::initialize(): the solver was "
595 <<
"initialized with a kernel of " <<
lclDim << endl
596 <<
"The block size however is only " << blockSize_ << endl
597 <<
"The last " <<
lclDim - blockSize_
598 <<
" vectors will be discarded." << endl;
600 std::vector<int>
nevind (curDim_ + blockSize_);
601 for (
int i = 0;
i < curDim_ + blockSize_; ++
i) {
609 lclV = Teuchos::null;
615 SDM
newZ (Teuchos::View, *
newstate.z, curDim_ + blockSize_, blockSize_);
617 lclz =
rcp (
new SDM (Teuchos::View, *z_, curDim_ + blockSize_, blockSize_));
619 lclz = Teuchos::null;
624 newstate.V == Teuchos::null,std::invalid_argument,
625 "Belos::BlockFGmresIter::initialize(): BlockFGmresStateIterState does not have initial kernel V_0.");
628 newstate.z == Teuchos::null,std::invalid_argument,
629 "Belos::BlockFGmresIter::initialize(): BlockFGmresStateIterState does not have initial norms z_0.");
637 template <
class ScalarType,
class MV,
class OP>
640 using Teuchos::Array;
645 typedef Teuchos::SerialDenseMatrix<int, ScalarType> SDM;
648 if (initialized_ ==
false) {
653 const int searchDim = blockSize_ * numBlocks_;
657 while (stest_->checkStatus (
this) !=
Passed && curDim_+blockSize_ <=
searchDim) {
661 const int lclDim = curDim_ + blockSize_;
664 std::vector<int>
curind (blockSize_);
665 for (
int i = 0;
i < blockSize_; ++
i) {
672 for (
int i = 0;
i < blockSize_; ++
i) {
706 "Belos::BlockFGmresIter::iterate(): After orthogonalization, the new "
707 "basis block does not have full rank. It contains " << blockSize_
708 <<
" vector" << (blockSize_ != 1 ?
"s" :
"")
709 <<
", but its rank is " <<
rank <<
".");
713 if (keepHessenberg_) {
730 curDim_ += blockSize_;
735 template<
class ScalarType,
class MV,
class OP>
738 typedef Teuchos::ScalarTraits<ScalarType> STS;
739 typedef Teuchos::ScalarTraits<MagnitudeType> STM;
744 Teuchos::BLAS<int, ScalarType>
blas;
750 int curDim = curDim_;
751 if (
dim >= curDim_ &&
dim < getMaxSubspaceDim ()) {
760 if (blockSize_ == 1) {
762 for (
int i = 0;
i < curDim; ++
i) {
764 blas.ROT (1, &(*R_)(
i, curDim), 1, &(*R_)(
i+1, curDim), 1, &cs[
i], &sn[
i]);
767 blas.ROTG (&(*R_)(curDim, curDim), &(*R_)(curDim+1, curDim), &cs[curDim], &sn[curDim]);
768 (*R_)(curDim+1, curDim) =
zero;
771 blas.ROT (1, &(*z_)(curDim,0), 1, &(*z_)(curDim+1,0), 1, &cs[curDim], &sn[curDim]);
775 for (
int j = 0;
j < blockSize_; ++
j) {
777 for (
int i = 0;
i < curDim +
j; ++
i) {
778 sigma =
blas.DOT (blockSize_, &(*R_)(
i+1,
i), 1, &(*R_)(
i+1,curDim+
j), 1);
786 const int maxidx =
blas.IAMAX (blockSize_+1, &(*R_)(curDim+
j,curDim+
j), 1);
788 for (
int i = 0;
i < blockSize_ + 1; ++
i) {
791 sigma =
blas.DOT (blockSize_, &(*R_)(curDim +
j + 1, curDim +
j), 1,
792 &(*R_)(curDim +
j + 1, curDim +
j), 1);
794 beta[curDim +
j] =
zero;
796 mu = STS::squareroot ((*R_)(curDim+
j,curDim+
j)*(*R_)(curDim+
j,curDim+
j)+
sigma);
797 if (STS::real ((*R_)(curDim +
j, curDim +
j)) < STM::zero ()) {
804 for (
int i = 0;
i < blockSize_; ++
i) {
810 for (
int i = 0;
i < blockSize_; ++
i) {
811 sigma =
blas.DOT (blockSize_, &(*R_)(curDim+
j+1,curDim+
j),
812 1, &(*z_)(curDim+
j+1,
i), 1);
816 1, &(*z_)(curDim+
j+1,
i), 1);
823 if (
dim >= curDim_ &&
dim < getMaxSubspaceDim ()) {
824 curDim_ =
dim + blockSize_;
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.
This class implements the block flexible GMRES iteration, where a block Krylov subspace is constructe...
int getNumIters() const
Get the current iteration count.
int getBlockSize() const
Get the blocksize to be used by the iterative solver in solving this linear problem.
OperatorTraits< ScalarType, MV, OP > OPT
Teuchos::RCP< const MV > getNativeResiduals(std::vector< MagnitudeType > *norms) const
Get the norms of the residuals native to the solver.
void iterate()
This method performs block FGmres iterations until the status test indicates the need to stop or an e...
MultiVecTraits< ScalarType, MV > MVT
GmresIterationState< ScalarType, MV > getState() const
Get the current state of the linear solver.
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.
SCT::magnitudeType MagnitudeType
const LinearProblem< ScalarType, MV, OP > & getProblem() const
Get a constant reference to the linear problem.
Teuchos::RCP< MV > getCurrentUpdate() const
Get the current update to the linear system.
void resetNumIters(int iter=0)
Reset the iteration count.
virtual ~BlockFGmresIter()
Destructor.
void updateLSQR(int dim=-1)
Method for updating QR factorization of upper Hessenberg matrix.
void setBlockSize(int blockSize)
Set the blocksize.
void initialize()
Initialize the solver with the initial vectors from the linear problem or random data.
BlockFGmresIter(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)
BlockFGmresIter constructor with linear problem, solver utilities, and parameter list of solver optio...
int getCurSubspaceDim() const
Get the dimension of the search subspace used to generate the current solution to the linear problem.
int getMaxSubspaceDim() const
Get the maximum dimension allocated for the search subspace.
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...
void setNumBlocks(int numBlocks)
Set the maximum number of blocks used by the iterative solver.
Teuchos::ScalarTraits< ScalarType > SCT
bool isInitialized()
States whether the solver has been initialized or not.
GmresIterationOrthoFailure is thrown when the GmresIteration object is unable to compute independent ...
Alternative run-time polymorphic interface for operators.
Operator()
Default constructor (does nothing).