18#ifndef BELOS_PSEUDO_BLOCK_TFQMR_ITER_HPP
19#define BELOS_PSEUDO_BLOCK_TFQMR_ITER_HPP
38#include "Teuchos_ScalarTraits.hpp"
39#include "Teuchos_ParameterList.hpp"
40#include "Teuchos_TimeMonitor.hpp"
59 template <
class ScalarType,
class MV>
62 typedef Teuchos::ScalarTraits<ScalarType>
SCT;
66 Teuchos::RCP<const MV>
W;
67 Teuchos::RCP<const MV>
U;
68 Teuchos::RCP<const MV>
AU;
70 Teuchos::RCP<const MV>
D;
71 Teuchos::RCP<const MV>
V;
77 Rtilde(Teuchos::null),
D(Teuchos::null),
V(Teuchos::null)
81 template <
class ScalarType,
class MV,
class OP>
89 typedef Teuchos::ScalarTraits<ScalarType>
SCT;
99 Teuchos::ParameterList ¶ms );
172 state.
alpha = alpha_;
176 state.
theta = theta_;
217 TEUCHOS_TEST_FOR_EXCEPTION(blockSize!=1,std::invalid_argument,
218 "Belos::PseudoBlockTFQMRIter::setBlockSize(): Cannot use a block size that is not one.");
232 const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > lp_;
233 const Teuchos::RCP<OutputManager<ScalarType> > om_;
234 const Teuchos::RCP<StatusTest<ScalarType,MV,OP> > stest_;
244 std::vector<ScalarType> alpha_, eta_, rho_, rho_old_;
245 std::vector<MagnitudeType> tau_, theta_;
262 Teuchos::RCP<MV> U_, AU_;
263 Teuchos::RCP<MV> Rtilde_;
266 Teuchos::RCP<MV> solnUpdate_;
277 template <
class ScalarType,
class MV,
class OP>
281 Teuchos::ParameterList &
294 template <
class ScalarType,
class MV,
class OP>
295 Teuchos::RCP<const MV>
298 MagnitudeType one = Teuchos::ScalarTraits<MagnitudeType>::one();
301 if ((
int) normvec->size() < numRHS_)
302 normvec->resize( numRHS_ );
305 for (
int i=0; i<numRHS_; i++) {
306 (*normvec)[i] = Teuchos::ScalarTraits<MagnitudeType>::squareroot( 2*iter_ + one )*tau_[i];
310 return Teuchos::null;
315 template <
class ScalarType,
class MV,
class OP>
319 const ScalarType STone = Teuchos::ScalarTraits<ScalarType>::one();
320 const ScalarType STzero = Teuchos::ScalarTraits<ScalarType>::zero();
321 const MagnitudeType MTzero = Teuchos::ScalarTraits<MagnitudeType>::zero();
324 TEUCHOS_TEST_FOR_EXCEPTION(newstate.
Rtilde == Teuchos::null,std::invalid_argument,
325 "Belos::PseudoBlockTFQMRIter::initialize(): PseudoBlockTFQMRIterState does not have initial residual.");
336 if ( Teuchos::is_null(D_) )
341 if ( Teuchos::is_null(solnUpdate_) )
346 if (newstate.
Rtilde != Rtilde_)
354 lp_->apply( *U_, *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 );
380 solnUpdate_ =
MVT::Clone( *Rtilde_, numRHS_ );
384 alpha_ = newstate.
alpha;
388 theta_ = newstate.
theta;
398 template <
class ScalarType,
class MV,
class OP>
404 if (initialized_ ==
false) {
409 const ScalarType STone = Teuchos::ScalarTraits<ScalarType>::one();
410 const ScalarType STzero = Teuchos::ScalarTraits<ScalarType>::zero();
411 const MagnitudeType MTone = Teuchos::ScalarTraits<MagnitudeType>::one();
412 const MagnitudeType MTzero = Teuchos::ScalarTraits<MagnitudeType>::zero();
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++)
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) {
463 MVT::MvAddMv( STone, *U_i, (theta_[i]*theta_[i]/alpha_[i])*eta_[i], *D_i, *D_i );
485 lp_->apply( *U_, *AU_ );
494 bool breakdown=
false;
495 for (
int i=0; i<numRHS_; ++i) {
496 theta_[i] /= tau_[i];
498 MagnitudeType cs = MTone / Teuchos::ScalarTraits<MagnitudeType>::squareroot(MTone + theta_[i]*theta_[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 );
533 for (
int i=0; i<numRHS_; ++i) {
534 beta[i] = rho_[i]/rho_old_[i];
554 lp_->apply( *U_, *AU_ );
557 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