Belos Package Browser (Single Doxygen Collection) Development
Private Types | Private Attributes
Belos::LinearProblem< ScalarType, MV, OP > Class Template Reference

The Belos::LinearProblem class is a wrapper that encapsulates the general information needed for solving a linear system of equations. More...

#include <BelosLinearProblem.hpp>

List of all members.

Private Types

typedef MultiVecTraits
< ScalarType, MV > 
MVT
typedef OperatorTraits
< ScalarType, MV, OP > 
OPT

Private Attributes

RCP< const OP > A_
 Operator of linear system.
RCP< MV > X_
 Solution std::vector of linear system.
RCP< MV > curX_
 Current solution std::vector of the linear system.
RCP< const MV > B_
 Right-hand side of linear system.
RCP< MV > curB_
 Current right-hand side of the linear system.
RCP< MV > R0_
 Initial residual of the linear system.
RCP< MV > PR0_
 Preconditioned initial residual of the linear system.
RCP< const OP > LP_
 Left preconditioning operator of linear system.
RCP< const OP > RP_
 Right preconditioning operator of linear system.
Teuchos::RCP< Teuchos::TimetimerOp_
 Timers.
Teuchos::RCP< Teuchos::TimetimerPrec_
int blocksize_
 Current block size of linear system.
int num2Solve_
 Number of linear systems that are currently being solver for ( <= blocksize_ )
std::vector< int > rhsIndex_
 Indices of current linear systems being solver for.
int lsNum_
 Number of linear systems that have been loaded in this linear problem object.
bool Left_Scale_
 Booleans to keep track of linear problem attributes/status.
bool Right_Scale_
bool isSet_
bool isHermitian_
bool solutionUpdated_
std::string label_
 Linear problem label that prefixes the timer labels.

Constructors/Destructor

 LinearProblem (void)
 Default Constructor.
 LinearProblem (const RCP< const OP > &A, const RCP< MV > &X, const RCP< const MV > &B)
 Unpreconditioned linear system constructor.
 LinearProblem (const LinearProblem< ScalarType, MV, OP > &Problem)
 Copy Constructor.
virtual ~LinearProblem (void)
 Destructor.

Set methods

void setOperator (const RCP< const OP > &A)
 Set Operator A of linear problem AX = B.
void setLHS (const RCP< MV > &X)
 Set left-hand-side X of linear problem AX = B.
void setRHS (const RCP< const MV > &B)
 Set right-hand-side B of linear problem AX = B.
void setLeftPrec (const RCP< const OP > &LP)
 Set left preconditioning operator (LP) of linear problem AX = B.
void setRightPrec (const RCP< const OP > &RP)
 Set right preconditioning operator (RP) of linear problem AX = B.
void setCurrLS ()
 Inform the linear problem that the solver is finished with the current linear system.
void setLSIndex (std::vector< int > &index)
 Inform the linear problem of the linear systems that need to be solved next.
void setHermitian ()
 Inform the linear problem that the operator is Hermitian.
void setLabel (const std::string &label)
 Set the label prefix used by the timers in this object. The default is "Belos".
RCP< MV > updateSolution (const RCP< MV > &update=null, bool updateLP=false, ScalarType scale=Teuchos::ScalarTraits< ScalarType >::one())
 Compute the new solution to the linear system given the /c update.
RCP< MV > updateSolution (const RCP< MV > &update=null, ScalarType scale=Teuchos::ScalarTraits< ScalarType >::one()) const
 Compute the new solution to the linear system given the /c update without updating the linear problem.

Set / Reset method

bool setProblem (const RCP< MV > &newX=null, const RCP< const MV > &newB=null)
 Setup the linear problem manager.

Accessor methods

RCP< const OP > getOperator () const
 Get a pointer to the operator A.
RCP< MV > getLHS () const
 Get a pointer to the left-hand side X.
RCP< const MV > getRHS () const
 Get a pointer to the right-hand side B.
RCP< const MV > getInitResVec () const
 Get a pointer to the initial residual vector.
RCP< const MV > getInitPrecResVec () const
 Get a pointer to the preconditioned initial residual vector.
RCP< MV > getCurrLHSVec ()
 Get a pointer to the current left-hand side (solution) of the linear system.
RCP< MV > getCurrRHSVec ()
 Get a pointer to the current right-hand side of the linear system.
RCP< const OP > getLeftPrec () const
 Get a pointer to the left preconditioning operator.
RCP< const OP > getRightPrec () const
 Get a pointer to the right preconditioning operator.
const std::vector< int > getLSIndex () const
 Get the 0-based index std::vector indicating the current linear systems being solved for.
int getLSNumber () const
 Get the number of linear systems that have been set with this LinearProblem object.
Teuchos::Array< Teuchos::RCP
< Teuchos::Time > > 
getTimers () const
 Return the timers for this object.

State methods

bool isSolutionUpdated () const
 Get the current status of the solution.
bool isProblemSet () const
 If the problem has been set, this will return true.
bool isHermitian () const
 Get the current symmetry of the operator.
bool isLeftPrec () const
 Get information on whether the linear system is being preconditioned on the left.
bool isRightPrec () const
 Get information on whether the linear system is being preconditioned on the right.

Apply / Compute methods

void apply (const MV &x, MV &y) const
 Apply the composite operator of this linear problem to x, returning y.
void applyOp (const MV &x, MV &y) const
 Apply ONLY the operator to x, returning y.
void applyLeftPrec (const MV &x, MV &y) const
 Apply ONLY the left preconditioner to x, returning y.
void applyRightPrec (const MV &x, MV &y) const
 Apply ONLY the right preconditioner to x, returning y.
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.
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.

Detailed Description

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

The Belos::LinearProblem class is a wrapper that encapsulates the general information needed for solving a linear system of equations.

Definition at line 80 of file BelosLinearProblem.hpp.


Member Typedef Documentation

template<class ScalarType, class MV, class OP>
typedef MultiVecTraits<ScalarType,MV> Belos::LinearProblem< ScalarType, MV, OP >::MVT [private]

Definition at line 415 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
typedef OperatorTraits<ScalarType,MV,OP> Belos::LinearProblem< ScalarType, MV, OP >::OPT [private]

Definition at line 416 of file BelosLinearProblem.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 424 of file BelosLinearProblem.hpp.

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

Unpreconditioned linear system constructor.

Creates an unpreconditioned LinearProblem instance with the Belos::Operator (A), initial guess (X), and right hand side (B). Preconditioners can be set using the setLeftPrec() and setRightPrec() methods, and scaling can also be set using the setLeftScale() and setRightScale() methods.

Definition at line 439 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 copy of an existing LinearProblem instance.

Definition at line 460 of file BelosLinearProblem.hpp.

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

Destructor.

Completely deletes a LinearProblem object.

Definition at line 486 of file BelosLinearProblem.hpp.


Member Function Documentation

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

Set Operator A of linear problem AX = B.

Sets a pointer to an Operator. No copy of the operator is made.

Definition at line 121 of file BelosLinearProblem.hpp.

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

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

Sets a pointer to a MultiVec. No copy of the object is made.

Definition at line 126 of file BelosLinearProblem.hpp.

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

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

Sets a pointer to a MultiVec. No copy of the object is made.

Definition at line 131 of file BelosLinearProblem.hpp.

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

Set left preconditioning operator (LP) of linear problem AX = B.

Sets a pointer to an Operator. No copy of the operator is made.

Definition at line 136 of file BelosLinearProblem.hpp.

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

Set right preconditioning operator (RP) of linear problem AX = B.

Sets a pointer to an Operator. No copy of the operator is made.

Definition at line 141 of file BelosLinearProblem.hpp.

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

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

Note:
This method is to be only used by the solver to inform the linear problem that it's finished with this block of linear systems. The next time the 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 blocksize is changed.

Definition at line 556 of file BelosLinearProblem.hpp.

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

Inform the linear problem of the linear systems that need to be solved next.

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

Definition at line 490 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
void Belos::LinearProblem< ScalarType, MV, OP >::setHermitian ( ) [inline]

Inform the linear problem that the operator is Hermitian.

This knowledge may allow the operator to take advantage of the linear problem symmetry. However, this should not be set to true if the preconditioner is not Hermitian, or symmetrically applied.

Definition at line 165 of file BelosLinearProblem.hpp.

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

Set the label prefix used by the timers in this object. The default is "Belos".

Note:
The timers are created during the first call to setProblem(). Any calls to this method to change the label after that will not change the label used in the timer.

Definition at line 171 of file BelosLinearProblem.hpp.

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

Compute the new solution to the linear system given the /c update.

Note:
If updateLP is true, then the next time GetCurrResVecs is called, a new residual will be computed. This keeps the linear problem from having to recompute the residual vector everytime it's asked for if the solution hasn't been updated. If updateLP is false, the new solution is computed without actually updating the linear problem.

Definition at line 590 of file BelosLinearProblem.hpp.

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

Compute the new solution to the linear system given the /c update without updating the linear problem.

Definition at line 184 of file BelosLinearProblem.hpp.

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

Setup the linear problem manager.

This is useful for solving the linear system with another right-hand side or getting the linear problem prepared to solve the linear system that was already passed in. The internal flags will be set as if the linear system manager was just initialized and the initial residual will be computed.

Definition at line 638 of file BelosLinearProblem.hpp.

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

Get a pointer to the operator A.

Definition at line 207 of file BelosLinearProblem.hpp.

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

Get a pointer to the left-hand side X.

Definition at line 210 of file BelosLinearProblem.hpp.

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

Get a pointer to the right-hand side B.

Definition at line 213 of file BelosLinearProblem.hpp.

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

Get a pointer to the initial residual vector.

Note:
This is the unpreconditioned residual.

Definition at line 218 of file BelosLinearProblem.hpp.

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

Get a pointer to the preconditioned initial residual vector.

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

Definition at line 223 of file BelosLinearProblem.hpp.

template<class ScalarType , class MV , class OP >
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 for.

  1. If the solution has been updated by the solver, then this vector is current ( see SolutionUpdated() ).
  2. If there is no linear system to solve, this method will return a NULL pointer

Definition at line 698 of file BelosLinearProblem.hpp.

template<class ScalarType , class MV , class OP >
RCP< 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 of any method that is interested in the current linear system being solved for.

  1. If the solution has been updated by the solver, then this vector is current ( see SolutionUpdated() ).
  2. If there is no linear system to solve, this method will return a NULL pointer

Definition at line 709 of file BelosLinearProblem.hpp.

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

Get a pointer to the left preconditioning operator.

Definition at line 246 of file BelosLinearProblem.hpp.

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

Get a pointer to the right preconditioning operator.

Definition at line 249 of file BelosLinearProblem.hpp.

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

Get the 0-based index std::vector indicating the current linear systems being solved for.

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 for. That may mean you need to update the information about the norms of your initial residual std::vector for weighting purposes. This information can keep you from querying the solver for information that rarely changes.

Note:
The length of the index std::vector is the number of right-hand sides being solved for. If an entry of this std::vector is -1 then that linear system is an augmented linear system and doesn't need to be considered for convergence.
The std::vector returned from this method is valid if isProblemSet() returns true.

Definition at line 263 of file BelosLinearProblem.hpp.

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

Get the number of linear systems that have been set with this LinearProblem object.

Definition at line 269 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]

Return the timers for this object.

The timers are ordered as follows:

  • time spent applying operator
  • time spent applying preconditioner

Definition at line 277 of file BelosLinearProblem.hpp.

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

Get the current status of the solution.

This only means that the current linear system being solved for ( 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 until restarts for GMRES.

Definition at line 292 of file BelosLinearProblem.hpp.

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

If the problem has been set, this will return true.

Definition at line 295 of file BelosLinearProblem.hpp.

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

Get the current symmetry of the operator.

Definition at line 298 of file BelosLinearProblem.hpp.

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

Get information on whether the linear system is being preconditioned on the left.

Definition at line 301 of file BelosLinearProblem.hpp.

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

Get information on whether the linear system is being preconditioned on the right.

Definition at line 304 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 720 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 application is only of the linear problem operator, no preconditioners are applied. Flexible variants of Krylov methods will use this application method within their code.

Precondition:

Definition at line 781 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 application is only of the left preconditioner, which may be required for flexible variants of Krylov methods.

Note:
This will return Undefined if the left preconditioner is not defined for this operator.

Definition at line 793 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 application is only of the right preconditioner, which may be required for flexible variants of Krylov methods.

Note:
This will return Undefined if the right preconditioner is not defined for this operator.

Definition at line 805 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 884 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 817 of file BelosLinearProblem.hpp.


Member Data Documentation

template<class ScalarType, class MV, class OP>
RCP<const OP> Belos::LinearProblem< ScalarType, MV, OP >::A_ [private]

Operator of linear system.

Definition at line 364 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
RCP<MV> Belos::LinearProblem< ScalarType, MV, OP >::X_ [private]

Solution std::vector of linear system.

Definition at line 367 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
RCP<MV> Belos::LinearProblem< ScalarType, MV, OP >::curX_ [private]

Current solution std::vector of the linear system.

Definition at line 370 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
RCP<const MV> Belos::LinearProblem< ScalarType, MV, OP >::B_ [private]

Right-hand side of linear system.

Definition at line 373 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
RCP<MV> Belos::LinearProblem< ScalarType, MV, OP >::curB_ [private]

Current right-hand side of the linear system.

Definition at line 376 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
RCP<MV> Belos::LinearProblem< ScalarType, MV, OP >::R0_ [private]

Initial residual of the linear system.

Definition at line 379 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
RCP<MV> Belos::LinearProblem< ScalarType, MV, OP >::PR0_ [private]

Preconditioned initial residual of the linear system.

Definition at line 382 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
RCP<const OP> Belos::LinearProblem< ScalarType, MV, OP >::LP_ [private]

Left preconditioning operator of linear system.

Definition at line 385 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
RCP<const OP> Belos::LinearProblem< ScalarType, MV, OP >::RP_ [private]

Right preconditioning operator of linear system.

Definition at line 388 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
Teuchos::RCP<Teuchos::Time> Belos::LinearProblem< ScalarType, MV, OP >::timerOp_ [mutable, private]

Timers.

Definition at line 391 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
Teuchos::RCP<Teuchos::Time> Belos::LinearProblem< ScalarType, MV, OP >::timerPrec_ [mutable, private]

Definition at line 391 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
int Belos::LinearProblem< ScalarType, MV, OP >::blocksize_ [private]

Current block size of linear system.

Definition at line 394 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
int Belos::LinearProblem< ScalarType, MV, OP >::num2Solve_ [private]

Number of linear systems that are currently being solver for ( <= blocksize_ )

Definition at line 397 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
std::vector<int> Belos::LinearProblem< ScalarType, MV, OP >::rhsIndex_ [private]

Indices of current linear systems being solver for.

Definition at line 400 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
int Belos::LinearProblem< ScalarType, MV, OP >::lsNum_ [private]

Number of linear systems that have been loaded in this linear problem object.

Definition at line 403 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
bool Belos::LinearProblem< ScalarType, MV, OP >::Left_Scale_ [private]

Booleans to keep track of linear problem attributes/status.

Definition at line 406 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
bool Belos::LinearProblem< ScalarType, MV, OP >::Right_Scale_ [private]

Definition at line 407 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
bool Belos::LinearProblem< ScalarType, MV, OP >::isSet_ [private]

Definition at line 408 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
bool Belos::LinearProblem< ScalarType, MV, OP >::isHermitian_ [private]

Definition at line 409 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
bool Belos::LinearProblem< ScalarType, MV, OP >::solutionUpdated_ [private]

Definition at line 410 of file BelosLinearProblem.hpp.

template<class ScalarType, class MV, class OP>
std::string Belos::LinearProblem< ScalarType, MV, OP >::label_ [private]

Linear problem label that prefixes the timer labels.

Definition at line 413 of file BelosLinearProblem.hpp.


The documentation for this class was generated from the following file:
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines