43 #ifndef BELOS_LINEAR_PROBLEM_HPP 44 #define BELOS_LINEAR_PROBLEM_HPP 51 #include "Teuchos_ParameterList.hpp" 52 #include "Teuchos_TimeMonitor.hpp" 81 template <
class ScalarType,
class MV,
class OP>
104 const Teuchos::RCP<MV> &X,
105 const Teuchos::RCP<const MV> &B);
134 void setLHS (
const Teuchos::RCP<MV> &X) {
143 void setRHS (
const Teuchos::RCP<const MV> &B) {
193 void setLSIndex (
const std::vector<int>& index);
213 void setLabel (
const std::string& label);
252 updateSolution (
const Teuchos::RCP<MV>& update = Teuchos::null,
253 bool updateLP =
false,
254 ScalarType scale = Teuchos::ScalarTraits<ScalarType>::one());
274 ScalarType scale = Teuchos::ScalarTraits<ScalarType>::one() )
const 308 setProblem (
const Teuchos::RCP<MV> &newX = Teuchos::null,
309 const Teuchos::RCP<const MV> &newB = Teuchos::null);
320 Teuchos::RCP<MV>
getLHS()
const {
return(X_); }
323 Teuchos::RCP<const MV>
getRHS()
const {
return(B_); }
326 Teuchos::RCP<const MV> getInitResVec()
const;
332 Teuchos::RCP<const MV> getInitPrecResVec()
const;
348 Teuchos::RCP<MV> getCurrLHSVec();
364 Teuchos::RCP<const MV> getCurrRHSVec();
393 const std::vector<int>
getLSIndex()
const {
return(rhsIndex_); }
409 Teuchos::Array<Teuchos::RCP<Teuchos::Time> >
getTimers()
const {
410 return Teuchos::tuple(timerOp_,timerPrec_);
457 void apply(
const MV& x, MV& y )
const;
466 void applyOp(
const MV& x, MV& y )
const;
474 void applyLeftPrec(
const MV& x, MV& y )
const;
482 void applyRightPrec(
const MV& x, MV& y )
const;
489 void computeCurrResVec( MV* R ,
const MV* X = 0,
const MV* B = 0 )
const;
496 void computeCurrPrecResVec( MV* R,
const MV* X = 0,
const MV* B = 0 )
const;
503 Teuchos::RCP<const OP> A_;
509 Teuchos::RCP<MV> curX_;
512 Teuchos::RCP<const MV> B_;
515 Teuchos::RCP<const MV> curB_;
518 Teuchos::RCP<MV> R0_;
521 Teuchos::RCP<MV> PR0_;
524 Teuchos::RCP<const MV> R0_user_;
527 Teuchos::RCP<const MV> PR0_user_;
530 Teuchos::RCP<const OP> LP_;
533 Teuchos::RCP<const OP> RP_;
536 mutable Teuchos::RCP<Teuchos::Time> timerOp_, timerPrec_;
545 std::vector<int> rhsIndex_;
567 bool solutionUpdated_;
582 template <
class ScalarType,
class MV,
class OP>
592 solutionUpdated_(false),
597 template <
class ScalarType,
class MV,
class OP>
599 const Teuchos::RCP<MV> &X,
600 const Teuchos::RCP<const MV> &B
613 solutionUpdated_(false),
618 template <
class ScalarType,
class MV,
class OP>
622 curX_(Problem.curX_),
624 curB_(Problem.curB_),
627 R0_user_(Problem.R0_user_),
628 PR0_user_(Problem.PR0_user_),
631 timerOp_(Problem.timerOp_),
632 timerPrec_(Problem.timerPrec_),
633 blocksize_(Problem.blocksize_),
634 num2Solve_(Problem.num2Solve_),
635 rhsIndex_(Problem.rhsIndex_),
636 lsNum_(Problem.lsNum_),
637 Left_Scale_(Problem.Left_Scale_),
638 Right_Scale_(Problem.Right_Scale_),
639 isSet_(Problem.isSet_),
640 isHermitian_(Problem.isHermitian_),
641 solutionUpdated_(Problem.solutionUpdated_),
642 label_(Problem.label_)
646 template <
class ScalarType,
class MV,
class OP>
650 template <
class ScalarType,
class MV,
class OP>
658 curB_ = Teuchos::null;
659 curX_ = Teuchos::null;
662 int validIdx = 0, ivalidIdx = 0;
663 blocksize_ = rhsIndex_.size();
664 std::vector<int> vldIndex( blocksize_ );
665 std::vector<int> newIndex( blocksize_ );
666 std::vector<int> iIndex( blocksize_ );
667 for (
int i=0; i<blocksize_; ++i) {
668 if (rhsIndex_[i] > -1) {
669 vldIndex[validIdx] = rhsIndex_[i];
670 newIndex[validIdx] = i;
674 iIndex[ivalidIdx] = i;
678 vldIndex.resize(validIdx);
679 newIndex.resize(validIdx);
680 iIndex.resize(ivalidIdx);
681 num2Solve_ = validIdx;
684 if (num2Solve_ != blocksize_) {
685 newIndex.resize(num2Solve_);
686 vldIndex.resize(num2Solve_);
692 Teuchos::RCP<MV> tmpCurB =
MVT::Clone( *B_, blocksize_ );
704 solutionUpdated_ =
false;
717 template <
class ScalarType,
class MV,
class OP>
724 if (num2Solve_ < blocksize_) {
729 std::vector<int> newIndex( num2Solve_ );
730 std::vector<int> vldIndex( num2Solve_ );
731 for (
int i=0; i<blocksize_; ++i) {
732 if ( rhsIndex_[i] > -1 ) {
733 vldIndex[validIdx] = rhsIndex_[i];
734 newIndex[validIdx] = i;
745 curX_ = Teuchos::null;
746 curB_ = Teuchos::null;
751 template <
class ScalarType,
class MV,
class OP>
762 if (update.is_null())
780 RCP<MV> rightPrecUpdate =
783 #ifdef BELOS_TEUCHOS_TIME_MONITOR 784 Teuchos::TimeMonitor PrecTimer (*timerPrec_);
786 OPT::Apply( *RP_, *update, *rightPrecUpdate );
789 MVT::MvAddMv( 1.0, *curX_, scale, *rightPrecUpdate, *curX_ );
791 solutionUpdated_ =
true;
807 RCP<MV> rightPrecUpdate =
810 #ifdef BELOS_TEUCHOS_TIME_MONITOR 811 Teuchos::TimeMonitor PrecTimer(*timerPrec_);
813 OPT::Apply( *RP_, *update, *rightPrecUpdate );
816 MVT::MvAddMv( 1.0, *curX_, scale, *rightPrecUpdate, *newSoln );
823 template <
class ScalarType,
class MV,
class OP>
826 if (label != label_) {
829 if (timerOp_ != Teuchos::null) {
830 std::string opLabel = label_ +
": Operation Op*x";
831 #ifdef BELOS_TEUCHOS_TIME_MONITOR 832 timerOp_ = Teuchos::TimeMonitor::getNewCounter( opLabel );
835 if (timerPrec_ != Teuchos::null) {
836 std::string precLabel = label_ +
": Operation Prec*x";
837 #ifdef BELOS_TEUCHOS_TIME_MONITOR 838 timerPrec_ = Teuchos::TimeMonitor::getNewCounter( precLabel );
844 template <
class ScalarType,
class MV,
class OP>
848 const Teuchos::RCP<const MV> &newB)
851 if (timerOp_ == Teuchos::null) {
852 std::string opLabel = label_ +
": Operation Op*x";
853 #ifdef BELOS_TEUCHOS_TIME_MONITOR 854 timerOp_ = Teuchos::TimeMonitor::getNewCounter( opLabel );
857 if (timerPrec_ == Teuchos::null) {
858 std::string precLabel = label_ +
": Operation Prec*x";
859 #ifdef BELOS_TEUCHOS_TIME_MONITOR 860 timerPrec_ = Teuchos::TimeMonitor::getNewCounter( precLabel );
865 if (newX != Teuchos::null)
867 if (newB != Teuchos::null)
872 curX_ = Teuchos::null;
873 curB_ = Teuchos::null;
877 if (A_ == Teuchos::null || X_ == Teuchos::null || B_ == Teuchos::null) {
885 solutionUpdated_ =
false;
888 if(Teuchos::is_null(R0_user_)) {
894 if (LP_!=Teuchos::null) {
899 #ifdef BELOS_TEUCHOS_TIME_MONITOR 900 Teuchos::TimeMonitor PrecTimer(*timerPrec_);
916 if (LP_!=Teuchos::null) {
920 #ifdef BELOS_TEUCHOS_TIME_MONITOR 921 Teuchos::TimeMonitor PrecTimer(*timerPrec_);
929 PR0_user_ = R0_user_;
940 template <
class ScalarType,
class MV,
class OP>
943 if(Teuchos::nonnull(R0_user_)) {
949 template <
class ScalarType,
class MV,
class OP>
952 if(Teuchos::nonnull(PR0_user_)) {
958 template <
class ScalarType,
class MV,
class OP>
965 return Teuchos::null;
969 template <
class ScalarType,
class MV,
class OP>
976 return Teuchos::null;
980 template <
class ScalarType,
class MV,
class OP>
986 const bool leftPrec = LP_ != null;
987 const bool rightPrec = RP_ != null;
998 if (!leftPrec && !rightPrec){
999 #ifdef BELOS_TEUCHOS_TIME_MONITOR 1000 Teuchos::TimeMonitor OpTimer(*timerOp_);
1007 else if( leftPrec && rightPrec )
1010 #ifdef BELOS_TEUCHOS_TIME_MONITOR 1011 Teuchos::TimeMonitor PrecTimer(*timerPrec_);
1016 #ifdef BELOS_TEUCHOS_TIME_MONITOR 1017 Teuchos::TimeMonitor OpTimer(*timerOp_);
1022 #ifdef BELOS_TEUCHOS_TIME_MONITOR 1023 Teuchos::TimeMonitor PrecTimer(*timerPrec_);
1034 #ifdef BELOS_TEUCHOS_TIME_MONITOR 1035 Teuchos::TimeMonitor OpTimer(*timerOp_);
1040 #ifdef BELOS_TEUCHOS_TIME_MONITOR 1041 Teuchos::TimeMonitor PrecTimer(*timerPrec_);
1052 #ifdef BELOS_TEUCHOS_TIME_MONITOR 1053 Teuchos::TimeMonitor PrecTimer(*timerPrec_);
1058 #ifdef BELOS_TEUCHOS_TIME_MONITOR 1059 Teuchos::TimeMonitor OpTimer(*timerOp_);
1066 template <
class ScalarType,
class MV,
class OP>
1069 #ifdef BELOS_TEUCHOS_TIME_MONITOR 1070 Teuchos::TimeMonitor OpTimer(*timerOp_);
1075 MVT::MvAddMv( Teuchos::ScalarTraits<ScalarType>::one(), x,
1076 Teuchos::ScalarTraits<ScalarType>::zero(), x, y );
1080 template <
class ScalarType,
class MV,
class OP>
1082 if (LP_!=Teuchos::null) {
1083 #ifdef BELOS_TEUCHOS_TIME_MONITOR 1084 Teuchos::TimeMonitor PrecTimer(*timerPrec_);
1089 MVT::MvAddMv( Teuchos::ScalarTraits<ScalarType>::one(), x,
1090 Teuchos::ScalarTraits<ScalarType>::zero(), x, y );
1094 template <
class ScalarType,
class MV,
class OP>
1096 if (RP_!=Teuchos::null) {
1097 #ifdef BELOS_TEUCHOS_TIME_MONITOR 1098 Teuchos::TimeMonitor PrecTimer(*timerPrec_);
1103 MVT::MvAddMv( Teuchos::ScalarTraits<ScalarType>::one(), x,
1104 Teuchos::ScalarTraits<ScalarType>::zero(), x, y );
1108 template <
class ScalarType,
class MV,
class OP>
1114 if (LP_!=Teuchos::null)
1118 #ifdef BELOS_TEUCHOS_TIME_MONITOR 1119 Teuchos::TimeMonitor OpTimer(*timerOp_);
1125 #ifdef BELOS_TEUCHOS_TIME_MONITOR 1126 Teuchos::TimeMonitor PrecTimer(*timerPrec_);
1134 #ifdef BELOS_TEUCHOS_TIME_MONITOR 1135 Teuchos::TimeMonitor OpTimer(*timerOp_);
1144 Teuchos::RCP<const MV> localB, localX;
1146 localB = Teuchos::rcp( B,
false );
1151 localX = Teuchos::rcp( X,
false );
1155 if (LP_!=Teuchos::null)
1159 #ifdef BELOS_TEUCHOS_TIME_MONITOR 1160 Teuchos::TimeMonitor OpTimer(*timerOp_);
1166 #ifdef BELOS_TEUCHOS_TIME_MONITOR 1167 Teuchos::TimeMonitor PrecTimer(*timerPrec_);
1175 #ifdef BELOS_TEUCHOS_TIME_MONITOR 1176 Teuchos::TimeMonitor OpTimer(*timerOp_);
1187 template <
class ScalarType,
class MV,
class OP>
1194 #ifdef BELOS_TEUCHOS_TIME_MONITOR 1195 Teuchos::TimeMonitor OpTimer(*timerOp_);
1203 Teuchos::RCP<const MV> localB, localX;
1205 localB = Teuchos::rcp( B,
false );
1210 localX = Teuchos::rcp( X,
false );
1215 #ifdef BELOS_TEUCHOS_TIME_MONITOR 1216 Teuchos::TimeMonitor OpTimer(*timerOp_);
Teuchos::RCP< const MV > getRHS() const
A pointer to the right-hand side B.
Exception thrown to signal error with the Belos::LinearProblem object.
void setHermitian(bool isSym=true)
Tell the linear problem that the (preconditioned) operator is Hermitian.
const std::vector< int > getLSIndex() const
(Zero-based) indices of the linear system(s) currently being solved.
Teuchos::RCP< const MV > getCurrRHSVec()
Get a pointer to the current right-hand side of the linear system.
static void MvRandom(MV &mv)
Replace the vectors in mv with random vectors.
Teuchos::RCP< MV > getCurrLHSVec()
Get a pointer to the current left-hand side (solution) of the linear system.
Teuchos::RCP< const MV > getInitResVec() const
A pointer to the initial unpreconditioned residual vector.
Teuchos::RCP< MV > updateSolution(const Teuchos::RCP< MV > &update=Teuchos::null, bool updateLP=false, ScalarType scale=Teuchos::ScalarTraits< ScalarType >::one())
Compute the new solution to the linear system using the given update vector.
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).
bool isRightPrec() const
Whether the linear system is being preconditioned on the right.
Declaration of basic traits for the multivector type.
virtual ~LinearProblem(void)
Destructor (declared virtual for memory safety of derived classes).
bool setProblem(const Teuchos::RCP< MV > &newX=Teuchos::null, const Teuchos::RCP< const MV > &newB=Teuchos::null)
Set up the linear problem manager.
Teuchos::Array< Teuchos::RCP< Teuchos::Time > > getTimers() const
The timers for this object.
bool isHermitian() const
Whether the (preconditioned) operator is Hermitian.
void applyLeftPrec(const MV &x, MV &y) const
Apply ONLY the left preconditioner to x, returning y.
LinearProblem(void)
Default constructor.
static void MvInit(MV &mv, const ScalarType alpha=Teuchos::ScalarTraits< ScalarType >::zero())
Replace each element of the vectors in mv with alpha.
void setOperator(const Teuchos::RCP< const OP > &A)
Set the operator A of the linear problem .
Class which defines basic traits for the operator type.
void setLSIndex(const std::vector< int > &index)
Tell the linear problem which linear system(s) need to be solved next.
static int GetNumberVecs(const MV &mv)
Obtain the number of vectors in mv.
static void MvAddMv(const ScalarType alpha, const MV &A, const ScalarType beta, const MV &B, MV &mv)
Replace mv with .
Teuchos::RCP< MV > updateSolution(const Teuchos::RCP< MV > &update=Teuchos::null, ScalarType scale=Teuchos::ScalarTraits< ScalarType >::one()) const
Compute the new solution to the linear system using the given update vector.
void setRHS(const Teuchos::RCP< const MV > &B)
Set right-hand-side B of linear problem .
Traits class which defines basic operations on multivectors.
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 > CloneViewNonConst(MV &mv, const std::vector< int > &index)
Creates a new MV that shares the selected contents of mv (shallow copy).
static Teuchos::RCP< MV > Clone(const MV &mv, const int numvecs)
Creates a new empty MV containing numvecs columns.
static void Apply(const OP &Op, const MV &x, MV &y, ETrans trans=NOTRANS)
Apply Op to x, putting the result into y.
void setLeftPrec(const Teuchos::RCP< const OP > &LP)
Set left preconditioner (LP) of linear problem .
Teuchos::RCP< MV > getLHS() const
A pointer to the left-hand side X.
A linear system to solve, and its associated information.
void applyRightPrec(const MV &x, MV &y) const
Apply ONLY the right preconditioner to x, returning y.
void applyOp(const MV &x, MV &y) const
Apply ONLY the operator to x, returning y.
Teuchos::RCP< const OP > getRightPrec() const
Get a pointer to the right preconditioner.
void setLabel(const std::string &label)
Set the label prefix used by the timers in this object.
int getLSNumber() const
The number of linear systems that have been set.
Teuchos::RCP< const OP > getLeftPrec() const
Get a pointer to the left preconditioner.
bool isProblemSet() const
Whether the problem has been set.
LinearProblemError(const std::string &what_arg)
bool isSolutionUpdated() const
Has the current approximate solution been updated?
void setCurrLS()
Tell the linear problem that the solver is finished with the current linear system.
bool isLeftPrec() const
Whether the linear system is being preconditioned on the left.
void computeCurrResVec(MV *R, const MV *X=0, const MV *B=0) const
Compute a residual R for this operator given a solution X, and right-hand side B. ...
Class which defines basic traits for the operator type.
Parent class to all Belos exceptions.
Teuchos::RCP< const OP > getOperator() const
A pointer to the (unpreconditioned) operator A.
void setRightPrec(const Teuchos::RCP< const OP > &RP)
Set right preconditioner (RP) of linear problem .
void setInitResVec(const Teuchos::RCP< const MV > &R0)
Set the user-defined residual of linear problem .
void setLHS(const Teuchos::RCP< MV > &X)
Set left-hand-side X of linear problem .
void apply(const MV &x, MV &y) const
Apply the composite operator of this linear problem to x, returning y.
Teuchos::RCP< const MV > getInitPrecResVec() const
A pointer to the preconditioned initial residual vector.
void setInitPrecResVec(const Teuchos::RCP< const MV > &PR0)
Set the user-defined preconditioned residual of linear problem .
void computeCurrPrecResVec(MV *R, const MV *X=0, const MV *B=0) const
Compute a residual R for this operator given a solution X, and right-hand side B. ...