42 #ifndef BELOS_PSEUDO_BLOCK_GMRES_ITER_HPP
43 #define BELOS_PSEUDO_BLOCK_GMRES_ITER_HPP
89 template <
class ScalarType,
class MV>
101 std::vector<Teuchos::RCP<const MV> >
V;
107 std::vector<Teuchos::RCP<const Teuchos::SerialDenseMatrix<int,ScalarType> > >
H;
109 std::vector<Teuchos::RCP<const Teuchos::SerialDenseMatrix<int,ScalarType> > >
R;
111 std::vector<Teuchos::RCP<const Teuchos::SerialDenseVector<int,ScalarType> > >
Z;
113 std::vector<Teuchos::RCP<const Teuchos::SerialDenseVector<int,ScalarType> > >
sn;
114 std::vector<Teuchos::RCP<const Teuchos::SerialDenseVector<int,MagnitudeType> > >
cs;
154 template<
class ScalarType,
class MV,
class OP>
262 for (
int i=0; i<
numRHS_; ++i) {
266 state.
sn[i] =
sn_[i];
267 state.
cs[i] =
cs_[i];
341 "Belos::PseudoBlockGmresIter::setBlockSize(): Cannot use a block size that is not one.");
374 std::vector<Teuchos::RCP<Teuchos::SerialDenseVector<int,ScalarType> > >
sn_;
375 std::vector<Teuchos::RCP<Teuchos::SerialDenseVector<int,MagnitudeType> > >
cs_;
397 std::vector<Teuchos::RCP<MV> >
V_;
402 std::vector<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > >
H_;
407 std::vector<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > >
R_;
408 std::vector<Teuchos::RCP<Teuchos::SerialDenseVector<int,ScalarType> > >
Z_;
413 template<
class ScalarType,
class MV,
class OP>
431 "Belos::PseudoBlockGmresIter::constructor: mandatory parameter 'Num Blocks' is not specified.");
432 int nb = Teuchos::getParameter<int>(params,
"Num Blocks");
439 template <
class ScalarType,
class MV,
class OP>
445 TEUCHOS_TEST_FOR_EXCEPTION(numBlocks <= 0, std::invalid_argument,
"Belos::PseudoBlockGmresIter::setNumBlocks was passed a non-positive argument.");
447 numBlocks_ = numBlocks;
450 initialized_ =
false;
455 template <
class ScalarType,
class MV,
class OP>
464 return currentUpdate;
466 currentUpdate = MVT::Clone(*(V_[0]), numRHS_);
467 std::vector<int> index(1), index2(curDim_);
468 for (
int i=0; i<curDim_; ++i) {
475 for (
int i=0; i<numRHS_; ++i) {
477 Teuchos::RCP<MV> cur_block_copy_vec = MVT::CloneViewNonConst( *currentUpdate, index );
487 H_[i]->values(), H_[i]->stride(), y.
values(), y.
stride() );
490 MVT::MvTimesMatAddMv( one, *Vjp1, y, zero, *cur_block_copy_vec );
493 return currentUpdate;
500 template <
class ScalarType,
class MV,
class OP>
511 if (static_cast<int> (norms->size()) < numRHS_)
512 norms->resize (numRHS_);
515 for (
int j = 0; j < numRHS_; ++j)
517 const ScalarType curNativeResid = (*Z_[j])(curDim_);
518 (*norms)[j] = STS::magnitude (curNativeResid);
525 template <
class ScalarType,
class MV,
class OP>
534 this->numRHS_ = MVT::GetNumberVecs (*(lp_->getCurrLHSVec()));
540 std::string errstr (
"Belos::PseudoBlockGmresIter::initialize(): "
541 "Specified multivectors must have a consistent "
542 "length and width.");
546 std::invalid_argument,
547 "Belos::PseudoBlockGmresIter::initialize(): "
548 "V and/or Z was not specified in the input state; "
549 "the V and/or Z arrays have length zero.");
560 RCP<const MV> lhsMV = lp_->getLHS();
561 RCP<const MV> rhsMV = lp_->getRHS();
565 RCP<const MV> vectorInBasisSpace = rhsMV.is_null() ? lhsMV : rhsMV;
569 std::invalid_argument,
570 "Belos::PseudoBlockGmresIter::initialize(): "
571 "The linear problem to solve does not specify multi"
572 "vectors from which we can clone basis vectors. The "
573 "right-hand side(s), left-hand side(s), or both should "
579 std::invalid_argument,
581 curDim_ = newstate.
curDim;
587 for (
int i=0; i<numRHS_; ++i) {
591 if (V_[i].
is_null() || MVT::GetNumberVecs(*V_[i]) < numBlocks_ + 1) {
592 V_[i] = MVT::Clone (*vectorInBasisSpace, numBlocks_ + 1);
596 TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetGlobalLength(*newstate.
V[i]) != MVT::GetGlobalLength(*V_[i]),
597 std::invalid_argument, errstr );
598 TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetNumberVecs(*newstate.
V[i]) < newstate.
curDim,
599 std::invalid_argument, errstr );
604 int lclDim = MVT::GetNumberVecs(*newstate.
V[i]);
605 if (newstate.
V[i] != V_[i]) {
607 if (curDim_ == 0 && lclDim > 1) {
609 <<
"Belos::PseudoBlockGmresIter::initialize(): the solver was "
610 <<
"initialized with a kernel of " << lclDim
612 <<
"The block size however is only " << 1
614 <<
"The last " << lclDim - 1 <<
" vectors will be discarded."
617 std::vector<int> nevind (curDim_ + 1);
618 for (
int j = 0; j < curDim_ + 1; ++j)
621 RCP<const MV> newV = MVT::CloneView (*newstate.
V[i], nevind);
622 RCP<MV> lclV = MVT::CloneViewNonConst( *V_[i], nevind );
625 MVT::MvAddMv (one, *newV, zero, *newV, *lclV);
635 for (
int i=0; i<numRHS_; ++i) {
640 if (Z_[i]->length() < numBlocks_+1) {
641 Z_[i]->shapeUninitialized(numBlocks_+1, 1);
645 TEUCHOS_TEST_FOR_EXCEPTION(newstate.
Z[i]->numRows() < curDim_, std::invalid_argument, errstr);
648 if (newstate.
Z[i] != Z_[i]) {
665 for (
int i=0; i<numRHS_; ++i) {
670 if (H_[i]->numRows() < numBlocks_+1 || H_[i]->numCols() < numBlocks_) {
671 H_[i]->shapeUninitialized(numBlocks_+1, numBlocks_);
675 if ((
int)newstate.
H.size() == numRHS_) {
678 TEUCHOS_TEST_FOR_EXCEPTION((newstate.
H[i]->numRows() < curDim_ || newstate.
H[i]->numCols() < curDim_), std::invalid_argument,
679 "Belos::PseudoBlockGmresIter::initialize(): Specified Hessenberg matrices must have a consistent size to the current subspace dimension");
681 if (newstate.
H[i] != H_[i]) {
702 if ((
int)newstate.
cs.size() == numRHS_ && (int)newstate.
sn.size() == numRHS_) {
703 for (
int i=0; i<numRHS_; ++i) {
704 if (cs_[i] != newstate.
cs[i])
706 if (sn_[i] != newstate.
sn[i])
712 for (
int i=0; i<numRHS_; ++i) {
716 cs_[i]->resize(numBlocks_+1);
720 sn_[i]->resize(numBlocks_+1);
742 template <
class ScalarType,
class MV,
class OP>
748 if (initialized_ ==
false) {
756 int searchDim = numBlocks_;
761 std::vector<int> index(1);
762 std::vector<int> index2(1);
769 for (
int i=0; i<numRHS_; ++i) {
773 MVT::MvAddMv( one, *tmp_vec, zero, *tmp_vec, *U_vec_view );
781 while (stest_->checkStatus(
this) !=
Passed && curDim_ < searchDim) {
787 lp_->apply( *U_vec, *AU_vec );
792 int num_prev = curDim_+1;
793 index.resize( num_prev );
794 for (
int i=0; i<num_prev; ++i) {
800 for (
int i=0; i<numRHS_; ++i) {
826 ortho_->projectAndNormalize( *V_new, h_array, r_new, V_array );
831 index2[0] = curDim_+1;
833 MVT::MvAddMv( one, *V_new, zero, *V_new, *tmp_vec );
875 template<
class ScalarType,
class MV,
class OP>
881 int curDim = curDim_;
882 if (dim >= curDim_ && dim < getMaxSubspaceDim()) {
891 for (i=0; i<numRHS_; ++i) {
897 for (j=0; j<curDim; j++) {
901 blas.
ROT( 1, &(*H_[i])(j,curDim), 1, &(*H_[i])(j+1, curDim), 1, &(*cs_[i])[j], &(*sn_[i])[j] );
906 blas.
ROTG( &(*H_[i])(curDim,curDim), &(*H_[i])(curDim+1,curDim), &(*cs_[i])[curDim], &(*sn_[i])[curDim] );
907 (*H_[i])(curDim+1,curDim) = zero;
911 blas.
ROT( 1, &(*Z_[i])(curDim), 1, &(*Z_[i])(curDim+1), 1, &(*cs_[i])[curDim], &(*sn_[i])[curDim] );
ScalarType * values() const
std::vector< Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > > H
The current Hessenberg matrix.
int getCurSubspaceDim() const
Get the dimension of the search subspace used to generate the current solution to the linear problem...
Collection of types and exceptions used within the Belos solvers.
void TRSM(ESide side, EUplo uplo, ETransp transa, EDiag diag, const OrdinalType &m, const OrdinalType &n, const alpha_type alpha, const A_type *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb) const
Teuchos::RCP< MV > U_vec_
Belos's basic output manager for sending information of select verbosity levels to the appropriate ou...
Class which manages the output and verbosity of the Belos solvers.
bool is_null(const boost::shared_ptr< T > &p)
PseudoBlockGmresIterState< ScalarType, MV > getState() const
Get the current state of the linear solver.
Teuchos::RCP< MV > cur_block_sol_
std::vector< Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > > R
The current upper-triangular matrix from the QR reduction of H.
Teuchos::ScalarTraits< ScalarType > SCT
void ROT(const OrdinalType &n, ScalarType *dx, const OrdinalType &incx, ScalarType *dy, const OrdinalType &incy, MagnitudeType *c, ScalarType *s) const
Pure virtual base class for defining the status testing capabilities of Belos.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Declaration of basic traits for the multivector type.
int getNumBlocks() const
Get the maximum number of blocks used by the iterative solver in solving this linear problem...
virtual ~PseudoBlockGmresIter()
Destructor.
std::vector< Teuchos::RCP< const MV > > V
The current Krylov basis.
Teuchos::RCP< const MV > getNativeResiduals(std::vector< MagnitudeType > *norms) const
Get the norms of the "native" residual vectors.
std::vector< Teuchos::RCP< Teuchos::SerialDenseVector< int, MagnitudeType > > > cs_
Teuchos::RCP< MV > AU_vec_
Pure virtual base class which describes the basic interface to the linear solver iteration.
const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > lp_
A pure virtual class for defining the status tests for the Belos iterative solvers.
Class which defines basic traits for the operator type.
int getNumIters() const
Get the current iteration count.
Structure to contain pointers to PseudoBlockGmresIter state variables.
const LinearProblem< ScalarType, MV, OP > & getProblem() const
Get a constant reference to the linear problem.
void initialize()
Initialize the solver with the initial vectors from the linear problem or random data.
bool isInitialized()
States whether the solver has been initialized or not.
int getBlockSize() const
Get the blocksize to be used by the iterative solver in solving this linear problem.
Traits class which defines basic operations on multivectors.
Teuchos::RCP< MV > getCurrentUpdate() const
Get the current update to the linear system.
bool isParameter(const std::string &name) const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::RCP< MV > cur_block_rhs_
void updateLSQR(int dim=-1)
Method for updating QR factorization of upper Hessenberg matrix.
PseudoBlockGmresIter(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)
PseudoBlockGmresIter constructor with linear problem, solver utilities, and parameter list of solver ...
std::vector< Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > > R_
PseudoBlockGmresIterState()
PseudoBlockGmresIterOrthoFailure(const std::string &what_arg)
std::vector< Teuchos::RCP< const Teuchos::SerialDenseVector< int, MagnitudeType > > > cs
A linear system to solve, and its associated information.
Class which describes the linear problem to be solved by the iterative solver.
const Teuchos::RCP< OrthoManager< ScalarType, MV > > ortho_
void setNumBlocks(int numBlocks)
Set the maximum number of blocks used by the iterative solver.
PseudoBlockGmresIterInitFailure(const std::string &what_arg)
SCT::magnitudeType MagnitudeType
Teuchos::ScalarTraits< ScalarType > SCT
std::vector< Teuchos::RCP< Teuchos::SerialDenseVector< int, ScalarType > > > Z_
std::vector< Teuchos::RCP< MV > > V_
int getMaxSubspaceDim() const
Get the maximum dimension allocated for the search subspace.
void resetNumIters(int iter=0)
Reset the iteration count.
MultiVecTraits< ScalarType, MV > MVT
std::vector< Teuchos::RCP< const Teuchos::SerialDenseVector< int, ScalarType > > > Z
The current right-hand side of the least squares system RY = Z.
OperatorTraits< ScalarType, MV, OP > OPT
void ROTG(ScalarType *da, ScalarType *db, rotg_c_type *c, ScalarType *s) const
void iterate()
This method performs block Gmres iterations until the status test indicates the need to stop or an er...
Class which defines basic traits for the operator type.
SCT::magnitudeType MagnitudeType
const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > stest_
int curDim
The current dimension of the reduction.
Parent class to all Belos exceptions.
std::vector< Teuchos::RCP< const Teuchos::SerialDenseVector< int, ScalarType > > > sn
The current Given's rotation coefficients.
Belos's templated virtual class for providing routines for orthogonalization and orthonormzalition of...
const Teuchos::RCP< OutputManager< ScalarType > > om_
Belos header file which uses auto-configuration information to include necessary C++ headers...
std::vector< Teuchos::RCP< Teuchos::SerialDenseVector< int, ScalarType > > > sn_
This class implements the pseudo-block GMRES iteration, where a block Krylov subspace is constructed ...
SerialDenseMatrix< OrdinalType, ScalarType > & assign(const SerialDenseMatrix< OrdinalType, ScalarType > &Source)
std::vector< Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > > H_
void setBlockSize(int blockSize)
Set the blocksize.
PseudoBlockGmresIterOrthoFailure is thrown when the orthogonalization manager is unable to generate o...
OrdinalType stride() const
Templated virtual class for providing orthogonalization/orthonormalization methods with matrix-based ...
PseudoBlockGmresIterInitFailure is thrown when the PseudoBlockGmresIter object is unable to generate ...