42 #ifndef BELOS_GMRESPOLYOP_HPP
43 #define BELOS_GMRESPOLYOP_HPP
73 #include "Teuchos_BLAS.hpp"
74 #include "Teuchos_LAPACK.hpp"
75 #include "Teuchos_as.hpp"
76 #include "Teuchos_RCP.hpp"
77 #include "Teuchos_SerialDenseMatrix.hpp"
78 #include "Teuchos_SerialDenseVector.hpp"
79 #include "Teuchos_SerialDenseSolver.hpp"
80 #include "Teuchos_ParameterList.hpp"
82 #ifdef BELOS_TEUCHOS_TIME_MONITOR
83 #include "Teuchos_TimeMonitor.hpp"
99 template <
class ScalarType,
class MV>
109 mv_ = Teuchos::rcp_const_cast<MV>( mv_in );
111 Teuchos::RCP<MV>
getMV() {
return mv_; }
143 const Teuchos::SerialDenseMatrix<int,ScalarType>& B,
const ScalarType beta)
166 void MvNorm ( std::vector<
typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>& normvec,
NormType type =
TwoNorm )
const
183 Teuchos::RCP<MV> mv_;
197 template <
class ScalarType,
class MV,
class OP>
206 const Teuchos::RCP<Teuchos::ParameterList>& params_in
208 : problem_(problem_in),
210 LP_(problem_in->getLeftPrec()),
211 RP_(problem_in->getRightPrec())
215 polyUpdateLabel_ = label_ +
": Hybrid Gmres: Vector Update";
216 #ifdef BELOS_TEUCHOS_TIME_MONITOR
217 timerPolyUpdate_ = Teuchos::TimeMonitor::getNewCounter(polyUpdateLabel_);
220 if (polyType_ ==
"Arnoldi" || polyType_==
"Roots")
222 else if (polyType_ ==
"Gmres")
225 TEUCHOS_TEST_FOR_EXCEPTION(polyType_!=
"Arnoldi"&&polyType_!=
"Gmres"&&polyType_!=
"Roots",std::invalid_argument,
226 "Belos::GmresPolyOp: \"Polynomial Type\" must be either \"Arnoldi\", \"Gmres\", or \"Roots\".");
231 : problem_(problem_in)
245 void setParameters(
const Teuchos::RCP<Teuchos::ParameterList>& params_in );
271 void ApplyPoly (
const MV& x, MV& y )
const;
290 #ifdef BELOS_TEUCHOS_TIME_MONITOR
291 Teuchos::RCP<Teuchos::Time> timerPolyUpdate_;
293 std::string polyUpdateLabel_;
297 typedef Teuchos::ScalarTraits<ScalarType> SCT ;
298 typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType MagnitudeType;
299 typedef Teuchos::ScalarTraits<MagnitudeType> MCT ;
302 static constexpr
int maxDegree_default_ = 25;
304 static constexpr
bool randomRHS_default_ =
true;
305 static constexpr
const char * label_default_ =
"Belos";
306 static constexpr
const char * polyType_default_ =
"Roots";
307 static constexpr
const char * orthoType_default_ =
"DGKS";
308 static constexpr std::ostream * outputStream_default_ = &std::cout;
309 static constexpr
bool damp_default_ =
false;
310 static constexpr
bool addRoots_default_ =
true;
313 Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > problem_;
314 Teuchos::RCP<Teuchos::ParameterList> params_;
315 Teuchos::RCP<const OP> LP_, RP_;
318 Teuchos::RCP<OutputManager<ScalarType> > printer_;
319 Teuchos::RCP<std::ostream> outputStream_ = Teuchos::rcp(outputStream_default_,
false);
322 Teuchos::RCP<MatOrthoManager<ScalarType,MV,OP> > ortho_;
326 int maxDegree_ = maxDegree_default_;
327 int verbosity_ = verbosity_default_;
328 bool randomRHS_ = randomRHS_default_;
329 std::string label_ = label_default_;
330 std::string polyType_ = polyType_default_;
331 std::string orthoType_ = orthoType_default_;
333 bool damp_ = damp_default_;
334 bool addRoots_ = addRoots_default_;
337 mutable Teuchos::RCP<MV> V_, wL_, wR_;
338 Teuchos::SerialDenseMatrix<OT,ScalarType> H_, y_;
339 Teuchos::SerialDenseVector<OT,ScalarType> r0_;
342 bool autoDeg =
false;
343 Teuchos::SerialDenseMatrix< OT, ScalarType > pCoeff_;
346 Teuchos::SerialDenseMatrix< OT, MagnitudeType > theta_;
350 void SortModLeja(Teuchos::SerialDenseMatrix< OT, MagnitudeType > &thetaN, std::vector<int> &index)
const ;
353 void ComputeAddedRoots();
356 template <
class ScalarType,
class MV,
class OP>
360 if (params_in->isParameter(
"Polynomial Type")) {
361 polyType_ = params_in->get(
"Polynomial Type", polyType_default_);
365 if (params_in->isParameter(
"Polynomial Tolerance")) {
366 if (params_in->isType<MagnitudeType> (
"Polynomial Tolerance")) {
367 polyTol_ = params_in->get (
"Polynomial Tolerance",
376 if (params_in->isParameter(
"Maximum Degree")) {
377 maxDegree_ = params_in->get(
"Maximum Degree", maxDegree_default_);
381 if (params_in->isParameter(
"Random RHS")) {
382 randomRHS_ = params_in->get(
"Random RHS", randomRHS_default_);
386 if (params_in->isParameter(
"Verbosity")) {
387 if (Teuchos::isParameterType<int>(*params_in,
"Verbosity")) {
388 verbosity_ = params_in->get(
"Verbosity", verbosity_default_);
391 verbosity_ = (int)Teuchos::getParameter<Belos::MsgType>(*params_in,
"Verbosity");
395 if (params_in->isParameter(
"Orthogonalization")) {
396 orthoType_ = params_in->get(
"Orthogonalization",orthoType_default_);
397 TEUCHOS_TEST_FOR_EXCEPTION( orthoType_ !=
"DGKS" && orthoType_ !=
"ICGS" && orthoType_ !=
"IMGS",
398 std::invalid_argument,
399 "Belos::GmresPolyOp: \"Orthogonalization\" must be either \"DGKS\", \"ICGS\", or \"IMGS\".");
403 if (params_in->isParameter(
"Timer Label")) {
404 label_ = params_in->get(
"Timer Label", label_default_);
408 if (params_in->isParameter(
"Output Stream")) {
409 outputStream_ = Teuchos::getParameter<Teuchos::RCP<std::ostream> >(*params_in,
"Output Stream");
413 if (params_in->isParameter(
"Damped Poly")) {
414 damp_ = params_in->get(
"Damped Poly", damp_default_);
418 if (params_in->isParameter(
"Add Roots")) {
419 addRoots_ = params_in->get(
"Add Roots", addRoots_default_);
423 template <
class ScalarType,
class MV,
class OP>
426 Teuchos::RCP< MV > V = MVT::Clone( *problem_->getRHS(), maxDegree_+2 );
429 std::vector<int> index(1,0);
430 Teuchos::RCP< MV > V0 = MVT::CloneViewNonConst(*V, index);
432 MVT::MvRandom( *V0 );
434 MVT::Assign( *problem_->getRHS(), *V0 );
436 if ( !LP_.is_null() ) {
437 Teuchos::RCP< MV > Vtemp = MVT::CloneCopy(*V0);
438 problem_->applyLeftPrec( *Vtemp, *V0);
441 Teuchos::RCP< MV > Vtemp = MVT::CloneCopy(*V0);
442 problem_->apply( *Vtemp, *V0);
445 for(
int i=0; i< maxDegree_+1; i++)
448 Teuchos::RCP< const MV > Vi = MVT::CloneView(*V, index);
450 Teuchos::RCP< MV > Vip1 = MVT::CloneViewNonConst(*V, index);
451 problem_->apply( *Vi, *Vip1);
455 Teuchos::Range1D range( 1, maxDegree_+1);
456 Teuchos::RCP< const MV > AV = MVT::CloneView( *V, range);
459 Teuchos::SerialDenseMatrix< OT, ScalarType > AVtransAV( maxDegree_+1, maxDegree_+1);
460 MVT::MvTransMv( SCT::one(), *AV, *AV, AVtransAV);
463 Teuchos::LAPACK< OT, ScalarType > lapack;
468 Teuchos::SerialDenseMatrix< OT, ScalarType > lhs;
469 while( status && dim_ >= 1)
471 Teuchos::SerialDenseMatrix< OT, ScalarType > lhstemp(Teuchos::Copy, AVtransAV, dim_+1, dim_+1);
472 lapack.POTRF(
'U', dim_+1, lhstemp.values(), lhstemp.stride(), &infoInt);
479 std::cout <<
"BelosGmresPolyOp.hpp: LAPACK POTRF was not successful!!" << std::endl;
480 std::cout <<
"Error code: " << infoInt << std::endl;
501 pCoeff_.shape( 1, 1);
502 pCoeff_(0,0) = SCT::one();
503 std::cout <<
"Poly Degree is zero. No preconditioner created." << std::endl;
507 pCoeff_.shape( dim_+1, 1);
509 Teuchos::Range1D rangeSub( 1, dim_+1);
510 Teuchos::RCP< const MV > AVsub = MVT::CloneView( *V, rangeSub);
513 MVT::MvTransMv( SCT::one(), *AVsub, *V0, pCoeff_);
514 lapack.POTRS(
'U', dim_+1, 1, lhs.values(), lhs.stride(), pCoeff_.values(), pCoeff_.stride(), &infoInt);
517 std::cout <<
"BelosGmresPolyOp.hpp: LAPACK POTRS was not successful!!" << std::endl;
518 std::cout <<
"Error code: " << infoInt << std::endl;
523 template <
class ScalarType,
class MV,
class OP>
526 std::string polyLabel = label_ +
": GmresPolyOp creation";
529 std::vector<int> idx(1,0);
530 Teuchos::RCP<MV> newX = MVT::Clone( *(problem_->getLHS()), 1 );
531 Teuchos::RCP<MV> newB = MVT::Clone( *(problem_->getRHS()), 1 );
532 MVT::MvInit( *newX, SCT::zero() );
534 MVT::MvRandom( *newB );
537 MVT::Assign( *(MVT::CloneView(*(problem_->getRHS()), idx)), *newB );
539 Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > newProblem =
541 newProblem->setInitResVec( newB );
542 newProblem->setLeftPrec( problem_->getLeftPrec() );
543 newProblem->setRightPrec( problem_->getRightPrec() );
544 newProblem->setLabel(polyLabel);
545 newProblem->setProblem();
546 newProblem->setLSIndex( idx );
549 Teuchos::ParameterList polyList;
552 polyList.set(
"Num Blocks",maxDegree_);
553 polyList.set(
"Block Size",1);
554 polyList.set(
"Keep Hessenberg",
true);
557 if (ortho_.is_null()) {
558 params_->set(
"Orthogonalization", orthoType_);
559 if (orthoType_==
"DGKS") {
562 else if (orthoType_==
"ICGS") {
565 else if (orthoType_==
"IMGS") {
569 TEUCHOS_TEST_FOR_EXCEPTION(orthoType_!=
"ICGS"&&orthoType_!=
"DGKS"&&orthoType_!=
"IMGS",std::invalid_argument,
570 "Belos::GmresPolyOp(): Invalid orthogonalization type.");
578 Teuchos::RCP<StatusTestMaxIters<ScalarType,MV,OP> > maxItrTst =
582 Teuchos::RCP<StatusTestGenResNorm<ScalarType,MV,OP> > convTst =
587 Teuchos::RCP<StatusTestCombo<ScalarType,MV,OP> > polyTest =
591 Teuchos::RCP<BlockGmresIter<ScalarType,MV,OP> > gmres_iter;
595 Teuchos::RCP<MV> V_0 = MVT::CloneCopy( *newB );
596 if ( !LP_.is_null() )
597 newProblem->applyLeftPrec( *newB, *V_0 );
600 Teuchos::RCP< MV > Vtemp = MVT::CloneCopy(*V_0);
601 newProblem->apply( *Vtemp, *V_0 );
608 int rank = ortho_->normalize( *V_0, Teuchos::rcpFromRef(r0_) );
610 "Belos::GmresPolyOp::generateArnoldiPoly(): Failed to compute initial block of orthonormal vectors for polynomial generation.");
615 newstate.
z = Teuchos::rcpFromRef( r0_);
617 gmres_iter->initializeGmres(newstate);
621 gmres_iter->iterate();
625 gmres_iter->updateLSQR( gmres_iter->getCurSubspaceDim() );
627 catch (std::exception e) {
629 printer_->stream(
Errors) <<
"Error! Caught exception in BlockGmresIter::iterate() at iteration "
630 << gmres_iter->getNumIters() << endl << e.what () << endl;
635 Teuchos::RCP<MV> currX = gmres_iter->getCurrentUpdate();
645 if(polyType_ ==
"Arnoldi"){
648 y_ = Teuchos::SerialDenseMatrix<OT,ScalarType>( Teuchos::Copy, *gmresState.
z, dim_, 1 );
654 Teuchos::BLAS<OT,ScalarType> blas;
655 blas.TRSM( Teuchos::LEFT_SIDE, Teuchos::UPPER_TRI, Teuchos::NO_TRANS,
656 Teuchos::NON_UNIT_DIAG, dim_, 1, SCT::one(),
657 gmresState.
R->values(), gmresState.
R->stride(),
658 y_.values(), y_.stride() );
664 H_ = Teuchos::SerialDenseMatrix<OT,ScalarType>(Teuchos::Copy, *gmresState.
H, dim_, dim_);
666 for(
int i=0; i <= dim_-3; i++) {
667 for(
int k=i+2; k <= dim_-1; k++) {
668 H_(k,i) = SCT::zero();
672 Teuchos::SerialDenseMatrix<OT,ScalarType> Htemp (Teuchos::Copy, H_, dim_, dim_);
675 ScalarType Hlast = (*gmresState.
H)(dim_,dim_-1);
676 Teuchos::SerialDenseMatrix<OT,ScalarType> HlastCol (Teuchos::View, H_, dim_, 1, 0, dim_-1);
679 Teuchos::SerialDenseMatrix< OT, ScalarType > F(dim_,1), E(dim_,1);
680 E.putScalar(SCT::zero());
681 E(dim_-1,0) = SCT::one();
683 Teuchos::SerialDenseSolver< OT, ScalarType > HSolver;
684 HSolver.setMatrix( Teuchos::rcpFromRef(Htemp));
685 HSolver.solveWithTransposeFlag( Teuchos::CONJ_TRANS );
686 HSolver.setVectors( Teuchos::rcpFromRef(F), Teuchos::rcpFromRef(E));
687 HSolver.factorWithEquilibration(
true );
691 info = HSolver.factor();
693 std::cout <<
"Hsolver factor: info = " << info << std::endl;
695 info = HSolver.solve();
697 std::cout <<
"Hsolver solve : info = " << info << std::endl;
701 F.scale(Hlast*Hlast);
705 Teuchos::LAPACK< OT, ScalarType > lapack;
706 theta_.shape(dim_,2);
713 std::vector<ScalarType> work(1);
714 std::vector<MagnitudeType> rwork(2*dim_);
717 lapack.GEEV(
'N',
'N',dim_,H_.values(),H_.stride(),theta_[0],theta_[1],vlr, ldv, vlr, ldv, &work[0], lwork, &rwork[0], &info);
718 lwork = std::abs (
static_cast<int> (Teuchos::ScalarTraits<ScalarType>::real (work[0])));
719 work.resize( lwork );
721 lapack.GEEV(
'N',
'N',dim_,H_.values(),H_.stride(),theta_[0],theta_[1],vlr, ldv, vlr, ldv, &work[0], lwork, &rwork[0], &info);
724 std::cout <<
"GEEV solve : info = " << info << std::endl;
728 std::vector<int> index(dim_);
729 for(
int i=0; i<dim_; ++i){
732 SortModLeja(theta_,index);
741 template <
class ScalarType,
class MV,
class OP>
746 std::vector<std::complex<MagnitudeType>> cmplxHRitz (dim_);
747 for(
unsigned int i=0; i<cmplxHRitz.size(); ++i){
748 cmplxHRitz[i] = std::complex<MagnitudeType>( theta_(i,0), theta_(i,1) );
752 const MagnitudeType one(1.0);
753 std::vector<MagnitudeType> pof (dim_,one);
754 for(
int j=0; j<dim_; ++j) {
755 for(
int i=0; i<dim_; ++i) {
757 pof[j] = std::abs(pof[j]*(one-(cmplxHRitz[j]/cmplxHRitz[i])));
763 std::vector<int> extra (dim_);
765 for(
int i=0; i<dim_; ++i){
766 extra[i] = ceil((log10(pof[i])-MagnitudeType(4.0))/MagnitudeType(14.0));
768 totalExtra += extra[i];
772 printer_->stream(
Warnings) <<
"Warning: Need to add " << totalExtra <<
" extra roots." << std::endl;}
775 if(addRoots_ && totalExtra>0)
777 theta_.reshape(dim_+totalExtra,2);
779 Teuchos::SerialDenseMatrix<OT,MagnitudeType> thetaPert (Teuchos::Copy, theta_, dim_+totalExtra, 2);
783 for(
int i=0; i<dim_; ++i){
784 for(
int j=0; j< extra[i]; ++j){
785 theta_(count,0) = theta_(i,0);
786 theta_(count,1) = theta_(i,1);
787 thetaPert(count,0) = theta_(i,0)+(j+MCT::one())*MagnitudeType(5e-8);
788 thetaPert(count,1) = theta_(i,1);
796 printer_->stream(
Warnings) <<
"New poly degree is: " << dim_ << std::endl;}
799 std::vector<int> index2(dim_);
800 for(
int i=0; i<dim_; ++i){
803 SortModLeja(thetaPert,index2);
805 for(
int i=0; i<dim_; ++i)
807 thetaPert(i,0) = theta_(index2[i],0);
808 thetaPert(i,1) = theta_(index2[i],1);
817 template <
class ScalarType,
class MV,
class OP>
818 void GmresPolyOp<ScalarType, MV, OP>::SortModLeja(Teuchos::SerialDenseMatrix< OT, MagnitudeType > &thetaN, std::vector<int> &index)
const
823 int dimN = index.size();
824 std::vector<int> newIndex(dimN);
825 Teuchos::SerialDenseMatrix< OT, MagnitudeType > sorted (thetaN.numRows(), thetaN.numCols());
826 Teuchos::SerialDenseVector< OT, MagnitudeType > absVal (thetaN.numRows());
827 Teuchos::SerialDenseVector< OT, MagnitudeType > prod (thetaN.numRows());
830 for(
int i = 0; i < dimN; i++){
831 absVal(i) = hypot(thetaN(i,0), thetaN(i,1));
833 MagnitudeType * maxPointer = std::max_element(absVal.values(), (absVal.values()+dimN));
834 int maxIndex = int (maxPointer- absVal.values());
837 sorted(0,0) = thetaN(maxIndex,0);
838 sorted(0,1) = thetaN(maxIndex,1);
839 newIndex[0] = index[maxIndex];
843 if(sorted(0,1)!= SCT::zero() && !SCT::isComplex)
845 sorted(1,0) = thetaN(maxIndex,0);
846 sorted(1,1) = -thetaN(maxIndex,1);
847 newIndex[1] = index[maxIndex+1];
860 for(
int i = 0; i < dimN; i++)
862 prod(i) = MCT::one();
863 for(
int k = 0; k < j; k++)
865 a = thetaN(i,0) - sorted(k,0);
866 b = thetaN(i,1) - sorted(k,1);
867 prod(i) = prod(i) + log10(sqrt(a*a + b*b));
872 MagnitudeType * maxPointer = std::max_element(prod.values(), (prod.values()+dimN));
873 int maxIndex = int (maxPointer- prod.values());
874 sorted(j,0) = thetaN(maxIndex,0);
875 sorted(j,1) = thetaN(maxIndex,1);
876 newIndex[j] = index[maxIndex];
879 if(sorted(j,1)!= SCT::zero() && !SCT::isComplex)
882 sorted(j,0) = thetaN(maxIndex,0);
883 sorted(j,1) = -thetaN(maxIndex,1);
884 newIndex[j] = index[maxIndex+1];
894 template <
class ScalarType,
class MV,
class OP>
898 if (polyType_ ==
"Arnoldi")
899 ApplyArnoldiPoly(x, y);
900 else if (polyType_ ==
"Gmres")
901 ApplyGmresPoly(x, y);
902 else if (polyType_ ==
"Roots")
903 ApplyRootsPoly(x, y);
907 problem_->applyOp( x, y );
911 template <
class ScalarType,
class MV,
class OP>
914 Teuchos::RCP<MV> AX = MVT::CloneCopy(x);
915 Teuchos::RCP<MV> AX2 = MVT::Clone( x, MVT::GetNumberVecs(x) );
918 if (!LP_.is_null()) {
919 Teuchos::RCP<MV> Xtmp = MVT::Clone( x, MVT::GetNumberVecs(x) );
920 problem_->applyLeftPrec( *AX, *Xtmp );
925 #ifdef BELOS_TEUCHOS_TIME_MONITOR
926 Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
928 MVT::MvAddMv(pCoeff_(0,0), *AX, SCT::zero(), y, y);
930 for(
int i=1; i < dim_+1; i++)
932 Teuchos::RCP<MV> X, Y;
943 problem_->apply(*X, *Y);
945 #ifdef BELOS_TEUCHOS_TIME_MONITOR
946 Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
948 MVT::MvAddMv(pCoeff_(i,0), *Y, SCT::one(), y, y);
953 if (!RP_.is_null()) {
954 Teuchos::RCP<MV> Ytmp = MVT::CloneCopy(y);
955 problem_->applyRightPrec( *Ytmp, y );
959 template <
class ScalarType,
class MV,
class OP>
962 MVT::MvInit( y, SCT::zero() );
963 Teuchos::RCP<MV> prod = MVT::CloneCopy(x);
964 Teuchos::RCP<MV> Xtmp = MVT::Clone( x, MVT::GetNumberVecs(x) );
965 Teuchos::RCP<MV> Xtmp2 = MVT::Clone( x, MVT::GetNumberVecs(x) );
968 if (!LP_.is_null()) {
969 problem_->applyLeftPrec( *prod, *Xtmp );
976 if(theta_(i,1)== SCT::zero() || SCT::isComplex)
979 #ifdef BELOS_TEUCHOS_TIME_MONITOR
980 Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
982 MVT::MvAddMv(SCT::one(), y, SCT::one()/theta_(i,0), *prod, y);
984 problem_->apply(*prod, *Xtmp);
986 #ifdef BELOS_TEUCHOS_TIME_MONITOR
987 Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
989 MVT::MvAddMv(SCT::one(), *prod, -SCT::one()/theta_(i,0), *Xtmp, *prod);
995 MagnitudeType mod = theta_(i,0)*theta_(i,0) + theta_(i,1)*theta_(i,1);
996 problem_->apply(*prod, *Xtmp);
998 #ifdef BELOS_TEUCHOS_TIME_MONITOR
999 Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
1001 MVT::MvAddMv(2*theta_(i,0), *prod, -SCT::one(), *Xtmp, *Xtmp);
1002 MVT::MvAddMv(SCT::one(), y, SCT::one()/mod, *Xtmp, y);
1006 problem_->apply(*Xtmp, *Xtmp2);
1008 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1009 Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
1011 MVT::MvAddMv(SCT::one(), *prod, -SCT::one()/mod, *Xtmp2, *prod);
1017 if(theta_(dim_-1,1)== SCT::zero() || SCT::isComplex)
1019 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1020 Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
1022 MVT::MvAddMv(SCT::one(), y, SCT::one()/theta_(dim_-1,0), *prod, y);
1026 if (!RP_.is_null()) {
1027 Teuchos::RCP<MV> Ytmp = MVT::CloneCopy(y);
1028 problem_->applyRightPrec( *Ytmp, y );
1032 template <
class ScalarType,
class MV,
class OP>
1037 V_ = MVT::Clone( x, dim_ );
1038 if (!LP_.is_null()) {
1039 wL_ = MVT::Clone( y, 1 );
1041 if (!RP_.is_null()) {
1042 wR_ = MVT::Clone( y, 1 );
1048 int n = MVT::GetNumberVecs( x );
1049 std::vector<int> idxi(1), idxi2, idxj(1);
1052 for (
int j=0; j<n; ++j) {
1056 Teuchos::RCP<const MV> x_view = MVT::CloneView( x, idxj );
1057 Teuchos::RCP<MV> y_view = MVT::CloneViewNonConst( y, idxj );
1058 if (!LP_.is_null()) {
1059 Teuchos::RCP<MV> v_curr = MVT::CloneViewNonConst( *V_, idxi );
1060 problem_->applyLeftPrec( *x_view, *v_curr );
1062 MVT::SetBlock( *x_view, idxi, *V_ );
1065 for (
int i=0; i<dim_-1; ++i) {
1069 for (
int ii=0; ii<i+1; ++ii) { idxi2[ii] = ii; }
1070 Teuchos::RCP<const MV> v_prev = MVT::CloneView( *V_, idxi2 );
1073 Teuchos::RCP<MV> v_curr = MVT::CloneViewNonConst( *V_, idxi );
1075 Teuchos::RCP<MV> v_next = MVT::CloneViewNonConst( *V_, idxi );
1081 if (!RP_.is_null()) {
1082 problem_->applyRightPrec( *v_curr, *wR_ );
1087 if (LP_.is_null()) {
1091 problem_->applyOp( *wR_, *wL_ );
1093 if (!LP_.is_null()) {
1094 problem_->applyLeftPrec( *wL_, *v_next );
1098 Teuchos::SerialDenseMatrix<OT,ScalarType> h(Teuchos::View,H_,i+1,1,0,i);
1100 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1101 Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
1103 MVT::MvTimesMatAddMv( -SCT::one(), *v_prev, h, SCT::one(), *v_next );
1107 MVT::MvScale( *v_next, SCT::one()/H_(i+1,i) );
1111 if (!RP_.is_null()) {
1113 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1114 Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
1116 MVT::MvTimesMatAddMv( SCT::one()/r0_(0), *V_, y_, SCT::zero(), *wR_ );
1118 problem_->applyRightPrec( *wR_, *y_view );
1121 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1122 Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
1124 MVT::MvTimesMatAddMv( SCT::one()/r0_(0), *V_, y_, SCT::zero(), *y_view );
Belos concrete class for performing the block GMRES iteration.
Belos header file which uses auto-configuration information to include necessary C++ headers.
Classical Gram-Schmidt (with DGKS correction) implementation of the Belos::OrthoManager class.
Pure virtual base class which augments the basic interface for a Gmres linear solver iteration.
Iterated Classical Gram-Schmidt (ICGS) implementation of the Belos::OrthoManager class.
Iterated Modified Gram-Schmidt (IMGS) implementation of the Belos::OrthoManager class.
Class which describes the linear problem to be solved by the iterative solver.
Declaration of basic traits for the multivector type.
Interface for multivectors used by Belos' linear solvers.
Class which defines basic traits for the operator type.
Alternative run-time polymorphic interface for operators.
Class which manages the output and verbosity of the Belos solvers.
Belos::StatusTest for logically combining several status tests.
Belos::StatusTestResNorm for specifying general residual norm stopping criteria.
Belos::StatusTest for specifying an implicit residual norm stopping criteria that checks for loss of ...
Belos::StatusTest class for specifying a maximum number of iterations.
A factory class for generating StatusTestOutput objects.
Collection of types and exceptions used within the Belos solvers.
Parent class to all Belos exceptions.
This class implements the block GMRES iteration, where a block Krylov subspace is constructed....
An implementation of the Belos::MatOrthoManager that performs orthogonalization using (potentially) m...
GmresIterationOrthoFailure is thrown when the GmresIteration object is unable to compute independent ...
int GetNumberVecs() const
The number of vectors (i.e., columns) in the multivector.
void MvDot(const MultiVec< ScalarType > &A, std::vector< ScalarType > &b) const
Compute the dot product of each column of *this with the corresponding column of A.
Teuchos::RCP< const MV > getConstMV() const
void MvScale(const ScalarType alpha)
Scale each element of the vectors in *this with alpha.
void MvAddMv(const ScalarType alpha, const MultiVec< ScalarType > &A, const ScalarType beta, const MultiVec< ScalarType > &B)
Replace *this with alpha * A + beta * B.
void MvPrint(std::ostream &os) const
Print *this multivector to the os output stream.
void MvScale(const std::vector< ScalarType > &alpha)
Scale each element of the i-th vector in *this with alpha[i].
void MvRandom()
Fill all the vectors in *this with random numbers.
GmresPolyMv * CloneCopy() const
Create a new MultiVec and copy contents of *this into it (deep copy).
void SetBlock(const MultiVec< ScalarType > &A, const std::vector< int > &index)
Copy the vectors in A to a set of vectors in *this.
GmresPolyMv(const Teuchos::RCP< const MV > &mv_in)
Teuchos::RCP< MV > getMV()
ptrdiff_t GetGlobalLength() const
The number of rows in the multivector.
GmresPolyMv(const Teuchos::RCP< MV > &mv_in)
void MvNorm(std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > &normvec, NormType type=TwoNorm) const
Compute the norm of each vector in *this.
void MvTransMv(const ScalarType alpha, const MultiVec< ScalarType > &A, Teuchos::SerialDenseMatrix< int, ScalarType > &B) const
Compute a dense matrix B through the matrix-matrix multiply alpha * A^T * (*this).
const GmresPolyMv * CloneView(const std::vector< int > &index) const
Creates a new Belos::MultiVec that shares the selected contents of *this. The index of the numvecs ve...
void MvTimesMatAddMv(const ScalarType alpha, const MultiVec< ScalarType > &A, const Teuchos::SerialDenseMatrix< int, ScalarType > &B, const ScalarType beta)
Update *this with alpha * A * B + beta * (*this).
GmresPolyMv * CloneViewNonConst(const std::vector< int > &index)
Creates a new Belos::MultiVec that shares the selected contents of *this. The index of the numvecs ve...
GmresPolyMv * Clone(const int numvecs) const
Create a new MultiVec with numvecs columns.
void MvInit(const ScalarType alpha)
Replace each element of the vectors in *this with alpha.
GmresPolyMv * CloneCopy(const std::vector< int > &index) const
Creates a new Belos::MultiVec and copies the selected contents of *this into the new multivector (dee...
GmresPolyOpOrthoFailure is thrown when the orthogonalization manager is unable to generate orthonorma...
GmresPolyOpOrthoFailure(const std::string &what_arg)
Belos's class for applying the GMRES polynomial operator that is used by the hybrid-GMRES linear solv...
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > ¶ms_in)
Process the passed in parameters.
void ApplyRootsPoly(const MV &x, MV &y) const
void generateGmresPoly()
This routine takes the matrix, preconditioner, and vectors from the linear problem as well as the par...
GmresPolyOp(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem_in)
Given no ParameterList, constructor creates no polynomial and only applies the given operator.
void ApplyPoly(const MV &x, MV &y) const
This routine takes the MV x and applies the polynomial operator phi(OP) to it resulting in the MV y,...
void ApplyGmresPoly(const MV &x, MV &y) const
void generateArnoldiPoly()
This routine takes the matrix, preconditioner, and vectors from the linear problem as well as the par...
void ApplyArnoldiPoly(const MV &x, MV &y) const
virtual ~GmresPolyOp()
Destructor.
void Apply(const MultiVec< ScalarType > &x, MultiVec< ScalarType > &y, ETrans=NOTRANS) const
This routine casts the MultiVec to GmresPolyMv to retrieve the MV. Then the above apply method is cal...
GmresPolyOp(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem_in, const Teuchos::RCP< Teuchos::ParameterList > ¶ms_in)
Basic contstructor.
An implementation of the Belos::MatOrthoManager that performs orthogonalization using multiple steps ...
An implementation of the Belos::MatOrthoManager that performs orthogonalization using multiple steps ...
A linear system to solve, and its associated information.
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 MvRandom(MV &mv)
Replace the vectors in mv with random vectors.
static void MvScale(MV &mv, const ScalarType alpha)
Scale each element of the vectors in mv with alpha.
static Teuchos::RCP< MV > Clone(const MV &mv, const int numvecs)
Creates a new empty MV containing numvecs columns.
static void MvTimesMatAddMv(const ScalarType alpha, const MV &A, const Teuchos::SerialDenseMatrix< int, ScalarType > &B, const ScalarType beta, MV &mv)
Update mv with .
static ptrdiff_t GetGlobalLength(const MV &mv)
Return the number of rows in the given multivector mv.
static void MvPrint(const MV &mv, std::ostream &os)
Print the mv multi-vector to the os output stream.
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 void MvDot(const MV &mv, const MV &A, std::vector< ScalarType > &b)
Compute a vector b where the components are the individual dot-products of the i-th columns of A and ...
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 MvInit(MV &mv, const ScalarType alpha=Teuchos::ScalarTraits< ScalarType >::zero())
Replace each element of the vectors in mv with alpha.
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 Teuchos::RCP< MV > CloneCopy(const MV &mv)
Creates a new MV and copies contents of mv into the new vector (deep copy).
static int GetNumberVecs(const MV &mv)
Obtain the number of vectors in mv.
Interface for multivectors used by Belos' linear solvers.
Alternative run-time polymorphic interface for operators.
Belos's basic output manager for sending information of select verbosity levels to the appropriate ou...
A class for extending the status testing capabilities of Belos via logical combinations.
An implementation of StatusTestResNorm using a family of residual norms.
A Belos::StatusTest class for specifying a maximum number of iterations.
ScaleType convertStringToScaleType(const std::string &scaleType)
Convert the given string to its ScaleType enum value.
NormType
The type of vector norm to compute.
ETrans
Whether to apply the (conjugate) transpose of an operator.
static const double polyTol
Relative residual tolerance for matrix polynomial construction.
Structure to contain pointers to GmresIteration state variables.
Teuchos::RCP< const MV > V
The current Krylov basis.
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > z
The current right-hand side of the least squares system RY = Z.
int curDim
The current dimension of the reduction.
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > H
The current Hessenberg matrix.
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > R
The current upper-triangular matrix from the QR reduction of H.