143 typedef Teuchos::ScalarTraits<ScalarType>
SCT;
145 typedef Teuchos::SerialDenseMatrix<int,ScalarType>
SDM;
146 typedef Teuchos::SerialDenseVector<int,ScalarType>
SDV;
163 Teuchos::ParameterList ¶ms );
289 if (!initialized_)
return 0;
317 void setSize(
int recycledBlocks,
int numBlocks ) {
319 if ( (recycledBlocks_ != recycledBlocks) || (numBlocks_ != numBlocks) ) {
320 recycledBlocks_ = recycledBlocks;
321 numBlocks_ = numBlocks;
322 cs_.sizeUninitialized( numBlocks_ );
323 sn_.sizeUninitialized( numBlocks_ );
324 Z_.shapeUninitialized( numBlocks_*blockSize_, blockSize_ );
340 const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > lp_;
341 const Teuchos::RCP<OutputManager<ScalarType> > om_;
342 const Teuchos::RCP<StatusTest<ScalarType,MV,OP> > stest_;
343 const Teuchos::RCP<OrthoManager<ScalarType,MV> > ortho_;
351 int numBlocks_, blockSize_;
356 std::vector<bool> trueRHSIndices_;
363 Teuchos::SerialDenseVector<int,MagnitudeType> cs_;
370 std::vector< SDM >House_;
382 int curDim_, iter_, lclIter_;
393 Teuchos::RCP<MV> U_, C_;
402 Teuchos::RCP<SDM > H_;
407 Teuchos::RCP<SDM > B_;
415 Teuchos::RCP<SDM> R_;
435 Teuchos::ParameterList ¶ms ):lp_(problem),
436 om_(printer), stest_(tester), ortho_(ortho) {
440 initialized_ =
false;
451 TEUCHOS_TEST_FOR_EXCEPTION(!params.isParameter(
"Num Blocks"), std::invalid_argument,
"Belos::BlockGCRODRIter::constructor: mandatory parameter \"Num Blocks\" is not specified.");
452 int nb = Teuchos::getParameter<int>(params,
"Num Blocks");
454 TEUCHOS_TEST_FOR_EXCEPTION(!params.isParameter(
"Recycled Blocks"), std::invalid_argument,
"Belos::BlockGCRODRIter::constructor: mandatory parameter \"Recycled Blocks\" is not specified.");
455 int rb = Teuchos::getParameter<int>(params,
"Recycled Blocks");
457 TEUCHOS_TEST_FOR_EXCEPTION(nb <= 0, std::invalid_argument,
"Belos::BlockGCRODRIter() was passed a non-positive argument for \"Num Blocks\".");
458 TEUCHOS_TEST_FOR_EXCEPTION(rb >= nb, std::invalid_argument,
"Belos::BlockGCRODRIter() the number of recycled blocks is larger than the allowable subspace.");
461 int bs = Teuchos::getParameter<int>(params,
"Block Size");
463 TEUCHOS_TEST_FOR_EXCEPTION(bs <= 0, std::invalid_argument,
"Belos::BlockGCRODRIter() the block size was passed a non-postitive argument.");
467 recycledBlocks_ = rb;
472 trueRHSIndices_.resize(blockSize_);
474 for(i=0; i<blockSize_; i++){
475 trueRHSIndices_[i] =
true;
480 cs_.sizeUninitialized( numBlocks_+1 );
481 sn_.sizeUninitialized( numBlocks_+1 );
482 Z_.shapeUninitialized( (numBlocks_+1)*blockSize_,blockSize_ );
484 House_.resize(numBlocks_);
486 for(i=0; i<numBlocks_;i++){
487 House_[i].shapeUninitialized(2*blockSize_, 2*blockSize_);
495 TEUCHOS_TEST_FOR_EXCEPTION( initialized_ ==
false,
BlockGCRODRIterInitFailure,
"Belos::BlockGCRODRIter::iterate(): GCRODRIter class not initialized." );
501 setSize( recycledBlocks_, numBlocks_ );
503 Teuchos::RCP<MV> Vnext;
504 Teuchos::RCP<const MV> Vprev;
505 std::vector<int> curind(blockSize_);
511 for(
int i = 0; i<blockSize_; i++){curind[i] = i;};
519 Teuchos::RCP<SDM > Z0 =
520 Teuchos::rcp(
new SDM(blockSize_,blockSize_) );
521 int rank = ortho_->normalize(*Vnext,Z0);
526 TEUCHOS_TEST_FOR_EXCEPTION(rank != blockSize_,
BlockGCRODRIterOrthoFailure,
"Belos::BlockGCRODRIter::iterate(): couldn't generate basis of full rank at the initial step.");
528 Teuchos::RCP<SDM > Z_block = Teuchos::rcp(
new SDM(Teuchos::View, Z_, blockSize_,blockSize_) );
529 Z_block->assign(*Z0);
531 std::vector<int> prevind(blockSize_*(numBlocks_ + 1));
538 while( (stest_->checkStatus(
this) !=
Passed) && (curDim_+blockSize_-1) < (numBlocks_*blockSize_)) {
544 int HFirstCol = curDim_-blockSize_;
545 int HLastCol = HFirstCol + blockSize_-1 ;
546 int HLastOrthRow = HLastCol;
547 int HFirstNormRow = HLastOrthRow + 1;
551 for(
int i = 0; i< blockSize_; i++){
552 curind[i] = curDim_ + i;
559 for(
int i = 0; i< blockSize_; i++){
560 curind[blockSize_ - 1 - i] = curDim_ - i - 1;
564 lp_->apply(*Vprev,*Vnext);
565 Vprev = Teuchos::null;
571 Teuchos::Array<Teuchos::RCP<const MV> > C(1, C_);
574 subB = Teuchos::rcp(
new SDM ( Teuchos::View,*B_,recycledBlocks_,blockSize_,0, HFirstCol ) );
576 Teuchos::Array<Teuchos::RCP<SDM > > AsubB;
577 AsubB.append( subB );
579 ortho_->project( *Vnext, AsubB, C );
583 prevind.resize(curDim_);
584 for (
int i=0; i<curDim_; i++) { prevind[i] = i; }
586 Teuchos::Array<Teuchos::RCP<const MV> > AVprev(1, Vprev);
589 Teuchos::RCP<SDM> subH = Teuchos::rcp(
new SDM ( Teuchos::View,*H_,curDim_,blockSize_,0,HFirstCol ) );
590 Teuchos::Array<Teuchos::RCP<SDM > > AsubH;
591 AsubH.append( subH );
593 Teuchos::RCP<SDM > subR = Teuchos::rcp(
new SDM ( Teuchos::View,*H_,blockSize_,blockSize_,HFirstNormRow,HFirstCol ) );
595 int rank = ortho_->projectAndNormalize(*Vnext,AsubH,subR,AVprev);
598 SDM subR2( Teuchos::View,*R_,(lclIter_+1)*blockSize_,blockSize_,0,HFirstCol);
599 SDM subH2( Teuchos::View,*H_,(lclIter_+1)*blockSize_,blockSize_,0,HFirstCol);
602 TEUCHOS_TEST_FOR_EXCEPTION(rank != blockSize_,
BlockGCRODRIterOrthoFailure,
"Belos::BlockGCRODRIter::iterate(): couldn't generate basis of full rank.");
608 curDim_ = curDim_ + blockSize_;
689 Teuchos::RCP<MV> currentUpdate = Teuchos::null;
691 if(curDim_<=blockSize_) {
692 return currentUpdate;
695 const ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
696 const ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
697 Teuchos::BLAS<int,ScalarType> blas;
698 currentUpdate =
MVT::Clone( *V_, blockSize_ );
702 SDM Y( Teuchos::Copy, Z_, curDim_-blockSize_, blockSize_ );
703 Teuchos::RCP<SDM> Rtmp = Teuchos::rcp(
new SDM(Teuchos::View, *R_, curDim_, curDim_-blockSize_));
720 blas.TRSM( Teuchos::LEFT_SIDE, Teuchos::UPPER_TRI, Teuchos::NO_TRANS,
721 Teuchos::NON_UNIT_DIAG, curDim_-blockSize_, blockSize_, one,
722 Rtmp->values(), Rtmp->stride(), Y.values(), Y.stride() );
729 std::vector<int> index(curDim_-blockSize_);
730 for (
int i=0; i<curDim_-blockSize_; i++ ) index[i] = i;
739 if (U_ != Teuchos::null) {
740 SDM z(recycledBlocks_,blockSize_);
741 SDM subB( Teuchos::View, *B_, recycledBlocks_, curDim_-blockSize_ );
742 z.multiply( Teuchos::NO_TRANS, Teuchos::NO_TRANS, one, subB, Y, zero );
751 return currentUpdate;
758 const ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
759 const ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
766 int curDim = curDim_;
771 Teuchos::BLAS<int, ScalarType> blas;
780 for (i=0; i<curDim-1; i++) {
784 blas.ROT( 1, &(*R_)(i,curDim-1), 1, &(*R_)(i+1, curDim-1), 1, &cs_[i], &sn_[i] );
789 blas.ROTG( &(*R_)(curDim-1,curDim-1), &(*R_)(curDim,curDim-1), &cs_[curDim-1], &sn_[curDim-1] );
790 (*R_)(curDim,curDim-1) = zero;
794 blas.ROT( 1, &Z_(curDim-1,0), 1, &Z_(curDim,0), 1, &cs_[curDim-1], &sn_[curDim-1] );
809 Teuchos::RCP< SDM > workmatrix = Teuchos::null;
810 Teuchos::RCP< SDV > workvec = Teuchos::null;
811 Teuchos::RCP<SDV> v_refl = Teuchos::null;
812 int R_colStart = curDim_-blockSize_;
813 Teuchos::RCP< SDM >Rblock = Teuchos::null;
818 for(i=0; i<lclIter_-1; i++){
819 int R_rowStart = i*blockSize_;
821 Teuchos::RCP< SDM > RblockCopy = rcp(
new SDM (Teuchos::Copy, *R_, 2*blockSize_,blockSize_, R_rowStart, R_colStart));
822 Teuchos::RCP< SDM > RblockView = rcp(
new SDM (Teuchos::View, *R_, 2*blockSize_,blockSize_, R_rowStart, R_colStart));
823 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());
830 Rblock = rcp(
new SDM (Teuchos::View, *R_, 2*blockSize_,blockSize_, curDim_-blockSize_, curDim_-blockSize_));
833 for(i=0; i<blockSize_; i++){
837 int curcol = (lclIter_ - 1)*blockSize_ + i;
839 ScalarType signDiag = (*R_)(curcol,curcol) / Teuchos::ScalarTraits<ScalarType>::magnitude((*R_)(curcol,curcol));
843 ScalarType nvs = blas.NRM2(blockSize_+1,&((*R_)[curcol][curcol]),1);
844 ScalarType alpha = -signDiag*nvs;
862 v_refl = rcp(
new SDV(Teuchos::Copy, &((*R_)(curcol,curcol)), blockSize_ + 1 ));
863 (*v_refl)[0] -= alpha;
864 nvs = blas.NRM2(blockSize_+1,v_refl -> values(),1);
880 if(i < blockSize_ - 1){
881 workvec = Teuchos::rcp(
new SDV(blockSize_ - i -1));
883 workmatrix = Teuchos::rcp(
new SDM (Teuchos::View, *Rblock, blockSize_+1, blockSize_ - i -1, lclCurcol, lclCurcol +1 ) );
884 blas.GEMV(Teuchos::TRANS, workmatrix->numRows(), workmatrix->numCols(), one, workmatrix->values(), workmatrix->stride(), v_refl->values(), 1, zero, workvec->values(), 1);
885 blas.GER(workmatrix->numRows(),workmatrix->numCols(), -2.0*one/nvs, v_refl->values(),1,workvec->values(),1,workmatrix->values(),workmatrix->stride());
892 workvec = Teuchos::rcp(
new SDV(2*blockSize_));
893 workmatrix = Teuchos::rcp(
new SDM (Teuchos::View, House_[lclIter_ -1], blockSize_+1, 2*blockSize_, i, 0 ) );
894 blas.GEMV(Teuchos::TRANS,workmatrix->numRows(),workmatrix->numCols(),one,workmatrix->values(),workmatrix->stride(), v_refl->values(), 1,zero,workvec->values(),1);
895 blas.GER(workmatrix->numRows(),workmatrix->numCols(), -2.0*one/nvs, v_refl -> values(),1,workvec->values(),1,workmatrix->values(),(*workmatrix).stride());
900 workvec = Teuchos::rcp(
new SDV(blockSize_));
901 workmatrix = Teuchos::rcp(
new SDM (Teuchos::View, Z_, blockSize_+1, blockSize_, curcol, 0 ) );
902 blas.GEMV(Teuchos::TRANS, workmatrix->numRows(), workmatrix->numCols(), one, workmatrix-> values(), workmatrix->stride(), v_refl -> values(), 1, zero, workvec->values(), 1);
903 blas.GER((*workmatrix).numRows(),(*workmatrix).numCols(), -2.0*one/nvs,v_refl -> values(), 1,&((*workvec)[0]),1,(*workmatrix)[0],(*workmatrix).stride());
908 (*R_)[curcol][curcol] = alpha;
909 for(
int ii=1; ii<= blockSize_; ii++){
910 (*R_)[curcol][curcol+ii] = 0;