42#ifndef BELOS_BLOCK_FGMRES_ITER_HPP
43#define BELOS_BLOCK_FGMRES_ITER_HPP
60#include "Teuchos_BLAS.hpp"
61#include "Teuchos_SerialDenseMatrix.hpp"
62#include "Teuchos_SerialDenseVector.hpp"
63#include "Teuchos_ScalarTraits.hpp"
64#include "Teuchos_ParameterList.hpp"
65#include "Teuchos_TimeMonitor.hpp"
82template<
class ScalarType,
class MV,
class OP>
92 typedef Teuchos::ScalarTraits<ScalarType>
SCT;
111 Teuchos::ParameterList ¶ms );
260 void setSize(
int blockSize,
int numBlocks);
278 const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> >
lp_;
279 const Teuchos::RCP<OutputManager<ScalarType> >
om_;
280 const Teuchos::RCP<StatusTest<ScalarType,MV,OP> >
stest_;
281 const Teuchos::RCP<OrthoManager<ScalarType,MV> >
ortho_;
293 Teuchos::SerialDenseVector<int,ScalarType>
beta,
sn;
294 Teuchos::SerialDenseVector<int,MagnitudeType>
cs;
321 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> >
H_;
326 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> >
R_;
327 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> >
z_;
332 template<
class ScalarType,
class MV,
class OP>
338 Teuchos::ParameterList ¶ms ):
351 TEUCHOS_TEST_FOR_EXCEPTION(
352 ! params.isParameter (
"Num Blocks"), std::invalid_argument,
353 "Belos::BlockFGmresIter::constructor: mandatory parameter 'Num Blocks' is not specified.");
354 const int nb = params.get<
int> (
"Num Blocks");
357 const int bs = params.get (
"Block Size", 1);
363 template <
class ScalarType,
class MV,
class OP>
369 TEUCHOS_TEST_FOR_EXCEPTION(numBlocks <= 0 || blockSize <= 0, std::invalid_argument,
"Belos::BlockFGmresIter::setSize was passed a non-positive argument.");
391 template <
class ScalarType,
class MV,
class OP>
396 typedef Teuchos::SerialDenseMatrix<int, ScalarType> SDM;
400 RCP<const MV> lhsMV =
lp_->getLHS();
401 RCP<const MV> rhsMV =
lp_->getRHS();
402 if (lhsMV == Teuchos::null && rhsMV == Teuchos::null) {
421 TEUCHOS_TEST_FOR_EXCEPTION(
423 std::invalid_argument,
"Belos::BlockFGmresIter::setStateSize(): "
424 "Cannot generate a Krylov basis with dimension larger the operator!");
427 if (
V_ == Teuchos::null) {
429 RCP<const MV> tmp = (rhsMV != Teuchos::null) ? rhsMV : lhsMV;
430 TEUCHOS_TEST_FOR_EXCEPTION(
431 tmp == Teuchos::null, std::invalid_argument,
432 "Belos::BlockFGmresIter::setStateSize(): "
433 "linear problem does not specify multivectors to clone from.");
439 RCP<const MV> tmp =
V_;
444 if (
Z_ == Teuchos::null) {
446 RCP<const MV> tmp = (rhsMV != Teuchos::null) ? rhsMV : lhsMV;
447 TEUCHOS_TEST_FOR_EXCEPTION(
448 tmp == Teuchos::null, std::invalid_argument,
449 "Belos::BlockFGmresIter::setStateSize(): "
450 "linear problem does not specify multivectors to clone from.");
456 RCP<const MV> tmp =
Z_;
462 if (
H_ == Teuchos::null) {
472 if (
z_ == Teuchos::null) {
486 template <
class ScalarType,
class MV,
class OP>
490 typedef Teuchos::SerialDenseMatrix<int, ScalarType> SDM;
492 Teuchos::RCP<MV> currentUpdate = Teuchos::null;
496 return currentUpdate;
499 const ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero ();
500 const ScalarType one = Teuchos::ScalarTraits<ScalarType>::one ();
501 Teuchos::BLAS<int,ScalarType> blas;
509 blas.TRSM (Teuchos::LEFT_SIDE, Teuchos::UPPER_TRI, Teuchos::NO_TRANS,
511 H_->values (),
H_->stride (), y.values (), y.stride ());
514 std::vector<int> index (
curDim_);
515 for (
int i = 0; i <
curDim_; ++i) {
521 return currentUpdate;
525 template <
class ScalarType,
class MV,
class OP>
526 Teuchos::RCP<const MV>
531 if (norms != NULL && (
int)norms->size() <
blockSize_) {
536 Teuchos::BLAS<int, ScalarType> blas;
543 return Teuchos::null;
547 template <
class ScalarType,
class MV,
class OP>
554 typedef Teuchos::ScalarTraits<ScalarType> STS;
555 typedef Teuchos::SerialDenseMatrix<int, ScalarType> SDM;
556 const ScalarType ZERO = STS::zero ();
557 const ScalarType ONE = STS::one ();
564 TEUCHOS_TEST_FOR_EXCEPTION(
566 "Belos::BlockFGmresIter::initialize(): Cannot initialize state storage!");
571 const char errstr[] =
"Belos::BlockFGmresIter::initialize(): The given "
572 "multivectors must have a consistent length and width.";
574 if (! newstate.
V.is_null () && ! newstate.
z.is_null ()) {
578 TEUCHOS_TEST_FOR_EXCEPTION(
580 std::invalid_argument, errstr );
581 TEUCHOS_TEST_FOR_EXCEPTION(
583 std::invalid_argument, errstr );
584 TEUCHOS_TEST_FOR_EXCEPTION(
586 std::invalid_argument, errstr );
592 TEUCHOS_TEST_FOR_EXCEPTION(
594 std::invalid_argument, errstr);
597 if (newstate.
V !=
V_) {
601 warn <<
"Belos::BlockFGmresIter::initialize(): the solver was "
602 <<
"initialized with a kernel of " << lclDim << endl
603 <<
"The block size however is only " <<
blockSize_ << endl
605 <<
" vectors will be discarded." << endl;
616 lclV = Teuchos::null;
620 if (newstate.
z !=
z_) {
626 lclz = Teuchos::null;
630 TEUCHOS_TEST_FOR_EXCEPTION(
631 newstate.
V == Teuchos::null,std::invalid_argument,
632 "Belos::BlockFGmresIter::initialize(): BlockFGmresStateIterState does not have initial kernel V_0.");
634 TEUCHOS_TEST_FOR_EXCEPTION(
635 newstate.
z == Teuchos::null,std::invalid_argument,
636 "Belos::BlockFGmresIter::initialize(): BlockFGmresStateIterState does not have initial norms z_0.");
644 template <
class ScalarType,
class MV,
class OP>
647 using Teuchos::Array;
652 typedef Teuchos::SerialDenseMatrix<int, ScalarType> SDM;
673 curind[i] = lclDim + i;
686 lp_->applyRightPrec (*Vprev, *Znext);
690 lp_->applyOp (*Znext, *Vnext);
695 std::vector<int> prevind (lclDim);
696 for (
int i = 0; i < lclDim; ++i) {
700 Array<RCP<const MV> > AVprev (1, Vprev);
704 Array<RCP<SDM> > AsubH;
709 const int rank =
ortho_->projectAndNormalize (*Vnext, AsubH, subR, AVprev);
710 TEUCHOS_TEST_FOR_EXCEPTION(
712 "Belos::BlockFGmresIter::iterate(): After orthogonalization, the new "
713 "basis block does not have full rank. It contains " <<
blockSize_
715 <<
", but its rank is " << rank <<
".");
732 template<
class ScalarType,
class MV,
class OP>
735 typedef Teuchos::ScalarTraits<ScalarType> STS;
736 typedef Teuchos::ScalarTraits<MagnitudeType> STM;
738 const ScalarType zero = STS::zero ();
739 const ScalarType two = (STS::one () + STS::one());
740 ScalarType sigma, mu, vscale, maxelem;
741 Teuchos::BLAS<int, ScalarType> blas;
759 for (
int i = 0; i < curDim; ++i) {
761 blas.ROT (1, &(*
H_)(i, curDim), 1, &(*
H_)(i+1, curDim), 1, &
cs[i], &
sn[i]);
764 blas.ROTG (&(*
H_)(curDim, curDim), &(*
H_)(curDim+1, curDim), &
cs[curDim], &
sn[curDim]);
765 (*H_)(curDim+1, curDim) = zero;
768 blas.ROT (1, &(*
z_)(curDim,0), 1, &(*
z_)(curDim+1,0), 1, &
cs[curDim], &
sn[curDim]);
774 for (
int i = 0; i < curDim + j; ++i) {
775 sigma = blas.DOT (
blockSize_, &(*
H_)(i+1,i), 1, &(*
H_)(i+1,curDim+j), 1);
776 sigma += (*H_)(i,curDim+j);
778 blas.AXPY (
blockSize_, ScalarType(-sigma), &(*
H_)(i+1,i), 1, &(*
H_)(i+1,curDim+j), 1);
779 (*H_)(i,curDim+j) -= sigma;
783 const int maxidx = blas.IAMAX (
blockSize_+1, &(*
H_)(curDim+j,curDim+j), 1);
784 maxelem = (*H_)(curDim + j + maxidx - 1, curDim + j);
786 (*H_)(curDim+j+i,curDim+j) /= maxelem;
788 sigma = blas.DOT (
blockSize_, &(*
H_)(curDim + j + 1, curDim + j), 1,
789 &(*
H_)(curDim + j + 1, curDim + j), 1);
791 beta[curDim + j] = zero;
793 mu = STS::squareroot ((*
H_)(curDim+j,curDim+j)*(*
H_)(curDim+j,curDim+j)+sigma);
794 if (STS::real ((*
H_)(curDim + j, curDim + j)) < STM::zero ()) {
795 vscale = (*H_)(curDim+j,curDim+j) - mu;
797 vscale = -sigma / ((*H_)(curDim+j, curDim+j) + mu);
799 beta[curDim+j] = two * vscale * vscale / (sigma + vscale*vscale);
800 (*H_)(curDim+j, curDim+j) = maxelem*mu;
802 (*H_)(curDim+j+1+i,curDim+j) /= vscale;
808 sigma = blas.DOT (
blockSize_, &(*
H_)(curDim+j+1,curDim+j),
809 1, &(*
z_)(curDim+j+1,i), 1);
810 sigma += (*z_)(curDim+j,i);
811 sigma *=
beta[curDim+j];
812 blas.AXPY (
blockSize_, ScalarType(-sigma), &(*
H_)(curDim+j+1,curDim+j),
813 1, &(*
z_)(curDim+j+1,i), 1);
814 (*z_)(curDim+j,i) -= sigma;
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.
bool stateStorageInitialized_
int getNumIters() const
Get the current iteration count.
void setStateSize()
Method for initalizing the state storage needed by block flexible GMRES.
int getBlockSize() const
Get the blocksize to be used by the iterative solver in solving this linear problem.
Teuchos::SerialDenseVector< int, ScalarType > beta
OperatorTraits< ScalarType, MV, OP > OPT
Teuchos::SerialDenseVector< int, MagnitudeType > cs
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...
Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > z_
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.
Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > R_
Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > H_
Teuchos::SerialDenseVector< int, ScalarType > sn
const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > stest_
SCT::magnitudeType MagnitudeType
const Teuchos::RCP< OutputManager< ScalarType > > om_
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.
const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > lp_
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...
const Teuchos::RCP< OrthoManager< ScalarType, MV > > ortho_
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.