Belos Package Browser (Single Doxygen Collection)
Development
|
This class implements the pseudo-block GMRES iteration, where a block Krylov subspace is constructed for all of the linear systems simultaneously. The QR decomposition of each block, upper Hessenberg matrix is performed each iteration to update the least squares system and give the current linear system residuals. More...
#include <BelosPseudoBlockGmresIter.hpp>
Public Types | |
typedef MultiVecTraits < ScalarType, MV > | MVT |
typedef OperatorTraits < ScalarType, MV, OP > | OPT |
typedef Teuchos::ScalarTraits < ScalarType > | SCT |
typedef SCT::magnitudeType | MagnitudeType |
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 | numRHS_ |
int | numBlocks_ |
std::vector< Teuchos::RCP < Teuchos::SerialDenseVector < int, ScalarType > > > | sn_ |
std::vector< Teuchos::RCP < Teuchos::SerialDenseVector < int, MagnitudeType > > > | cs_ |
Teuchos::RCP< MV > | U_vec_ |
Teuchos::RCP< MV > | AU_vec_ |
Teuchos::RCP< MV > | cur_block_rhs_ |
Teuchos::RCP< MV > | cur_block_sol_ |
bool | initialized_ |
int | curDim_ |
int | iter_ |
std::vector< Teuchos::RCP< MV > > | V_ |
std::vector< Teuchos::RCP < Teuchos::SerialDenseMatrix < int, ScalarType > > > | H_ |
std::vector< Teuchos::RCP < Teuchos::SerialDenseMatrix < int, ScalarType > > > | R_ |
std::vector< Teuchos::RCP < Teuchos::SerialDenseVector < int, ScalarType > > > | Z_ |
Constructors/Destructor | |
PseudoBlockGmresIter (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) | |
PseudoBlockGmresIter constructor with linear problem, solver utilities, and parameter list of solver options. More... | |
virtual | ~PseudoBlockGmresIter () |
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 (const PseudoBlockGmresIterState< 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... | |
PseudoBlockGmresIterState < 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 "native" residual vectors. 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... | |
bool | isInitialized () |
States whether the solver has been initialized or not. More... | |
This class implements the pseudo-block GMRES iteration, where a block Krylov subspace is constructed for all of the linear systems simultaneously. The QR decomposition of each block, upper Hessenberg matrix is performed each iteration to update the least squares system and give the current linear system residuals.
Definition at line 155 of file BelosPseudoBlockGmresIter.hpp.
typedef MultiVecTraits<ScalarType,MV> Belos::PseudoBlockGmresIter< ScalarType, MV, OP >::MVT |
Definition at line 162 of file BelosPseudoBlockGmresIter.hpp.
typedef OperatorTraits<ScalarType,MV,OP> Belos::PseudoBlockGmresIter< ScalarType, MV, OP >::OPT |
Definition at line 163 of file BelosPseudoBlockGmresIter.hpp.
typedef Teuchos::ScalarTraits<ScalarType> Belos::PseudoBlockGmresIter< ScalarType, MV, OP >::SCT |
Definition at line 164 of file BelosPseudoBlockGmresIter.hpp.
typedef SCT::magnitudeType Belos::PseudoBlockGmresIter< ScalarType, MV, OP >::MagnitudeType |
Definition at line 165 of file BelosPseudoBlockGmresIter.hpp.
Belos::PseudoBlockGmresIter< ScalarType, MV, OP >::PseudoBlockGmresIter | ( | 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 | ||
) |
PseudoBlockGmresIter 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: false Definition at line 414 of file BelosPseudoBlockGmresIter.hpp.
|
inlinevirtual |
Destructor.
Definition at line 185 of file BelosPseudoBlockGmresIter.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 PseudoBlockGmresIterOrthoFailure.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 743 of file BelosPseudoBlockGmresIter.hpp.
void Belos::PseudoBlockGmresIter< ScalarType, MV, OP >::initialize | ( | const PseudoBlockGmresIterState< ScalarType, MV > & | newstate | ) |
Initialize the solver to an iterate, providing a complete state.
The PseudoBlockGmresIter 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 528 of file BelosPseudoBlockGmresIter.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 241 of file BelosPseudoBlockGmresIter.hpp.
|
inline |
Get the current state of the linear solver.
The data is only valid if isInitialized() == true
.
Definition at line 254 of file BelosPseudoBlockGmresIter.hpp.
|
inlinevirtual |
Get the current iteration count.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 279 of file BelosPseudoBlockGmresIter.hpp.
|
inlinevirtual |
Reset the iteration count.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 282 of file BelosPseudoBlockGmresIter.hpp.
Teuchos::RCP< const MV > Belos::PseudoBlockGmresIter< ScalarType, MV, OP >::getNativeResiduals | ( | std::vector< MagnitudeType > * | norms | ) | const |
Get the norms of the "native" residual vectors.
If norms != NULL, fill *norms with the native residual norms. There are numRHS_ of them. *norms will be resized if it has too few entries to hold the data.
For an explanation of "native" vs. "exact" (also known as "implicit" vs. "explicit") residuals, see the documentation of PseudoBlockGmresSolMgr::isLOADetected()
. In brief: "Native" residuals are cheaper to compute than "exact" residuals, but the two may differ, especially when using a left preconditioner.
PseudoBlockGmresSolMgr
knows that this method always returns null. Definition at line 503 of file BelosPseudoBlockGmresIter.hpp.
|
virtual |
Get the current update to the linear system.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 456 of file BelosPseudoBlockGmresIter.hpp.
void Belos::PseudoBlockGmresIter< 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 876 of file BelosPseudoBlockGmresIter.hpp.
|
inline |
Get the dimension of the search subspace used to generate the current solution to the linear problem.
Definition at line 318 of file BelosPseudoBlockGmresIter.hpp.
|
inline |
Get the maximum dimension allocated for the search subspace.
Definition at line 324 of file BelosPseudoBlockGmresIter.hpp.
|
inlinevirtual |
Get a constant reference to the linear problem.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 333 of file BelosPseudoBlockGmresIter.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 336 of file BelosPseudoBlockGmresIter.hpp.
|
inlinevirtual |
Set the blocksize.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 339 of file BelosPseudoBlockGmresIter.hpp.
|
inline |
Get the maximum number of blocks used by the iterative solver in solving this linear problem.
Definition at line 345 of file BelosPseudoBlockGmresIter.hpp.
void Belos::PseudoBlockGmresIter< ScalarType, MV, OP >::setNumBlocks | ( | int | numBlocks | ) |
Set the maximum number of blocks used by the iterative solver.
Definition at line 440 of file BelosPseudoBlockGmresIter.hpp.
|
inlinevirtual |
States whether the solver has been initialized or not.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 351 of file BelosPseudoBlockGmresIter.hpp.
|
private |
Definition at line 360 of file BelosPseudoBlockGmresIter.hpp.
|
private |
Definition at line 361 of file BelosPseudoBlockGmresIter.hpp.
|
private |
Definition at line 362 of file BelosPseudoBlockGmresIter.hpp.
|
private |
Definition at line 363 of file BelosPseudoBlockGmresIter.hpp.
|
private |
Definition at line 369 of file BelosPseudoBlockGmresIter.hpp.
|
private |
Definition at line 371 of file BelosPseudoBlockGmresIter.hpp.
|
private |
Definition at line 374 of file BelosPseudoBlockGmresIter.hpp.
|
private |
Definition at line 375 of file BelosPseudoBlockGmresIter.hpp.
|
private |
Definition at line 378 of file BelosPseudoBlockGmresIter.hpp.
|
private |
Definition at line 378 of file BelosPseudoBlockGmresIter.hpp.
|
private |
Definition at line 381 of file BelosPseudoBlockGmresIter.hpp.
|
private |
Definition at line 381 of file BelosPseudoBlockGmresIter.hpp.
|
private |
Definition at line 389 of file BelosPseudoBlockGmresIter.hpp.
|
private |
Definition at line 392 of file BelosPseudoBlockGmresIter.hpp.
|
private |
Definition at line 392 of file BelosPseudoBlockGmresIter.hpp.
|
private |
Definition at line 397 of file BelosPseudoBlockGmresIter.hpp.
|
private |
Definition at line 402 of file BelosPseudoBlockGmresIter.hpp.
|
private |
Definition at line 407 of file BelosPseudoBlockGmresIter.hpp.
|
private |
Definition at line 408 of file BelosPseudoBlockGmresIter.hpp.