Belos
Version of the Day
|
This class implements the preconditioned transpose-free QMR algorithm for solving non-Hermitian linear systems of equations Ax = b, where b is the right-hand side vector and x is the corresponding solution. More...
#include <BelosPseudoBlockTFQMRIter.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 | |
Constructor/Destructor. | |
PseudoBlockTFQMRIter (const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem, const Teuchos::RCP< OutputManager< ScalarType > > &printer, const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &tester, Teuchos::ParameterList ¶ms) | |
Belos::PseudoBlockTFQMRIter constructor. More... | |
virtual | ~PseudoBlockTFQMRIter () |
Belos::PseudoBlockTFQMRIter destructor. More... | |
Solver methods | |
void | iterate () |
This method performs pseudo-block TFQMR iterations until the status test indicates the need to stop or an error occurs (in which case, an std::exception is thrown). More... | |
void | initializeTFQMR (const PseudoBlockTFQMRIterState< 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... | |
PseudoBlockTFQMRIterState < 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... | |
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... | |
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 preconditioned transpose-free QMR algorithm for solving non-Hermitian linear systems of equations Ax = b, where b is the right-hand side vector and x is the corresponding solution.
Definition at line 147 of file BelosPseudoBlockTFQMRIter.hpp.
typedef MultiVecTraits<ScalarType,MV> Belos::PseudoBlockTFQMRIter< ScalarType, MV, OP >::MVT |
Definition at line 152 of file BelosPseudoBlockTFQMRIter.hpp.
typedef OperatorTraits<ScalarType,MV,OP> Belos::PseudoBlockTFQMRIter< ScalarType, MV, OP >::OPT |
Definition at line 153 of file BelosPseudoBlockTFQMRIter.hpp.
typedef Teuchos::ScalarTraits<ScalarType> Belos::PseudoBlockTFQMRIter< ScalarType, MV, OP >::SCT |
Definition at line 154 of file BelosPseudoBlockTFQMRIter.hpp.
typedef SCT::magnitudeType Belos::PseudoBlockTFQMRIter< ScalarType, MV, OP >::MagnitudeType |
Definition at line 155 of file BelosPseudoBlockTFQMRIter.hpp.
Belos::PseudoBlockTFQMRIter< ScalarType, MV, OP >::PseudoBlockTFQMRIter | ( | const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > & | problem, |
const Teuchos::RCP< OutputManager< ScalarType > > & | printer, | ||
const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > & | tester, | ||
Teuchos::ParameterList & | params | ||
) |
Belos::PseudoBlockTFQMRIter constructor.
Definition at line 343 of file BelosPseudoBlockTFQMRIter.hpp.
|
inlinevirtual |
Belos::PseudoBlockTFQMRIter destructor.
Definition at line 167 of file BelosPseudoBlockTFQMRIter.hpp.
|
virtual |
This method performs pseudo-block TFQMR 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 pseudo-block TFQMR iterations until the status test evaluates as Passed, at which point the method returns to the caller.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 464 of file BelosPseudoBlockTFQMRIter.hpp.
void Belos::PseudoBlockTFQMRIter< ScalarType, MV, OP >::initializeTFQMR | ( | const PseudoBlockTFQMRIterState< ScalarType, MV > & | newstate | ) |
Initialize the solver to an iterate, providing a complete state.
The PseudoBlockTFQMRIter 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 381 of file BelosPseudoBlockTFQMRIter.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 212 of file BelosPseudoBlockTFQMRIter.hpp.
|
inline |
Get the current state of the linear solver.
The data is only valid if isInitialized() == true
.
Definition at line 225 of file BelosPseudoBlockTFQMRIter.hpp.
|
inlinevirtual |
Get the current iteration count.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 253 of file BelosPseudoBlockTFQMRIter.hpp.
|
inlinevirtual |
Reset the iteration count.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 256 of file BelosPseudoBlockTFQMRIter.hpp.
Teuchos::RCP< const MV > Belos::PseudoBlockTFQMRIter< ScalarType, MV, OP >::getNativeResiduals | ( | std::vector< MagnitudeType > * | norms | ) | const |
Get the norms of the residuals native to the solver.
Definition at line 361 of file BelosPseudoBlockTFQMRIter.hpp.
|
inlinevirtual |
Get the current update to the linear system.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 266 of file BelosPseudoBlockTFQMRIter.hpp.
|
inlinevirtual |
Get a constant reference to the linear problem.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 275 of file BelosPseudoBlockTFQMRIter.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 278 of file BelosPseudoBlockTFQMRIter.hpp.
|
inlinevirtual |
Set the blocksize.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 281 of file BelosPseudoBlockTFQMRIter.hpp.
|
inlinevirtual |
States whether the solver has been initialized or not.
Implements Belos::Iteration< ScalarType, MV, OP >.
Definition at line 287 of file BelosPseudoBlockTFQMRIter.hpp.