10#ifndef BELOS_PCPG_ITER_HPP
11#define BELOS_PCPG_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"
54 template <
class ScalarType,
class MV>
87 Teuchos::RCP<const Teuchos::SerialDenseMatrix<int,ScalarType> >
D;
91 R(Teuchos::null),
Z(Teuchos::null),
92 P(Teuchos::null),
AP(Teuchos::null),
93 U(Teuchos::null),
C(Teuchos::null),
100 template<
class ScalarType,
class MV,
class OP>
110 typedef Teuchos::ScalarTraits<ScalarType>
SCT;
127 Teuchos::ParameterList ¶ms );
231 if (!initialized_)
return 0;
237 if (!initialized_)
return 0;
260 TEUCHOS_TEST_FOR_EXCEPTION(blockSize!=1,std::invalid_argument,
261 "Belos::PCPGIter::setBlockSize(): Cannot use a block size that is not one.");
265 void setSize(
int savedBlocks );
286 const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > lp_;
287 const Teuchos::RCP<OutputManager<ScalarType> > om_;
288 const Teuchos::RCP<StatusTest<ScalarType,MV,OP> > stest_;
289 const Teuchos::RCP<OrthoManager<ScalarType,MV> > ortho_;
307 bool stateStorageInitialized_;
337 Teuchos::RCP<MV> AP_;
347 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > D_;
352 template<
class ScalarType,
class MV,
class OP>
357 Teuchos::ParameterList ¶ms ):
364 stateStorageInitialized_(false),
365 keepDiagonal_(false),
366 initDiagonal_(false),
373 TEUCHOS_TEST_FOR_EXCEPTION(!params.isParameter(
"Saved Blocks"), std::invalid_argument,
374 "Belos::PCPGIter::constructor: mandatory parameter \"Saved Blocks\" is not specified.");
375 int rb = Teuchos::getParameter<int>(params,
"Saved Blocks");
378 keepDiagonal_ = params.get(
"Keep Diagonal",
false);
381 initDiagonal_ = params.get(
"Initialize Diagonal",
false);
389 template <
class ScalarType,
class MV,
class OP>
395 TEUCHOS_TEST_FOR_EXCEPTION(savedBlocks <= 0, std::invalid_argument,
"Belos::PCPGIter::setSize() was passed a non-positive argument for \"Num Saved Blocks\".");
397 if ( savedBlocks_ != savedBlocks) {
398 stateStorageInitialized_ =
false;
399 savedBlocks_ = savedBlocks;
400 initialized_ =
false;
409 template <
class ScalarType,
class MV,
class OP>
412 stateStorageInitialized_ =
false;
413 initialized_ =
false;
421 template <
class ScalarType,
class MV,
class OP>
422 void PCPGIter<ScalarType,MV,OP>::setStateSize ()
424 if (!stateStorageInitialized_) {
427 Teuchos::RCP<const MV> lhsMV = lp_->getLHS();
428 Teuchos::RCP<const MV> rhsMV = lp_->getRHS();
429 if (lhsMV == Teuchos::null && rhsMV == Teuchos::null) {
436 int newsd = savedBlocks_ ;
441 if (Z_ == Teuchos::null) {
442 Teuchos::RCP<const MV> tmp = ( (rhsMV!=Teuchos::null)? rhsMV: lhsMV );
443 Z_ = MVT::Clone( *tmp, 1 );
445 if (P_ == Teuchos::null) {
446 Teuchos::RCP<const MV> tmp = ( (rhsMV!=Teuchos::null)? rhsMV: lhsMV );
447 P_ = MVT::Clone( *tmp, 1 );
449 if (AP_ == Teuchos::null) {
450 Teuchos::RCP<const MV> tmp = ( (rhsMV!=Teuchos::null)? rhsMV: lhsMV );
451 AP_ = MVT::Clone( *tmp, 1 );
454 if (C_ == Teuchos::null) {
457 Teuchos::RCP<const MV> tmp = ( (rhsMV!=Teuchos::null)? rhsMV: lhsMV );
458 TEUCHOS_TEST_FOR_EXCEPTION(tmp == Teuchos::null,std::invalid_argument,
459 "Belos::PCPGIter::setStateSize(): linear problem does not specify multivectors to clone from.");
460 TEUCHOS_TEST_FOR_EXCEPTION( 0 != prevUdim_,std::invalid_argument,
461 "Belos::PCPGIter::setStateSize(): prevUdim not zero and C is null.");
462 C_ = MVT::Clone( *tmp, savedBlocks_ );
466 if (MVT::GetNumberVecs(*C_) < savedBlocks_ ) {
467 Teuchos::RCP<const MV> tmp = C_;
468 C_ = MVT::Clone( *tmp, savedBlocks_ );
471 if (U_ == Teuchos::null) {
472 Teuchos::RCP<const MV> tmp = ( (rhsMV!=Teuchos::null)? rhsMV: lhsMV );
473 TEUCHOS_TEST_FOR_EXCEPTION( 0 != prevUdim_,std::invalid_argument,
474 "Belos::PCPGIter::setStateSize(): prevUdim not zero and U is null.");
475 U_ = MVT::Clone( *tmp, savedBlocks_ );
479 if (MVT::GetNumberVecs(*U_) < savedBlocks_ ) {
480 Teuchos::RCP<const MV> tmp = U_;
481 U_ = MVT::Clone( *tmp, savedBlocks_ );
485 if (D_ == Teuchos::null) {
486 D_ = Teuchos::rcp(
new Teuchos::SerialDenseMatrix<int,ScalarType>() );
489 D_->shape( newsd, newsd );
492 if (D_->numRows() < newsd || D_->numCols() < newsd) {
493 D_->shapeUninitialized( newsd, newsd );
498 stateStorageInitialized_ =
true;
505 template <
class ScalarType,
class MV,
class OP>
509 TEUCHOS_TEST_FOR_EXCEPTION(!stateStorageInitialized_,std::invalid_argument,
510 "Belos::PCPGIter::initialize(): Cannot initialize state storage!");
514 std::string errstr(
"Belos::PCPGIter::initialize(): Specified multivectors must have a consistent length and width.");
516 if (newstate.
R != Teuchos::null){
519 if (newstate.
U == Teuchos::null){
525 prevUdim_ = newstate.
curDim;
526 if (newstate.
C == Teuchos::null){
527 std::vector<int> index(prevUdim_);
528 for (
int i=0; i< prevUdim_; ++i)
533 lp_->apply( *Ukeff, *Ckeff );
535 curDim_ = prevUdim_ ;
539 if (!stateStorageInitialized_)
546 newstate.
curDim = curDim_;
550 std::vector<int> zero_index(1);
552 if ( lp_->getLeftPrec() != Teuchos::null ) {
553 lp_->applyLeftPrec( *R_, *Z_ );
560 std::vector<int> nextind(1);
561 nextind[0] = curDim_;
566 newstate.
curDim = curDim_;
569 std::invalid_argument, errstr );
570 if (newstate.
U != U_) {
575 std::invalid_argument, errstr );
576 if (newstate.
C != C_) {
582 TEUCHOS_TEST_FOR_EXCEPTION(newstate.
R == Teuchos::null,std::invalid_argument,
583 "Belos::PCPGIter::initialize(): PCPGStateIterState does not have initial kernel R_0.");
593 template <
class ScalarType,
class MV,
class OP>
599 if (initialized_ ==
false) {
602 const bool debug =
false;
605 Teuchos::SerialDenseMatrix<int,ScalarType> alpha( 1, 1 );
606 Teuchos::SerialDenseMatrix<int,ScalarType> pAp( 1, 1 );
607 Teuchos::SerialDenseMatrix<int,ScalarType> beta( 1, 1 );
608 Teuchos::SerialDenseMatrix<int,ScalarType> rHz( 1, 1 ), rHz_old( 1, 1 );
611 std::cout <<
" Iterate Warning: begin from nonzero iter_ ?" << std::endl;
614 std::vector<int> prevInd;
615 Teuchos::RCP<const MV> Uprev;
616 Teuchos::RCP<const MV> Cprev;
617 Teuchos::SerialDenseMatrix<int,ScalarType> CZ;
620 prevInd.resize( prevUdim_ );
621 for(
int i=0; i<prevUdim_ ; i++) prevInd[i] = i;
622 CZ.reshape( prevUdim_ , 1 );
628 Teuchos::RCP<MV> cur_soln_vec = lp_->getCurrLHSVec();
632 "Belos::PCPGIter::iterate(): current linear system has more than one std::vector!" );
636 "Belos::PCPGIter::iterate(): mistake in initialization !" );
639 const ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
640 const ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
643 std::vector<int> curind(1);
647 curind[0] = curDim_ - 1;
654 std::cout <<
" Input CZ post ortho " << std::endl;
655 CZ.print( std::cout );
657 if( curDim_ == savedBlocks_ ){
658 std::vector<int> zero_index(1);
671 while (stest_->checkStatus(
this) !=
Passed ) {
672 Teuchos::RCP<const MV> P;
676 curind[0] = curDim_ - 1;
679 std::cout << iter_ <<
" " << curDim_ <<
" " << rnorm[0] << std::endl;
681 if( prevUdim_ + iter_ < savedBlocks_ ){
684 lp_->applyOp( *P, *AP );
687 if( prevUdim_ + iter_ == savedBlocks_ ){
689 lp_->applyOp( *P_, *AP );
692 lp_->applyOp( *P_, *AP_ );
697 if( keepDiagonal_ && prevUdim_ + iter_ <= savedBlocks_ )
698 (*D_)(iter_ -1 ,iter_ -1 ) = pAp(0,0);
702 "Belos::PCPGIter::iterate(): non-positive value for p^H*A*p encountered!" );
705 alpha(0,0) = rHz(0,0) / pAp(0,0);
709 "Belos::PCPGIter::iterate(): non-positive value for alpha encountered!" );
712 if( curDim_ < savedBlocks_ ){
713 MVT::MvAddMv( one, *cur_soln_vec, alpha(0,0), *P, *cur_soln_vec );
715 MVT::MvAddMv( one, *cur_soln_vec, alpha(0,0), *P_, *cur_soln_vec );
721 rHz_old(0,0) = rHz(0,0);
725 if( prevUdim_ + iter_ <= savedBlocks_ ){
734 if ( lp_->getLeftPrec() != Teuchos::null ) {
735 lp_->applyLeftPrec( *R_, *Z_ );
742 beta(0,0) = rHz(0,0) / rHz_old(0,0);
744 if( curDim_ < savedBlocks_ ){
746 curind[0] = curDim_ - 1;
753 std::cout <<
" Check CZ " << std::endl;
755 CZ.print( std::cout );
759 if( curDim_ == savedBlocks_ ){
760 std::vector<int> zero_index(1);
764 Pnext = Teuchos::null;
772 std::cout <<
" Check CZ " << std::endl;
774 CZ.print( std::cout );
783 TEUCHOS_TEST_FOR_EXCEPTION( AP != Teuchos::null || P != Teuchos::null, std::logic_error,
"Loop recurrence violated. Please contact Belos team.");
785 if( prevUdim_ + iter_ < savedBlocks_ ) --curDim_;
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.
CGIterationInitFailure is thrown when the CGIteration object is unable to generate an initial iterate...
CGPositiveDefiniteFailure is thrown when the the CG 'alpha = p^H*A*P' value is less than zero,...
Iteration()
Default Constructor.
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 void MvNorm(const MV &mv, std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > &normvec, NormType type=TwoNorm)
Compute the norm of each individual vector of mv. Upon return, normvec[i] holds the value of ,...
static Teuchos::RCP< MV > Clone(const MV &mv, const int numvecs)
Creates a new empty MV containing numvecs columns.
static Teuchos::RCP< const MV > CloneView(const MV &mv, const std::vector< int > &index)
Creates a new const MV that shares the selected contents of mv (shallow copy).
static Teuchos::RCP< MV > CloneViewNonConst(MV &mv, const std::vector< int > &index)
Creates a new MV that shares the selected contents of mv (shallow copy).
static void MvAddMv(const ScalarType alpha, const MV &A, const ScalarType beta, const MV &B, MV &mv)
Replace mv with .
static void SetBlock(const MV &A, const std::vector< int > &index, MV &mv)
Copy the vectors in A to a set of vectors in mv indicated by the indices given in index.
static int GetNumberVecs(const MV &mv)
Obtain the number of vectors in mv.
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...
void initialize(PCPGIterState< ScalarType, MV > &newstate)
Initialize the solver to an iterate, providing a complete state.
PCPGIterState< ScalarType, MV > getState() const
Get the current state of the linear solver.
MultiVecTraits< ScalarType, MV > MVT
int getNumIters() const
Get the current iteration count.
Teuchos::RCP< MV > getCurrentUpdate() const
Get the current update to the linear system solution?.
void initialize()
Initialize the solver with the initial vectors from the linear problem. An exception is thrown if ini...
void setBlockSize(int blockSize)
Get the blocksize to be used by the iterative solver in solving this linear problem.
void resetNumIters(int iter=0)
Reset the iteration count.
virtual ~PCPGIter()
Destructor.
int getNumRecycledBlocks() const
Get the maximum number of recycled blocks used by the iterative solver in solving this linear problem...
void resetState()
tell the Iterator to "reset" itself; delete and rebuild the seed space.
const LinearProblem< ScalarType, MV, OP > & getProblem() const
Get a constant reference to the linear problem.
Teuchos::RCP< const MV > getNativeResiduals(std::vector< MagnitudeType > *) const
int getCurSubspaceDim() const
Get the current dimension of the whole seed subspace.
void setSize(int savedBlocks)
Set the maximum number of saved or recycled blocks used by the iterative solver.
SCT::magnitudeType MagnitudeType
OperatorTraits< ScalarType, MV, OP > OPT
int getBlockSize() const
Get the maximum number of blocks used by the iterative solver in solving this linear problem.
int getPrevSubspaceDim() const
Get the dimension of the search subspace used to solve the current solution to the linear problem.
Teuchos::ScalarTraits< ScalarType > SCT
PCPGIter(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)
PCPGIter constructor with linear problem, solver utilities, and parameter list of solver options.
void iterate()
PCPGIter iterates CG until the status test either requests a stop or detects an error....
bool isInitialized()
States whether the solver has been initialized or not.
A pure virtual class for defining the status tests for the Belos iterative solvers.
Structure to contain pointers to PCPGIter state variables.
Teuchos::RCP< MV > AP
The matrix A applied to current decent direction std::vector.
int prevUdim
Number of block columns in matrices C and U before current iteration.
Teuchos::RCP< MV > Z
The current preconditioned residual.
Teuchos::RCP< MV > U
The recycled subspace.
Teuchos::RCP< MV > C
C = AU, U spans recycled subspace.
int curDim
The current dimension of the reduction.
Teuchos::RCP< MV > P
The current decent direction std::vector.
Teuchos::RCP< MV > R
The current residual.
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > D
The current Hessenberg matrix.