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 ¶ms );
194 if (!initialized_)
return 0;
228 void setSize(
int blockSize,
int numBlocks);
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_;
289 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > H_;
294 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > R_;
295 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > z_;
300 template<
class ScalarType,
class MV,
class OP>
306 Teuchos::ParameterList ¶ms ):
314 stateStorageInitialized_(false),
319 TEUCHOS_TEST_FOR_EXCEPTION(
320 ! params.isParameter (
"Num Blocks"), std::invalid_argument,
321 "Belos::BlockFGmresIter::constructor: mandatory parameter 'Num Blocks' is not specified.");
322 const int nb = params.get<
int> (
"Num Blocks");
325 const int bs = params.get (
"Block Size", 1);
331 template <
class ScalarType,
class MV,
class OP>
337 TEUCHOS_TEST_FOR_EXCEPTION(numBlocks <= 0 || blockSize <= 0, std::invalid_argument,
"Belos::BlockFGmresIter::setSize was passed a non-positive argument.");
338 if (blockSize == blockSize_ && numBlocks == numBlocks_) {
343 if (blockSize!=blockSize_ || numBlocks!=numBlocks_)
344 stateStorageInitialized_ =
false;
346 blockSize_ = blockSize;
347 numBlocks_ = numBlocks;
349 initialized_ =
false;
359 template <
class ScalarType,
class MV,
class OP>
360 void BlockFGmresIter<ScalarType,MV,OP>::setStateSize ()
364 typedef Teuchos::SerialDenseMatrix<int, ScalarType> SDM;
366 if (! stateStorageInitialized_) {
368 RCP<const MV> lhsMV = lp_->getLHS();
369 RCP<const MV> rhsMV = lp_->getRHS();
370 if (lhsMV == Teuchos::null && rhsMV == Teuchos::null) {
371 stateStorageInitialized_ =
false;
378 int newsd = blockSize_*(numBlocks_+1);
389 TEUCHOS_TEST_FOR_EXCEPTION(
390 blockSize_ *
static_cast<ptrdiff_t
> (numBlocks_) > MVT::GetGlobalLength (*rhsMV),
391 std::invalid_argument,
"Belos::BlockFGmresIter::setStateSize(): "
392 "Cannot generate a Krylov basis with dimension larger the operator!");
395 if (V_ == Teuchos::null) {
397 RCP<const MV> tmp = (rhsMV != Teuchos::null) ? rhsMV : lhsMV;
398 TEUCHOS_TEST_FOR_EXCEPTION(
399 tmp == Teuchos::null, std::invalid_argument,
400 "Belos::BlockFGmresIter::setStateSize(): "
401 "linear problem does not specify multivectors to clone from.");
402 V_ = MVT::Clone (*tmp, newsd);
406 if (MVT::GetNumberVecs (*V_) < newsd) {
407 RCP<const MV> tmp = V_;
408 V_ = MVT::Clone (*tmp, newsd);
412 if (Z_ == Teuchos::null) {
414 RCP<const MV> tmp = (rhsMV != Teuchos::null) ? rhsMV : lhsMV;
415 TEUCHOS_TEST_FOR_EXCEPTION(
416 tmp == Teuchos::null, std::invalid_argument,
417 "Belos::BlockFGmresIter::setStateSize(): "
418 "linear problem does not specify multivectors to clone from.");
419 Z_ = MVT::Clone (*tmp, newsd);
423 if (MVT::GetNumberVecs (*Z_) < newsd) {
424 RCP<const MV> tmp = Z_;
425 Z_ = MVT::Clone (*tmp, newsd);
430 if (H_ == Teuchos::null) {
431 H_ = rcp (
new SDM (newsd, newsd-blockSize_));
434 H_->shapeUninitialized (newsd, newsd - blockSize_);
440 if (z_ == Teuchos::null) {
441 z_ = rcp (
new SDM (newsd, blockSize_));
444 z_->shapeUninitialized (newsd, blockSize_);
448 stateStorageInitialized_ =
true;
454 template <
class ScalarType,
class MV,
class OP>
458 typedef Teuchos::SerialDenseMatrix<int, ScalarType> SDM;
460 Teuchos::RCP<MV> currentUpdate = Teuchos::null;
464 return currentUpdate;
467 const ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero ();
468 const ScalarType one = Teuchos::ScalarTraits<ScalarType>::one ();
469 Teuchos::BLAS<int,ScalarType> blas;
474 SDM y (Teuchos::Copy, *z_, curDim_, blockSize_);
477 blas.TRSM (Teuchos::LEFT_SIDE, Teuchos::UPPER_TRI, Teuchos::NO_TRANS,
478 Teuchos::NON_UNIT_DIAG, curDim_, blockSize_, one,
479 H_->values (), H_->stride (), y.values (), y.stride ());
482 std::vector<int> index (curDim_);
483 for (
int i = 0; i < curDim_; ++i) {
489 return currentUpdate;
493 template <
class ScalarType,
class MV,
class OP>
494 Teuchos::RCP<const MV>
499 if (norms != NULL && (
int)norms->size() < blockSize_) {
500 norms->resize (blockSize_);
504 Teuchos::BLAS<int, ScalarType> blas;
505 for (
int j = 0; j < blockSize_; ++j) {
506 (*norms)[j] = blas.NRM2 (blockSize_, &(*z_)(curDim_, j), 1);
511 return Teuchos::null;
515 template <
class ScalarType,
class MV,
class OP>
522 typedef Teuchos::ScalarTraits<ScalarType> STS;
523 typedef Teuchos::SerialDenseMatrix<int, ScalarType> SDM;
524 const ScalarType ZERO = STS::zero ();
525 const ScalarType ONE = STS::one ();
528 if (! stateStorageInitialized_) {
532 TEUCHOS_TEST_FOR_EXCEPTION(
533 ! stateStorageInitialized_, std::invalid_argument,
534 "Belos::BlockFGmresIter::initialize(): Cannot initialize state storage!");
539 const char errstr[] =
"Belos::BlockFGmresIter::initialize(): The given "
540 "multivectors must have a consistent length and width.";
542 if (! newstate.
V.is_null () && ! newstate.
z.is_null ()) {
546 TEUCHOS_TEST_FOR_EXCEPTION(
548 std::invalid_argument, errstr );
549 TEUCHOS_TEST_FOR_EXCEPTION(
551 std::invalid_argument, errstr );
552 TEUCHOS_TEST_FOR_EXCEPTION(
553 newstate.
curDim > blockSize_*(numBlocks_+1),
554 std::invalid_argument, errstr );
556 curDim_ = newstate.
curDim;
560 TEUCHOS_TEST_FOR_EXCEPTION(
561 newstate.
z->numRows() < curDim_ || newstate.
z->numCols() < blockSize_,
562 std::invalid_argument, errstr);
565 if (newstate.
V != V_) {
567 if (curDim_ == 0 && lclDim > blockSize_) {
568 std::ostream& warn = om_->stream (
Warnings);
569 warn <<
"Belos::BlockFGmresIter::initialize(): the solver was "
570 <<
"initialized with a kernel of " << lclDim << endl
571 <<
"The block size however is only " << blockSize_ << endl
572 <<
"The last " << lclDim - blockSize_
573 <<
" vectors will be discarded." << endl;
575 std::vector<int> nevind (curDim_ + blockSize_);
576 for (
int i = 0; i < curDim_ + blockSize_; ++i) {
584 lclV = Teuchos::null;
588 if (newstate.
z != z_) {
590 SDM newZ (Teuchos::View, *newstate.
z, curDim_ + blockSize_, blockSize_);
592 lclz = rcp (
new SDM (Teuchos::View, *z_, curDim_ + blockSize_, blockSize_));
594 lclz = Teuchos::null;
598 TEUCHOS_TEST_FOR_EXCEPTION(
599 newstate.
V == Teuchos::null,std::invalid_argument,
600 "Belos::BlockFGmresIter::initialize(): BlockFGmresStateIterState does not have initial kernel V_0.");
602 TEUCHOS_TEST_FOR_EXCEPTION(
603 newstate.
z == Teuchos::null,std::invalid_argument,
604 "Belos::BlockFGmresIter::initialize(): BlockFGmresStateIterState does not have initial norms z_0.");
612 template <
class ScalarType,
class MV,
class OP>
615 using Teuchos::Array;
620 typedef Teuchos::SerialDenseMatrix<int, ScalarType> SDM;
623 if (initialized_ ==
false) {
628 const int searchDim = blockSize_ * numBlocks_;
632 while (stest_->checkStatus (
this) !=
Passed && curDim_+blockSize_ <= searchDim) {
636 const int lclDim = curDim_ + blockSize_;
639 std::vector<int> curind (blockSize_);
640 for (
int i = 0; i < blockSize_; ++i) {
641 curind[i] = lclDim + i;
647 for (
int i = 0; i < blockSize_; ++i) {
648 curind[i] = curDim_ + i;
654 lp_->applyRightPrec (*Vprev, *Znext);
658 lp_->applyOp (*Znext, *Vnext);
663 std::vector<int> prevind (lclDim);
664 for (
int i = 0; i < lclDim; ++i) {
668 Array<RCP<const MV> > AVprev (1, Vprev);
671 RCP<SDM> subH = rcp (
new SDM (View, *H_, lclDim, blockSize_, 0, curDim_));
672 Array<RCP<SDM> > AsubH;
676 RCP<SDM> subR = rcp (
new SDM (View, *H_, blockSize_, blockSize_, lclDim, curDim_));
677 const int rank = ortho_->projectAndNormalize (*Vnext, AsubH, subR, AVprev);
678 TEUCHOS_TEST_FOR_EXCEPTION(
680 "Belos::BlockFGmresIter::iterate(): After orthogonalization, the new "
681 "basis block does not have full rank. It contains " << blockSize_
682 <<
" vector" << (blockSize_ != 1 ?
"s" :
"")
683 <<
", but its rank is " << rank <<
".");
695 curDim_ += blockSize_;
700 template<
class ScalarType,
class MV,
class OP>
703 typedef Teuchos::ScalarTraits<ScalarType> STS;
704 typedef Teuchos::ScalarTraits<MagnitudeType> STM;
706 const ScalarType zero = STS::zero ();
707 const ScalarType two = (STS::one () + STS::one());
708 ScalarType sigma, mu, vscale, maxelem;
709 Teuchos::BLAS<int, ScalarType> blas;
715 int curDim = curDim_;
725 if (blockSize_ == 1) {
727 for (
int i = 0; i < curDim; ++i) {
729 blas.ROT (1, &(*H_)(i, curDim), 1, &(*H_)(i+1, curDim), 1, &cs[i], &sn[i]);
732 blas.ROTG (&(*H_)(curDim, curDim), &(*H_)(curDim+1, curDim), &cs[curDim], &sn[curDim]);
733 (*H_)(curDim+1, curDim) = zero;
736 blas.ROT (1, &(*z_)(curDim,0), 1, &(*z_)(curDim+1,0), 1, &cs[curDim], &sn[curDim]);
740 for (
int j = 0; j < blockSize_; ++j) {
742 for (
int i = 0; i < curDim + j; ++i) {
743 sigma = blas.DOT (blockSize_, &(*H_)(i+1,i), 1, &(*H_)(i+1,curDim+j), 1);
744 sigma += (*H_)(i,curDim+j);
746 blas.AXPY (blockSize_, ScalarType(-sigma), &(*H_)(i+1,i), 1, &(*H_)(i+1,curDim+j), 1);
747 (*H_)(i,curDim+j) -= sigma;
751 const int maxidx = blas.IAMAX (blockSize_+1, &(*H_)(curDim+j,curDim+j), 1);
752 maxelem = (*H_)(curDim + j + maxidx - 1, curDim + j);
753 for (
int i = 0; i < blockSize_ + 1; ++i) {
754 (*H_)(curDim+j+i,curDim+j) /= maxelem;
756 sigma = blas.DOT (blockSize_, &(*H_)(curDim + j + 1, curDim + j), 1,
757 &(*H_)(curDim + j + 1, curDim + j), 1);
759 beta[curDim + j] = zero;
761 mu = STS::squareroot ((*H_)(curDim+j,curDim+j)*(*H_)(curDim+j,curDim+j)+sigma);
762 if (STS::real ((*H_)(curDim + j, curDim + j)) < STM::zero ()) {
763 vscale = (*H_)(curDim+j,curDim+j) - mu;
765 vscale = -sigma / ((*H_)(curDim+j, curDim+j) + mu);
767 beta[curDim+j] = two * vscale * vscale / (sigma + vscale*vscale);
768 (*H_)(curDim+j, curDim+j) = maxelem*mu;
769 for (
int i = 0; i < blockSize_; ++i) {
770 (*H_)(curDim+j+1+i,curDim+j) /= vscale;
775 for (
int i = 0; i < blockSize_; ++i) {
776 sigma = blas.DOT (blockSize_, &(*H_)(curDim+j+1,curDim+j),
777 1, &(*z_)(curDim+j+1,i), 1);
778 sigma += (*z_)(curDim+j,i);
779 sigma *= beta[curDim+j];
780 blas.AXPY (blockSize_, ScalarType(-sigma), &(*H_)(curDim+j+1,curDim+j),
781 1, &(*z_)(curDim+j+1,i), 1);
782 (*z_)(curDim+j,i) -= sigma;
789 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.
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
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 ...
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.
Teuchos::RCP< const MV > Z
The current preconditioned Krylov basis (only used in flexible GMRES).
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.