50#ifndef BELOS_PSEUDO_BLOCK_TFQMR_ITER_HPP
51#define BELOS_PSEUDO_BLOCK_TFQMR_ITER_HPP
70#include "Teuchos_ScalarTraits.hpp"
71#include "Teuchos_ParameterList.hpp"
72#include "Teuchos_TimeMonitor.hpp"
91 template <
class ScalarType,
class MV>
94 typedef Teuchos::ScalarTraits<ScalarType>
SCT;
98 Teuchos::RCP<const MV>
W;
99 Teuchos::RCP<const MV>
U;
100 Teuchos::RCP<const MV>
AU;
102 Teuchos::RCP<const MV>
D;
103 Teuchos::RCP<const MV>
V;
109 Rtilde(Teuchos::null),
D(Teuchos::null),
V(Teuchos::null)
113 template <
class ScalarType,
class MV,
class OP>
121 typedef Teuchos::ScalarTraits<ScalarType>
SCT;
131 Teuchos::ParameterList ¶ms );
204 state.
alpha = alpha_;
208 state.
theta = theta_;
249 TEUCHOS_TEST_FOR_EXCEPTION(blockSize!=1,std::invalid_argument,
250 "Belos::PseudoBlockTFQMRIter::setBlockSize(): Cannot use a block size that is not one.");
264 const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > lp_;
265 const Teuchos::RCP<OutputManager<ScalarType> > om_;
266 const Teuchos::RCP<StatusTest<ScalarType,MV,OP> > stest_;
276 std::vector<ScalarType> alpha_, eta_, rho_, rho_old_;
277 std::vector<MagnitudeType> tau_, theta_;
294 Teuchos::RCP<MV> U_, AU_;
295 Teuchos::RCP<MV> Rtilde_;
298 Teuchos::RCP<MV> solnUpdate_;
309 template <
class ScalarType,
class MV,
class OP>
313 Teuchos::ParameterList &
326 template <
class ScalarType,
class MV,
class OP>
327 Teuchos::RCP<const MV>
330 MagnitudeType one = Teuchos::ScalarTraits<MagnitudeType>::one();
333 if ((
int) normvec->size() < numRHS_)
334 normvec->resize( numRHS_ );
337 for (
int i=0; i<numRHS_; i++) {
338 (*normvec)[i] = Teuchos::ScalarTraits<MagnitudeType>::squareroot( 2*iter_ + one )*tau_[i];
342 return Teuchos::null;
347 template <
class ScalarType,
class MV,
class OP>
351 const ScalarType STone = Teuchos::ScalarTraits<ScalarType>::one();
352 const ScalarType STzero = Teuchos::ScalarTraits<ScalarType>::zero();
353 const MagnitudeType MTzero = Teuchos::ScalarTraits<MagnitudeType>::zero();
356 TEUCHOS_TEST_FOR_EXCEPTION(newstate.
Rtilde == Teuchos::null,std::invalid_argument,
357 "Belos::PseudoBlockTFQMRIter::initialize(): PseudoBlockTFQMRIterState does not have initial residual.");
368 if ( Teuchos::is_null(D_) )
373 if ( Teuchos::is_null(solnUpdate_) )
378 if (newstate.
Rtilde != Rtilde_)
386 lp_->apply( *U_, *V_ );
390 alpha_.resize( numRHS_, STone );
391 eta_.resize( numRHS_, STzero );
392 rho_.resize( numRHS_ );
393 rho_old_.resize( numRHS_ );
394 tau_.resize( numRHS_ );
395 theta_.resize( numRHS_, MTzero );
412 solnUpdate_ =
MVT::Clone( *Rtilde_, numRHS_ );
416 alpha_ = newstate.
alpha;
420 theta_ = newstate.
theta;
430 template <
class ScalarType,
class MV,
class OP>
436 if (initialized_ ==
false) {
441 const ScalarType STone = Teuchos::ScalarTraits<ScalarType>::one();
442 const ScalarType STzero = Teuchos::ScalarTraits<ScalarType>::zero();
443 const MagnitudeType MTone = Teuchos::ScalarTraits<MagnitudeType>::one();
444 const MagnitudeType MTzero = Teuchos::ScalarTraits<MagnitudeType>::zero();
445 std::vector< ScalarType > beta(numRHS_,STzero);
446 std::vector<int> index(1);
454 while (stest_->checkStatus(
this) !=
Passed) {
456 for (
int iIter=0; iIter<2; iIter++)
465 for (
int i=0; i<numRHS_; i++) {
466 rho_old_[i] = rho_[i];
467 alpha_[i] = rho_old_[i]/alpha_[i];
475 for (
int i=0; i<numRHS_; ++i) {
495 MVT::MvAddMv( STone, *U_i, (theta_[i]*theta_[i]/alpha_[i])*eta_[i], *D_i, *D_i );
517 lp_->apply( *U_, *AU_ );
526 bool breakdown=
false;
527 for (
int i=0; i<numRHS_; ++i) {
528 theta_[i] /= tau_[i];
530 MagnitudeType cs = MTone / Teuchos::ScalarTraits<MagnitudeType>::squareroot(MTone + theta_[i]*theta_[i]);
531 tau_[i] *= theta_[i]*cs;
532 if ( tau_[i] == MTzero ) {
535 eta_[i] = cs*cs*alpha_[i];
542 for (
int i=0; i<numRHS_; ++i) {
546 MVT::MvAddMv( STone, *update_i, eta_[i], *D_i, *update_i );
565 for (
int i=0; i<numRHS_; ++i) {
566 beta[i] = rho_[i]/rho_old_[i];
586 lp_->apply( *U_, *AU_ );
589 for (
int i=0; i<numRHS_; ++i) {
Belos header file which uses auto-configuration information to include necessary C++ headers.
Pure virtual base class which describes the basic interface to the linear solver iteration.
Class which describes the linear problem to be solved by the iterative solver.
Declaration of basic traits for the multivector type.
Class which defines basic traits for the operator type.
Class which manages the output and verbosity of the Belos solvers.
Pure virtual base class for defining the status testing capabilities of Belos.
Collection of types and exceptions used within the Belos solvers.
Iteration()
Default Constructor.
Traits class which defines basic operations on multivectors.
static Teuchos::RCP< MV > CloneCopy(const MV &mv)
Creates a new MV and copies contents of mv into the new vector (deep copy).
static void MvNorm(const MV &mv, std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > &normvec, NormType type=TwoNorm)
Compute the norm of each individual vector of mv. Upon return, normvec[i] holds the value of ,...
static void MvDot(const MV &mv, const MV &A, std::vector< ScalarType > &b)
Compute a vector b where the components are the individual dot-products of the i-th columns of A and ...
static Teuchos::RCP< MV > Clone(const MV &mv, const int numvecs)
Creates a new empty MV containing numvecs columns.
static Teuchos::RCP< const MV > CloneView(const MV &mv, const std::vector< int > &index)
Creates a new const MV that shares the selected contents of mv (shallow copy).
static Teuchos::RCP< MV > CloneViewNonConst(MV &mv, const std::vector< int > &index)
Creates a new MV that shares the selected contents of mv (shallow copy).
static void MvAddMv(const ScalarType alpha, const MV &A, const ScalarType beta, const MV &B, MV &mv)
Replace mv with .
static void MvInit(MV &mv, const ScalarType alpha=Teuchos::ScalarTraits< ScalarType >::zero())
Replace each element of the vectors in mv with alpha.
static int GetNumberVecs(const MV &mv)
Obtain the number of vectors in mv.
Class which defines basic traits for the operator type.
Belos's basic output manager for sending information of select verbosity levels to the appropriate ou...
int getNumIters() const
Get the current iteration count.
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.
int getBlockSize() const
Get the blocksize to be used by the iterative solver in solving this linear problem.
bool isInitialized()
States whether the solver has been initialized or not.
virtual ~PseudoBlockTFQMRIter()
Belos::PseudoBlockTFQMRIter destructor.
PseudoBlockTFQMRIterState< ScalarType, MV > getState() const
Get the current state of the linear solver.
MultiVecTraits< ScalarType, MV > MVT
void initialize()
Initialize the solver with the initial vectors from the linear problem or random data.
Teuchos::RCP< const MV > getNativeResiduals(std::vector< MagnitudeType > *norms) const
Teuchos::RCP< MV > getCurrentUpdate() const
Get the current update to the linear system.
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.
Teuchos::ScalarTraits< ScalarType > SCT
void initializeTFQMR(const PseudoBlockTFQMRIterState< ScalarType, MV > &newstate)
Initialize the solver to an iterate, providing a complete state.
SCT::magnitudeType MagnitudeType
void setBlockSize(int blockSize)
Set the blocksize.
const LinearProblem< ScalarType, MV, OP > & getProblem() const
Get a constant reference to the linear problem.
OperatorTraits< ScalarType, MV, OP > OPT
A pure virtual class for defining the status tests for the Belos iterative solvers.
Structure to contain pointers to PseudoBlockTFQMRIter state variables.
std::vector< MagnitudeType > theta
std::vector< MagnitudeType > tau
SCT::magnitudeType MagnitudeType
Teuchos::RCP< const MV > AU
std::vector< ScalarType > rho
Teuchos::RCP< const MV > D
std::vector< ScalarType > alpha
Teuchos::RCP< const MV > Rtilde
Teuchos::RCP< const MV > W
The current residual basis.
Teuchos::RCP< const MV > U
Teuchos::RCP< const MV > V
PseudoBlockTFQMRIterState()
Teuchos::ScalarTraits< ScalarType > SCT
std::vector< ScalarType > eta