42 #ifndef BELOS_BLOCK_CG_ITER_HPP
43 #define BELOS_BLOCK_CG_ITER_HPP
80 template<
class ScalarType,
class MV,
class OP,
81 const bool lapackSupportsScalarType =
165 template<
class ScalarType,
class MV,
class OP>
363 template<
class ScalarType,
class MV,
class OP>
376 stateStorageInitialized_(false),
380 int bs = params.
get(
"Block Size", 1);
384 template <
class ScalarType,
class MV,
class OP>
387 if (! stateStorageInitialized_) {
392 stateStorageInitialized_ =
false;
404 "Belos::BlockCGIter::setStateSize: LinearProblem lacks "
405 "multivectors from which to clone.");
413 stateStorageInitialized_ =
true;
418 template <
class ScalarType,
class MV,
class OP>
424 (blockSize <= 0, std::invalid_argument,
"Belos::BlockGmresIter::"
425 "setBlockSize: blockSize = " << blockSize <<
" <= 0.");
426 if (blockSize == blockSize_) {
429 if (blockSize!=blockSize_) {
430 stateStorageInitialized_ =
false;
432 blockSize_ = blockSize;
433 initialized_ =
false;
438 template <
class ScalarType,
class MV,
class OP>
442 const char prefix[] =
"Belos::BlockCGIter::initialize: ";
445 if (! stateStorageInitialized_) {
450 (! stateStorageInitialized_, std::invalid_argument,
451 prefix <<
"Cannot initialize state storage!");
454 const char errstr[] =
"Specified multivectors must have a consistent "
464 std::invalid_argument, prefix << errstr );
467 std::invalid_argument, prefix << errstr );
470 if (newstate.
R != R_) {
478 lp_->applyLeftPrec( *R_, *Z_ );
481 lp_->applyRightPrec( *Z_, *tmp );
486 lp_->applyRightPrec( *R_, *Z_ );
496 prefix <<
"BlockCGStateIterState does not have initial residual.");
503 template <
class ScalarType,
class MV,
class OP>
506 const char prefix[] =
"Belos::BlockCGIter::iterate: ";
511 if (initialized_ ==
false) {
525 rHz_old( blockSize_, blockSize_ ), pAp( blockSize_, blockSize_ );
540 prefix <<
"Current linear system does not have the right number of vectors!" );
544 prefix <<
"Failed to compute initial block of orthonormal direction vectors.");
549 while (stest_->checkStatus(
this) !=
Passed) {
554 lp_->applyOp( *P_, *AP_ );
566 lltSolver.factorWithEquilibration(
true );
567 info = lltSolver.factor();
570 prefix <<
"Failed to compute Cholesky factorization using LAPACK routine POTRF.");
574 lltSolver.setVectors (Teuchos::rcpFromRef (alpha), Teuchos::rcpFromRef (alpha));
575 info = lltSolver.solve();
578 prefix <<
"Failed to compute alpha using Cholesky factorization (POTRS).");
582 lp_->updateSolution();
589 lp_->applyLeftPrec( *R_, *Z_ );
592 lp_->applyRightPrec( *Z_, *tmp );
597 lp_->applyRightPrec( *R_, *Z_ );
612 info = lltSolver.solve();
615 prefix <<
"Failed to compute beta using Cholesky factorization (POTRS).");
626 prefix <<
"Failed to compute block of orthonormal direction vectors.");
ScalarType * values() const
Teuchos::RCP< const MV > R
The current residual.
Collection of types and exceptions used within the Belos solvers.
Belos's basic output manager for sending information of select verbosity levels to the appropriate ou...
CGIterationState< ScalarType, MV > getState() const
Get the current state of the linear solver.
void setDoCondEst(bool val)
Sets whether or not to store the diagonal for condition estimation.
Teuchos::ArrayView< MagnitudeType > getOffDiag()
Gets the off-diagonal for condition estimation (NOT_IMPLEMENTED)
int getBlockSize() const
Get the blocksize to be used by the iterative solver in solving this linear problem.
Class which manages the output and verbosity of the Belos solvers.
Teuchos::ScalarTraits< ScalarType > SCT
void initialize()
Initialize the solver with the initial vectors from the linear problem or random data.
int getNumIters() const
Get the current iteration count.
const Teuchos::RCP< OutputManager< ScalarType > > om_
BlockCGIter(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &, const Teuchos::RCP< OutputManager< ScalarType > > &, const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &, const Teuchos::RCP< MatOrthoManager< ScalarType, MV, OP > > &, Teuchos::ParameterList &)
const LinearProblem< ScalarType, MV, OP > & getProblem() const
Get a constant reference to the linear problem.
T & get(ParameterList &l, const std::string &name)
Structure to contain pointers to CGIteration state variables.
Pure virtual base class which augments the basic interface for a conjugate gradient linear solver ite...
void resetNumIters(int iter=0)
Reset the iteration count.
Pure virtual base class for defining the status testing capabilities of Belos.
void setDoCondEst(bool)
Sets whether or not to store the diagonal for condition estimation.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
MultiVecTraits< ScalarType, MV > MVT
MultiVecTraits< ScalarType, MV > MVT
Declaration of basic traits for the multivector type.
Teuchos::ScalarTraits< ScalarType > SCT
static void MvTransMv(const ScalarType alpha, const MV &A, const MV &mv, Teuchos::SerialDenseMatrix< int, ScalarType > &B)
Compute a dense matrix B through the matrix-matrix multiply .
A pure virtual class for defining the status tests for the Belos iterative solvers.
Class which defines basic traits for the operator type.
CGIterateFailure is thrown when the CGIteration object is unable to compute the next iterate in the C...
Teuchos::RCP< const MV > AP
The matrix A applied to current decent direction vector.
const LinearProblem< ScalarType, MV, OP > & getProblem() const
Get a constant reference to the linear problem.
static int GetNumberVecs(const MV &mv)
Obtain the number of vectors in mv.
static void Assign(const MV &A, MV &mv)
mv := A
Teuchos::RCP< MV > getCurrentUpdate() const
Get the current update to the linear system.
Traits class which defines basic operations on multivectors.
OperatorTraits< ScalarType, MV, OP > OPT
bool stateStorageInitialized_
void initializeCG(CGIterationState< ScalarType, MV > &)
Initialize the solver to an iterate, providing a complete state.
CGIterationLAPACKFailure is thrown when a nonzero return value is passed back from an LAPACK routine...
static Teuchos::RCP< MV > Clone(const MV &mv, const int numvecs)
Creates a new empty MV containing numvecs columns.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::RCP< const MV > getNativeResiduals(std::vector< MagnitudeType > *) const
Get the norms of the residuals native to the solver.
void iterate()
This method performs linear solver iterations until the status test indicates the need to stop or an ...
A linear system to solve, and its associated information.
int getNumIters() const
Get the current iteration count.
Class which describes the linear problem to be solved by the iterative solver.
SCT::magnitudeType MagnitudeType
SCT::magnitudeType MagnitudeType
Type traits class that says whether Teuchos::LAPACK has a valid implementation for the given ScalarTy...
static Teuchos::RCP< MV > CloneCopy(const MV &mv)
Creates a new MV and copies contents of mv into the new vector (deep copy).
static void MvTimesMatAddMv(const ScalarType alpha, const MV &A, const Teuchos::SerialDenseMatrix< int, ScalarType > &B, const ScalarType beta, MV &mv)
Update mv with .
const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > lp_
Teuchos::RCP< MV > getCurrentUpdate() const
Get the current update to the linear system.
bool isInitialized()
States whether the solver has been initialized or not.
bool isInitialized()
States whether the solver has been initialized or not.
void resetNumIters(int iter=0)
Reset the iteration count to iter.
const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > stest_
const Teuchos::RCP< OrthoManager< ScalarType, MV > > ortho_
Teuchos::RCP< const MV > Z
The current preconditioned residual.
Teuchos::ArrayView< MagnitudeType > getDiag()
Gets the diagonal for condition estimation (NOT_IMPLEMENTED)
OperatorTraits< ScalarType, MV, OP > OPT
static ptrdiff_t GetGlobalLength(const MV &mv)
Return the number of rows in the given multivector mv.
CGIterationOrthoFailure is thrown when the CGIteration object is unable to compute independent direct...
void setBlockSize(int blockSize)
Set the blocksize to be used by the iterative solver in solving this linear problem.
Teuchos::RCP< const MV > getNativeResiduals(std::vector< MagnitudeType > *) const
int getBlockSize() const
Get the block size to be used by the iterative solver in solving this linear problem.
void initialize()
Initialize the solver with the initial vectors from the linear problem or random data.
Class which defines basic traits for the operator type.
virtual ~BlockCGIter()
Destructor.
Belos's templated virtual class for providing routines for orthogonalization and orthonormzalition of...
Belos header file which uses auto-configuration information to include necessary C++ headers...
CGIterationState< ScalarType, MV > getState() const
Get the current state of the linear solver.
Teuchos::RCP< const MV > P
The current decent direction vector.
Templated virtual class for providing orthogonalization/orthonormalization methods with matrix-based ...
Stub implementation of BlockCGIter, for ScalarType types for which Teuchos::LAPACK does NOT have a va...