576 TEUCHOS_TEST_FOR_EXCEPTION(!stateStorageInitialized_,std::invalid_argument,
577 "Belos::PCPGIter::initialize(): Cannot initialize state storage!");
581 std::string errstr(
"Belos::PCPGIter::initialize(): Specified multivectors must have a consistent length and width.");
583 if (newstate.
R != Teuchos::null){
586 if (newstate.
U == Teuchos::null){
592 prevUdim_ = newstate.
curDim;
593 if (newstate.
C == Teuchos::null){
594 std::vector<int> index(prevUdim_);
595 for (
int i=0; i< prevUdim_; ++i)
597 Teuchos::RCP<const MV> Ukeff = MVT::CloneView( *newstate.
U, index );
598 newstate.
C = MVT::Clone( *newstate.
U, prevUdim_ );
599 Teuchos::RCP<MV> Ckeff = MVT::CloneViewNonConst( *newstate.
C, index );
600 lp_->apply( *Ukeff, *Ckeff );
602 curDim_ = prevUdim_ ;
606 if (!stateStorageInitialized_)
613 newstate.
curDim = curDim_;
617 std::vector<int> zero_index(1);
619 if ( lp_->getLeftPrec() != Teuchos::null ) {
620 lp_->applyLeftPrec( *R_, *Z_ );
621 MVT::SetBlock( *Z_, zero_index , *P_ );
624 MVT::SetBlock( *R_, zero_index, *P_ );
627 std::vector<int> nextind(1);
628 nextind[0] = curDim_;
630 MVT::SetBlock( *P_, nextind, *newstate.
U );
633 newstate.
curDim = curDim_;
635 TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetNumberVecs(*newstate.
U) != savedBlocks_ ,
636 std::invalid_argument, errstr );
637 if (newstate.
U != U_) {
641 TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetNumberVecs(*newstate.
C) != savedBlocks_ ,
642 std::invalid_argument, errstr );
643 if (newstate.
C != C_) {
649 TEUCHOS_TEST_FOR_EXCEPTION(newstate.
R == Teuchos::null,std::invalid_argument,
650 "Belos::PCPGIter::initialize(): PCPGStateIterState does not have initial kernel R_0.");
666 if (initialized_ ==
false) {
669 const bool debug =
false;
672 Teuchos::SerialDenseMatrix<int,ScalarType> alpha( 1, 1 );
673 Teuchos::SerialDenseMatrix<int,ScalarType> pAp( 1, 1 );
674 Teuchos::SerialDenseMatrix<int,ScalarType> beta( 1, 1 );
675 Teuchos::SerialDenseMatrix<int,ScalarType> rHz( 1, 1 ), rHz_old( 1, 1 );
678 std::cout <<
" Iterate Warning: begin from nonzero iter_ ?" << std::endl;
681 std::vector<int> prevInd;
682 Teuchos::RCP<const MV> Uprev;
683 Teuchos::RCP<const MV> Cprev;
684 Teuchos::SerialDenseMatrix<int,ScalarType> CZ;
687 prevInd.resize( prevUdim_ );
688 for(
int i=0; i<prevUdim_ ; i++) prevInd[i] = i;
689 CZ.reshape( prevUdim_ , 1 );
690 Uprev = MVT::CloneView(*U_, prevInd);
691 Cprev = MVT::CloneView(*C_, prevInd);
695 Teuchos::RCP<MV> cur_soln_vec = lp_->getCurrLHSVec();
699 "Belos::CGIter::iterate(): current linear system has more than one std::vector!" );
703 "Belos::CGIter::iterate(): mistake in initialization !" );
706 const ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
707 const ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
710 std::vector<int> curind(1);
711 std::vector<ScalarType> rnorm(MVT::GetNumberVecs(*cur_soln_vec));
714 curind[0] = curDim_ - 1;
715 P = MVT::CloneViewNonConst(*U_,curind);
716 MVT::MvTransMv( one, *Cprev, *P, CZ );
717 MVT::MvTimesMatAddMv( -one, *Uprev, CZ, one, *P );
720 MVT::MvTransMv( one, *Cprev, *P, CZ );
721 std::cout <<
" Input CZ post ortho " << std::endl;
722 CZ.print( std::cout );
724 if( curDim_ == savedBlocks_ ){
725 std::vector<int> zero_index(1);
727 MVT::SetBlock( *P, zero_index, *P_ );
733 MVT::MvTransMv( one, *R_, *Z_, rHz );
738 while (stest_->checkStatus(
this) !=
Passed ) {
739 Teuchos::RCP<const MV> P;
743 curind[0] = curDim_ - 1;
745 MVT::MvNorm(*R_, rnorm);
746 std::cout << iter_ <<
" " << curDim_ <<
" " << rnorm[0] << std::endl;
748 if( prevUdim_ + iter_ < savedBlocks_ ){
749 P = MVT::CloneView(*U_,curind);
750 AP = MVT::CloneViewNonConst(*C_,curind);
751 lp_->applyOp( *P, *AP );
752 MVT::MvTransMv( one, *P, *AP, pAp );
754 if( prevUdim_ + iter_ == savedBlocks_ ){
755 AP = MVT::CloneViewNonConst(*C_,curind);
756 lp_->applyOp( *P_, *AP );
757 MVT::MvTransMv( one, *P_, *AP, pAp );
759 lp_->applyOp( *P_, *AP_ );
760 MVT::MvTransMv( one, *P_, *AP_, pAp );
764 if( keepDiagonal_ && prevUdim_ + iter_ <= savedBlocks_ )
765 (*D_)(iter_ -1 ,iter_ -1 ) = pAp(0,0);
769 "Belos::CGIter::iterate(): non-positive value for p^H*A*p encountered!" );
772 alpha(0,0) = rHz(0,0) / pAp(0,0);
776 "Belos::CGIter::iterate(): non-positive value for alpha encountered!" );
779 if( curDim_ < savedBlocks_ ){
780 MVT::MvAddMv( one, *cur_soln_vec, alpha(0,0), *P, *cur_soln_vec );
782 MVT::MvAddMv( one, *cur_soln_vec, alpha(0,0), *P_, *cur_soln_vec );
788 rHz_old(0,0) = rHz(0,0);
792 if( prevUdim_ + iter_ <= savedBlocks_ ){
793 MVT::MvAddMv( one, *R_, -alpha(0,0), *AP, *R_ );
796 MVT::MvAddMv( one, *R_, -alpha(0,0), *AP_, *R_ );
801 if ( lp_->getLeftPrec() != Teuchos::null ) {
802 lp_->applyLeftPrec( *R_, *Z_ );
807 MVT::MvTransMv( one, *R_, *Z_, rHz );
809 beta(0,0) = rHz(0,0) / rHz_old(0,0);
811 if( curDim_ < savedBlocks_ ){
813 curind[0] = curDim_ - 1;
814 Teuchos::RCP<MV> Pnext = MVT::CloneViewNonConst(*U_,curind);
815 MVT::MvAddMv( one, *Z_, beta(0,0), *P, *Pnext );
817 MVT::MvTransMv( one, *Cprev, *Z_, CZ );
818 MVT::MvTimesMatAddMv( -one, *Uprev, CZ, one, *Pnext );
820 std::cout <<
" Check CZ " << std::endl;
821 MVT::MvTransMv( one, *Cprev, *Pnext, CZ );
822 CZ.print( std::cout );
826 if( curDim_ == savedBlocks_ ){
827 std::vector<int> zero_index(1);
829 MVT::SetBlock( *Pnext, zero_index, *P_ );
831 Pnext = Teuchos::null;
833 MVT::MvAddMv( one, *Z_, beta(0,0), *P_, *P_ );
835 MVT::MvTransMv( one, *Cprev, *Z_, CZ );
836 MVT::MvTimesMatAddMv( -one, *Uprev, CZ, one, *P_ );
839 std::cout <<
" Check CZ " << std::endl;
840 MVT::MvTransMv( one, *Cprev, *P_, CZ );
841 CZ.print( std::cout );
850 TEUCHOS_TEST_FOR_EXCEPTION( AP != Teuchos::null || P != Teuchos::null, std::logic_error,
"Loop recurrence violated. Please contact Belos team.");
852 if( prevUdim_ + iter_ < savedBlocks_ ) --curDim_;