|
Belos
Version of the Day
|
This class implements the block flexible GMRES iteration, where a block 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 <BelosBlockFGmresIter.hpp>

Public Types | |
| typedef MultiVecTraits < ScalarType, MV > | MVT |
| typedef MultiVecTraitsExt < ScalarType, MV > | MVText |
| typedef OperatorTraits < ScalarType, MV, OP > | OPT |
| typedef Teuchos::ScalarTraits < ScalarType > | SCT |
| typedef SCT::magnitudeType | MagnitudeType |
Public Member Functions | |
Constructors/Destructor | |
| BlockFGmresIter (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) | |
| BlockFGmresIter constructor with linear problem, solver utilities, and parameter list of solver options. | |
| virtual | ~BlockFGmresIter () |
| Destructor. | |
Solver methods | |
| void | iterate () |
| This method performs block FGmres iterations until the status test indicates the need to stop or an error occurs (in which case, an std::exception is thrown). | |
| void | initializeGmres (GmresIterationState< ScalarType, MV > newstate) |
| Initialize the solver to an iterate, providing a complete state. | |
| void | initialize () |
| Initialize the solver with the initial vectors from the linear problem or random data. | |
| GmresIterationState < ScalarType, MV > | getState () const |
| Get the current state of the linear solver. | |
Status methods | |
| int | getNumIters () const |
| Get the current iteration count. | |
| void | resetNumIters (int iter=0) |
| Reset the iteration count. | |
| Teuchos::RCP< const MV > | getNativeResiduals (std::vector< MagnitudeType > *norms) const |
| Get the norms of the residuals native to the solver. | |
| Teuchos::RCP< MV > | getCurrentUpdate () const |
| Get the current update to the linear system. | |
| void | updateLSQR (int dim=-1) |
| Method for updating QR factorization of upper Hessenberg matrix. | |
| int | getCurSubspaceDim () const |
| Get the dimension of the search subspace used to generate the current solution to the linear problem. | |
| int | getMaxSubspaceDim () const |
| Get the maximum dimension allocated for the search subspace. | |
Accessor methods | |
| const LinearProblem < ScalarType, MV, OP > & | getProblem () const |
| Get a constant reference to the linear problem. | |
| int | getBlockSize () const |
| Get the blocksize to be used by the iterative solver in solving this linear problem. | |
| void | setBlockSize (int blockSize) |
| Set the blocksize. | |
| int | getNumBlocks () const |
| Get the maximum number of blocks used by the iterative solver in solving this linear problem. | |
| void | setNumBlocks (int numBlocks) |
| Set the maximum number of blocks used by the iterative solver. | |
| void | setSize (int blockSize, int numBlocks) |
| Set the blocksize and number of blocks to be used by the iterative solver in solving this linear problem. | |
| bool | isInitialized () |
| States whether the solver has been initialized or not. | |
This class implements the block flexible GMRES iteration, where a block 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 83 of file BelosBlockFGmresIter.hpp.
| typedef MultiVecTraits<ScalarType,MV> Belos::BlockFGmresIter< ScalarType, MV, OP >::MVT |
Definition at line 90 of file BelosBlockFGmresIter.hpp.
| typedef MultiVecTraitsExt<ScalarType,MV> Belos::BlockFGmresIter< ScalarType, MV, OP >::MVText |
Definition at line 91 of file BelosBlockFGmresIter.hpp.
| typedef OperatorTraits<ScalarType,MV,OP> Belos::BlockFGmresIter< ScalarType, MV, OP >::OPT |
Definition at line 92 of file BelosBlockFGmresIter.hpp.
| typedef Teuchos::ScalarTraits<ScalarType> Belos::BlockFGmresIter< ScalarType, MV, OP >::SCT |
Definition at line 93 of file BelosBlockFGmresIter.hpp.
| typedef SCT::magnitudeType Belos::BlockFGmresIter< ScalarType, MV, OP >::MagnitudeType |
Definition at line 94 of file BelosBlockFGmresIter.hpp.
| Belos::BlockFGmresIter< ScalarType, MV, OP >::BlockFGmresIter | ( | 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 | ||
| ) |
BlockFGmresIter 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 block size used by the algorithm. This can also be specified using the setBlockSize() method. Default: 1int 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 335 of file BelosBlockFGmresIter.hpp.
| virtual Belos::BlockFGmresIter< ScalarType, MV, OP >::~BlockFGmresIter | ( | ) | [inline, virtual] |
Destructor.
Definition at line 115 of file BelosBlockFGmresIter.hpp.
| void Belos::BlockFGmresIter< ScalarType, MV, OP >::iterate | ( | ) | [virtual] |
This method performs block FGmres 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 FGmres iterations until the status test evaluates as Passed, at which point the method returns to the caller.
The block FGmres 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 GmresIterationOrthoFailure.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 646 of file BelosBlockFGmresIter.hpp.
| void Belos::BlockFGmresIter< ScalarType, MV, OP >::initializeGmres | ( | GmresIterationState< ScalarType, MV > | newstate | ) | [virtual] |
Initialize the solver to an iterate, providing a complete state.
The BlockFGmresIter 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. Implements Belos::GmresIteration< ScalarType, MV, OP >.
Definition at line 550 of file BelosBlockFGmresIter.hpp.
| void Belos::BlockFGmresIter< ScalarType, MV, OP >::initialize | ( | ) | [inline, virtual] |
Initialize the solver with the initial vectors from the linear problem or random data.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 171 of file BelosBlockFGmresIter.hpp.
| GmresIterationState<ScalarType,MV> Belos::BlockFGmresIter< ScalarType, MV, OP >::getState | ( | ) | const [inline, virtual] |
Get the current state of the linear solver.
The data is only valid if isInitialized() == true.
Implements Belos::GmresIteration< ScalarType, MV, OP >.
Definition at line 184 of file BelosBlockFGmresIter.hpp.
| int Belos::BlockFGmresIter< ScalarType, MV, OP >::getNumIters | ( | ) | const [inline, virtual] |
Get the current iteration count.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 202 of file BelosBlockFGmresIter.hpp.
| void Belos::BlockFGmresIter< ScalarType, MV, OP >::resetNumIters | ( | int | iter = 0 | ) | [inline, virtual] |
Reset the iteration count.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 205 of file BelosBlockFGmresIter.hpp.
| Teuchos::RCP< const MV > Belos::BlockFGmresIter< ScalarType, MV, OP >::getNativeResiduals | ( | std::vector< MagnitudeType > * | norms | ) | const |
Get the norms of the residuals native to the solver.
Definition at line 529 of file BelosBlockFGmresIter.hpp.
| Teuchos::RCP< MV > Belos::BlockFGmresIter< ScalarType, MV, OP >::getCurrentUpdate | ( | ) | const [virtual] |
Get the current update to the linear system.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 489 of file BelosBlockFGmresIter.hpp.
| void Belos::BlockFGmresIter< ScalarType, MV, OP >::updateLSQR | ( | int | dim = -1 | ) | [virtual] |
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. Implements Belos::GmresIteration< ScalarType, MV, OP >.
Definition at line 734 of file BelosBlockFGmresIter.hpp.
| int Belos::BlockFGmresIter< ScalarType, MV, OP >::getCurSubspaceDim | ( | ) | const [inline, virtual] |
Get the dimension of the search subspace used to generate the current solution to the linear problem.
Implements Belos::GmresIteration< ScalarType, MV, OP >.
Definition at line 226 of file BelosBlockFGmresIter.hpp.
| int Belos::BlockFGmresIter< ScalarType, MV, OP >::getMaxSubspaceDim | ( | ) | const [inline, virtual] |
Get the maximum dimension allocated for the search subspace.
Implements Belos::GmresIteration< ScalarType, MV, OP >.
Definition at line 232 of file BelosBlockFGmresIter.hpp.
| const LinearProblem<ScalarType,MV,OP>& Belos::BlockFGmresIter< ScalarType, MV, OP >::getProblem | ( | ) | const [inline, virtual] |
Get a constant reference to the linear problem.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 241 of file BelosBlockFGmresIter.hpp.
| int Belos::BlockFGmresIter< ScalarType, MV, OP >::getBlockSize | ( | ) | const [inline, virtual] |
Get the blocksize to be used by the iterative solver in solving this linear problem.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 244 of file BelosBlockFGmresIter.hpp.
| void Belos::BlockFGmresIter< ScalarType, MV, OP >::setBlockSize | ( | int | blockSize | ) | [inline, virtual] |
Set the blocksize.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 247 of file BelosBlockFGmresIter.hpp.
| int Belos::BlockFGmresIter< ScalarType, MV, OP >::getNumBlocks | ( | ) | const [inline] |
Get the maximum number of blocks used by the iterative solver in solving this linear problem.
Definition at line 250 of file BelosBlockFGmresIter.hpp.
| void Belos::BlockFGmresIter< ScalarType, MV, OP >::setNumBlocks | ( | int | numBlocks | ) | [inline] |
Set the maximum number of blocks used by the iterative solver.
Definition at line 253 of file BelosBlockFGmresIter.hpp.
| void Belos::BlockFGmresIter< ScalarType, MV, OP >::setSize | ( | int | blockSize, |
| int | numBlocks | ||
| ) | [virtual] |
Set the blocksize and number of blocks to be used by the iterative solver in solving this linear problem.
Changing either the block size or the number of blocks will reset the solver to an uninitialized state.
Implements Belos::GmresIteration< ScalarType, MV, OP >.
Definition at line 365 of file BelosBlockFGmresIter.hpp.
| bool Belos::BlockFGmresIter< ScalarType, MV, OP >::isInitialized | ( | ) | [inline, virtual] |
States whether the solver has been initialized or not.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 264 of file BelosBlockFGmresIter.hpp.
1.7.6.1