Belos Package Browser (Single Doxygen Collection)
Development
|
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 OperatorTraits < ScalarType, MV, OP > | OPT |
typedef Teuchos::ScalarTraits < ScalarType > | SCT |
typedef SCT::magnitudeType | MagnitudeType |
Private Member Functions | |
void | setStateSize () |
Method for initalizing the state storage needed by block flexible GMRES. More... | |
Private Attributes | |
const Teuchos::RCP < LinearProblem< ScalarType, MV, OP > > | lp_ |
const Teuchos::RCP < OutputManager< ScalarType > > | om_ |
const Teuchos::RCP< StatusTest < ScalarType, MV, OP > > | stest_ |
const Teuchos::RCP < OrthoManager< ScalarType, MV > > | ortho_ |
int | blockSize_ |
int | numBlocks_ |
Teuchos::SerialDenseVector < int, ScalarType > | beta |
Teuchos::SerialDenseVector < int, ScalarType > | sn |
Teuchos::SerialDenseVector < int, MagnitudeType > | cs |
bool | initialized_ |
bool | stateStorageInitialized_ |
int | curDim_ |
int | iter_ |
Teuchos::RCP< MV > | V_ |
Teuchos::RCP< MV > | Z_ |
Teuchos::RCP < Teuchos::SerialDenseMatrix < int, ScalarType > > | H_ |
Teuchos::RCP < Teuchos::SerialDenseMatrix < int, ScalarType > > | R_ |
Teuchos::RCP < Teuchos::SerialDenseMatrix < int, ScalarType > > | z_ |
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. More... | |
virtual | ~BlockFGmresIter () |
Destructor. More... | |
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). More... | |
void | initializeGmres (GmresIterationState< ScalarType, MV > &newstate) |
Initialize the solver to an iterate, providing a complete state. More... | |
void | initialize () |
Initialize the solver with the initial vectors from the linear problem or random data. More... | |
GmresIterationState < 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 | 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... | |
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... | |
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. More... | |
bool | isInitialized () |
States whether the solver has been initialized or not. More... | |
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 51 of file BelosBlockFGmresIter.hpp.
typedef MultiVecTraits<ScalarType,MV> Belos::BlockFGmresIter< ScalarType, MV, OP >::MVT |
Definition at line 58 of file BelosBlockFGmresIter.hpp.
typedef OperatorTraits<ScalarType,MV,OP> Belos::BlockFGmresIter< ScalarType, MV, OP >::OPT |
Definition at line 59 of file BelosBlockFGmresIter.hpp.
typedef Teuchos::ScalarTraits<ScalarType> Belos::BlockFGmresIter< ScalarType, MV, OP >::SCT |
Definition at line 60 of file BelosBlockFGmresIter.hpp.
typedef SCT::magnitudeType Belos::BlockFGmresIter< ScalarType, MV, OP >::MagnitudeType |
Definition at line 61 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 302 of file BelosBlockFGmresIter.hpp.
|
inlinevirtual |
Destructor.
Definition at line 82 of file BelosBlockFGmresIter.hpp.
|
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 613 of file BelosBlockFGmresIter.hpp.
|
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 517 of file BelosBlockFGmresIter.hpp.
|
inlinevirtual |
Initialize the solver with the initial vectors from the linear problem or random data.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 138 of file BelosBlockFGmresIter.hpp.
|
inlinevirtual |
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 151 of file BelosBlockFGmresIter.hpp.
|
inlinevirtual |
Get the current iteration count.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 169 of file BelosBlockFGmresIter.hpp.
|
inlinevirtual |
Reset the iteration count.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 172 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 496 of file BelosBlockFGmresIter.hpp.
|
virtual |
Get the current update to the linear system.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 456 of file BelosBlockFGmresIter.hpp.
|
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 701 of file BelosBlockFGmresIter.hpp.
|
inlinevirtual |
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 193 of file BelosBlockFGmresIter.hpp.
|
inlinevirtual |
Get the maximum dimension allocated for the search subspace.
Implements Belos::GmresIteration< ScalarType, MV, OP >.
Definition at line 199 of file BelosBlockFGmresIter.hpp.
|
inlinevirtual |
Get a constant reference to the linear problem.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 208 of file BelosBlockFGmresIter.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 211 of file BelosBlockFGmresIter.hpp.
|
inlinevirtual |
Set the blocksize.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 214 of file BelosBlockFGmresIter.hpp.
|
inline |
Get the maximum number of blocks used by the iterative solver in solving this linear problem.
Definition at line 217 of file BelosBlockFGmresIter.hpp.
|
inline |
Set the maximum number of blocks used by the iterative solver.
Definition at line 220 of file BelosBlockFGmresIter.hpp.
|
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 332 of file BelosBlockFGmresIter.hpp.
|
inlinevirtual |
States whether the solver has been initialized or not.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 231 of file BelosBlockFGmresIter.hpp.
|
private |
Method for initalizing the state storage needed by block flexible GMRES.
Definition at line 360 of file BelosBlockFGmresIter.hpp.
|
private |
Definition at line 246 of file BelosBlockFGmresIter.hpp.
|
private |
Definition at line 247 of file BelosBlockFGmresIter.hpp.
|
private |
Definition at line 248 of file BelosBlockFGmresIter.hpp.
|
private |
Definition at line 249 of file BelosBlockFGmresIter.hpp.
|
private |
Definition at line 256 of file BelosBlockFGmresIter.hpp.
|
private |
Definition at line 258 of file BelosBlockFGmresIter.hpp.
|
private |
Definition at line 261 of file BelosBlockFGmresIter.hpp.
|
private |
Definition at line 261 of file BelosBlockFGmresIter.hpp.
|
private |
Definition at line 262 of file BelosBlockFGmresIter.hpp.
|
private |
Definition at line 270 of file BelosBlockFGmresIter.hpp.
|
private |
Definition at line 275 of file BelosBlockFGmresIter.hpp.
|
private |
Definition at line 278 of file BelosBlockFGmresIter.hpp.
|
private |
Definition at line 278 of file BelosBlockFGmresIter.hpp.
|
private |
Definition at line 283 of file BelosBlockFGmresIter.hpp.
|
private |
Definition at line 284 of file BelosBlockFGmresIter.hpp.
|
private |
Definition at line 289 of file BelosBlockFGmresIter.hpp.
|
private |
Definition at line 294 of file BelosBlockFGmresIter.hpp.
|
private |
Definition at line 295 of file BelosBlockFGmresIter.hpp.