10#ifndef BELOS_GMRESPOLYOP_HPP
11#define BELOS_GMRESPOLYOP_HPP
38#include "Teuchos_BLAS.hpp"
39#include "Teuchos_LAPACK.hpp"
40#include "Teuchos_as.hpp"
41#include "Teuchos_RCP.hpp"
42#include "Teuchos_SerialDenseMatrix.hpp"
43#include "Teuchos_SerialDenseVector.hpp"
44#include "Teuchos_SerialDenseSolver.hpp"
45#include "Teuchos_ParameterList.hpp"
47#ifdef BELOS_TEUCHOS_TIME_MONITOR
48 #include "Teuchos_TimeMonitor.hpp"
64 template <
class ScalarType,
class MV>
74 mv_ = Teuchos::rcp_const_cast<MV>( mv_in );
76 Teuchos::RCP<MV>
getMV() {
return mv_; }
77 Teuchos::RCP<const MV>
getConstMV()
const {
return mv_; }
108 const Teuchos::SerialDenseMatrix<int,ScalarType>& B,
const ScalarType beta)
131 void MvNorm ( std::vector<
typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>& normvec,
NormType type =
TwoNorm )
const
148 Teuchos::RCP<MV> mv_;
162 template <
class ScalarType,
class MV,
class OP>
171 const Teuchos::RCP<Teuchos::ParameterList>& params_in
173 : problem_(problem_in),
175 LP_(problem_in->getLeftPrec()),
176 RP_(problem_in->getRightPrec())
180 polyUpdateLabel_ = label_ +
": Hybrid Gmres: Vector Update";
181#ifdef BELOS_TEUCHOS_TIME_MONITOR
182 timerPolyUpdate_ = Teuchos::TimeMonitor::getNewCounter(polyUpdateLabel_);
185 if (polyType_ ==
"Arnoldi" || polyType_==
"Roots")
187 else if (polyType_ ==
"Gmres")
190 TEUCHOS_TEST_FOR_EXCEPTION(polyType_!=
"Arnoldi"&&polyType_!=
"Gmres"&&polyType_!=
"Roots",std::invalid_argument,
191 "Belos::GmresPolyOp: \"Polynomial Type\" must be either \"Arnoldi\", \"Gmres\", or \"Roots\".");
196 : problem_(problem_in)
210 void setParameters(
const Teuchos::RCP<Teuchos::ParameterList>& params_in );
236 void ApplyPoly (
const MV& x, MV& y )
const;
255#ifdef BELOS_TEUCHOS_TIME_MONITOR
256 Teuchos::RCP<Teuchos::Time> timerPolyUpdate_;
258 std::string polyUpdateLabel_;
262 typedef Teuchos::ScalarTraits<ScalarType> SCT ;
263 typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType MagnitudeType;
264 typedef Teuchos::ScalarTraits<MagnitudeType> MCT ;
267 static constexpr int maxDegree_default_ = 25;
269 static constexpr bool randomRHS_default_ =
true;
270 static constexpr const char * label_default_ =
"Belos";
271 static constexpr const char * polyType_default_ =
"Roots";
272 static constexpr const char * orthoType_default_ =
"DGKS";
273 static constexpr bool damp_default_ =
false;
274 static constexpr bool addRoots_default_ =
true;
277 Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > problem_;
278 Teuchos::RCP<Teuchos::ParameterList> params_;
279 Teuchos::RCP<const OP> LP_, RP_;
282 Teuchos::RCP<OutputManager<ScalarType> > printer_;
283 Teuchos::RCP<std::ostream> outputStream_ = Teuchos::rcpFromRef(std::cout);
286 Teuchos::RCP<MatOrthoManager<ScalarType,MV,OP> > ortho_;
290 int maxDegree_ = maxDegree_default_;
291 int verbosity_ = verbosity_default_;
292 bool randomRHS_ = randomRHS_default_;
293 std::string label_ = label_default_;
294 std::string polyType_ = polyType_default_;
295 std::string orthoType_ = orthoType_default_;
297 bool damp_ = damp_default_;
298 bool addRoots_ = addRoots_default_;
301 mutable Teuchos::RCP<MV> V_, wL_, wR_;
302 Teuchos::SerialDenseMatrix<OT,ScalarType> H_, y_;
303 Teuchos::SerialDenseVector<OT,ScalarType> r0_;
306 bool autoDeg =
false;
307 Teuchos::SerialDenseMatrix< OT, ScalarType > pCoeff_;
310 Teuchos::SerialDenseMatrix< OT, MagnitudeType > theta_;
314 void SortModLeja(Teuchos::SerialDenseMatrix< OT, MagnitudeType > &thetaN, std::vector<int> &index)
const ;
317 void ComputeAddedRoots();
320 template <
class ScalarType,
class MV,
class OP>
324 if (params_in->isParameter(
"Polynomial Type")) {
325 polyType_ = params_in->get(
"Polynomial Type", polyType_default_);
329 if (params_in->isParameter(
"Polynomial Tolerance")) {
330 if (params_in->isType<MagnitudeType> (
"Polynomial Tolerance")) {
331 polyTol_ = params_in->get (
"Polynomial Tolerance",
340 if (params_in->isParameter(
"Maximum Degree")) {
341 maxDegree_ = params_in->get(
"Maximum Degree", maxDegree_default_);
345 if (params_in->isParameter(
"Random RHS")) {
346 randomRHS_ = params_in->get(
"Random RHS", randomRHS_default_);
350 if (params_in->isParameter(
"Verbosity")) {
351 if (Teuchos::isParameterType<int>(*params_in,
"Verbosity")) {
352 verbosity_ = params_in->get(
"Verbosity", verbosity_default_);
355 verbosity_ = (int)Teuchos::getParameter<Belos::MsgType>(*params_in,
"Verbosity");
359 if (params_in->isParameter(
"Orthogonalization")) {
360 orthoType_ = params_in->get(
"Orthogonalization",orthoType_default_);
364 if (params_in->isParameter(
"Timer Label")) {
365 label_ = params_in->get(
"Timer Label", label_default_);
369 if (params_in->isParameter(
"Output Stream")) {
370 outputStream_ = Teuchos::getParameter<Teuchos::RCP<std::ostream> >(*params_in,
"Output Stream");
374 if (params_in->isParameter(
"Damped Poly")) {
375 damp_ = params_in->get(
"Damped Poly", damp_default_);
379 if (params_in->isParameter(
"Add Roots")) {
380 addRoots_ = params_in->get(
"Add Roots", addRoots_default_);
384 template <
class ScalarType,
class MV,
class OP>
387 Teuchos::RCP< MV > V =
MVT::Clone( *problem_->getRHS(), maxDegree_+1 );
390 std::vector<int> index(1,0);
397 if ( !LP_.is_null() ) {
399 problem_->applyLeftPrec( *Vtemp, *V0);
403 problem_->apply( *Vtemp, *V0);
406 for(
int i=0; i< maxDegree_; i++)
412 problem_->apply( *Vi, *Vip1);
416 Teuchos::Range1D range( 1, maxDegree_);
420 Teuchos::SerialDenseMatrix< OT, ScalarType > AVtransAV( maxDegree_, maxDegree_);
424 Teuchos::LAPACK< OT, ScalarType > lapack;
429 Teuchos::SerialDenseMatrix< OT, ScalarType > lhs;
430 while( status && dim_ >= 1)
432 Teuchos::SerialDenseMatrix< OT, ScalarType > lhstemp(Teuchos::Copy, AVtransAV, dim_, dim_);
433 lapack.POTRF(
'U', dim_, lhstemp.values(), lhstemp.stride(), &infoInt);
440 std::cout <<
"BelosGmresPolyOp.hpp: LAPACK POTRF was not successful!!" << std::endl;
441 std::cout <<
"Error code: " << infoInt << std::endl;
462 pCoeff_.shape( 1, 1);
463 pCoeff_(0,0) = SCT::one();
464 std::cout <<
"Poly Degree is zero. No preconditioner created." << std::endl;
468 pCoeff_.shape( dim_, 1);
470 Teuchos::Range1D rangeSub( 1, dim_);
475 lapack.POTRS(
'U', dim_, 1, lhs.values(), lhs.stride(), pCoeff_.values(), pCoeff_.stride(), &infoInt);
478 std::cout <<
"BelosGmresPolyOp.hpp: LAPACK POTRS was not successful!!" << std::endl;
479 std::cout <<
"Error code: " << infoInt << std::endl;
484 template <
class ScalarType,
class MV,
class OP>
487 std::string polyLabel = label_ +
": GmresPolyOp creation";
490 std::vector<int> idx(1,0);
491 Teuchos::RCP<MV> newX =
MVT::Clone( *(problem_->getLHS()), 1 );
492 Teuchos::RCP<MV> newB =
MVT::Clone( *(problem_->getRHS()), 1 );
500 Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > newProblem =
502 newProblem->setInitResVec( newB );
503 newProblem->setLeftPrec( problem_->getLeftPrec() );
504 newProblem->setRightPrec( problem_->getRightPrec() );
505 newProblem->setLabel(polyLabel);
506 newProblem->setProblem();
507 newProblem->setLSIndex( idx );
510 Teuchos::ParameterList polyList;
513 polyList.set(
"Num Blocks",maxDegree_);
514 polyList.set(
"Block Size",1);
515 polyList.set(
"Keep Hessenberg",
true);
521 if (ortho_.is_null()) {
522 params_->set(
"Orthogonalization", orthoType_);
524 Teuchos::RCP<Teuchos::ParameterList> paramsOrtho;
526 ortho_ = factory.
makeMatOrthoManager (orthoType_, Teuchos::null, printer_, polyLabel, paramsOrtho);
530 Teuchos::RCP<StatusTestMaxIters<ScalarType,MV,OP> > maxItrTst =
534 Teuchos::RCP<StatusTestGenResNorm<ScalarType,MV,OP> > convTst =
539 Teuchos::RCP<StatusTestCombo<ScalarType,MV,OP> > polyTest =
543 Teuchos::RCP<BlockGmresIter<ScalarType,MV,OP> > gmres_iter;
548 if ( !LP_.is_null() )
549 newProblem->applyLeftPrec( *newB, *V_0 );
553 newProblem->apply( *Vtemp, *V_0 );
560 int rank = ortho_->normalize( *V_0, Teuchos::rcpFromRef(r0_) );
562 "Belos::GmresPolyOp::generateArnoldiPoly(): Failed to compute initial block of orthonormal vectors for polynomial generation.");
567 newstate.
z = Teuchos::rcpFromRef( r0_);
569 gmres_iter->initializeGmres(newstate);
573 gmres_iter->iterate();
577 gmres_iter->updateLSQR( gmres_iter->getCurSubspaceDim() );
579 catch (std::exception& e) {
581 printer_->stream(
Errors) <<
"Error! Caught exception in BlockGmresIter::iterate() at iteration "
582 << gmres_iter->getNumIters() << endl << e.what () << endl;
587 Teuchos::RCP<MV> currX = gmres_iter->getCurrentUpdate();
597 if(polyType_ ==
"Arnoldi"){
600 y_ = Teuchos::SerialDenseMatrix<OT,ScalarType>( Teuchos::Copy, *gmresState.
z, dim_, 1 );
606 Teuchos::BLAS<OT,ScalarType> blas;
607 blas.TRSM( Teuchos::LEFT_SIDE, Teuchos::UPPER_TRI, Teuchos::NO_TRANS,
608 Teuchos::NON_UNIT_DIAG, dim_, 1, SCT::one(),
609 gmresState.
R->values(), gmresState.
R->stride(),
610 y_.values(), y_.stride() );
616 H_ = Teuchos::SerialDenseMatrix<OT,ScalarType>(Teuchos::Copy, *gmresState.
H, dim_, dim_);
618 for(
int i=0; i <= dim_-3; i++) {
619 for(
int k=i+2; k <= dim_-1; k++) {
620 H_(k,i) = SCT::zero();
624 Teuchos::SerialDenseMatrix<OT,ScalarType> Htemp (Teuchos::Copy, H_, dim_, dim_);
627 ScalarType Hlast = (*gmresState.
H)(dim_,dim_-1);
628 Teuchos::SerialDenseMatrix<OT,ScalarType> HlastCol (Teuchos::View, H_, dim_, 1, 0, dim_-1);
631 Teuchos::SerialDenseMatrix< OT, ScalarType > F(dim_,1), E(dim_,1);
632 E.putScalar(SCT::zero());
633 E(dim_-1,0) = SCT::one();
635 Teuchos::SerialDenseSolver< OT, ScalarType > HSolver;
636 HSolver.setMatrix( Teuchos::rcpFromRef(Htemp));
637 HSolver.solveWithTransposeFlag( Teuchos::CONJ_TRANS );
638 HSolver.setVectors( Teuchos::rcpFromRef(F), Teuchos::rcpFromRef(E));
639 HSolver.factorWithEquilibration(
true );
643 info = HSolver.factor();
645 std::cout <<
"Hsolver factor: info = " << info << std::endl;
647 info = HSolver.solve();
649 std::cout <<
"Hsolver solve : info = " << info << std::endl;
653 F.scale(Hlast*Hlast);
657 Teuchos::LAPACK< OT, ScalarType > lapack;
658 theta_.shape(dim_,2);
665 std::vector<ScalarType> work(1);
666 std::vector<MagnitudeType> rwork(2*dim_);
669 lapack.GEEV(
'N',
'N',dim_,H_.values(),H_.stride(),theta_[0],theta_[1],vlr, ldv, vlr, ldv, &work[0], lwork, &rwork[0], &info);
670 lwork = std::abs (
static_cast<int> (Teuchos::ScalarTraits<ScalarType>::real (work[0])));
671 work.resize( lwork );
673 lapack.GEEV(
'N',
'N',dim_,H_.values(),H_.stride(),theta_[0],theta_[1],vlr, ldv, vlr, ldv, &work[0], lwork, &rwork[0], &info);
676 std::cout <<
"GEEV solve : info = " << info << std::endl;
681 const MagnitudeType tol = 10.0 * Teuchos::ScalarTraits<MagnitudeType>::eps();
682 std::vector<int> index(dim_);
683 for(
int i=0; i<dim_; ++i){
686 TEUCHOS_TEST_FOR_EXCEPTION(hypot(theta_(i,0),theta_(i,1)) < tol, std::runtime_error,
"BelosGmresPolyOp Error: One of the computed polynomial roots is approximately zero. This will cause a divide by zero error! Your matrix may be close to singular. Please select a lower polynomial degree or give a shifted matrix.");
688 SortModLeja(theta_,index);
697 template <
class ScalarType,
class MV,
class OP>
698 void GmresPolyOp<ScalarType, MV, OP>::ComputeAddedRoots()
702 std::vector<std::complex<MagnitudeType>> cmplxHRitz (dim_);
703 for(
unsigned int i=0; i<cmplxHRitz.size(); ++i){
704 cmplxHRitz[i] = std::complex<MagnitudeType>( theta_(i,0), theta_(i,1) );
708 const MagnitudeType one(1.0);
709 std::vector<MagnitudeType> pof (dim_,one);
710 for(
int j=0; j<dim_; ++j) {
711 for(
int i=0; i<dim_; ++i) {
713 pof[j] = std::abs(pof[j]*(one-(cmplxHRitz[j]/cmplxHRitz[i])));
719 std::vector<int> extra (dim_);
721 for(
int i=0; i<dim_; ++i){
722 if (pof[i] > MCT::zero())
723 extra[i] = ceil((log10(pof[i])-MagnitudeType(4.0))/MagnitudeType(14.0));
727 totalExtra += extra[i];
731 printer_->stream(
Warnings) <<
"Warning: Need to add " << totalExtra <<
" extra roots." << std::endl;}
734 if(addRoots_ && totalExtra>0)
736 theta_.reshape(dim_+totalExtra,2);
738 Teuchos::SerialDenseMatrix<OT,MagnitudeType> thetaPert (Teuchos::Copy, theta_, dim_+totalExtra, 2);
742 for(
int i=0; i<dim_; ++i){
743 for(
int j=0; j< extra[i]; ++j){
744 theta_(count,0) = theta_(i,0);
745 theta_(count,1) = theta_(i,1);
746 thetaPert(count,0) = theta_(i,0)+(j+MCT::one())*MagnitudeType(5e-8);
747 thetaPert(count,1) = theta_(i,1);
755 printer_->stream(
Warnings) <<
"New poly degree is: " << dim_ << std::endl;}
758 std::vector<int> index2(dim_);
759 for(
int i=0; i<dim_; ++i){
762 SortModLeja(thetaPert,index2);
764 for(
int i=0; i<dim_; ++i)
766 thetaPert(i,0) = theta_(index2[i],0);
767 thetaPert(i,1) = theta_(index2[i],1);
776 template <
class ScalarType,
class MV,
class OP>
777 void GmresPolyOp<ScalarType, MV, OP>::SortModLeja(Teuchos::SerialDenseMatrix< OT, MagnitudeType > &thetaN, std::vector<int> &index)
const
782 int dimN = index.size();
783 std::vector<int> newIndex(dimN);
784 Teuchos::SerialDenseMatrix< OT, MagnitudeType > sorted (thetaN.numRows(), thetaN.numCols());
785 Teuchos::SerialDenseVector< OT, MagnitudeType > absVal (thetaN.numRows());
786 Teuchos::SerialDenseVector< OT, MagnitudeType > prod (thetaN.numRows());
789 for(
int i = 0; i < dimN; i++){
790 absVal(i) = hypot(thetaN(i,0), thetaN(i,1));
792 MagnitudeType * maxPointer = std::max_element(absVal.values(), (absVal.values()+dimN));
793 int maxIndex = int (maxPointer- absVal.values());
796 sorted(0,0) = thetaN(maxIndex,0);
797 sorted(0,1) = thetaN(maxIndex,1);
798 newIndex[0] = index[maxIndex];
802 if(sorted(0,1)!= SCT::zero() && !SCT::isComplex)
804 sorted(1,0) = thetaN(maxIndex,0);
805 sorted(1,1) = -thetaN(maxIndex,1);
806 newIndex[1] = index[maxIndex+1];
819 for(
int i = 0; i < dimN; i++)
821 prod(i) = MCT::one();
822 for(
int k = 0; k < j; k++)
824 a = thetaN(i,0) - sorted(k,0);
825 b = thetaN(i,1) - sorted(k,1);
826 if (a*a + b*b > MCT::zero())
827 prod(i) = prod(i) + log10(hypot(a,b));
829 prod(i) = -std::numeric_limits<MagnitudeType>::infinity();
836 maxPointer = std::max_element(prod.values(), (prod.values()+dimN));
837 maxIndex = int (maxPointer- prod.values());
838 sorted(j,0) = thetaN(maxIndex,0);
839 sorted(j,1) = thetaN(maxIndex,1);
840 newIndex[j] = index[maxIndex];
843 if(sorted(j,1)!= SCT::zero() && !SCT::isComplex)
846 sorted(j,0) = thetaN(maxIndex,0);
847 sorted(j,1) = -thetaN(maxIndex,1);
848 newIndex[j] = index[maxIndex+1];
858 template <
class ScalarType,
class MV,
class OP>
862 if (polyType_ ==
"Arnoldi")
864 else if (polyType_ ==
"Gmres")
866 else if (polyType_ ==
"Roots")
871 problem_->applyOp( x, y );
875 template <
class ScalarType,
class MV,
class OP>
882 if (!LP_.is_null()) {
884 problem_->applyLeftPrec( *AX, *Xtmp );
889#ifdef BELOS_TEUCHOS_TIME_MONITOR
890 Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
894 for(
int i=1; i < dim_; i++)
896 Teuchos::RCP<MV> X, Y;
907 problem_->apply(*X, *Y);
909#ifdef BELOS_TEUCHOS_TIME_MONITOR
910 Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
917 if (!RP_.is_null()) {
919 problem_->applyRightPrec( *Ytmp, y );
923 template <
class ScalarType,
class MV,
class OP>
932 if (!LP_.is_null()) {
933 problem_->applyLeftPrec( *prod, *Xtmp );
940 if(theta_(i,1)== SCT::zero() || SCT::isComplex)
943#ifdef BELOS_TEUCHOS_TIME_MONITOR
944 Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
946 MVT::MvAddMv(SCT::one(), y, SCT::one()/theta_(i,0), *prod, y);
948 problem_->apply(*prod, *Xtmp);
950#ifdef BELOS_TEUCHOS_TIME_MONITOR
951 Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
953 MVT::MvAddMv(SCT::one(), *prod, -SCT::one()/theta_(i,0), *Xtmp, *prod);
959 MagnitudeType mod = theta_(i,0)*theta_(i,0) + theta_(i,1)*theta_(i,1);
960 problem_->apply(*prod, *Xtmp);
962#ifdef BELOS_TEUCHOS_TIME_MONITOR
963 Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
965 MVT::MvAddMv(2*theta_(i,0), *prod, -SCT::one(), *Xtmp, *Xtmp);
970 problem_->apply(*Xtmp, *Xtmp2);
972#ifdef BELOS_TEUCHOS_TIME_MONITOR
973 Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
975 MVT::MvAddMv(SCT::one(), *prod, -SCT::one()/mod, *Xtmp2, *prod);
981 if(theta_(dim_-1,1)== SCT::zero() || SCT::isComplex)
983#ifdef BELOS_TEUCHOS_TIME_MONITOR
984 Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
986 MVT::MvAddMv(SCT::one(), y, SCT::one()/theta_(dim_-1,0), *prod, y);
990 if (!RP_.is_null()) {
992 problem_->applyRightPrec( *Ytmp, y );
996 template <
class ScalarType,
class MV,
class OP>
1002 if (!LP_.is_null()) {
1005 if (!RP_.is_null()) {
1013 std::vector<int> idxi(1), idxi2, idxj(1);
1016 for (
int j=0; j<n; ++j) {
1022 if (!LP_.is_null()) {
1024 problem_->applyLeftPrec( *x_view, *v_curr );
1029 for (
int i=0; i<dim_-1; ++i) {
1033 for (
int ii=0; ii<i+1; ++ii) { idxi2[ii] = ii; }
1045 if (!RP_.is_null()) {
1046 problem_->applyRightPrec( *v_curr, *wR_ );
1051 if (LP_.is_null()) {
1055 problem_->applyOp( *wR_, *wL_ );
1057 if (!LP_.is_null()) {
1058 problem_->applyLeftPrec( *wL_, *v_next );
1062 Teuchos::SerialDenseMatrix<OT,ScalarType> h(Teuchos::View,H_,i+1,1,0,i);
1064#ifdef BELOS_TEUCHOS_TIME_MONITOR
1065 Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
1075 if (!RP_.is_null()) {
1077#ifdef BELOS_TEUCHOS_TIME_MONITOR
1078 Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
1082 problem_->applyRightPrec( *wR_, *y_view );
1085#ifdef BELOS_TEUCHOS_TIME_MONITOR
1086 Teuchos::TimeMonitor updateTimer( *timerPolyUpdate_ );
Belos concrete class for performing the block GMRES iteration.
Belos header file which uses auto-configuration information to include necessary C++ headers.
Pure virtual base class which augments the basic interface for a Gmres linear solver iteration.
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.
BelosError(const std::string &what_arg)
This class implements the block GMRES iteration, where a block Krylov subspace is constructed....
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.
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.
Teuchos::RCP< const MV > getConstMV() const
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.
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)
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...
ptrdiff_t GetGlobalLength() const
The number of rows in the multivector.
GmresPolyMv(const Teuchos::RCP< MV > &mv_in)
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 MvNorm(std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > &normvec, NormType type=TwoNorm) const
Compute the norm of each vector in *this.
Teuchos::RCP< MV > getMV()
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).
GmresPolyMv * Clone(const int numvecs) const
Create a new MultiVec with numvecs columns.
GmresPolyMv * CloneCopy() const
Create a new MultiVec and copy contents of *this into it (deep copy).
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 * CloneCopy(const std::vector< int > &index) const
Creates a new Belos::MultiVec and copies the selected contents of *this into the new multivector (dee...
void MvInit(const ScalarType alpha)
Replace each element of the vectors in *this with alpha.
GmresPolyOpOrthoFailure is thrown when the orthogonalization manager is unable to generate orthonorma...
GmresPolyOpOrthoFailure(const std::string &what_arg)
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.
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 void MvTimesMatAddMv(const ScalarType alpha, const MV &A, const Teuchos::SerialDenseMatrix< int, ScalarType > &B, const ScalarType beta, MV &mv)
Update mv with .
static Teuchos::RCP< MV > CloneCopy(const MV &mv)
Creates a new MV and copies contents of mv into the new vector (deep copy).
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< MV > Clone(const MV &mv, const int numvecs)
Creates a new empty MV containing numvecs columns.
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 int GetNumberVecs(const MV &mv)
Obtain the number of vectors in mv.
static void Assign(const MV &A, MV &mv)
mv := A
Interface for multivectors used by Belos' linear solvers.
MultiVec()
Default constructor.
Operator()
Default constructor (does nothing).
Enumeration of all valid Belos (Mat)OrthoManager classes.
Teuchos::RCP< Belos::MatOrthoManager< Scalar, MV, OP > > makeMatOrthoManager(const std::string &ortho, const Teuchos::RCP< const OP > &M, const Teuchos::RCP< OutputManager< Scalar > > &, const std::string &label, const Teuchos::RCP< Teuchos::ParameterList > ¶ms)
Return an instance of the specified MatOrthoManager subclass.
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.