|
Belos Package Browser (Single Doxygen Collection)
Development
|
Go to the documentation of this file.
42 #ifndef BELOS_GCRODR_ITER_HPP
43 #define BELOS_GCRODR_ITER_HPP
87 template <
class ScalarType,
class MV>
113 U(Teuchos::null),
C(Teuchos::null),
114 H(Teuchos::null),
B(Teuchos::null)
153 template<
class ScalarType,
class MV,
class OP>
328 TEUCHOS_TEST_FOR_EXCEPTION(blockSize!=1,std::invalid_argument,
"Belos::GCRODRIter::setBlockSize(): Cannot use a block size that is not one.");
332 void setSize(
int recycledBlocks,
int numBlocks ) {
412 template<
class ScalarType,
class MV,
class OP>
418 lp_(problem), om_(printer), stest_(tester), ortho_(ortho) {
432 int nb = Teuchos::getParameter<int>(params,
"Num Blocks");
434 TEUCHOS_TEST_FOR_EXCEPTION(!params.
isParameter(
"Recycled Blocks"), std::invalid_argument,
"Belos::GCRODRIter::constructor: mandatory parameter \"Recycled Blocks\" is not specified.");
435 int rb = Teuchos::getParameter<int>(params,
"Recycled Blocks");
437 TEUCHOS_TEST_FOR_EXCEPTION(nb <= 0, std::invalid_argument,
"Belos::GCRODRIter() was passed a non-positive argument for \"Num Blocks\".");
438 TEUCHOS_TEST_FOR_EXCEPTION(rb >= nb, std::invalid_argument,
"Belos::GCRODRIter() the number of recycled blocks is larger than the allowable subspace.");
452 template <
class ScalarType,
class MV,
class OP>
460 return currentUpdate;
465 currentUpdate = MVT::Clone( *V_, 1 );
479 std::vector<int> index(curDim_);
480 for (
int i=0; i<curDim_; i++ ) index[i] = i;
482 MVT::MvTimesMatAddMv( one, *Vjp1, y, zero, *currentUpdate );
486 if (U_ != Teuchos::null) {
490 MVT::MvTimesMatAddMv( -one, *U_, z, one, *currentUpdate );
493 return currentUpdate;
500 template <
class ScalarType,
class MV,
class OP>
505 if ( norms && (
int)norms->size()==0 )
510 (*norms)[0] = blas.
NRM2( 1, &z_(curDim_), 1);
512 return Teuchos::null;
519 template <
class ScalarType,
class MV,
class OP>
522 if (newstate.
V != Teuchos::null && newstate.
H != Teuchos::null) {
523 curDim_ = newstate.
curDim;
531 TEUCHOS_TEST_FOR_EXCEPTION(newstate.
V == Teuchos::null,std::invalid_argument,
"Belos::GCRODRIter::initialize(): GCRODRIterState does not have V initialized.");
532 TEUCHOS_TEST_FOR_EXCEPTION(newstate.
H == Teuchos::null,std::invalid_argument,
"Belos::GCRODRIter::initialize(): GCRODRIterState does not have H initialized.");
543 template <
class ScalarType,
class MV,
class OP>
549 setSize( recycledBlocks_, numBlocks_ );
553 std::vector<int> curind(1);
560 Vnext = MVT::CloneViewNonConst(*V_,curind);
565 int rank = ortho_->normalize( *Vnext, z0 );
570 std::vector<int> prevind(numBlocks_+1);
577 if (U_ == Teuchos::null) {
578 while (stest_->checkStatus(
this) !=
Passed && curDim_+1 <= numBlocks_) {
580 int lclDim = curDim_ + 1;
584 Vnext = MVT::CloneViewNonConst(*V_,curind);
588 Vprev = MVT::CloneView(*V_,curind);
591 lp_->apply(*Vprev,*Vnext);
596 prevind.resize(lclDim);
597 for (
int i=0; i<lclDim; i++) { prevind[i] = i; }
598 Vprev = MVT::CloneView(*V_,prevind);
612 rank = ortho_->projectAndNormalize(*Vnext,AsubH,subR,AVprev);
629 while (stest_->checkStatus(
this) !=
Passed && curDim_+1 <= numBlocks_) {
631 int lclDim = curDim_ + 1;
635 Vnext = MVT::CloneViewNonConst(*V_,curind);
640 Vprev = MVT::CloneView(*V_,curind);
643 lp_->apply(*Vprev,*Vnext);
644 Vprev = Teuchos::null;
654 ortho_->project( *Vnext, AsubB,
C );
658 prevind.resize(lclDim);
659 for (
int i=0; i<lclDim; i++) { prevind[i] = i; }
660 Vprev = MVT::CloneView(*V_,prevind);
672 rank = ortho_->projectAndNormalize(*Vnext,AsubH,subR,AVprev);
693 template<
class ScalarType,
class MV,
class OP>
702 int curDim = curDim_;
703 if ( (dim >= curDim_) && (dim < getMaxSubspaceDim()) )
713 for (i=0; i<curDim; i++) {
717 blas.
ROT( 1, &R_(i,curDim), 1, &R_(i+1, curDim), 1, &cs_[i], &sn_[i] );
722 blas.
ROTG( &R_(curDim,curDim), &R_(curDim+1,curDim), &cs_[curDim], &sn_[curDim] );
723 R_(curDim+1,curDim) = zero;
727 blas.
ROT( 1, &z_(curDim), 1, &z_(curDim+1), 1, &cs_[curDim], &sn_[curDim] );
GCRODRIterOrthoFailure is thrown when the GCRODRIter object is unable to compute independent directio...
const Teuchos::RCP< OutputManager< ScalarType > > om_
Belos header file which uses auto-configuration information to include necessary C++ headers.
const LinearProblem< ScalarType, MV, OP > & getProblem() const
Get a constant reference to the linear problem.
const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > stest_
Teuchos::ScalarTraits< ScalarType > SCT
int getNumIters() const
Get the current iteration count.
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...
Class which describes the linear problem to be solved by the iterative solver.
Pure virtual base class for defining the status testing capabilities of Belos.
This class implements the GCRODR iteration, where a single-std::vector Krylov subspace is constructed...
int getNumBlocks() const
Get the maximum number of blocks used by the iterative solver in solving this linear problem.
Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > B
The projection of the Krylov subspace against the recycled subspace.
OperatorTraits< ScalarType, MV, OP > OPT
void resetNumIters(int iter=0)
Reset the iteration count.
Teuchos::SerialDenseVector< int, ScalarType > z_
Teuchos::SerialDenseVector< int, ScalarType > sn_
Structure to contain pointers to GCRODRIter state variables.
const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > lp_
Teuchos::SerialDenseMatrix< int, ScalarType > R_
ScalarType * values() const
bool isInitialized()
States whether the solver has been initialized or not.
void initialize()
Initialize the solver with empty data. Calling this method will result in error, as GCRODRIter must b...
int getRecycledBlocks() const
Get the maximum number of recycled blocks used by the iterative solver in solving this linear problem...
GCRODRIterInitFailure(const std::string &what_arg)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
ScalarTraits< ScalarType >::magnitudeType NRM2(const OrdinalType &n, const ScalarType *x, const OrdinalType &incx) const
Teuchos::RCP< MV > V
The current Krylov basis.
void setBlockSize(int blockSize)
Set the blocksize.
A linear system to solve, and its associated information.
SCT::magnitudeType MagnitudeType
OrdinalType stride() const
Teuchos::RCP< MV > getCurrentUpdate() const
Get the current update to the linear system.
Class which manages the output and verbosity of the Belos solvers.
GCRODRIterInitFailure is thrown when the GCRODRIter object is unable to generate an initial iterate i...
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
void setSize(int recycledBlocks, int numBlocks)
Set the maximum number of blocks used by the iterative solver and the number of recycled vectors.
MultiVecTraits< ScalarType, MV > MVT
void setNumBlocks(int numBlocks)
Set the maximum number of blocks used by the iterative solver.
Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > B_
bool isParameter(const std::string &name) const
void setRecycledBlocks(int recycledBlocks)
Set the maximum number of recycled blocks used by the iterative solver.
Class which defines basic traits for the operator type.
Array< T > & append(const T &x)
Templated virtual class for providing orthogonalization/orthonormalization methods with matrix-based ...
Teuchos::RCP< MV > U
The recycled subspace and its projection.
int curDim
The current dimension of the reduction.
void ROTG(ScalarType *da, ScalarType *db, rotg_c_type *c, ScalarType *s) const
Belos's templated virtual class for providing routines for orthogonalization and orthonormzalition of...
SerialDenseMatrix< OrdinalType, ScalarType > & assign(const SerialDenseMatrix< OrdinalType, ScalarType > &Source)
A pure virtual class for defining the status tests for the Belos iterative solvers.
Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > H
The current Hessenberg matrix.
int sizeUninitialized(OrdinalType length_in)
Parent class to all Belos exceptions.
int shapeUninitialized(OrdinalType numRows, OrdinalType numCols)
Declaration of basic traits for the multivector type.
GCRODRIterState< ScalarType, MV > getState() const
Get the current state of the linear solver.
GCRODRIter(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)
GCRODRIter constructor with linear problem, solver utilities, and parameter list of solver options.
Teuchos::RCP< const MV > getNativeResiduals(std::vector< MagnitudeType > *norms) const
Get the norms of the residuals native to the solver.
Collection of types and exceptions used within the Belos solvers.
GCRODRIterOrthoFailure(const std::string &what_arg)
Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > H_
Teuchos::SerialDenseVector< int, MagnitudeType > cs_
int getBlockSize() const
Get the blocksize to be used by the iterative solver in solving this linear problem.
void updateLSQR(int dim=-1)
Method for updating QR factorization of upper Hessenberg matrix.
int multiply(ETransp transa, ETransp transb, ScalarType alpha, const SerialDenseMatrix< OrdinalType, ScalarType > &A, const SerialDenseMatrix< OrdinalType, ScalarType > &B, ScalarType beta)
const Teuchos::RCP< OrthoManager< ScalarType, MV > > ortho_
int getMaxSubspaceDim() const
Get the maximum dimension allocated for the search subspace.
void iterate()
This method performs block Gmres iterations until the status test indicates the need to stop or an er...
virtual ~GCRODRIter()
Destructor.
Traits class which defines basic operations on multivectors.
int getCurSubspaceDim() const
Get the dimension of the search subspace used to generate the current solution to the linear problem.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
void ROT(const OrdinalType &n, ScalarType *dx, const OrdinalType &incx, ScalarType *dy, const OrdinalType &incy, MagnitudeType *c, ScalarType *s) const