42 #ifndef BELOS_BLOCK_GCRODR_ITER_HPP
43 #define BELOS_BLOCK_GCRODR_ITER_HPP
94 template <
class ScalarType,
class MV>
121 U(Teuchos::null),
C(Teuchos::null),
122 H(Teuchos::null),
B(Teuchos::null)
166 template<
class ScalarType,
class MV,
class OP>
321 if (!initialized_)
return 0;
349 void setSize(
int recycledBlocks,
int numBlocks ) {
351 if ( (recycledBlocks_ != recycledBlocks) || (numBlocks_ != numBlocks) ) {
352 recycledBlocks_ = recycledBlocks;
353 numBlocks_ = numBlocks;
383 int numBlocks_, blockSize_;
388 std::vector<bool> trueRHSIndices_;
402 std::vector< SDM >House_;
414 int curDim_, iter_, lclIter_;
462 template<
class ScalarType,
class MV,
class OP>
468 om_(printer), stest_(tester), ortho_(ortho) {
472 initialized_ =
false;
484 int nb = Teuchos::getParameter<int>(params,
"Num Blocks");
486 TEUCHOS_TEST_FOR_EXCEPTION(!params.
isParameter(
"Recycled Blocks"), std::invalid_argument,
"Belos::BlockGCRODRIter::constructor: mandatory parameter \"Recycled Blocks\" is not specified.");
487 int rb = Teuchos::getParameter<int>(params,
"Recycled Blocks");
489 TEUCHOS_TEST_FOR_EXCEPTION(nb <= 0, std::invalid_argument,
"Belos::BlockGCRODRIter() was passed a non-positive argument for \"Num Blocks\".");
490 TEUCHOS_TEST_FOR_EXCEPTION(rb >= nb, std::invalid_argument,
"Belos::BlockGCRODRIter() the number of recycled blocks is larger than the allowable subspace.");
493 int bs = Teuchos::getParameter<int>(params,
"Block Size");
495 TEUCHOS_TEST_FOR_EXCEPTION(bs <= 0, std::invalid_argument,
"Belos::BlockGCRODRIter() the block size was passed a non-postitive argument.");
499 recycledBlocks_ = rb;
504 trueRHSIndices_.resize(blockSize_);
506 for(i=0; i<blockSize_; i++){
507 trueRHSIndices_[i] =
true;
516 House_.resize(numBlocks_);
518 for(i=0; i<numBlocks_;i++){
519 House_[i].shapeUninitialized(2*blockSize_, 2*blockSize_);
525 template <
class ScalarType,
class MV,
class OP>
533 setSize( recycledBlocks_, numBlocks_ );
537 std::vector<int> curind(blockSize_);
543 for(
int i = 0; i<blockSize_; i++){curind[i] = i;};
548 Vnext = MVT::CloneViewNonConst(*V_,curind);
553 int rank = ortho_->normalize(*Vnext,Z0);
561 Z_block->assign(*Z0);
563 std::vector<int> prevind(blockSize_*(numBlocks_ + 1));
570 while( (stest_->checkStatus(
this) !=
Passed) && (curDim_+blockSize_-1) < (numBlocks_*blockSize_)) {
576 int HFirstCol = curDim_-blockSize_;
577 int HLastCol = HFirstCol + blockSize_-1 ;
578 int HLastOrthRow = HLastCol;
579 int HFirstNormRow = HLastOrthRow + 1;
583 for(
int i = 0; i< blockSize_; i++){
584 curind[i] = curDim_ + i;
586 Vnext = MVT::CloneViewNonConst(*V_,curind);
591 for(
int i = 0; i< blockSize_; i++){
592 curind[blockSize_ - 1 - i] = curDim_ - i - 1;
594 Vprev = MVT::CloneView(*V_,curind);
596 lp_->apply(*Vprev,*Vnext);
597 Vprev = Teuchos::null;
611 ortho_->project( *Vnext, AsubB, C );
615 prevind.resize(curDim_);
616 for (
int i=0; i<curDim_; i++) { prevind[i] = i; }
617 Vprev = MVT::CloneView(*V_,prevind);
627 int rank = ortho_->projectAndNormalize(*Vnext,AsubH,subR,AVprev);
630 SDM subR2(
Teuchos::View,*R_,(lclIter_+1)*blockSize_,blockSize_,0,HFirstCol);
631 SDM subH2(
Teuchos::View,*H_,(lclIter_+1)*blockSize_,blockSize_,0,HFirstCol);
640 curDim_ = curDim_ + blockSize_;
650 template <
class ScalarType,
class MV,
class OP>
652 if (newstate.
V != Teuchos::null && newstate.
H != Teuchos::null) {
653 curDim_ = newstate.
curDim;
664 SDM Identity(2*blockSize_, 2*blockSize_);
665 for(
int i=0;i<2*blockSize_; i++){
668 for(
int i=0; i<numBlocks_;i++){
669 House_[i].
assign(Identity);
673 TEUCHOS_TEST_FOR_EXCEPTION(newstate.
V == Teuchos::null,std::invalid_argument,
"Belos::GCRODRIter::initialize(): BlockGCRODRIterState does not have V initialized.");
674 TEUCHOS_TEST_FOR_EXCEPTION(newstate.
H == Teuchos::null,std::invalid_argument,
"Belos::GCRODRIter::initialize(): BlockGCRODRIterState does not have H initialized.");
686 template <
class ScalarType,
class MV,
class OP>
694 if (static_cast<int> (norms->size()) < blockSize_) {
695 norms->resize( blockSize_ );
698 for (
int j=0; j<blockSize_; j++) {
699 if(trueRHSIndices_[j]){
700 (*norms)[j] = blas.
NRM2( blockSize_, &Z_(curDim_-blockSize_+j, j), 1);
706 return Teuchos::null;
709 return Teuchos::null;
715 template <
class ScalarType,
class MV,
class OP>
723 if(curDim_<=blockSize_) {
724 return currentUpdate;
730 currentUpdate = MVT::Clone( *V_, blockSize_ );
754 Rtmp->values(), Rtmp->stride(), Y.
values(), Y.
stride() );
761 std::vector<int> index(curDim_-blockSize_);
762 for (
int i=0; i<curDim_-blockSize_; i++ ) index[i] = i;
764 MVT::MvTimesMatAddMv( one, *Vjp1, Y, zero, *currentUpdate );
771 if (U_ != Teuchos::null) {
772 SDM z(recycledBlocks_,blockSize_);
777 MVT::MvTimesMatAddMv( -one, *U_, z, one, *currentUpdate );
783 return currentUpdate;
786 template<
class ScalarType,
class MV,
class OP>
798 int curDim = curDim_;
799 if ( (dim >= curDim_) && (dim < getMaxSubspaceDim()) ){
812 for (i=0; i<curDim-1; i++) {
816 blas.
ROT( 1, &(*R_)(i,curDim-1), 1, &(*R_)(i+1, curDim-1), 1, &cs_[i], &sn_[i] );
821 blas.
ROTG( &(*R_)(curDim-1,curDim-1), &(*R_)(curDim,curDim-1), &cs_[curDim-1], &sn_[curDim-1] );
822 (*R_)(curDim,curDim-1) = zero;
826 blas.
ROT( 1, &Z_(curDim-1,0), 1, &Z_(curDim,0), 1, &cs_[curDim-1], &sn_[curDim-1] );
844 int R_colStart = curDim_-blockSize_;
850 for(i=0; i<lclIter_-1; i++){
851 int R_rowStart = i*blockSize_;
855 blas.
GEMM(
Teuchos::NO_TRANS,
Teuchos::NO_TRANS, 2*blockSize_,blockSize_,2*blockSize_,one,House_[i].values(),House_[i].stride(), RblockCopy->values(),RblockCopy -> stride(), zero, RblockView->values(),RblockView -> stride());
862 Rblock =
rcp(
new SDM (
Teuchos::View, *R_, 2*blockSize_,blockSize_, curDim_-blockSize_, curDim_-blockSize_));
865 for(i=0; i<blockSize_; i++){
869 int curcol = (lclIter_ - 1)*blockSize_ + i;
875 ScalarType nvs = blas.
NRM2(blockSize_+1,&((*R_)[curcol][curcol]),1);
876 ScalarType alpha = -signDiag*nvs;
895 (*v_refl)[0] -= alpha;
896 nvs = blas.
NRM2(blockSize_+1,v_refl -> values(),1);
912 if(i < blockSize_ - 1){
916 blas.
GEMV(
Teuchos::TRANS, workmatrix->numRows(), workmatrix->numCols(), one, workmatrix->values(), workmatrix->stride(), v_refl->values(), 1, zero, workvec->values(), 1);
917 blas.
GER(workmatrix->numRows(),workmatrix->numCols(), -2.0*one/nvs, v_refl->values(),1,workvec->values(),1,workmatrix->values(),workmatrix->stride());
926 blas.
GEMV(
Teuchos::TRANS,workmatrix->numRows(),workmatrix->numCols(),one,workmatrix->values(),workmatrix->stride(), v_refl->values(), 1,zero,workvec->values(),1);
927 blas.
GER(workmatrix->numRows(),workmatrix->numCols(), -2.0*one/nvs, v_refl -> values(),1,workvec->values(),1,workmatrix->values(),(*workmatrix).stride());
934 blas.
GEMV(
Teuchos::TRANS, workmatrix->numRows(), workmatrix->numCols(), one, workmatrix-> values(), workmatrix->stride(), v_refl -> values(), 1, zero, workvec->values(), 1);
935 blas.
GER((*workmatrix).numRows(),(*workmatrix).numCols(), -2.0*one/nvs,v_refl -> values(), 1,&((*workvec)[0]),1,(*workmatrix)[0],(*workmatrix).stride());
940 (*R_)[curcol][curcol] = alpha;
941 for(
int ii=1; ii<= blockSize_; ii++){
942 (*R_)[curcol][curcol+ii] = 0;