10 #ifndef BELOS_PCPG_ITER_HPP
11 #define BELOS_PCPG_ITER_HPP
54 template <
class ScalarType,
class MV>
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>
231 if (!initialized_)
return 0;
237 if (!initialized_)
return 0;
261 "Belos::PCPGIter::setBlockSize(): Cannot use a block size that is not one.");
265 void setSize(
int savedBlocks );
307 bool stateStorageInitialized_;
352 template<
class ScalarType,
class MV,
class OP>
364 stateStorageInitialized_(false),
365 keepDiagonal_(false),
366 initDiagonal_(false),
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>
424 if (!stateStorageInitialized_) {
429 if (lhsMV == Teuchos::null && rhsMV == Teuchos::null) {
436 int newsd = savedBlocks_ ;
441 if (Z_ == Teuchos::null) {
443 Z_ = MVT::Clone( *tmp, 1 );
445 if (P_ == Teuchos::null) {
447 P_ = MVT::Clone( *tmp, 1 );
449 if (AP_ == Teuchos::null) {
451 AP_ = MVT::Clone( *tmp, 1 );
454 if (C_ == Teuchos::null) {
459 "Belos::PCPGIter::setStateSize(): linear problem does not specify multivectors to clone from.");
461 "Belos::PCPGIter::setStateSize(): prevUdim not zero and C is null.");
462 C_ = MVT::Clone( *tmp, savedBlocks_ );
466 if (MVT::GetNumberVecs(*C_) < savedBlocks_ ) {
468 C_ = MVT::Clone( *tmp, savedBlocks_ );
471 if (U_ == Teuchos::null) {
474 "Belos::PCPGIter::setStateSize(): prevUdim not zero and U is null.");
475 U_ = MVT::Clone( *tmp, savedBlocks_ );
479 if (MVT::GetNumberVecs(*U_) < savedBlocks_ ) {
481 U_ = MVT::Clone( *tmp, savedBlocks_ );
485 if (D_ == Teuchos::null) {
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>
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)
531 newstate.
C = MVT::Clone( *newstate.
U, prevUdim_ );
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_ );
554 MVT::SetBlock( *Z_, zero_index , *P_ );
557 MVT::SetBlock( *R_, zero_index, *P_ );
560 std::vector<int> nextind(1);
561 nextind[0] = curDim_;
563 MVT::SetBlock( *P_, nextind, *newstate.
U );
566 newstate.
curDim = curDim_;
569 std::invalid_argument, errstr );
570 if (newstate.
U != U_) {
575 std::invalid_argument, errstr );
576 if (newstate.
C != C_) {
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;
611 std::cout <<
" Iterate Warning: begin from nonzero iter_ ?" << std::endl;
614 std::vector<int> prevInd;
620 prevInd.resize( prevUdim_ );
621 for(
int i=0; i<prevUdim_ ; i++) prevInd[i] = i;
623 Uprev = MVT::CloneView(*U_, prevInd);
624 Cprev = MVT::CloneView(*C_, prevInd);
632 "Belos::PCPGIter::iterate(): current linear system has more than one std::vector!" );
636 "Belos::PCPGIter::iterate(): mistake in initialization !" );
643 std::vector<int> curind(1);
644 std::vector<ScalarType> rnorm(MVT::GetNumberVecs(*cur_soln_vec));
647 curind[0] = curDim_ - 1;
648 P = MVT::CloneViewNonConst(*U_,curind);
649 MVT::MvTransMv( one, *Cprev, *P, CZ );
650 MVT::MvTimesMatAddMv( -one, *Uprev, CZ, one, *P );
653 MVT::MvTransMv( one, *Cprev, *P, CZ );
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);
660 MVT::SetBlock( *P, zero_index, *P_ );
666 MVT::MvTransMv( one, *R_, *Z_, rHz );
671 while (stest_->checkStatus(
this) !=
Passed ) {
676 curind[0] = curDim_ - 1;
678 MVT::MvNorm(*R_, rnorm);
679 std::cout << iter_ <<
" " << curDim_ <<
" " << rnorm[0] << std::endl;
681 if( prevUdim_ + iter_ < savedBlocks_ ){
682 P = MVT::CloneView(*U_,curind);
683 AP = MVT::CloneViewNonConst(*C_,curind);
684 lp_->applyOp( *P, *AP );
685 MVT::MvTransMv( one, *P, *AP, pAp );
687 if( prevUdim_ + iter_ == savedBlocks_ ){
688 AP = MVT::CloneViewNonConst(*C_,curind);
689 lp_->applyOp( *P_, *AP );
690 MVT::MvTransMv( one, *P_, *AP, pAp );
692 lp_->applyOp( *P_, *AP_ );
693 MVT::MvTransMv( one, *P_, *AP_, pAp );
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_ ){
726 MVT::MvAddMv( one, *R_, -alpha(0,0), *AP, *R_ );
729 MVT::MvAddMv( one, *R_, -alpha(0,0), *AP_, *R_ );
734 if ( lp_->getLeftPrec() != Teuchos::null ) {
735 lp_->applyLeftPrec( *R_, *Z_ );
740 MVT::MvTransMv( one, *R_, *Z_, rHz );
742 beta(0,0) = rHz(0,0) / rHz_old(0,0);
744 if( curDim_ < savedBlocks_ ){
746 curind[0] = curDim_ - 1;
748 MVT::MvAddMv( one, *Z_, beta(0,0), *P, *Pnext );
750 MVT::MvTransMv( one, *Cprev, *Z_, CZ );
751 MVT::MvTimesMatAddMv( -one, *Uprev, CZ, one, *Pnext );
753 std::cout <<
" Check CZ " << std::endl;
754 MVT::MvTransMv( one, *Cprev, *Pnext, CZ );
755 CZ.
print( std::cout );
759 if( curDim_ == savedBlocks_ ){
760 std::vector<int> zero_index(1);
762 MVT::SetBlock( *Pnext, zero_index, *P_ );
764 Pnext = Teuchos::null;
766 MVT::MvAddMv( one, *Z_, beta(0,0), *P_, *P_ );
768 MVT::MvTransMv( one, *Cprev, *Z_, CZ );
769 MVT::MvTimesMatAddMv( -one, *Uprev, CZ, one, *P_ );
772 std::cout <<
" Check CZ " << std::endl;
773 MVT::MvTransMv( one, *Cprev, *P_, CZ );
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_;
Collection of types and exceptions used within the Belos solvers.
int prevUdim
Number of block columns in matrices C and U before current iteration.
Belos's basic output manager for sending information of select verbosity levels to the appropriate ou...
bool isInitialized()
States whether the solver has been initialized or not.
Class which manages the output and verbosity of the Belos solvers.
Teuchos::RCP< MV > R
The current residual.
Pure virtual base class which augments the basic interface for a conjugate gradient linear solver ite...
PCPGIterState< ScalarType, MV > getState() const
Get the current state of the linear solver.
Pure virtual base class for defining the status testing capabilities of Belos.
Teuchos::RCP< MV > C
C = AU, U spans recycled subspace.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
int getCurSubspaceDim() const
Get the current dimension of the whole seed subspace.
Declaration of basic traits for the multivector type.
virtual std::ostream & print(std::ostream &os) const
void setSize(int savedBlocks)
Set the maximum number of saved or recycled blocks used by the iterative solver.
int curDim
The current dimension of the reduction.
void resetNumIters(int iter=0)
Reset the iteration count.
Teuchos::RCP< MV > P
The current decent direction std::vector.
Teuchos::RCP< MV > AP
The matrix A applied to current decent direction std::vector.
Structure to contain pointers to PCPGIter state variables.
OperatorTraits< ScalarType, MV, OP > OPT
int getPrevSubspaceDim() const
Get the dimension of the search subspace used to solve the current solution to the linear problem...
A pure virtual class for defining the status tests for the Belos iterative solvers.
int getBlockSize() const
Get the maximum number of blocks used by the iterative solver in solving this linear problem...
Class which defines basic traits for the operator type.
Traits class which defines basic operations on multivectors.
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > D
The current Hessenberg matrix.
CGPositiveDefiniteFailure is thrown when the the CG 'alpha = p^H*A*P' value is less than zero...
bool isParameter(const std::string &name) const
virtual ~PCPGIter()
Destructor.
void setBlockSize(int blockSize)
Get the blocksize to be used by the iterative solver in solving this linear problem.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
int getNumIters() const
Get the current iteration count.
void initialize()
Initialize the solver with the initial vectors from the linear problem. An exception is thrown if ini...
A linear system to solve, and its associated information.
Class which describes the linear problem to be solved by the iterative solver.
int getNumRecycledBlocks() const
Get the maximum number of recycled blocks used by the iterative solver in solving this linear problem...
const LinearProblem< ScalarType, MV, OP > & getProblem() const
Get a constant reference to the linear problem.
SCT::magnitudeType MagnitudeType
Teuchos::RCP< MV > U
The recycled subspace.
void resetState()
tell the Iterator to "reset" itself; delete and rebuild the seed space.
int reshape(OrdinalType numRows, OrdinalType numCols)
void iterate()
PCPGIter iterates CG until the status test either requests a stop or detects an error. In the latter case, std::exception is thrown.
Teuchos::RCP< MV > Z
The current preconditioned residual.
Teuchos::RCP< const MV > getNativeResiduals(std::vector< MagnitudeType > *) const
Get the norms of the residuals native to the solver.
Class which defines basic traits for the operator type.
CGIterationInitFailure is thrown when the CGIteration object is unable to generate an initial iterate...
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...
Belos's templated virtual class for providing routines for orthogonalization and orthonormzalition of...
This class implements the PCPG iteration, where a single-std::vector Krylov subspace is constructed...
Belos header file which uses auto-configuration information to include necessary C++ headers...
MultiVecTraits< ScalarType, MV > MVT
Teuchos::RCP< MV > getCurrentUpdate() const
Get the current update to the linear system solution?.
Templated virtual class for providing orthogonalization/orthonormalization methods with matrix-based ...
Teuchos::ScalarTraits< ScalarType > SCT