11#ifndef BELOS_LINEAR_PROBLEM_HPP
12#define BELOS_LINEAR_PROBLEM_HPP
19#include "Teuchos_ParameterList.hpp"
20#include "Teuchos_TimeMonitor.hpp"
49 template <
class ScalarType,
class MV,
class OP>
71 const Teuchos::RCP<MV> &
X,
72 const Teuchos::RCP<const MV> &B);
110 void setRHS (
const Teuchos::RCP<const MV> &B) {
164 virtual void setLSIndex (
const std::vector<int>& index);
284 const Teuchos::RCP<const MV> &
newB = Teuchos::null);
298 Teuchos::RCP<const MV>
getRHS()
const {
return(
B_); }
384 Teuchos::Array<Teuchos::RCP<Teuchos::Time> >
getTimers()
const {
432 virtual void apply(
const MV&
x, MV& y )
const;
441 virtual void applyOp(
const MV&
x, MV& y )
const;
478 Teuchos::RCP<const OP>
A_;
487 Teuchos::RCP<const MV>
B_;
505 Teuchos::RCP<const OP>
LP_;
508 Teuchos::RCP<const OP>
RP_;
556 template <
class ScalarType,
class MV,
class OP>
564 solutionUpdated_(
false),
569 template <
class ScalarType,
class MV,
class OP>
571 const Teuchos::RCP<MV> &
X,
572 const Teuchos::RCP<const MV> &B
583 solutionUpdated_(
false),
588 template <
class ScalarType,
class MV,
class OP>
602 timerPrec_(
Problem.timerPrec_),
603 blocksize_(
Problem.blocksize_),
604 num2Solve_(
Problem.num2Solve_),
608 isHermitian_(
Problem.isHermitian_),
609 solutionUpdated_(
Problem.solutionUpdated_),
614 template <
class ScalarType,
class MV,
class OP>
622 curB_ = Teuchos::null;
623 curX_ = Teuchos::null;
627 blocksize_ = rhsIndex_.size();
628 std::vector<int>
vldIndex( blocksize_ );
629 std::vector<int>
newIndex( blocksize_ );
630 std::vector<int>
iIndex( blocksize_ );
631 for (
int i=0;
i<blocksize_; ++
i) {
632 if (rhsIndex_[
i] > -1) {
648 if (num2Solve_ != blocksize_) {
654 curX_ = MVT::Clone( *X_, blocksize_ );
656 Teuchos::RCP<MV>
tmpCurB = MVT::Clone( *B_, blocksize_ );
660 Teuchos::RCP<const MV>
tptr = MVT::CloneView( *B_,
vldIndex );
668 solutionUpdated_ =
false;
671 curX_ = MVT::CloneViewNonConst( *X_, rhsIndex_ );
672 curB_ = MVT::CloneView( *B_, rhsIndex_ );
681 template <
class ScalarType,
class MV,
class OP>
688 if (num2Solve_ < blocksize_) {
693 std::vector<int>
newIndex( num2Solve_ );
694 std::vector<int>
vldIndex( num2Solve_ );
695 for (
int i=0;
i<blocksize_; ++
i) {
696 if ( rhsIndex_[
i] > -1 ) {
702 Teuchos::RCP<MV>
tptr = MVT::CloneViewNonConst( *curX_,
newIndex );
709 curX_ = Teuchos::null;
710 curB_ = Teuchos::null;
715 template <
class ScalarType,
class MV,
class OP>
739 MVT::MvAddMv( 1.0, *curX_,
scale, *
update, *curX_ );
750 solutionUpdated_ =
true;
777 template <
class ScalarType,
class MV,
class OP>
780 if (
label != label_) {
783 if (timerOp_ != Teuchos::null) {
784 std::string
opLabel = label_ +
": Operation Op*x";
785#ifdef BELOS_TEUCHOS_TIME_MONITOR
786 timerOp_ = Teuchos::TimeMonitor::getNewCounter(
opLabel );
789 if (timerPrec_ != Teuchos::null) {
790 std::string
precLabel = label_ +
": Operation Prec*x";
791#ifdef BELOS_TEUCHOS_TIME_MONITOR
792 timerPrec_ = Teuchos::TimeMonitor::getNewCounter(
precLabel );
798 template <
class ScalarType,
class MV,
class OP>
802 const Teuchos::RCP<const MV> &
newB)
805 if (timerOp_ == Teuchos::null) {
806 std::string
opLabel = label_ +
": Operation Op*x";
807#ifdef BELOS_TEUCHOS_TIME_MONITOR
808 timerOp_ = Teuchos::TimeMonitor::getNewCounter(
opLabel );
811 if (timerPrec_ == Teuchos::null) {
812 std::string
precLabel = label_ +
": Operation Prec*x";
813#ifdef BELOS_TEUCHOS_TIME_MONITOR
814 timerPrec_ = Teuchos::TimeMonitor::getNewCounter(
precLabel );
818 Y_temp_ = Teuchos::null;
821 if (
newX != Teuchos::null)
823 if (
newB != Teuchos::null)
828 curX_ = Teuchos::null;
829 curB_ = Teuchos::null;
833 if (A_ == Teuchos::null || X_ == Teuchos::null || B_ == Teuchos::null) {
841 solutionUpdated_ =
false;
844 if(Teuchos::is_null(R0_user_)) {
845 if (R0_==Teuchos::null || MVT::GetNumberVecs( *R0_ )!=MVT::GetNumberVecs( *B_ )) {
846 R0_ = MVT::Clone( *B_, MVT::GetNumberVecs( *B_ ) );
848 computeCurrResVec( &*R0_, &*X_, &*B_ );
850 if (LP_!=Teuchos::null) {
851 if (PR0_==Teuchos::null || (PR0_==R0_) || (MVT::GetNumberVecs(*PR0_)!=MVT::GetNumberVecs(*B_))) {
852 PR0_ = MVT::Clone( *B_, MVT::GetNumberVecs( *B_ ) );
854 applyLeftPrec( *R0_, *PR0_ );
862 if (MVT::GetNumberVecs( *R0_user_ )!=MVT::GetNumberVecs( *B_ )) {
863 Teuchos::RCP<MV>
helper = MVT::Clone( *B_, MVT::GetNumberVecs( *B_ ) );
864 computeCurrResVec( &*
helper, &*X_, &*B_ );
865 R0_user_ = Teuchos::null;
869 if (LP_!=Teuchos::null) {
872 if (PR0_user_==Teuchos::null || (PR0_user_==R0_) || (PR0_user_==R0_user_)
873 || (MVT::GetNumberVecs(*PR0_user_)!=MVT::GetNumberVecs(*B_))) {
874 Teuchos::RCP<MV>
helper = MVT::Clone( *B_, MVT::GetNumberVecs( *B_ ) );
877 applyLeftPrec( *getInitResVec(), *
helper );
878 PR0_user_ = Teuchos::null;
885 if (R0_user_!=Teuchos::null)
887 PR0_user_ = R0_user_;
891 PR0_user_ = Teuchos::null;
904 template <
class ScalarType,
class MV,
class OP>
907 if(Teuchos::nonnull(R0_user_)) {
913 template <
class ScalarType,
class MV,
class OP>
916 if(Teuchos::nonnull(PR0_user_)) {
922 template <
class ScalarType,
class MV,
class OP>
929 return Teuchos::null;
933 template <
class ScalarType,
class MV,
class OP>
940 return Teuchos::null;
944 template <
class ScalarType,
class MV,
class OP>
954 const auto num_vecs_y = MVT::GetNumberVecs(y);
955 if(!Y_temp_ || MVT::GetNumberVecs(*Y_temp_) !=
num_vecs_y) {
971 applyRightPrec(
x, y );
972 applyOp( y, *Y_temp_ );
973 applyLeftPrec( *Y_temp_, y );
980 applyOp(
x, *Y_temp_ );
981 applyLeftPrec( *Y_temp_, y );
988 applyRightPrec(
x, *Y_temp_ );
989 applyOp( *Y_temp_, y );
993 template <
class ScalarType,
class MV,
class OP>
996#ifdef BELOS_TEUCHOS_TIME_MONITOR
997 Teuchos::TimeMonitor
OpTimer(*timerOp_);
999 OPT::Apply( *A_,
x, y);
1002 MVT::MvAddMv( Teuchos::ScalarTraits<ScalarType>::one(),
x,
1003 Teuchos::ScalarTraits<ScalarType>::zero(),
x, y );
1007 template <
class ScalarType,
class MV,
class OP>
1009 if (LP_!=Teuchos::null) {
1010#ifdef BELOS_TEUCHOS_TIME_MONITOR
1011 Teuchos::TimeMonitor
PrecTimer(*timerPrec_);
1013 OPT::Apply( *LP_,
x, y);
1016 MVT::MvAddMv( Teuchos::ScalarTraits<ScalarType>::one(),
x,
1017 Teuchos::ScalarTraits<ScalarType>::zero(),
x, y );
1021 template <
class ScalarType,
class MV,
class OP>
1023 if (RP_!=Teuchos::null) {
1024#ifdef BELOS_TEUCHOS_TIME_MONITOR
1025 Teuchos::TimeMonitor
PrecTimer(*timerPrec_);
1027 OPT::Apply( *RP_,
x, y);
1030 MVT::MvAddMv( Teuchos::ScalarTraits<ScalarType>::one(),
x,
1031 Teuchos::ScalarTraits<ScalarType>::zero(),
x, y );
1035 template <
class ScalarType,
class MV,
class OP>
1041 if (LP_!=Teuchos::null)
1043 Teuchos::RCP<MV>
R_temp = MVT::Clone( *B, MVT::GetNumberVecs( *B ) );
1046 applyLeftPrec( *
R_temp, *R );
1051 MVT::MvAddMv( -1.0, *R, 1.0, *B, *R );
1058 localB = Teuchos::rcp( B,
false );
1063 localX = Teuchos::rcp(
X,
false );
1067 if (LP_!=Teuchos::null)
1072 applyLeftPrec( *
R_temp, *R );
1077 MVT::MvAddMv( -1.0, *R, 1.0, *
localB, *R );
1084 template <
class ScalarType,
class MV,
class OP>
1091 MVT::MvAddMv( -1.0, *R, 1.0, *B, *R );
1097 localB = Teuchos::rcp( B,
false );
1102 localX = Teuchos::rcp(
X,
false );
1107 MVT::MvAddMv( -1.0, *R, 1.0, *
localB, *R );
Declaration of basic traits for the multivector type.
Class which defines basic traits for the operator type.
Parent class to all Belos exceptions.
Exception thrown to signal error with the Belos::LinearProblem object.
LinearProblemError(const std::string &what_arg)
A linear system to solve, and its associated information.
Teuchos::RCP< const MV > PR0_user_
User-defined preconditioned initial residual of the linear system.
std::vector< int > rhsIndex_
Indices of current linear systems being solver for.
bool isSet_
Has the linear problem to solve been set?
void setLabel(const std::string &label)
Set the label prefix used by the timers in this object.
bool isRightPrec() const
Whether the linear system is being preconditioned on the right.
void setInitPrecResVec(const Teuchos::RCP< const MV > &PR0)
Set the user-defined preconditioned residual of linear problem .
int blocksize_
Current block size of linear system.
virtual void applyLeftPrec(const MV &x, MV &y) const
Apply ONLY the left preconditioner to x, returning y.
bool isSolutionUpdated() const
Has the current approximate solution been updated?
std::string label_
Linear problem label that prefixes the timer labels.
Teuchos::RCP< const MV > curB_
Current right-hand side of the linear system.
void setLeftPrec(const Teuchos::RCP< const OP > &LP)
Set left preconditioner (LP) of linear problem .
Teuchos::RCP< const MV > B_
Right-hand side of linear system.
virtual 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.
virtual ~LinearProblem(void)
Destructor (declared virtual for memory safety of derived classes).
virtual void apply(const MV &x, MV &y) const
Apply the composite operator of this linear problem to x, returning y.
void setLHS(const Teuchos::RCP< MV > &X)
Set left-hand-side X of linear problem .
virtual void setLSIndex(const std::vector< int > &index)
Tell the linear problem which linear system(s) need to be solved next.
bool solutionUpdated_
Has the current approximate solution been updated?
Teuchos::RCP< MV > PR0_
Preconditioned initial residual of the linear system.
virtual void applyRightPrec(const MV &x, MV &y) const
Apply ONLY the right preconditioner to x, returning y.
Teuchos::RCP< const OP > LP_
Left preconditioning operator of linear system.
Teuchos::RCP< const OP > getLeftPrec() const
Get a pointer to the left preconditioner.
Teuchos::RCP< const MV > R0_user_
User-defined initial residual of the linear system.
Teuchos::RCP< MV > curX_
Current solution vector of the linear system.
Teuchos::RCP< MV > R0_
Initial residual of the linear system.
void setRightPrec(const Teuchos::RCP< const OP > &RP)
Set right preconditioner (RP) of linear problem .
virtual void setCurrLS()
Tell the linear problem that the solver is finished with the current linear system.
int getLSNumber() const
The number of linear systems that have been set.
void setHermitian(bool isSym=true)
Tell the linear problem that the (preconditioned) operator is Hermitian.
void setOperator(const Teuchos::RCP< const OP > &A)
Set the operator A of the linear problem .
bool isProblemSet() const
Whether the problem has been set.
int num2Solve_
Number of linear systems that are currently being solver for ( <= blocksize_ )
Teuchos::RCP< Teuchos::Time > timerPrec_
Teuchos::RCP< const MV > getInitPrecResVec() const
A pointer to the preconditioned initial residual vector.
Teuchos::RCP< const OP > getRightPrec() const
Get a pointer to the right preconditioner.
Teuchos::RCP< const OP > getOperator() const
A pointer to the (unpreconditioned) operator A.
Teuchos::RCP< Teuchos::Time > timerOp_
Timers.
const std::vector< int > & getLSIndex() const
(Zero-based) indices of the linear system(s) currently being solved.
Teuchos::Array< Teuchos::RCP< Teuchos::Time > > getTimers() const
The timers for this object.
bool isHermitian_
Whether the operator A is symmetric (in real arithmetic, or Hermitian in complex arithmetic).
LinearProblem(void)
Default constructor.
void setRHS(const Teuchos::RCP< const MV > &B)
Set right-hand-side B of linear problem .
Teuchos::RCP< const MV > getRHS() const
A pointer to the right-hand side B.
bool isHermitian() const
Whether the (preconditioned) operator is Hermitian.
Teuchos::RCP< MV > getLHS() const
A pointer to the left-hand side X.
virtual void applyOp(const MV &x, MV &y) const
Apply ONLY the operator to x, returning y.
Teuchos::RCP< const MV > getInitResVec() const
A pointer to the initial unpreconditioned residual vector.
Teuchos::RCP< MV > getCurrLHSVec()
Get a pointer to the current left-hand side (solution) of the linear system.
bool isLeftPrec() const
Whether the linear system is being preconditioned on the left.
Teuchos::RCP< const MV > getCurrRHSVec()
Get a pointer to the current right-hand side of the linear system.
virtual 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.
int lsNum_
Number of linear systems that have been loaded in this linear problem object.
Teuchos::RCP< MV > Y_temp_
Scratch vector for use in apply()
Teuchos::RCP< MV > X_
Solution vector of linear system.
Teuchos::RCP< const OP > A_
Operator of linear system.
virtual 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.
void setInitResVec(const Teuchos::RCP< const MV > &R0)
Set the user-defined residual of linear problem .
virtual 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.
virtual bool setProblem(const Teuchos::RCP< MV > &newX=Teuchos::null, const Teuchos::RCP< const MV > &newB=Teuchos::null)
Set up the linear problem manager.
Teuchos::RCP< const OP > RP_
Right preconditioning operator of linear system.
Alternative run-time polymorphic interface for operators.
Operator()
Default constructor (does nothing).