Belos
Version of the Day
|
This class implements the GCRODR iteration, where a single-std::vector Krylov subspace is constructed. The QR decomposition of block, upper Hessenberg matrix is performed each iteration to update the least squares system and give the current linear system residuals. More...
#include <BelosGCRODRIter.hpp>
Public Types | |
typedef MultiVecTraits < ScalarType, MV > | MVT |
typedef OperatorTraits < ScalarType, MV, OP > | OPT |
typedef Teuchos::ScalarTraits < ScalarType > | SCT |
typedef SCT::magnitudeType | MagnitudeType |
Public Member Functions | |
Constructors/Destructor | |
GCRODRIter (const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem, const Teuchos::RCP< OutputManager< ScalarType > > &printer, const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &tester, const Teuchos::RCP< MatOrthoManager< ScalarType, MV, OP > > &ortho, Teuchos::ParameterList ¶ms) | |
GCRODRIter constructor with linear problem, solver utilities, and parameter list of solver options. More... | |
virtual | ~GCRODRIter () |
Destructor. More... | |
Solver methods | |
void | iterate () |
This method performs block Gmres iterations until the status test indicates the need to stop or an error occurs (in which case, an std::exception is thrown). More... | |
void | initialize (GCRODRIterState< ScalarType, MV > &newstate) |
Initialize the solver to an iterate, providing a complete state. More... | |
void | initialize () |
Initialize the solver with empty data. Calling this method will result in error, as GCRODRIter must be initialized with a valid state. More... | |
GCRODRIterState< ScalarType, MV > | getState () const |
Get the current state of the linear solver. More... | |
Status methods | |
int | getNumIters () const |
Get the current iteration count. More... | |
void | resetNumIters (int iter=0) |
Reset the iteration count. More... | |
Teuchos::RCP< const MV > | getNativeResiduals (std::vector< MagnitudeType > *norms) const |
Get the norms of the residuals native to the solver. More... | |
Teuchos::RCP< MV > | getCurrentUpdate () const |
Get the current update to the linear system. More... | |
void | updateLSQR (int dim=-1) |
Method for updating QR factorization of upper Hessenberg matrix. More... | |
int | getCurSubspaceDim () const |
Get the dimension of the search subspace used to generate the current solution to the linear problem. More... | |
int | getMaxSubspaceDim () const |
Get the maximum dimension allocated for the search subspace. More... | |
Accessor methods | |
const LinearProblem < ScalarType, MV, OP > & | getProblem () const |
Get a constant reference to the linear problem. More... | |
int | getNumBlocks () const |
Get the maximum number of blocks used by the iterative solver in solving this linear problem. More... | |
void | setNumBlocks (int numBlocks) |
Set the maximum number of blocks used by the iterative solver. More... | |
int | getRecycledBlocks () const |
Get the maximum number of recycled blocks used by the iterative solver in solving this linear problem. More... | |
void | setRecycledBlocks (int recycledBlocks) |
Set the maximum number of recycled blocks used by the iterative solver. More... | |
int | getBlockSize () const |
Get the blocksize to be used by the iterative solver in solving this linear problem. More... | |
void | setBlockSize (int blockSize) |
Set the blocksize. More... | |
void | setSize (int recycledBlocks, int numBlocks) |
Set the maximum number of blocks used by the iterative solver and the number of recycled vectors. More... | |
bool | isInitialized () |
States whether the solver has been initialized or not. More... | |
Public Member Functions inherited from Belos::Iteration< ScalarType, MV, OP > | |
Iteration () | |
Default Constructor. More... | |
virtual | ~Iteration () |
Destructor. More... | |
virtual Teuchos::RCP< const MV > | getNativeResiduals (std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > *norms) const =0 |
Get the residuals native to the solver. More... | |
This class implements the GCRODR iteration, where a single-std::vector Krylov subspace is constructed. The QR decomposition of block, upper Hessenberg matrix is performed each iteration to update the least squares system and give the current linear system residuals.
Definition at line 154 of file BelosGCRODRIter.hpp.
typedef MultiVecTraits<ScalarType,MV> Belos::GCRODRIter< ScalarType, MV, OP >::MVT |
Definition at line 161 of file BelosGCRODRIter.hpp.
typedef OperatorTraits<ScalarType,MV,OP> Belos::GCRODRIter< ScalarType, MV, OP >::OPT |
Definition at line 162 of file BelosGCRODRIter.hpp.
typedef Teuchos::ScalarTraits<ScalarType> Belos::GCRODRIter< ScalarType, MV, OP >::SCT |
Definition at line 163 of file BelosGCRODRIter.hpp.
typedef SCT::magnitudeType Belos::GCRODRIter< ScalarType, MV, OP >::MagnitudeType |
Definition at line 164 of file BelosGCRODRIter.hpp.
Belos::GCRODRIter< ScalarType, MV, OP >::GCRODRIter | ( | const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > & | problem, |
const Teuchos::RCP< OutputManager< ScalarType > > & | printer, | ||
const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > & | tester, | ||
const Teuchos::RCP< MatOrthoManager< ScalarType, MV, OP > > & | ortho, | ||
Teuchos::ParameterList & | params | ||
) |
GCRODRIter constructor with linear problem, solver utilities, and parameter list of solver options.
This constructor takes pointers required by the linear solver, in addition to a parameter list of options for the linear solver. These options include the following:
int
specifying the maximum number of blocks allocated for the solver basis. Default: 25bool
specifying whether the timers should be restarted each time iterate() is called. Default: falsebool
specifying whether the upper Hessenberg should be stored separately from the least squares system. Default: false Definition at line 413 of file BelosGCRODRIter.hpp.
|
inlinevirtual |
Destructor.
Definition at line 184 of file BelosGCRODRIter.hpp.
|
virtual |
This method performs block Gmres iterations until the status test indicates the need to stop or an error occurs (in which case, an std::exception is thrown).
iterate() will first determine whether the solver is inintialized; if not, it will call initialize() using default arguments. After initialization, the solver performs block Gmres iterations until the status test evaluates as Passed, at which point the method returns to the caller.
The block Gmres iteration proceeds as follows:
blockSize
vectors in the Krylov basis.The status test is queried at the beginning of the iteration.
Possible exceptions thrown include the GCRODRIterOrthoFailure.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 544 of file BelosGCRODRIter.hpp.
void Belos::GCRODRIter< ScalarType, MV, OP >::initialize | ( | GCRODRIterState< ScalarType, MV > & | newstate | ) |
Initialize the solver to an iterate, providing a complete state.
The GCRODRIter contains a certain amount of state, consisting of the current Krylov basis and the associated Hessenberg matrix.
initialize() gives the user the opportunity to manually set these, although this must be done with caution, abiding by the rules given below. All notions of orthogonality and orthonormality are derived from the inner product specified by the orthogonalization manager.
true
(see post-conditions of isInitialize())The user has the option of specifying any component of the state using initialize(). However, these arguments are assumed to match the post-conditions specified under isInitialized(). Any necessary component of the state not given to initialize() will be generated.
newstate
which directly points to the multivectors in the solver, the data is not copied. Definition at line 520 of file BelosGCRODRIter.hpp.
|
inlinevirtual |
Initialize the solver with empty data. Calling this method will result in error, as GCRODRIter must be initialized with a valid state.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 240 of file BelosGCRODRIter.hpp.
|
inline |
Get the current state of the linear solver.
The data is only valid if isInitialized() == true
.
Definition at line 252 of file BelosGCRODRIter.hpp.
|
inlinevirtual |
Get the current iteration count.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 270 of file BelosGCRODRIter.hpp.
|
inlinevirtual |
Reset the iteration count.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 273 of file BelosGCRODRIter.hpp.
Teuchos::RCP< const MV > Belos::GCRODRIter< ScalarType, MV, OP >::getNativeResiduals | ( | std::vector< MagnitudeType > * | norms | ) | const |
Get the norms of the residuals native to the solver.
Definition at line 501 of file BelosGCRODRIter.hpp.
|
virtual |
Get the current update to the linear system.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 453 of file BelosGCRODRIter.hpp.
void Belos::GCRODRIter< ScalarType, MV, OP >::updateLSQR | ( | int | dim = -1 | ) |
Method for updating QR factorization of upper Hessenberg matrix.
dim
>= getCurSubspaceDim()
and dim
< getMaxSubspaceDim()
, then the dim-th
equations of the least squares problem will be updated. Definition at line 694 of file BelosGCRODRIter.hpp.
|
inline |
Get the dimension of the search subspace used to generate the current solution to the linear problem.
Definition at line 294 of file BelosGCRODRIter.hpp.
|
inline |
Get the maximum dimension allocated for the search subspace.
Definition at line 300 of file BelosGCRODRIter.hpp.
|
inlinevirtual |
Get a constant reference to the linear problem.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 309 of file BelosGCRODRIter.hpp.
|
inline |
Get the maximum number of blocks used by the iterative solver in solving this linear problem.
Definition at line 312 of file BelosGCRODRIter.hpp.
|
inline |
Set the maximum number of blocks used by the iterative solver.
Definition at line 315 of file BelosGCRODRIter.hpp.
|
inline |
Get the maximum number of recycled blocks used by the iterative solver in solving this linear problem.
Definition at line 318 of file BelosGCRODRIter.hpp.
|
inline |
Set the maximum number of recycled blocks used by the iterative solver.
Definition at line 321 of file BelosGCRODRIter.hpp.
|
inlinevirtual |
Get the blocksize to be used by the iterative solver in solving this linear problem.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 324 of file BelosGCRODRIter.hpp.
|
inlinevirtual |
Set the blocksize.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 327 of file BelosGCRODRIter.hpp.
|
inline |
Set the maximum number of blocks used by the iterative solver and the number of recycled vectors.
Definition at line 332 of file BelosGCRODRIter.hpp.
|
inlinevirtual |
States whether the solver has been initialized or not.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 345 of file BelosGCRODRIter.hpp.