10 #ifndef BELOS_CG_ITER_HPP
11 #define BELOS_CG_ITER_HPP
46 template<
class ScalarType,
class MV,
class OP>
176 "Belos::CGIter::setBlockSize(): Cannot use a block size that is not one.");
287 template<
class ScalarType,
class MV,
class OP>
296 convTest_(convTester),
298 stateStorageInitialized_(false),
300 assertPositiveDefiniteness_( params.get(
"Assert Positive Definiteness", true) ),
301 numEntriesForCondEst_(params.get(
"Max Size For Condest",0) ),
309 template <
class ScalarType,
class MV,
class OP>
312 if (!stateStorageInitialized_) {
318 stateStorageInitialized_ =
false;
329 "Belos::CGIter::setStateSize(): linear problem does not specify multivectors to clone from.");
330 S_ = MVT::Clone( *tmp, 2 );
331 std::vector<int> index(1,0);
333 R_ = MVT::CloneViewNonConst( *S_, index );
335 Z_ = MVT::CloneViewNonConst( *S_, index );
336 P_ = MVT::Clone( *tmp, 1 );
337 AP_ = MVT::Clone( *tmp, 1 );
342 if(numEntriesForCondEst_ > 0) {
343 diag_.resize(numEntriesForCondEst_);
344 offdiag_.resize(numEntriesForCondEst_-1);
348 stateStorageInitialized_ =
true;
356 template <
class ScalarType,
class MV,
class OP>
360 if (!stateStorageInitialized_)
364 "Belos::CGIter::initialize(): Cannot initialize state storage!");
368 std::string errstr(
"Belos::CGIter::initialize(): Specified multivectors must have a consistent length and width.");
373 std::invalid_argument, errstr );
375 std::invalid_argument, errstr );
378 if (newstate.
R != R_) {
380 MVT::Assign( *newstate.
R, *R_ );
387 lp_->applyLeftPrec( *R_, *Z_ );
390 lp_->applyRightPrec( *tmp, *Z_ );
394 lp_->applyRightPrec( *R_, *Z_ );
397 MVT::Assign( *R_, *Z_ );
399 MVT::Assign( *Z_, *P_ );
404 "Belos::CGIter::initialize(): CGIterationState does not have initial residual.");
414 template <
class ScalarType,
class MV,
class OP>
420 if (initialized_ ==
false) {
425 std::vector<ScalarType> alpha(1);
426 std::vector<ScalarType> beta(1);
427 std::vector<ScalarType> rHz(1), rHz_old(1), pAp(1);
435 ScalarType pAp_old = one, beta_old = one, rHz_old2 = one;
442 "Belos::CGIter::iterate(): current linear system has more than one vector!" );
445 if (foldConvergenceDetectionIntoAllreduce_ && convTest_->getResNormType() ==
Belos::TwoNorm) {
446 MVT::MvTransMv( one, *R_, *S_, rHs );
450 MVT::MvDot( *R_, *Z_, rHz );
455 while (stest_->checkStatus(
this) !=
Passed) {
461 lp_->applyOp( *P_, *AP_ );
464 MVT::MvDot( *P_, *AP_, pAp );
465 alpha[0] = rHz[0] / pAp[0];
468 if(assertPositiveDefiniteness_) {
470 "Belos::CGIter::iterate(): non-positive value for p^H*A*p encountered!" );
476 MVT::MvAddMv( one, *cur_soln_vec, alpha[0], *P_, *cur_soln_vec );
477 lp_->updateSolution();
485 MVT::MvAddMv( one, *R_, -alpha[0], *AP_, *R_ );
491 lp_->applyLeftPrec( *R_, *Z_ );
494 lp_->applyRightPrec( *tmp, *Z_ );
498 lp_->applyRightPrec( *R_, *Z_ );
501 MVT::Assign( *R_, *Z_ );
504 if (foldConvergenceDetectionIntoAllreduce_ && convTest_->getResNormType() ==
Belos::TwoNorm) {
505 MVT::MvTransMv( one, *R_, *S_, rHs );
509 MVT::MvDot( *R_, *Z_, rHz );
511 beta[0] = rHz[0] / rHz_old[0];
513 MVT::MvAddMv( one, *Z_, beta[0], *P_, *P_ );
524 rHz_old2 = rHz_old[0];
Teuchos::RCP< const MV > R
The current residual.
Collection of types and exceptions used within the Belos solvers.
const Teuchos::RCP< OutputManager< ScalarType > > om_
Teuchos::ArrayRCP< MagnitudeType > diag_
Teuchos::ScalarTraits< ScalarType > SCT
Belos's basic output manager for sending information of select verbosity levels to the appropriate ou...
const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > stest_
Class which manages the output and verbosity of the Belos solvers.
bool assertPositiveDefiniteness_
void setStateSize()
Method for initalizing the state storage needed by CG.
void setDoCondEst(bool val)
Sets whether or not to store the diagonal for condition estimation.
Teuchos::ArrayRCP< MagnitudeType > offdiag_
bool foldConvergenceDetectionIntoAllreduce_
This class implements the preconditioned Conjugate Gradient (CG) iteration.
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...
Pure virtual base class for defining the status testing capabilities of Belos.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
static magnitudeType real(T a)
Declaration of basic traits for the multivector type.
void setBlockSize(int blockSize)
Set the blocksize to be used by the iterative solver in solving this linear problem.
An implementation of StatusTestResNorm using a family of residual norms.
MultiVecTraits< ScalarType, MV > MVT
A pure virtual class for defining the status tests for the Belos iterative solvers.
CGIter(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< StatusTestGenResNorm< ScalarType, MV, OP > > &convTester, Teuchos::ParameterList ¶ms)
CGIter constructor with linear problem, solver utilities, and parameter list of solver options...
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.
Teuchos::RCP< const MV > getNativeResiduals(std::vector< MagnitudeType > *norms) const
Get the norms of the residuals native to the solver.
CGIterationState< ScalarType, MV > getState() const
Get the current state of the linear solver.
bool isInitialized()
States whether the solver has been initialized or not.
OperatorTraits< ScalarType, MV, OP > OPT
int getNumIters() const
Get the current iteration count.
Traits class which defines basic operations on multivectors.
CGPositiveDefiniteFailure is thrown when the the CG 'alpha = p^H*A*P' value is less than zero...
Teuchos::ArrayView< MagnitudeType > getDiag()
Gets the diagonal for condition estimation.
int getBlockSize() const
Get the blocksize to be used by the iterative solver in solving this linear problem.
void iterate()
This method performs CG iterations until the status test indicates the need to stop or an error occur...
A linear system to solve, and its associated information.
Class which describes the linear problem to be solved by the iterative solver.
bool stateStorageInitialized_
Teuchos::ArrayView< MagnitudeType > getOffDiag()
Gets the off-diagonal for condition estimation.
Teuchos::RCP< const MV > Z
The current preconditioned residual.
virtual ~CGIter()
Destructor.
void initialize()
Initialize the solver with the initial vectors from the linear problem or random data.
Teuchos::RCP< MV > getCurrentUpdate() const
Get the current update to the linear system.
SCT::magnitudeType MagnitudeType
Belos::StatusTestResNorm for specifying general residual norm stopping criteria.
const Teuchos::RCP< StatusTestGenResNorm< ScalarType, MV, OP > > convTest_
int numEntriesForCondEst_
void resetNumIters(int iter=0)
Reset the iteration count.
const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > lp_
const LinearProblem< ScalarType, MV, OP > & getProblem() const
Get a constant reference to the linear problem.
Class which defines basic traits for the operator type.
Belos header file which uses auto-configuration information to include necessary C++ headers...
Teuchos::RCP< const MV > P
The current decent direction vector.
void initializeCG(CGIterationState< ScalarType, MV > &newstate)
Initialize the solver to an iterate, providing a complete state.