18 #ifndef BELOS_TFQMR_ITER_HPP
19 #define BELOS_TFQMR_ITER_HPP
61 template <
class ScalarType,
class MV>
73 Rtilde(Teuchos::null),
D(Teuchos::null),
V(Teuchos::null)
94 template <
class ScalarType,
class MV,
class OP>
181 state.solnUpdate = solnUpdate_;
222 "Belos::TFQMRIter::setBlockSize(): Cannot use a block size that is not one.");
252 std::vector<ScalarType> alpha_, rho_, rho_old_;
253 std::vector<MagnitudeType> tau_, cs_, theta_;
266 bool stateStorageInitialized_;
290 template <
class ScalarType,
class MV,
class OP>
306 stateStorageInitialized_(false),
313 template <
class ScalarType,
class MV,
class OP>
321 return Teuchos::null;
327 template <
class ScalarType,
class MV,
class OP>
330 if (!stateStorageInitialized_) {
335 if (lhsMV == Teuchos::null && rhsMV == Teuchos::null) {
336 stateStorageInitialized_ =
false;
343 if (R_ == Teuchos::null) {
347 "Belos::TFQMRIter::setStateSize(): linear problem does not specify multivectors to clone from.");
348 R_ = MVT::Clone( *tmp, 1 );
349 D_ = MVT::Clone( *tmp, 1 );
350 V_ = MVT::Clone( *tmp, 1 );
351 solnUpdate_ = MVT::Clone( *tmp, 1 );
355 stateStorageInitialized_ =
true;
362 template <
class ScalarType,
class MV,
class OP>
366 if (!stateStorageInitialized_)
370 "Belos::TFQMRIter::initialize(): Cannot initialize state storage!");
374 std::string errstr(
"Belos::TFQMRIter::initialize(): Specified multivectors must have a consistent length and width.");
379 if (newstate.
R != Teuchos::null) {
382 std::invalid_argument, errstr );
384 std::invalid_argument, errstr );
387 if (newstate.
R != R_) {
389 MVT::Assign( *newstate.
R, *R_ );
395 W_ = MVT::CloneCopy( *R_ );
396 U_ = MVT::CloneCopy( *R_ );
397 Rtilde_ = MVT::CloneCopy( *R_ );
399 MVT::MvInit( *solnUpdate_ );
403 lp_->apply( *U_, *V_ );
404 AU_ = MVT::CloneCopy( *V_ );
409 MVT::MvNorm( *R_, tau_ );
410 MVT::MvDot( *R_, *Rtilde_, rho_old_ );
415 "Belos::TFQMRIter::initialize(): TFQMRIterState does not have initial residual.");
425 template <
class ScalarType,
class MV,
class OP>
431 if (initialized_ ==
false) {
440 ScalarType eta = STzero, beta = STzero;
449 "Belos::TFQMRIter::iterate(): current linear system has more than one vector!" );
455 while (stest_->checkStatus(
this) !=
Passed) {
457 for (
int iIter=0; iIter<2; iIter++)
465 MVT::MvDot( *V_, *Rtilde_, alpha_ );
466 alpha_[0] = rho_old_[0]/alpha_[0];
474 MVT::MvAddMv( STone, *W_, -alpha_[0], *AU_, *W_ );
481 MVT::MvAddMv( STone, *U_, (theta_[0]*theta_[0]/alpha_[0])*eta, *D_, *D_ );
492 MVT::MvAddMv( STone, *U_, -alpha_[0], *V_, *U_ );
495 lp_->apply( *U_, *AU_ );
502 MVT::MvNorm( *W_, theta_ );
503 theta_[0] /= tau_[0];
506 tau_[0] *= theta_[0]*cs_[0];
507 eta = cs_[0]*cs_[0]*alpha_[0];
514 MVT::MvAddMv( STone, *solnUpdate_, eta, *D_, *solnUpdate_ );
519 if ( tau_[0] == MTzero ) {
529 MVT::MvDot( *W_, *Rtilde_, rho_ );
530 beta = rho_[0]/rho_old_[0];
531 rho_old_[0] = rho_[0];
538 MVT::MvAddMv( STone, *W_, beta, *U_, *U_ );
541 MVT::MvAddMv( STone, *AU_, beta, *V_, *V_ );
544 lp_->apply( *U_, *AU_ );
547 MVT::MvAddMv( STone, *AU_, beta, *V_, *V_ );
560 #endif // BELOS_TFQMR_ITER_HPP
Teuchos::RCP< const MV > V
Collection of types and exceptions used within the Belos solvers.
This class implements the preconditioned transpose-free QMR algorithm for solving non-Hermitian linea...
Belos's basic output manager for sending information of select verbosity levels to the appropriate ou...
Teuchos::ScalarTraits< ScalarType > SCT
SCT::magnitudeType MagnitudeType
OperatorTraits< ScalarType, MV, OP > OPT
TFQMRIterateFailure(const std::string &what_arg)
void setBlockSize(int blockSize)
Set the blocksize.
Class which manages the output and verbosity of the Belos solvers.
Pure virtual base class for defining the status testing capabilities of Belos.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Declaration of basic traits for the multivector type.
void initialize()
Initialize the solver with the initial vectors from the linear problem or random data.
TFQMRIter(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::TFQMRIter constructor.
virtual ~TFQMRIter()
Belos::TFQMRIter destructor.
Teuchos::RCP< const MV > R
The current residual basis.
Teuchos::RCP< const MV > getNativeResiduals(std::vector< MagnitudeType > *norms) const
Get the norms of the residuals native to the solver.
Teuchos::RCP< const MV > U
void resetNumIters(int iter=0)
Reset the iteration count.
Pure virtual base class which describes the basic interface to the linear solver iteration.
A pure virtual class for defining the status tests for the Belos iterative solvers.
Class which defines basic traits for the operator type.
Traits class which defines basic operations on multivectors.
Teuchos::RCP< MV > getCurrentUpdate() const
Get the current update to the linear system.
void initializeTFQMR(const TFQMRIterState< ScalarType, MV > &newstate)
Initialize the solver to an iterate, providing a complete state.
Teuchos::RCP< const MV > Rtilde
void iterate()
This method performs TFQMR iterations until the status test indicates the need to stop or an error oc...
bool isInitialized()
States whether the solver has been initialized or not.
Teuchos::RCP< const MV > W
A linear system to solve, and its associated information.
Structure to contain pointers to TFQMRIter state variables.
Class which describes the linear problem to be solved by the iterative solver.
int getNumIters() const
Get the current iteration count.
MultiVecTraits< ScalarType, MV > MVT
TFQMRIterateFailure is thrown when the TFQMRIter object is unable to compute the next iterate in the ...
int getBlockSize() const
Get the blocksize to be used by the iterative solver in solving this linear problem.
Class which defines basic traits for the operator type.
Parent class to all Belos exceptions.
Teuchos::RCP< const MV > D
Belos header file which uses auto-configuration information to include necessary C++ headers...
const LinearProblem< ScalarType, MV, OP > & getProblem() const
Get a constant reference to the linear problem.
TFQMRIterState< ScalarType, MV > getState() const
Get the current state of the linear solver.