50 #ifndef BELOS_PSEUDO_BLOCK_TFQMR_ITER_HPP
51 #define BELOS_PSEUDO_BLOCK_TFQMR_ITER_HPP
92 template <
class ScalarType,
class MV>
110 Rtilde(Teuchos::null),
D(Teuchos::null),
V(Teuchos::null)
146 template <
class ScalarType,
class MV,
class OP>
237 state.
alpha = alpha_;
241 state.
theta = theta_;
283 "Belos::PseudoBlockTFQMRIter::setBlockSize(): Cannot use a block size that is not one.");
309 std::vector<ScalarType> alpha_, eta_, rho_, rho_old_;
310 std::vector<MagnitudeType> tau_, theta_;
342 template <
class ScalarType,
class MV,
class OP>
359 template <
class ScalarType,
class MV,
class OP>
366 if ((
int) normvec->size() < numRHS_)
367 normvec->resize( numRHS_ );
370 for (
int i=0; i<numRHS_; i++) {
375 return Teuchos::null;
380 template <
class ScalarType,
class MV,
class OP>
390 "Belos::PseudoBlockTFQMRIter::initialize(): PseudoBlockTFQMRIterState does not have initial residual.");
393 int numRHS = MVT::GetNumberVecs(*newstate.
Rtilde);
398 if (
Teuchos::is_null(Rtilde_) || (MVT::GetNumberVecs(*Rtilde_) == numRHS_) )
402 D_ = MVT::Clone( *newstate.
Rtilde, numRHS_ );
403 MVT::MvInit( *D_, STzero );
407 solnUpdate_ = MVT::Clone( *newstate.
Rtilde, numRHS_ );
408 MVT::MvInit( *solnUpdate_, STzero );
411 if (newstate.
Rtilde != Rtilde_)
412 Rtilde_ = MVT::CloneCopy( *newstate.
Rtilde );
413 W_ = MVT::CloneCopy( *Rtilde_ );
414 U_ = MVT::CloneCopy( *Rtilde_ );
415 V_ = MVT::Clone( *Rtilde_, numRHS_ );
419 lp_->apply( *U_, *V_ );
420 AU_ = MVT::CloneCopy( *V_ );
423 alpha_.resize( numRHS_, STone );
424 eta_.resize( numRHS_, STzero );
425 rho_.resize( numRHS_ );
426 rho_old_.resize( numRHS_ );
427 tau_.resize( numRHS_ );
428 theta_.resize( numRHS_, MTzero );
430 MVT::MvNorm( *Rtilde_, tau_ );
431 MVT::MvDot( *Rtilde_, *Rtilde_, rho_ );
436 Rtilde_ = MVT::CloneCopy( *newstate.
Rtilde );
437 W_ = MVT::CloneCopy( *newstate.
W );
438 U_ = MVT::CloneCopy( *newstate.
U );
439 AU_ = MVT::CloneCopy( *newstate.
AU );
440 D_ = MVT::CloneCopy( *newstate.
D );
441 V_ = MVT::CloneCopy( *newstate.
V );
445 solnUpdate_ = MVT::Clone( *Rtilde_, numRHS_ );
446 MVT::MvInit( *solnUpdate_, STzero );
449 alpha_ = newstate.
alpha;
453 theta_ = newstate.
theta;
463 template <
class ScalarType,
class MV,
class OP>
469 if (initialized_ ==
false) {
478 std::vector< ScalarType > beta(numRHS_,STzero);
479 std::vector<int> index(1);
487 while (stest_->checkStatus(
this) !=
Passed) {
489 for (
int iIter=0; iIter<2; iIter++)
497 MVT::MvDot( *V_, *Rtilde_, alpha_ );
498 for (
int i=0; i<numRHS_; i++) {
499 rho_old_[i] = rho_[i];
500 alpha_[i] = rho_old_[i]/alpha_[i];
508 for (
int i=0; i<numRHS_; ++i) {
519 MVT::MvAddMv( STone, *W_i, -alpha_[i], *AU_i, *W_i );
528 MVT::MvAddMv( STone, *U_i, (theta_[i]*theta_[i]/alpha_[i])*eta_[i], *D_i, *D_i );
541 MVT::MvAddMv( STone, *U2_i, -alpha_[i], *V_i, *U2_i );
550 lp_->apply( *U_, *AU_ );
557 MVT::MvNorm( *W_, theta_ );
559 bool breakdown=
false;
560 for (
int i=0; i<numRHS_; ++i) {
561 theta_[i] /= tau_[i];
564 tau_[i] *= theta_[i]*cs;
565 if ( tau_[i] == MTzero ) {
568 eta_[i] = cs*cs*alpha_[i];
575 for (
int i=0; i<numRHS_; ++i) {
579 MVT::MvAddMv( STone, *update_i, eta_[i], *D_i, *update_i );
596 MVT::MvDot( *W_, *Rtilde_, rho_ );
598 for (
int i=0; i<numRHS_; ++i) {
599 beta[i] = rho_[i]/rho_old_[i];
610 MVT::MvAddMv( STone, *W_i, beta[i], *U_i, *U_i );
615 MVT::MvAddMv( STone, *AU_i, beta[i], *V_i, *V_i );
619 lp_->apply( *U_, *AU_ );
622 for (
int i=0; i<numRHS_; ++i) {
626 MVT::MvAddMv( STone, *AU_i, beta[i], *V_i, *V_i );
639 #endif // BELOS_PSEUDO_BLOCK_TFQMR_ITER_HPP
Collection of types and exceptions used within the Belos solvers.
Belos's basic output manager for sending information of select verbosity levels to the appropriate ou...
OperatorTraits< ScalarType, MV, OP > OPT
Class which manages the output and verbosity of the Belos solvers.
Teuchos::RCP< const MV > Rtilde
SCT::magnitudeType MagnitudeType
This class implements the preconditioned transpose-free QMR algorithm for solving non-Hermitian linea...
void initialize()
Initialize the solver with the initial vectors from the linear problem or random data.
PseudoBlockTFQMRIterateFailure is thrown when the PseudoBlockTFQMRIter object is unable to compute th...
Structure to contain pointers to PseudoBlockTFQMRIter state variables.
bool is_null(const std::shared_ptr< T > &p)
MultiVecTraits< ScalarType, MV > MVT
Pure virtual base class for defining the status testing capabilities of Belos.
PseudoBlockTFQMRIterState< ScalarType, MV > getState() const
Get the current state of the linear solver.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Declaration of basic traits for the multivector type.
const LinearProblem< ScalarType, MV, OP > & getProblem() const
Get a constant reference to the linear problem.
Teuchos::RCP< const MV > AU
int getNumIters() const
Get the current iteration count.
Pure virtual base class which describes the basic interface to the linear solver iteration.
std::vector< ScalarType > alpha
Teuchos::RCP< const MV > getNativeResiduals(std::vector< MagnitudeType > *norms) const
Get the norms of the residuals native to the solver.
A pure virtual class for defining the status tests for the Belos iterative solvers.
Class which defines basic traits for the operator type.
PseudoBlockTFQMRIterInitFailure(const std::string &what_arg)
int getBlockSize() const
Get the blocksize to be used by the iterative solver in solving this linear problem.
Traits class which defines basic operations on multivectors.
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.
PseudoBlockTFQMRIterState()
PseudoBlockTFQMRIterInitFailure is thrown when the PseudoBlockTFQMRIter object is unable to generate ...
Teuchos::RCP< MV > getCurrentUpdate() const
Get the current update to the linear system.
Teuchos::ScalarTraits< ScalarType > SCT
A linear system to solve, and its associated information.
std::vector< MagnitudeType > theta
Class which describes the linear problem to be solved by the iterative solver.
void setBlockSize(int blockSize)
Set the blocksize.
Teuchos::RCP< const MV > D
std::vector< ScalarType > eta
Teuchos::RCP< const MV > U
std::vector< MagnitudeType > tau
PseudoBlockTFQMRIterateFailure(const std::string &what_arg)
virtual ~PseudoBlockTFQMRIter()
Belos::PseudoBlockTFQMRIter destructor.
Teuchos::RCP< const MV > W
The current residual basis.
std::vector< ScalarType > rho
Class which defines basic traits for the operator type.
Teuchos::ScalarTraits< ScalarType > SCT
void iterate()
This method performs pseudo-block TFQMR iterations until the status test indicates the need to stop o...
Parent class to all Belos exceptions.
void resetNumIters(int iter=0)
Reset the iteration count.
Belos header file which uses auto-configuration information to include necessary C++ headers...
SCT::magnitudeType MagnitudeType
Teuchos::RCP< const MV > V
bool isInitialized()
States whether the solver has been initialized or not.
void initializeTFQMR(const PseudoBlockTFQMRIterState< ScalarType, MV > &newstate)
Initialize the solver to an iterate, providing a complete state.