10#ifndef BELOS_BLOCK_CG_ITER_HPP
11#define BELOS_BLOCK_CG_ITER_HPP
28#include "Teuchos_LAPACK.hpp"
29#include "Teuchos_SerialDenseMatrix.hpp"
30#include "Teuchos_SerialDenseVector.hpp"
31#include "Teuchos_SerialSymDenseMatrix.hpp"
32#include "Teuchos_SerialSpdDenseSolver.hpp"
33#include "Teuchos_ScalarTraits.hpp"
34#include "Teuchos_ParameterList.hpp"
35#include "Teuchos_TimeMonitor.hpp"
46 template <
class ScalarType,
class MV>
58 void initialize(Teuchos::RCP<const MV> tmp,
int _numVectors) {
60 this->
R = MVT::Clone( *tmp, _numVectors );
61 this->
Z = MVT::Clone( *tmp, _numVectors );
62 this->
P = MVT::Clone( *tmp, _numVectors );
63 this->
AP = MVT::Clone(*tmp, _numVectors );
68 bool matches(Teuchos::RCP<const MV> tmp,
int _numVectors=1)
const {
82template<
class ScalarType,
class MV,
class OP,
83 const bool lapackSupportsScalarType =
89 typedef Teuchos::ScalarTraits<ScalarType>
SCT;
96 Teuchos::ParameterList & )
98 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::logic_error,
"Stub");
104 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::logic_error,
"Stub");
108 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::logic_error,
"Stub");
112 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::logic_error,
"Stub");
115 Teuchos::RCP<CGIterationStateBase<ScalarType,MV> >
getState ()
const {
116 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::logic_error,
"Stub");
120 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::logic_error,
"Stub");
124 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::logic_error,
"Stub");
128 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::logic_error,
"Stub");
131 Teuchos::RCP<const MV>
133 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::logic_error,
"Stub");
137 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::logic_error,
"Stub");
141 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::logic_error,
"Stub");
145 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::logic_error,
"Stub");
149 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::logic_error,
"Stub");
153 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::logic_error,
"Stub");
157 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::logic_error,
"Stub");
166template<
class ScalarType,
class MV,
class OP>
176 using SCT = Teuchos::ScalarTraits<ScalarType>;
191 Teuchos::ParameterList ¶ms );
245 Teuchos::RCP<CGIterationStateBase<ScalarType,MV> >
getState()
const {
255 auto s = Teuchos::rcp_dynamic_cast<BlockCGIterationState<ScalarType,MV> >(state,
true);
305 Teuchos::ArrayView<MagnitudeType> temp;
311 Teuchos::ArrayView<MagnitudeType> temp;
324 const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > lp_;
325 const Teuchos::RCP<OutputManager<ScalarType> > om_;
326 const Teuchos::RCP<StatusTest<ScalarType,MV,OP> > stest_;
327 const Teuchos::RCP<OrthoManager<ScalarType,MV> > ortho_;
346 bool stateStorageInitialized_;
364 Teuchos::RCP<MV> AP_;
368 template<
class ScalarType,
class MV,
class OP>
374 Teuchos::ParameterList& params) :
381 stateStorageInitialized_(false),
385 int bs = params.get(
"Block Size", 1);
389 template <
class ScalarType,
class MV,
class OP>
394 TEUCHOS_TEST_FOR_EXCEPTION
395 (blockSize <= 0, std::invalid_argument,
"Belos::BlockGmresIter::"
396 "setBlockSize: blockSize = " << blockSize <<
" <= 0.");
397 if (blockSize == blockSize_) {
400 if (blockSize!=blockSize_) {
401 stateStorageInitialized_ =
false;
403 blockSize_ = blockSize;
404 initialized_ =
false;
407 template <
class ScalarType,
class MV,
class OP>
411 const char prefix[] =
"Belos::BlockCGIter::initialize: ";
414 Teuchos::RCP<const MV> lhsMV = lp_->getLHS();
415 Teuchos::RCP<const MV> rhsMV = lp_->getRHS();
416 Teuchos::RCP<const MV> tmp = ( (rhsMV!=Teuchos::null)? rhsMV: lhsMV );
417 TEUCHOS_ASSERT(!newstate.is_null());
419 newstate->initialize(tmp, blockSize_);
423 const char errstr[] =
"Specified multivectors must have a consistent "
428 TEUCHOS_TEST_FOR_EXCEPTION
430 std::invalid_argument, prefix << errstr );
431 TEUCHOS_TEST_FOR_EXCEPTION
433 std::invalid_argument, prefix << errstr );
443 if ( lp_->getLeftPrec() != Teuchos::null ) {
444 lp_->applyLeftPrec( *R_, *Z_ );
445 if ( lp_->getRightPrec() != Teuchos::null ) {
446 Teuchos::RCP<MV> tmp2 =
MVT::Clone( *Z_, blockSize_ );
447 lp_->applyRightPrec( *Z_, *tmp2 );
451 else if ( lp_->getRightPrec() != Teuchos::null ) {
452 lp_->applyRightPrec( *R_, *Z_ );
464 template <
class ScalarType,
class MV,
class OP>
467 const char prefix[] =
"Belos::BlockCGIter::iterate: ";
472 if (initialized_ ==
false) {
480 Teuchos::LAPACK<int,ScalarType> lapack;
483 Teuchos::SerialDenseMatrix<int,ScalarType> alpha( blockSize_, blockSize_ );
484 Teuchos::SerialDenseMatrix<int,ScalarType> beta( blockSize_, blockSize_ );
485 Teuchos::SerialDenseMatrix<int,ScalarType> rHz( blockSize_, blockSize_ ),
486 rHz_old( blockSize_, blockSize_ ), pAp( blockSize_, blockSize_ );
487 Teuchos::SerialSymDenseMatrix<int,ScalarType> pApHerm(Teuchos::View, uplo, pAp.values(), blockSize_, blockSize_);
490 Teuchos::SerialSpdDenseSolver<int,ScalarType> lltSolver;
493 const ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
496 Teuchos::RCP<MV> cur_soln_vec = lp_->getCurrLHSVec();
499 TEUCHOS_TEST_FOR_EXCEPTION
501 prefix <<
"Current linear system does not have the right number of vectors!" );
502 int rank = ortho_->normalize( *P_, Teuchos::null );
503 TEUCHOS_TEST_FOR_EXCEPTION
505 prefix <<
"Failed to compute initial block of orthonormal direction vectors.");
510 while (stest_->checkStatus(
this) !=
Passed) {
515 lp_->applyOp( *P_, *AP_ );
526 lltSolver.setMatrix( Teuchos::rcp(&pApHerm,
false) );
527 lltSolver.factorWithEquilibration(
true );
528 info = lltSolver.factor();
529 TEUCHOS_TEST_FOR_EXCEPTION
531 prefix <<
"Failed to compute Cholesky factorization using LAPACK routine POTRF.");
535 lltSolver.setVectors (Teuchos::rcpFromRef (alpha), Teuchos::rcpFromRef (alpha));
536 info = lltSolver.solve();
537 TEUCHOS_TEST_FOR_EXCEPTION
539 prefix <<
"Failed to compute alpha using Cholesky factorization (POTRS).");
543 lp_->updateSolution();
549 if ( lp_->getLeftPrec() != Teuchos::null ) {
550 lp_->applyLeftPrec( *R_, *Z_ );
551 if ( lp_->getRightPrec() != Teuchos::null ) {
552 Teuchos::RCP<MV> tmp =
MVT::Clone( *Z_, blockSize_ );
553 lp_->applyRightPrec( *Z_, *tmp );
557 else if ( lp_->getRightPrec() != Teuchos::null ) {
558 lp_->applyRightPrec( *R_, *Z_ );
572 lltSolver.setVectors( Teuchos::rcp( &beta,
false ), Teuchos::rcp( &beta,
false ) );
573 info = lltSolver.solve();
574 TEUCHOS_TEST_FOR_EXCEPTION
576 prefix <<
"Failed to compute beta using Cholesky factorization (POTRS).");
584 rank = ortho_->normalize( *P_, Teuchos::null );
585 TEUCHOS_TEST_FOR_EXCEPTION
587 prefix <<
"Failed to compute block of orthonormal direction vectors.");
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.
Templated virtual class for providing orthogonalization/orthonormalization methods with matrix-based ...
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.
Pure virtual base class for defining the status testing capabilities of Belos.
Collection of types and exceptions used within the Belos solvers.
void initialize()
Initialize the solver with the initial vectors from the linear problem or random data.
OperatorTraits< ScalarType, MV, OP > OPT
Teuchos::RCP< CGIterationStateBase< ScalarType, MV > > getState() const
Get the current state of the linear solver.
void setDoCondEst(bool)
Sets whether or not to store the diagonal for condition estimation.
const LinearProblem< ScalarType, MV, OP > & getProblem() const
Get a constant reference to the linear problem.
Teuchos::ArrayView< MagnitudeType > getDiag()
Gets the diagonal for condition estimation (NOT_IMPLEMENTED).
virtual ~BlockCGIter()=default
Destructor.
Teuchos::ScalarTraits< ScalarType > SCT
void resetNumIters(int iter=0)
Reset the iteration count.
typename SCT::magnitudeType MagnitudeType
void setState(Teuchos::RCP< CGIterationStateBase< ScalarType, MV > > state)
MultiVecTraits< ScalarType, MV > MVT
Teuchos::ArrayView< MagnitudeType > getOffDiag()
Gets the off-diagonal for condition estimation (NOT_IMPLEMENTED).
int getNumIters() const
Get the current iteration count.
bool isInitialized()
States whether the solver has been initialized or not.
int getBlockSize() const
Get the block size to be used by the iterative solver in solving this linear problem.
Teuchos::RCP< const MV > getNativeResiduals(std::vector< MagnitudeType > *) const
BlockCGIter(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)
BlockCGIter constructor with linear problem, solver utilities, and parameter list of solver options.
Teuchos::RCP< MV > getCurrentUpdate() const
Get the current update to the linear system.
OperatorTraits< ScalarType, MV, OP > OPT
Teuchos::RCP< const MV > getNativeResiduals(std::vector< MagnitudeType > *) const
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 &)
SCT::magnitudeType MagnitudeType
void resetNumIters(int iter=0)
Reset the iteration count to iter.
void setBlockSize(int blockSize)
Set the blocksize to be used by the iterative solver in solving this linear problem.
int getNumIters() const
Get the current iteration count.
void iterate()
This method performs linear solver iterations until the status test indicates the need to stop or an ...
void setState(Teuchos::RCP< CGIterationStateBase< ScalarType, MV > > state)
void initialize()
Initialize the solver with the initial vectors from the linear problem or random data.
void initializeCG(Teuchos::RCP< BlockCGIterationState< ScalarType, MV > >, Teuchos::RCP< MV >)
const LinearProblem< ScalarType, MV, OP > & getProblem() const
Get a constant reference to the linear problem.
int getBlockSize() const
Get the blocksize to be used by the iterative solver in solving this linear problem.
Teuchos::RCP< MV > getCurrentUpdate() const
Get the current update to the linear system.
Teuchos::RCP< CGIterationStateBase< ScalarType, MV > > getState() const
Get the current state of the linear solver.
Teuchos::ScalarTraits< ScalarType > SCT
bool isInitialized()
States whether the solver has been initialized or not.
void setDoCondEst(bool val)
Sets whether or not to store the diagonal for condition estimation.
MultiVecTraits< ScalarType, MV > MVT
Structure to contain pointers to BlockCGIteration state variables.
virtual ~BlockCGIterationState()=default
void initialize(Teuchos::RCP< const MV > tmp, int _numVectors)
bool matches(Teuchos::RCP< const MV > tmp, int _numVectors=1) const
BlockCGIterationState()=default
BlockCGIterationState(Teuchos::RCP< const MV > tmp)
CGIterateFailure is thrown when the CGIteration object is unable to compute the next iterate in the C...
CGIterationLAPACKFailure is thrown when a nonzero return value is passed back from an LAPACK routine.
CGIterationOrthoFailure is thrown when the CGIteration object is unable to compute independent direct...
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.
Belos's templated virtual class for providing routines for orthogonalization and orthonormzalition of...
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 void MvTimesMatAddMv(const ScalarType alpha, const MV &A, const Teuchos::SerialDenseMatrix< int, ScalarType > &B, const ScalarType beta, MV &mv)
Update mv with .
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 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...
A pure virtual class for defining the status tests for the Belos iterative solvers.