10#ifndef BELOS_CG_SINGLE_RED_ITER_HPP
11#define BELOS_CG_SINGLE_RED_ITER_HPP
28#include "Teuchos_SerialDenseMatrix.hpp"
29#include "Teuchos_SerialDenseVector.hpp"
30#include "Teuchos_ScalarTraits.hpp"
31#include "Teuchos_ParameterList.hpp"
32#include "Teuchos_TimeMonitor.hpp"
53 template <
class ScalarType,
class MV>
65 void initialize(Teuchos::RCP<const MV> tmp,
int _numVectors) {
68 TEUCHOS_ASSERT(_numVectors == 1);
71 W = MVT::Clone( *tmp, 3 );
72 std::vector<int> index2(2,0);
73 std::vector<int> index(1,0);
78 S = MVT::CloneViewNonConst( *
W, index2 );
83 U = MVT::CloneViewNonConst( *
W, index2 );
86 this->
R = MVT::CloneViewNonConst( *
W, index );
88 AZ = MVT::CloneViewNonConst( *
W, index );
90 this->
Z = MVT::CloneViewNonConst( *
W, index );
95 T = MVT::CloneViewNonConst( *
W, index2 );
98 V = MVT::Clone( *tmp, 2 );
100 this->
AP = MVT::CloneViewNonConst( *
V, index );
102 this->
P = MVT::CloneViewNonConst( *
V, index );
107 bool matches(Teuchos::RCP<const MV> tmp,
int _numVectors=1)
const {
126template<
class ScalarType,
class MV,
class OP>
136 using SCT = Teuchos::ScalarTraits<ScalarType>;
151 Teuchos::ParameterList ¶ms );
205 Teuchos::RCP<CGIterationStateBase<ScalarType,MV> >
getState()
const {
221 auto s = Teuchos::rcp_dynamic_cast<CGSingleRedIterationState<ScalarType,MV> >(state,
true);
268 TEUCHOS_TEST_FOR_EXCEPTION(blockSize!=1,std::invalid_argument,
269 "Belos::CGSingleRedIter::setBlockSize(): Cannot use a block size that is not one.");
280 Teuchos::ArrayView<MagnitudeType> temp;
286 Teuchos::ArrayView<MagnitudeType> temp;
300 const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > lp_;
301 const Teuchos::RCP<OutputManager<ScalarType> > om_;
302 const Teuchos::RCP<StatusTest<ScalarType,MV,OP> > stest_;
303 const Teuchos::RCP<StatusTestGenResNorm<ScalarType,MV,OP> > convTest_;
317 bool foldConvergenceDetectionIntoAllreduce_;
332 Teuchos::RCP<MV> AZ_;
336 Teuchos::RCP<MV> AP_;
353 template<
class ScalarType,
class MV,
class OP>
358 Teuchos::ParameterList ¶ms ):
362 convTest_(convTester),
366 foldConvergenceDetectionIntoAllreduce_ = params.get<
bool>(
"Fold Convergence Detection Into Allreduce",
false);
371 template <
class ScalarType,
class MV,
class OP>
375 Teuchos::RCP<const MV> lhsMV = lp_->getLHS();
376 Teuchos::RCP<const MV> rhsMV = lp_->getRHS();
377 Teuchos::RCP<const MV> tmp = ( (rhsMV!=Teuchos::null)? rhsMV: lhsMV );
378 TEUCHOS_ASSERT(!newstate.is_null());
380 newstate->initialize(tmp, 1);
383 std::string errstr(
"Belos::CGSingleRedIter::initialize(): Specified multivectors must have a consistent length and width.");
388 std::invalid_argument, errstr );
390 std::invalid_argument, errstr );
401 if ( lp_->getLeftPrec() != Teuchos::null ) {
402 lp_->applyLeftPrec( *R_, *Z_ );
403 if ( lp_->getRightPrec() != Teuchos::null ) {
405 lp_->applyRightPrec( *Z_, *tmp2 );
409 else if ( lp_->getRightPrec() != Teuchos::null ) {
410 lp_->applyRightPrec( *R_, *Z_ );
417 lp_->applyOp( *Z_, *AZ_ );
431 template <
class ScalarType,
class MV,
class OP>
432 Teuchos::RCP<const MV>
435 (*norms)[0] = std::sqrt(Teuchos::ScalarTraits<ScalarType>::magnitude(rHz_));
436 return Teuchos::null;
437 }
else if (foldConvergenceDetectionIntoAllreduce_ && convTest_->getResNormType() ==
Belos::TwoNorm) {
438 (*norms)[0] = std::sqrt(Teuchos::ScalarTraits<ScalarType>::magnitude(rHr_));
439 return Teuchos::null;
447 template <
class ScalarType,
class MV,
class OP>
458 Teuchos::SerialDenseMatrix<int,ScalarType> sHz( 2, 1 );
459 Teuchos::SerialDenseMatrix<int,ScalarType> sHt( 2, 2 );
466 const ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
467 const MagnitudeType zero = Teuchos::ScalarTraits<MagnitudeType>::zero();
470 Teuchos::RCP<MV> cur_soln_vec = lp_->getCurrLHSVec();
474 "Belos::CGSingleRedIter::iterate(): current linear system has more than one vector!" );
476 if (foldConvergenceDetectionIntoAllreduce_ && convTest_->getResNormType() ==
Belos::TwoNorm) {
488 if ((Teuchos::ScalarTraits<ScalarType>::magnitude(delta) < Teuchos::ScalarTraits<ScalarType>::eps()) &&
489 (stest_->checkStatus(
this) ==
Passed))
491 alpha = rHz_ / delta;
495 "Belos::CGSingleRedIter::iterate(): non-positive value for p^H*A*p encountered!" );
500 if (foldConvergenceDetectionIntoAllreduce_ && convTest_->getResNormType() ==
Belos::TwoNorm) {
508 MVT::MvAddMv( one, *cur_soln_vec, alpha, *P_, *cur_soln_vec );
516 if ( lp_->getLeftPrec() != Teuchos::null ) {
517 lp_->applyLeftPrec( *R_, *Z_ );
518 if ( lp_->getRightPrec() != Teuchos::null ) {
520 lp_->applyRightPrec( *tmp, *Z_ );
523 else if ( lp_->getRightPrec() != Teuchos::null ) {
524 lp_->applyRightPrec( *R_, *Z_ );
531 lp_->applyOp( *Z_, *AZ_ );
547 if (stest_->checkStatus(
this) ==
Passed) {
551 beta = rHz_ / rHz_old;
552 alpha = rHz_ / (delta - (beta*rHz_ / alpha));
556 "Belos::CGSingleRedIter::iterate(): non-positive value for p^H*A*p encountered!" );
574 MVT::MvAddMv( one, *cur_soln_vec, alpha, *P_, *cur_soln_vec );
585 if (stest_->checkStatus(
this) ==
Passed) {
591 if ( lp_->getLeftPrec() != Teuchos::null ) {
592 lp_->applyLeftPrec( *R_, *Z_ );
593 if ( lp_->getRightPrec() != Teuchos::null ) {
595 lp_->applyRightPrec( *tmp, *Z_ );
598 else if ( lp_->getRightPrec() != Teuchos::null ) {
599 lp_->applyRightPrec( *R_, *Z_ );
606 lp_->applyOp( *Z_, *AZ_ );
616 beta = rHz_ / rHz_old;
617 alpha = rHz_ / (delta - (beta*rHz_ / alpha));
621 "Belos::CGSingleRedIter::iterate(): non-positive value for p^H*A*p encountered!" );
Pure virtual base class which augments the basic interface for a conjugate gradient linear solver ite...
Belos header file which uses auto-configuration information to include necessary C++ headers.
Class which describes the linear problem to be solved by the iterative solver.
Declaration of basic traits for the multivector type.
Class which defines basic traits for the operator type.
Class which manages the output and verbosity of the Belos solvers.
Belos::StatusTestResNorm for specifying general residual norm stopping criteria.
Pure virtual base class for defining the status testing capabilities of Belos.
Collection of types and exceptions used within the Belos solvers.
CGIterateFailure is thrown when the CGIteration object is unable to compute the next iterate in the C...
Structure to contain pointers to CGIteration state variables.
Teuchos::RCP< MV > R
The current residual.
Teuchos::RCP< MV > P
The current decent direction vector.
virtual void initialize(Teuchos::RCP< const MV > tmp, int _numVectors)
virtual bool matches(Teuchos::RCP< const MV > tmp, int _numVectors=1) const
Teuchos::RCP< MV > Z
The current preconditioned residual.
Teuchos::RCP< MV > AP
The matrix A applied to current decent direction vector.
CGPositiveDefiniteFailure is thrown when the the CG 'alpha = p^H*A*P' value is less than zero,...
int getBlockSize() const
Get the blocksize to be used by the iterative solver in solving this linear problem.
typename SCT::magnitudeType MagnitudeType
OperatorTraits< ScalarType, MV, OP > OPT
const LinearProblem< ScalarType, MV, OP > & getProblem() const
Get a constant reference to the linear problem.
void setState(Teuchos::RCP< CGIterationStateBase< ScalarType, MV > > state)
Teuchos::RCP< MV > getCurrentUpdate() const
Get the current update to the linear system.
Teuchos::ArrayView< MagnitudeType > getOffDiag()
Gets the off-diagonal for condition estimation (NOT_IMPLEMENTED).
void resetNumIters(int iter=0)
Reset the iteration count.
Teuchos::RCP< CGIterationStateBase< ScalarType, MV > > getState() const
Get the current state of the linear solver.
void initializeCG(Teuchos::RCP< CGIterationStateBase< ScalarType, MV > > newstate, Teuchos::RCP< MV > R_0)
Initialize the solver to an iterate, providing a complete state.
void setDoCondEst(bool)
Sets whether or not to store the diagonal for condition estimation.
virtual ~CGSingleRedIter()=default
Destructor.
CGSingleRedIter(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)
CGSingleRedIter constructor with linear problem, solver utilities, and parameter list of solver optio...
Teuchos::RCP< const MV > getNativeResiduals(std::vector< MagnitudeType > *norms) const
MultiVecTraits< ScalarType, MV > MVT
void setBlockSize(int blockSize)
Set the blocksize to be used by the iterative solver in solving this linear problem.
Teuchos::ArrayView< MagnitudeType > getDiag()
Gets the diagonal for condition estimation (NOT_IMPLEMENTED).
void iterate()
This method performs CG iterations until the status test indicates the need to stop or an error occur...
Teuchos::ScalarTraits< ScalarType > SCT
bool isInitialized()
States whether the solver has been initialized or not.
int getNumIters() const
Get the current iteration count.
void initialize()
Initialize the solver with the initial vectors from the linear problem or random data.
Structure to contain pointers to CGSingleRedIteration state variables.
virtual ~CGSingleRedIterationState()=default
CGSingleRedIterationState(Teuchos::RCP< const MV > tmp)
bool matches(Teuchos::RCP< const MV > tmp, int _numVectors=1) const
CGSingleRedIterationState()=default
void initialize(Teuchos::RCP< const MV > tmp, int _numVectors)
Traits class which defines basic operations on multivectors.
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 .
static Teuchos::RCP< MV > CloneCopy(const MV &mv)
Creates a new MV and copies contents of mv into the new vector (deep copy).
static ptrdiff_t GetGlobalLength(const MV &mv)
Return the number of rows in the given multivector mv.
static Teuchos::RCP< MV > Clone(const MV &mv, const int numvecs)
Creates a new empty MV containing numvecs columns.
static void MvAddMv(const ScalarType alpha, const MV &A, const ScalarType beta, const MV &B, MV &mv)
Replace mv with .
static int GetNumberVecs(const MV &mv)
Obtain the number of vectors in mv.
static void Assign(const MV &A, MV &mv)
mv := A
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...
An implementation of StatusTestResNorm using a family of residual norms.
A pure virtual class for defining the status tests for the Belos iterative solvers.