Belos
Version of the Day
|
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 . More... | |
void | setLHS (const Teuchos::RCP< MV > &X) |
Set left-hand-side X of linear problem . More... | |
void | setRHS (const Teuchos::RCP< const MV > &B) |
Set right-hand-side B of linear problem . More... | |
void | setInitResVec (const Teuchos::RCP< const MV > &R0) |
Set the user-defined residual of linear problem . More... | |
void | setInitPrecResVec (const Teuchos::RCP< const MV > &PR0) |
Set the user-defined preconditioned residual of linear problem . More... | |
void | setLeftPrec (const Teuchos::RCP< const OP > &LP) |
Set left preconditioner (LP ) of linear problem . More... | |
void | setRightPrec (const Teuchos::RCP< const OP > &RP) |
Set right preconditioner (RP ) of linear problem . More... | |
virtual void | setCurrLS () |
Tell the linear problem that the solver is finished with the current linear system. More... | |
virtual 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... | |
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. More... | |
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. More... | |
Set / Reset method | |
virtual 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 | |
virtual void | apply (const MV &x, MV &y) const |
Apply the composite operator of this linear problem to x , returning y . More... | |
virtual void | applyOp (const MV &x, MV &y) const |
Apply ONLY the operator to x , returning y . More... | |
virtual void | applyLeftPrec (const MV &x, MV &y) const |
Apply ONLY the left preconditioner to x , returning y . More... | |
virtual void | applyRightPrec (const MV &x, MV &y) const |
Apply ONLY the right preconditioner to x , returning y . More... | |
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 . More... | |
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 . More... | |
Protected Attributes | |
Teuchos::RCP< const OP > | A_ |
Operator of linear system. More... | |
Teuchos::RCP< MV > | X_ |
Solution vector of linear system. More... | |
Teuchos::RCP< MV > | curX_ |
Current solution vector of the linear system. More... | |
Teuchos::RCP< const MV > | B_ |
Right-hand side of linear system. More... | |
Teuchos::RCP< const MV > | curB_ |
Current right-hand side of the linear system. More... | |
Teuchos::RCP< MV > | R0_ |
Initial residual of the linear system. More... | |
Teuchos::RCP< MV > | PR0_ |
Preconditioned initial residual of the linear system. More... | |
Teuchos::RCP< const MV > | R0_user_ |
User-defined initial residual of the linear system. More... | |
Teuchos::RCP< const MV > | PR0_user_ |
User-defined preconditioned initial residual of the linear system. More... | |
Teuchos::RCP< const OP > | LP_ |
Left preconditioning operator of linear system. More... | |
Teuchos::RCP< const OP > | RP_ |
Right preconditioning operator of linear system. More... | |
Teuchos::RCP< Teuchos::Time > | timerOp_ |
Timers. More... | |
Teuchos::RCP< Teuchos::Time > | timerPrec_ |
int | blocksize_ |
Current block size of linear system. More... | |
int | num2Solve_ |
Number of linear systems that are currently being solver for ( <= blocksize_ ) More... | |
std::vector< int > | rhsIndex_ |
Indices of current linear systems being solver for. More... | |
int | lsNum_ |
Number of linear systems that have been loaded in this linear problem object. More... | |
std::string | label_ |
Linear problem label that prefixes the timer labels. More... | |
Booleans to keep track of linear problem attributes and status. | |
bool | isSet_ |
Has the linear problem to solve been set? More... | |
bool | isHermitian_ |
Whether the operator A is symmetric (in real arithmetic, or Hermitian in complex arithmetic). More... | |
bool | solutionUpdated_ |
Has the current approximate solution been updated? More... | |
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.
ScalarType | The type of the entries in the matrix and vectors. |
MV | The (multi)vector type. |
OP | The operator type. (Operators are functions that take a multivector as input and compute a multivector as output.) |
Definition at line 61 of file BelosIteration.hpp.
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 586 of file BelosLinearProblem.hpp.
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 599 of file BelosLinearProblem.hpp.
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 618 of file BelosLinearProblem.hpp.
|
inlinevirtual |
Destructor (declared virtual for memory safety of derived classes).
Definition at line 112 of file BelosLinearProblem.hpp.
|
inline |
Set the operator A of the linear problem .
The operator is set by pointer; no copy of the operator is made.
Definition at line 122 of file BelosLinearProblem.hpp.
|
inline |
Set left-hand-side X of linear problem .
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.
|
inline |
Set right-hand-side B of linear problem .
The multivector is set by pointer; no copy of the object is made.
Definition at line 142 of file BelosLinearProblem.hpp.
|
inline |
Set the user-defined residual of linear problem .
The multivector is set by pointer; no copy of the object is made.
Definition at line 152 of file BelosLinearProblem.hpp.
|
inline |
Set the user-defined preconditioned residual of linear problem .
The multivector is set by pointer; no copy of the object is made.
Definition at line 162 of file BelosLinearProblem.hpp.
|
inline |
Set left preconditioner (LP
) of linear problem .
The operator is set by pointer; no copy of the operator is made.
Definition at line 170 of file BelosLinearProblem.hpp.
|
inline |
Set right preconditioner (RP
) of linear problem .
The operator is set by pointer; no copy of the operator is made.
Definition at line 175 of file BelosLinearProblem.hpp.
|
virtual |
Tell the linear problem that the solver is finished with the current linear system.
Definition at line 711 of file BelosLinearProblem.hpp.
|
virtual |
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 644 of file BelosLinearProblem.hpp.
|
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.
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 807 of file BelosLinearProblem.hpp.
|
virtual |
Compute the new solution to the linear system using the given update vector.
Let be the update vector, the scale factor, and the current solution. If there is a right preconditioner , then we compute the new solution as . Otherwise, if there is no right preconditioner, we compute the new solution as .
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.
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 by which to multiply the solution update vector when computing the update. This is ignored if the update vector is null. |
Definition at line 747 of file BelosLinearProblem.hpp.
|
inlinevirtual |
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.
update | [in/out] The solution update vector. If null, this method returns a pointer to the new solution. |
scale | [in] The factor by which to multiply the solution update vector when computing the update. This is ignored if the update vector is null. |
Definition at line 279 of file BelosLinearProblem.hpp.
|
virtual |
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.
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). |
Definition at line 830 of file BelosLinearProblem.hpp.
|
inline |
A pointer to the (unpreconditioned) operator A.
Definition at line 324 of file BelosLinearProblem.hpp.
|
inline |
A pointer to the left-hand side X.
Definition at line 327 of file BelosLinearProblem.hpp.
|
inline |
A pointer to the right-hand side B.
Definition at line 330 of file BelosLinearProblem.hpp.
Teuchos::RCP< const MV > Belos::LinearProblem< ScalarType, MV, OP >::getInitResVec | ( | ) | const |
A pointer to the initial unpreconditioned residual vector.
Definition at line 932 of file BelosLinearProblem.hpp.
Teuchos::RCP< const MV > Belos::LinearProblem< ScalarType, MV, OP >::getInitPrecResVec | ( | ) | const |
A pointer to the preconditioned initial residual vector.
Definition at line 941 of file BelosLinearProblem.hpp.
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.
isSolutionUpdated()
).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 950 of file BelosLinearProblem.hpp.
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.
isSolutionUpdated()
).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 961 of file BelosLinearProblem.hpp.
|
inline |
Get a pointer to the left preconditioner.
Definition at line 374 of file BelosLinearProblem.hpp.
|
inline |
Get a pointer to the right preconditioner.
Definition at line 377 of file BelosLinearProblem.hpp.
|
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.
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 400 of file BelosLinearProblem.hpp.
|
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 408 of file BelosLinearProblem.hpp.
|
inline |
The timers for this object.
The timers are ordered as follows:
Definition at line 416 of file BelosLinearProblem.hpp.
|
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 433 of file BelosLinearProblem.hpp.
|
inline |
Whether the problem has been set.
Definition at line 436 of file BelosLinearProblem.hpp.
|
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 443 of file BelosLinearProblem.hpp.
|
inline |
Whether the linear system is being preconditioned on the left.
Definition at line 446 of file BelosLinearProblem.hpp.
|
inline |
Whether the linear system is being preconditioned on the right.
Definition at line 449 of file BelosLinearProblem.hpp.
|
virtual |
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:
getOperator().get()!=NULL
Definition at line 972 of file BelosLinearProblem.hpp.
|
virtual |
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 1020 of file BelosLinearProblem.hpp.
|
virtual |
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 1034 of file BelosLinearProblem.hpp.
|
virtual |
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 1048 of file BelosLinearProblem.hpp.
|
virtual |
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.
Definition at line 1111 of file BelosLinearProblem.hpp.
|
virtual |
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.
Definition at line 1062 of file BelosLinearProblem.hpp.
|
protected |
Operator of linear system.
Definition at line 510 of file BelosLinearProblem.hpp.
|
protected |
Solution vector of linear system.
Definition at line 513 of file BelosLinearProblem.hpp.
|
protected |
Current solution vector of the linear system.
Definition at line 516 of file BelosLinearProblem.hpp.
|
protected |
Right-hand side of linear system.
Definition at line 519 of file BelosLinearProblem.hpp.
|
protected |
Current right-hand side of the linear system.
Definition at line 522 of file BelosLinearProblem.hpp.
|
protected |
Initial residual of the linear system.
Definition at line 525 of file BelosLinearProblem.hpp.
|
protected |
Preconditioned initial residual of the linear system.
Definition at line 528 of file BelosLinearProblem.hpp.
|
protected |
User-defined initial residual of the linear system.
Definition at line 531 of file BelosLinearProblem.hpp.
|
protected |
User-defined preconditioned initial residual of the linear system.
Definition at line 534 of file BelosLinearProblem.hpp.
|
protected |
Left preconditioning operator of linear system.
Definition at line 537 of file BelosLinearProblem.hpp.
|
protected |
Right preconditioning operator of linear system.
Definition at line 540 of file BelosLinearProblem.hpp.
|
mutableprotected |
Timers.
Definition at line 543 of file BelosLinearProblem.hpp.
|
mutableprotected |
Definition at line 543 of file BelosLinearProblem.hpp.
|
protected |
Current block size of linear system.
Definition at line 546 of file BelosLinearProblem.hpp.
|
protected |
Number of linear systems that are currently being solver for ( <= blocksize_ )
Definition at line 549 of file BelosLinearProblem.hpp.
|
protected |
Indices of current linear systems being solver for.
Definition at line 552 of file BelosLinearProblem.hpp.
|
protected |
Number of linear systems that have been loaded in this linear problem object.
Definition at line 555 of file BelosLinearProblem.hpp.
|
protected |
Has the linear problem to solve been set?
Definition at line 561 of file BelosLinearProblem.hpp.
|
protected |
Whether the operator A is symmetric (in real arithmetic, or Hermitian in complex arithmetic).
Definition at line 565 of file BelosLinearProblem.hpp.
|
protected |
Has the current approximate solution been updated?
Definition at line 568 of file BelosLinearProblem.hpp.
|
protected |
Linear problem label that prefixes the timer labels.
Definition at line 573 of file BelosLinearProblem.hpp.