10 #ifndef BELOS_BLOCK_CG_ITER_HPP
11 #define BELOS_BLOCK_CG_ITER_HPP
47 template<
class ScalarType,
class MV,
class OP,
48 const bool lapackSupportsScalarType =
132 template<
class ScalarType,
class MV,
class OP>
330 template<
class ScalarType,
class MV,
class OP>
343 stateStorageInitialized_(false),
347 int bs = params.
get(
"Block Size", 1);
351 template <
class ScalarType,
class MV,
class OP>
354 if (! stateStorageInitialized_) {
359 stateStorageInitialized_ =
false;
371 "Belos::BlockCGIter::setStateSize: LinearProblem lacks "
372 "multivectors from which to clone.");
380 stateStorageInitialized_ =
true;
385 template <
class ScalarType,
class MV,
class OP>
391 (blockSize <= 0, std::invalid_argument,
"Belos::BlockGmresIter::"
392 "setBlockSize: blockSize = " << blockSize <<
" <= 0.");
393 if (blockSize == blockSize_) {
396 if (blockSize!=blockSize_) {
397 stateStorageInitialized_ =
false;
399 blockSize_ = blockSize;
400 initialized_ =
false;
405 template <
class ScalarType,
class MV,
class OP>
409 const char prefix[] =
"Belos::BlockCGIter::initialize: ";
412 if (! stateStorageInitialized_) {
417 (! stateStorageInitialized_, std::invalid_argument,
418 prefix <<
"Cannot initialize state storage!");
421 const char errstr[] =
"Specified multivectors must have a consistent "
431 std::invalid_argument, prefix << errstr );
434 std::invalid_argument, prefix << errstr );
437 if (newstate.
R != R_) {
445 lp_->applyLeftPrec( *R_, *Z_ );
448 lp_->applyRightPrec( *Z_, *tmp );
453 lp_->applyRightPrec( *R_, *Z_ );
463 prefix <<
"BlockCGStateIterState does not have initial residual.");
470 template <
class ScalarType,
class MV,
class OP>
473 const char prefix[] =
"Belos::BlockCGIter::iterate: ";
478 if (initialized_ ==
false) {
492 rHz_old( blockSize_, blockSize_ ), pAp( blockSize_, blockSize_ );
507 prefix <<
"Current linear system does not have the right number of vectors!" );
511 prefix <<
"Failed to compute initial block of orthonormal direction vectors.");
516 while (stest_->checkStatus(
this) !=
Passed) {
521 lp_->applyOp( *P_, *AP_ );
533 lltSolver.factorWithEquilibration(
true );
534 info = lltSolver.factor();
537 prefix <<
"Failed to compute Cholesky factorization using LAPACK routine POTRF.");
541 lltSolver.setVectors (Teuchos::rcpFromRef (alpha), Teuchos::rcpFromRef (alpha));
542 info = lltSolver.solve();
545 prefix <<
"Failed to compute alpha using Cholesky factorization (POTRS).");
549 lp_->updateSolution();
556 lp_->applyLeftPrec( *R_, *Z_ );
559 lp_->applyRightPrec( *Z_, *tmp );
564 lp_->applyRightPrec( *R_, *Z_ );
579 info = lltSolver.solve();
582 prefix <<
"Failed to compute beta using Cholesky factorization (POTRS).");
593 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...