Belos  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
List of all members
Belos::LinearProblem< ScalarType, MV, OP > Class Template Reference

A linear system to solve, and its associated information. More...

#include <BelosIteration.hpp>

Public Member Functions

Constructors/Destructor
 LinearProblem (void)
 Default constructor. More...
 
 LinearProblem (const Teuchos::RCP< const OP > &A, const Teuchos::RCP< MV > &X, const Teuchos::RCP< const MV > &B)
 Unpreconditioned linear system constructor. More...
 
 LinearProblem (const LinearProblem< ScalarType, MV, OP > &Problem)
 Copy constructor. More...
 
virtual ~LinearProblem (void)
 Destructor (declared virtual for memory safety of derived classes). More...
 
Set methods
void setOperator (const Teuchos::RCP< const OP > &A)
 Set the operator A of the linear problem $AX = B$. More...
 
void setLHS (const Teuchos::RCP< MV > &X)
 Set left-hand-side X of linear problem $AX = B$. More...
 
void setRHS (const Teuchos::RCP< const MV > &B)
 Set right-hand-side B of linear problem $AX = B$. More...
 
void setInitResVec (const Teuchos::RCP< const MV > &R0)
 Set the user-defined residual of linear problem $AX = B$. More...
 
void setInitPrecResVec (const Teuchos::RCP< const MV > &PR0)
 Set the user-defined preconditioned residual of linear problem $AX = B$. More...
 
void setLeftPrec (const Teuchos::RCP< const OP > &LP)
 Set left preconditioner (LP) of linear problem $AX = B$. More...
 
void setRightPrec (const Teuchos::RCP< const OP > &RP)
 Set right preconditioner (RP) of linear problem $AX = B$. More...
 
void setCurrLS ()
 Tell the linear problem that the solver is finished with the current linear system. More...
 
void setLSIndex (const std::vector< int > &index)
 Tell the linear problem which linear system(s) need to be solved next. More...
 
void setHermitian (bool isSym=true)
 Tell the linear problem that the (preconditioned) operator is Hermitian. More...
 
void setLabel (const std::string &label)
 Set the label prefix used by the timers in this object. More...
 
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. More...
 
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. More...
 
Set / Reset method
bool setProblem (const Teuchos::RCP< MV > &newX=Teuchos::null, const Teuchos::RCP< const MV > &newB=Teuchos::null)
 Set up the linear problem manager. More...
 
Accessor methods
Teuchos::RCP< const OP > getOperator () const
 A pointer to the (unpreconditioned) operator A. More...
 
Teuchos::RCP< MV > getLHS () const
 A pointer to the left-hand side X. More...
 
Teuchos::RCP< const MV > getRHS () const
 A pointer to the right-hand side B. More...
 
Teuchos::RCP< const MV > getInitResVec () const
 A pointer to the initial unpreconditioned residual vector. More...
 
Teuchos::RCP< const MV > getInitPrecResVec () const
 A pointer to the preconditioned initial residual vector. More...
 
Teuchos::RCP< MV > getCurrLHSVec ()
 Get a pointer to the current left-hand side (solution) of the linear system. More...
 
Teuchos::RCP< const MV > getCurrRHSVec ()
 Get a pointer to the current right-hand side of the linear system. More...
 
Teuchos::RCP< const OP > getLeftPrec () const
 Get a pointer to the left preconditioner. More...
 
Teuchos::RCP< const OP > getRightPrec () const
 Get a pointer to the right preconditioner. More...
 
const std::vector< int > getLSIndex () const
 (Zero-based) indices of the linear system(s) currently being solved. More...
 
int getLSNumber () const
 The number of linear systems that have been set. More...
 
Teuchos::Array< Teuchos::RCP
< Teuchos::Time > > 
getTimers () const
 The timers for this object. More...
 
State methods
bool isSolutionUpdated () const
 Has the current approximate solution been updated? More...
 
bool isProblemSet () const
 Whether the problem has been set. More...
 
bool isHermitian () const
 Whether the (preconditioned) operator is Hermitian. More...
 
bool isLeftPrec () const
 Whether the linear system is being preconditioned on the left. More...
 
bool isRightPrec () const
 Whether the linear system is being preconditioned on the right. More...
 
Apply / Compute methods
void apply (const MV &x, MV &y) const
 Apply the composite operator of this linear problem to x, returning y. More...
 
void applyOp (const MV &x, MV &y) const
 Apply ONLY the operator to x, returning y. More...
 
void applyLeftPrec (const MV &x, MV &y) const
 Apply ONLY the left preconditioner to x, returning y. More...
 
void applyRightPrec (const MV &x, MV &y) const
 Apply ONLY the right preconditioner to x, returning y. More...
 
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. More...
 
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. More...
 

Detailed Description

template<class ScalarType, class MV, class OP>
class Belos::LinearProblem< ScalarType, MV, OP >

A linear system to solve, and its associated information.

This class encapsulates the general information needed for solving a linear system of equations using an iterative method.

Template Parameters
ScalarTypeThe type of the entries in the matrix and vectors.
MVThe (multi)vector type.
OPThe operator type. (Operators are functions that take a multivector as input and compute a multivector as output.)
Examples:
BlockCG/BlockCGEpetraExFile.cpp, BlockCG/BlockPrecCGEpetraExFile.cpp, BlockCG/PseudoBlockCGEpetraExFile.cpp, BlockCG/PseudoBlockPrecCGEpetraExFile.cpp, BlockGmres/BlockFlexGmresEpetraExFile.cpp, BlockGmres/BlockGmresEpetraExFile.cpp, BlockGmres/BlockGmresPolyEpetraExFile.cpp, BlockGmres/BlockPrecGmresEpetraExFile.cpp, BlockGmres/PseudoBlockGmresEpetraExFile.cpp, BlockGmres/PseudoBlockPrecGmresEpetraExFile.cpp, GCRODR/GCRODREpetraExFile.cpp, GCRODR/PrecGCRODREpetraExFile.cpp, PCPG/PCPGEpetraExFile.cpp, TFQMR/PseudoBlockTFQMREpetraExFile.cpp, and TFQMR/TFQMREpetraExFile.cpp.

Definition at line 61 of file BelosIteration.hpp.

Constructor & Destructor Documentation

template<class ScalarType , class MV , class OP >
Belos::LinearProblem< ScalarType, MV, OP >::LinearProblem ( void  )

Default constructor.

Creates an empty Belos::LinearProblem instance. The operator A, left-hand-side X and right-hand-side B must be set using the setOperator(), setLHS(), and setRHS() methods respectively.

Definition at line 580 of file BelosLinearProblem.hpp.

template<class ScalarType , class MV , class OP >
Belos::LinearProblem< ScalarType, MV, OP >::LinearProblem ( const Teuchos::RCP< const OP > &  A,
const Teuchos::RCP< MV > &  X,
const Teuchos::RCP< const MV > &  B 
)

Unpreconditioned linear system constructor.

Creates an unpreconditioned LinearProblem instance with the operator (A), initial guess (X), and right hand side (B). Preconditioners can be set using the setLeftPrec() and setRightPrec() methods.

Definition at line 593 of file BelosLinearProblem.hpp.

template<class ScalarType , class MV , class OP >
Belos::LinearProblem< ScalarType, MV, OP >::LinearProblem ( const LinearProblem< ScalarType, MV, OP > &  Problem)

Copy constructor.

Makes a copy of an existing LinearProblem instance.

Definition at line 612 of file BelosLinearProblem.hpp.

template<class ScalarType , class MV , class OP >
Belos::LinearProblem< ScalarType, MV, OP >::~LinearProblem ( void  )
virtual

Destructor (declared virtual for memory safety of derived classes).

Definition at line 638 of file BelosLinearProblem.hpp.

Member Function Documentation

template<class ScalarType, class MV, class OP>
void Belos::LinearProblem< ScalarType, MV, OP >::setOperator ( const Teuchos::RCP< const OP > &  A)
inline

Set the operator A of the linear problem $AX = B$.

The operator is set by pointer; no copy of the operator is made.

Definition at line 122 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
void Belos::LinearProblem< ScalarType, MV, OP >::setLHS ( const Teuchos::RCP< MV > &  X)
inline

Set left-hand-side X of linear problem $AX = B$.

Setting the "left-hand side" sets the starting vector (also called "initial guess") of an iterative method. The multivector is set by pointer; no copy of the object is made. Belos' solvers will modify this multivector in place.

Definition at line 133 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
void Belos::LinearProblem< ScalarType, MV, OP >::setRHS ( const Teuchos::RCP< const MV > &  B)
inline

Set right-hand-side B of linear problem $AX = B$.

The multivector is set by pointer; no copy of the object is made.

Definition at line 142 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
void Belos::LinearProblem< ScalarType, MV, OP >::setInitResVec ( const Teuchos::RCP< const MV > &  R0)
inline

Set the user-defined residual of linear problem $AX = B$.

The multivector is set by pointer; no copy of the object is made.

Note
If the passed in residual vector is not compatible with B, this vector will be ignored.

Definition at line 152 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
void Belos::LinearProblem< ScalarType, MV, OP >::setInitPrecResVec ( const Teuchos::RCP< const MV > &  PR0)
inline

Set the user-defined preconditioned residual of linear problem $AX = B$.

The multivector is set by pointer; no copy of the object is made.

Note
If the passed in residual vector is not compatible with B, this vector will be ignored.

Definition at line 162 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
void Belos::LinearProblem< ScalarType, MV, OP >::setLeftPrec ( const Teuchos::RCP< const OP > &  LP)
inline

Set left preconditioner (LP) of linear problem $AX = B$.

The operator is set by pointer; no copy of the operator is made.

Definition at line 170 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
void Belos::LinearProblem< ScalarType, MV, OP >::setRightPrec ( const Teuchos::RCP< const OP > &  RP)
inline

Set right preconditioner (RP) of linear problem $AX = B$.

The operator is set by pointer; no copy of the operator is made.

Definition at line 175 of file BelosLinearProblem.hpp.

template<class ScalarType , class MV , class OP >
void Belos::LinearProblem< ScalarType, MV, OP >::setCurrLS ( )

Tell the linear problem that the solver is finished with the current linear system.

Note
This method is only to be used by the solver to inform the linear problem that it is finished with the current block of linear systems. The next time that Curr{RHS, LHS}Vec() is called, the next linear system will be returned. Computing the next linear system isn't done in this method in case the block size is changed.

Definition at line 709 of file BelosLinearProblem.hpp.

template<class ScalarType , class MV , class OP >
void Belos::LinearProblem< ScalarType, MV, OP >::setLSIndex ( const std::vector< int > &  index)

Tell the linear problem which linear system(s) need to be solved next.

Any calls to get the current RHS/LHS vectors after this method is called will return the new linear system(s) indicated by index. The length of index is assumed to be the blocksize. Entries of index must be between 0 and the number of vectors in the RHS/LHS multivector. An entry of index may also be -1, which means this column of the linear system is augmented using a random vector.

Definition at line 642 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
void Belos::LinearProblem< ScalarType, MV, OP >::setHermitian ( bool  isSym = true)
inline

Tell the linear problem that the (preconditioned) operator is Hermitian.

This knowledge may allow the operator to take advantage of the linear problem's symmetry. However, this method should not be called if the preconditioned operator is not Hermitian (or symmetric in real arithmetic).

We make no attempt to detect the symmetry of the operators, so we cannot check whether this method has been called incorrectly.

Definition at line 208 of file BelosLinearProblem.hpp.

template<class ScalarType , class MV , class OP >
void Belos::LinearProblem< ScalarType, MV, OP >::setLabel ( const std::string &  label)

Set the label prefix used by the timers in this object.

The default label prefix is "Belos". The timers are created during the call to setProblem(). If they have already been created and this label is different than the current one, then this method will generate a new timer.

Definition at line 815 of file BelosLinearProblem.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP< MV > Belos::LinearProblem< ScalarType, MV, OP >::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.

Let $\delta$ be the update vector, $\alpha$ the scale factor, and $x$ the current solution. If there is a right preconditioner $M_R^{-1}$, then we compute the new solution as $x + \alpha M_R^{-1} \delta$. Otherwise, if there is no right preconditioner, we compute the new solution as $x + \alpha \delta$.

This method always returns the new solution. If updateLP is false, it computes the new solution as a deep copy, without modifying the internally stored current solution. If updateLP is true, it computes the new solution in place, and returns a pointer to the internally stored solution.

Parameters
update[in/out] The solution update vector. If null, this method returns a pointer to the new solution.
updateLP[in] This is ignored if the update vector is null. Otherwise, if updateLP is true, the following things happen: (a) this LinearProblem's stored solution is updated in place, and (b) the next time GetCurrResVecs() is called, a new residual will be computed. If updateLP is false, then the new solution is computed and returned as a copy, without modifying this LinearProblem's stored current solution.
scale[in] The factor $\alpha$ by which to multiply the solution update vector when computing the update. This is ignored if the update vector is null.
Returns
A pointer to the new solution. This is freshly allocated if updateLP is false, otherwise it is a view of the LinearProblem's stored current solution.

Definition at line 745 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
Teuchos::RCP<MV> Belos::LinearProblem< ScalarType, MV, OP >::updateSolution ( const Teuchos::RCP< MV > &  update = Teuchos::null,
ScalarType  scale = Teuchos::ScalarTraits<ScalarType>::one() 
) const
inline

Compute the new solution to the linear system using the given update vector.

This method does the same thing as calling the three-argument version of updateSolution() with updateLP = false. It does not update the linear problem or change the linear problem's state in any way.

Parameters
update[in/out] The solution update vector. If null, this method returns a pointer to the new solution.
scale[in] The factor $\alpha$ by which to multiply the solution update vector when computing the update. This is ignored if the update vector is null.
Returns
A pointer to the new solution.

Definition at line 276 of file BelosLinearProblem.hpp.

template<class ScalarType , class MV , class OP >
bool Belos::LinearProblem< ScalarType, MV, OP >::setProblem ( const Teuchos::RCP< MV > &  newX = Teuchos::null,
const Teuchos::RCP< const MV > &  newB = Teuchos::null 
)

Set up the linear problem manager.

Call this method if you want to solve the linear system with a different left- or right-hand side, or if you want to prepare the linear problem to solve the linear system that was already passed in. (In the latter case, call this method with the default arguments.) The internal flags will be set as if the linear system manager was just initialized, and the initial residual will be computed.

Many of Belos' solvers require that this method has been called on the linear problem, before they can solve it.

Parameters
newX[in/out] If you want to solve the linear system with a different left-hand side, pass it in here. Otherwise, set this to null (the default value).
newB[in] If you want to solve the linear system with a different right-hand side, pass it in here. Otherwise, set this to null (the default value).
Returns
Whether the linear problem was successfully set up. Successful setup requires at least that the matrix operator A, the left-hand side X, and the right-hand side B all be non-null.

Definition at line 838 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
Teuchos::RCP<const OP> Belos::LinearProblem< ScalarType, MV, OP >::getOperator ( ) const
inline

A pointer to the (unpreconditioned) operator A.

Definition at line 320 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
Teuchos::RCP<MV> Belos::LinearProblem< ScalarType, MV, OP >::getLHS ( ) const
inline

A pointer to the left-hand side X.

Definition at line 323 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
Teuchos::RCP<const MV> Belos::LinearProblem< ScalarType, MV, OP >::getRHS ( ) const
inline

A pointer to the right-hand side B.

Definition at line 326 of file BelosLinearProblem.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP< const MV > Belos::LinearProblem< ScalarType, MV, OP >::getInitResVec ( ) const

A pointer to the initial unpreconditioned residual vector.

Definition at line 949 of file BelosLinearProblem.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP< const MV > Belos::LinearProblem< ScalarType, MV, OP >::getInitPrecResVec ( ) const

A pointer to the preconditioned initial residual vector.

Note
This is the preconditioned residual if the linear system is left preconditioned.

Definition at line 958 of file BelosLinearProblem.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP< MV > Belos::LinearProblem< ScalarType, MV, OP >::getCurrLHSVec ( )

Get a pointer to the current left-hand side (solution) of the linear system.

This method is called by the solver or any method that is interested in the current linear system being solved.

  • If the solution has been updated by the solver, then this vector is current ( see isSolutionUpdated() ).
  • If there is no linear system to solve, this method returns a null pointer.

This method is not the same thing as getLHS(). The getLHS() method just returns a pointer to the original left-hand side vector. This method only returns a valid vector if the current subset of right-hand side(s) to solve has been set (via the setLSIndex() method).

Definition at line 967 of file BelosLinearProblem.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP< const MV > Belos::LinearProblem< ScalarType, MV, OP >::getCurrRHSVec ( )

Get a pointer to the current right-hand side of the linear system.

This method is called by the solver or any method that is interested in the current linear system being solved.

  • If the solution has been updated by the solver, then this vector is current ( see isSolutionUpdated() ).
  • If there is no linear system to solve, this method returns a null pointer.

This method is not the same thing as getRHS(). The getRHS() method just returns a pointer to the original right-hand side vector. This method only returns a valid vector if the current subset of right-hand side(s) to solve has been set (via the setLSIndex() method).

Definition at line 978 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
Teuchos::RCP<const OP> Belos::LinearProblem< ScalarType, MV, OP >::getLeftPrec ( ) const
inline

Get a pointer to the left preconditioner.

Definition at line 370 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
Teuchos::RCP<const OP> Belos::LinearProblem< ScalarType, MV, OP >::getRightPrec ( ) const
inline

Get a pointer to the right preconditioner.

Definition at line 373 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
const std::vector<int> Belos::LinearProblem< ScalarType, MV, OP >::getLSIndex ( ) const
inline

(Zero-based) indices of the linear system(s) currently being solved.

Since the block size is independent of the number of right-hand sides for some solvers (GMRES, CG, etc.), it is important to know which linear systems are being solved. That may mean you need to update the information about the norms of your initial residual vector for weighting purposes. This information can help you avoid querying the solver for information that rarely changes.

Note
The length of the returned index vector is the number of right-hand sides currently being solved. If an entry of the index vector is -1, then the corresponding linear system is an augmented linear system and doesn't need to be considered for convergence.
The index vector returned from this method can only be nonempty if setLSIndex() has been called with a nonempty index vector argument, or if this linear problem was constructed via the copy constructor of a linear problem with a nonempty index vector.

Definition at line 396 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
int Belos::LinearProblem< ScalarType, MV, OP >::getLSNumber ( ) const
inline

The number of linear systems that have been set.

This can be used by status test classes to determine if the solver manager has advanced and is solving another linear system. This is incremented by one every time that setLSIndex() completes successfully.

Definition at line 404 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
Teuchos::Array<Teuchos::RCP<Teuchos::Time> > Belos::LinearProblem< ScalarType, MV, OP >::getTimers ( ) const
inline

The timers for this object.

The timers are ordered as follows:

  • time spent applying operator
  • time spent applying preconditioner

Definition at line 412 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
bool Belos::LinearProblem< ScalarType, MV, OP >::isSolutionUpdated ( ) const
inline

Has the current approximate solution been updated?

This only means that the current linear system for which the solver is solving (as obtained by getCurr{LHS, RHS}Vec()) has been updated by the solver. This will be true every iteration for solvers like CG, but not true for solvers like GMRES until the solver restarts.

Definition at line 429 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
bool Belos::LinearProblem< ScalarType, MV, OP >::isProblemSet ( ) const
inline

Whether the problem has been set.

Definition at line 432 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
bool Belos::LinearProblem< ScalarType, MV, OP >::isHermitian ( ) const
inline

Whether the (preconditioned) operator is Hermitian.

If preconditioner(s) are defined and this method returns true, then the entire preconditioned operator is Hermitian (or symmetric in real arithmetic).

Definition at line 439 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
bool Belos::LinearProblem< ScalarType, MV, OP >::isLeftPrec ( ) const
inline

Whether the linear system is being preconditioned on the left.

Definition at line 442 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
bool Belos::LinearProblem< ScalarType, MV, OP >::isRightPrec ( ) const
inline

Whether the linear system is being preconditioned on the right.

Definition at line 445 of file BelosLinearProblem.hpp.

template<class ScalarType , class MV , class OP >
void Belos::LinearProblem< ScalarType, MV, OP >::apply ( const MV &  x,
MV &  y 
) const

Apply the composite operator of this linear problem to x, returning y.

This application is the composition of the left/right preconditioner and operator. Most Krylov methods will use this application method within their code.

Precondition:

Definition at line 989 of file BelosLinearProblem.hpp.

template<class ScalarType , class MV , class OP >
void Belos::LinearProblem< ScalarType, MV, OP >::applyOp ( const MV &  x,
MV &  y 
) const

Apply ONLY the operator to x, returning y.

This method only applies the linear problem's operator, without any preconditioners that may have been defined. Flexible variants of Krylov methods will use this method. If no operator has been defined, this method just copies x into y.

Definition at line 1075 of file BelosLinearProblem.hpp.

template<class ScalarType , class MV , class OP >
void Belos::LinearProblem< ScalarType, MV, OP >::applyLeftPrec ( const MV &  x,
MV &  y 
) const

Apply ONLY the left preconditioner to x, returning y.

This method only applies the left preconditioner. This may be required for flexible variants of Krylov methods. If no left preconditioner has been defined, this method just copies x into y.

Definition at line 1089 of file BelosLinearProblem.hpp.

template<class ScalarType , class MV , class OP >
void Belos::LinearProblem< ScalarType, MV, OP >::applyRightPrec ( const MV &  x,
MV &  y 
) const

Apply ONLY the right preconditioner to x, returning y.

This method only applies the right preconditioner. This may be required for flexible variants of Krylov methods. If no right preconditioner has been defined, this method just copies x into y.

Definition at line 1103 of file BelosLinearProblem.hpp.

template<class ScalarType , class MV , class OP >
void Belos::LinearProblem< ScalarType, MV, OP >::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.

This method will compute the residual for the current linear system if X and B are null pointers. The result will be returned into R. Otherwise R = OP(A)X - B will be computed and returned.

Note
This residual will not be preconditioned if the system has a left preconditioner.

Definition at line 1196 of file BelosLinearProblem.hpp.

template<class ScalarType , class MV , class OP >
void Belos::LinearProblem< ScalarType, MV, OP >::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.

This method will compute the residual for the current linear system if X and B are null pointers. The result will be returned into R. Otherwise R = OP(A)X - B will be computed and returned.

Note
This residual will be preconditioned if the system has a left preconditioner.

Definition at line 1117 of file BelosLinearProblem.hpp.


The documentation for this class was generated from the following files:

Generated on Fri Jun 5 2020 10:20:58 for Belos by doxygen 1.8.5