18 #ifndef BELOS_PSEUDO_BLOCK_TFQMR_ITER_HPP
19 #define BELOS_PSEUDO_BLOCK_TFQMR_ITER_HPP
59 template <
class ScalarType,
class MV>
77 Rtilde(Teuchos::null),
D(Teuchos::null),
V(Teuchos::null)
81 template <
class ScalarType,
class MV,
class OP>
218 "Belos::PseudoBlockTFQMRIter::setBlockSize(): Cannot use a block size that is not one.");
277 template <
class ScalarType,
class MV,
class OP>
294 template <
class ScalarType,
class MV,
class OP>
301 if ((
int) normvec->size() < numRHS_)
302 normvec->resize( numRHS_ );
305 for (
int i=0; i<numRHS_; i++) {
315 template <
class ScalarType,
class MV,
class OP>
325 "Belos::PseudoBlockTFQMRIter::initialize(): PseudoBlockTFQMRIterState does not have initial residual.");
328 int numRHS = MVT::GetNumberVecs(*newstate.
Rtilde);
333 if (
Teuchos::is_null(Rtilde_) || (MVT::GetNumberVecs(*Rtilde_) == numRHS_) )
337 D_ = MVT::Clone( *newstate.
Rtilde, numRHS_ );
338 MVT::MvInit( *D_, STzero );
342 solnUpdate_ = MVT::Clone( *newstate.
Rtilde, numRHS_ );
343 MVT::MvInit( *solnUpdate_, STzero );
346 if (newstate.
Rtilde != Rtilde_)
347 Rtilde_ = MVT::CloneCopy( *newstate.
Rtilde );
348 W_ = MVT::CloneCopy( *Rtilde_ );
349 U_ = MVT::CloneCopy( *Rtilde_ );
350 V_ = MVT::Clone( *Rtilde_, numRHS_ );
354 lp_->apply( *U_, *V_ );
355 AU_ = MVT::CloneCopy( *V_ );
358 alpha_.resize( numRHS_, STone );
359 eta_.resize( numRHS_, STzero );
360 rho_.resize( numRHS_ );
361 rho_old_.resize( numRHS_ );
362 tau_.resize( numRHS_ );
363 theta_.resize( numRHS_, MTzero );
365 MVT::MvNorm( *Rtilde_, tau_ );
366 MVT::MvDot( *Rtilde_, *Rtilde_, rho_ );
371 Rtilde_ = MVT::CloneCopy( *newstate.
Rtilde );
372 W_ = MVT::CloneCopy( *newstate.
W );
373 U_ = MVT::CloneCopy( *newstate.
U );
374 AU_ = MVT::CloneCopy( *newstate.
AU );
375 D_ = MVT::CloneCopy( *newstate.
D );
376 V_ = MVT::CloneCopy( *newstate.
V );
380 solnUpdate_ = MVT::Clone( *Rtilde_, numRHS_ );
381 MVT::MvInit( *solnUpdate_, STzero );
384 alpha_ = newstate.
alpha;
388 theta_ = newstate.
theta;
398 template <
class ScalarType,
class MV,
class OP>
404 if (initialized_ ==
false) {
413 std::vector< ScalarType > beta(numRHS_,STzero);
414 std::vector<int> index(1);
422 while (stest_->checkStatus(
this) !=
Passed) {
424 for (
int iIter=0; iIter<2; iIter++)
432 MVT::MvDot( *V_, *Rtilde_, alpha_ );
433 for (
int i=0; i<numRHS_; i++) {
434 rho_old_[i] = rho_[i];
435 alpha_[i] = rho_old_[i]/alpha_[i];
443 for (
int i=0; i<numRHS_; ++i) {
454 MVT::MvAddMv( STone, *W_i, -alpha_[i], *AU_i, *W_i );
463 MVT::MvAddMv( STone, *U_i, (theta_[i]*theta_[i]/alpha_[i])*eta_[i], *D_i, *D_i );
476 MVT::MvAddMv( STone, *U2_i, -alpha_[i], *V_i, *U2_i );
485 lp_->apply( *U_, *AU_ );
492 MVT::MvNorm( *W_, theta_ );
494 bool breakdown=
false;
495 for (
int i=0; i<numRHS_; ++i) {
496 theta_[i] /= tau_[i];
499 tau_[i] *= theta_[i]*cs;
500 if ( tau_[i] == MTzero ) {
503 eta_[i] = cs*cs*alpha_[i];
510 for (
int i=0; i<numRHS_; ++i) {
514 MVT::MvAddMv( STone, *update_i, eta_[i], *D_i, *update_i );
531 MVT::MvDot( *W_, *Rtilde_, rho_ );
533 for (
int i=0; i<numRHS_; ++i) {
534 beta[i] = rho_[i]/rho_old_[i];
545 MVT::MvAddMv( STone, *W_i, beta[i], *U_i, *U_i );
550 MVT::MvAddMv( STone, *AU_i, beta[i], *V_i, *V_i );
554 lp_->apply( *U_, *AU_ );
557 for (
int i=0; i<numRHS_; ++i) {
561 MVT::MvAddMv( STone, *AU_i, beta[i], *V_i, *V_i );
574 #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...
std::vector< ScalarType > rho_old_
OperatorTraits< ScalarType, MV, OP > OPT
const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > stest_
Class which manages the output and verbosity of the Belos solvers.
Teuchos::RCP< MV > solnUpdate_
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.
Structure to contain pointers to PseudoBlockTFQMRIter state variables.
bool is_null(const std::shared_ptr< T > &p)
std::vector< ScalarType > alpha_
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.
Teuchos::RCP< MV > Rtilde_
A pure virtual class for defining the status tests for the Belos iterative solvers.
Class which defines basic traits for the operator type.
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.
std::vector< MagnitudeType > theta_
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()
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.
std::vector< MagnitudeType > tau_
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
std::vector< ScalarType > rho_
const Teuchos::RCP< OutputManager< ScalarType > > om_
virtual ~PseudoBlockTFQMRIter()
Belos::PseudoBlockTFQMRIter destructor.
Teuchos::RCP< const MV > W
The current residual basis.
std::vector< ScalarType > eta_
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...
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
const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > lp_
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.