10#ifndef BELOS_BICGSTAB_ITER_HPP
11#define BELOS_BICGSTAB_ITER_HPP
28#include "Teuchos_SerialDenseMatrix.hpp"
29#include "Teuchos_SerialDenseVector.hpp"
30#include "Teuchos_ScalarTraits.hpp"
31#include "Teuchos_ParameterList.hpp"
32#include "Teuchos_TimeMonitor.hpp"
54 template <
class ScalarType,
class MV>
58 Teuchos::RCP<const MV>
R;
61 Teuchos::RCP<const MV>
Rhat;
64 Teuchos::RCP<const MV>
P;
67 Teuchos::RCP<const MV>
V;
72 P(Teuchos::null),
V(Teuchos::null)
80 template<
class ScalarType,
class MV,
class OP>
90 typedef Teuchos::ScalarTraits<ScalarType>
SCT;
92 typedef Teuchos::ScalarTraits<MagnitudeType>
MT;
105 Teuchos::ParameterList ¶ms );
175 state.
alpha = alpha_;
176 state.
omega = omega_;
219 TEUCHOS_TEST_FOR_EXCEPTION(blockSize!=1,std::invalid_argument,
220 "Belos::BiCGStabIter::setBlockSize(): Cannot use a block size that is not one.");
230 void axpy(
const ScalarType alpha,
const MV & A,
231 const std::vector<ScalarType> beta,
const MV& B, MV& mv,
bool minus=
false);
236 const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > lp_;
237 const Teuchos::RCP<OutputManager<ScalarType> > om_;
238 const Teuchos::RCP<StatusTest<ScalarType,MV,OP> > stest_;
264 Teuchos::RCP<MV> Rhat_;
275 std::vector<ScalarType> rho_old_, alpha_, omega_;
280 template<
class ScalarType,
class MV,
class OP>
284 Teuchos::ParameterList & ):
298 template <
class ScalarType,
class MV,
class OP>
302 Teuchos::RCP<const MV> lhsMV = lp_->getCurrLHSVec();
303 Teuchos::RCP<const MV> rhsMV = lp_->getCurrRHSVec();
304 TEUCHOS_TEST_FOR_EXCEPTION((lhsMV==Teuchos::null && rhsMV==Teuchos::null),std::invalid_argument,
305 "Belos::BiCGStabIter::initialize(): Cannot initialize state storage!");
308 Teuchos::RCP<const MV> tmp = ( (rhsMV!=Teuchos::null)? rhsMV: lhsMV );
322 rho_old_.resize(numRHS_);
323 alpha_.resize(numRHS_);
324 omega_.resize(numRHS_);
332 std::string errstr(
"Belos::BlockPseudoCGIter::initialize(): Specified multivectors must have a consistent length and width.");
335 const ScalarType one = SCT::one();
337 if (!Teuchos::is_null(newstate.
R)) {
340 std::invalid_argument, errstr );
342 std::invalid_argument, errstr );
345 if (newstate.
R != R_) {
351 lp_->computeCurrResVec(R_.get());
355 if (!Teuchos::is_null(newstate.
Rhat) && newstate.
Rhat != Rhat_) {
365 if (!Teuchos::is_null(newstate.
V) && newstate.
V != V_) {
375 if (!Teuchos::is_null(newstate.
P) && newstate.
P != P_) {
385 if (newstate.
rho_old.size () ==
static_cast<size_t> (numRHS_)) {
391 rho_old_.assign(numRHS_,one);
395 if (newstate.
alpha.size() ==
static_cast<size_t> (numRHS_)) {
397 alpha_ = newstate.
alpha;
401 alpha_.assign(numRHS_,one);
405 if (newstate.
omega.size() ==
static_cast<size_t> (numRHS_)) {
407 omega_ = newstate.
omega;
411 omega_.assign(numRHS_,one);
417 TEUCHOS_TEST_FOR_EXCEPTION(Teuchos::is_null(newstate.
R),std::invalid_argument,
418 "Belos::BiCGStabIter::initialize(): BiCGStabStateIterState does not have initial residual.");
428 template <
class ScalarType,
class MV,
class OP>
436 if (initialized_ ==
false) {
442 std::vector<ScalarType> rho_new( numRHS_ ), beta( numRHS_ );
443 std::vector<ScalarType> rhatV( numRHS_ ), tT( numRHS_ ), tS( numRHS_ );
446 const ScalarType one = SCT::one();
449 RCP<MV> leftPrecVec, leftPrecVec2;
454 if (lp_->isLeftPrec() || lp_->isRightPrec()) {
464 Teuchos::RCP<MV> X = lp_->getCurrLHSVec();
469 while (stest_->checkStatus(
this) !=
Passed && !breakdown_) {
479 for(i=0; i<numRHS_; i++) {
482 if (SCT::magnitude(rho_new[i]) < MT::sfmin())
485 beta[i] = (rho_new[i] / rho_old_[i]) * (alpha_[i] / omega_[i]);
491 axpy(one, *P_, omega_, *V_, *P_,
true);
492 axpy(one, *R_, beta, *P_, *P_);
496 if(lp_->isLeftPrec()) {
497 if(lp_->isRightPrec()) {
498 if(leftPrecVec == Teuchos::null) {
501 lp_->applyLeftPrec(*P_,*leftPrecVec);
502 lp_->applyRightPrec(*leftPrecVec,*Y);
505 lp_->applyLeftPrec(*P_,*Y);
508 else if(lp_->isRightPrec()) {
509 lp_->applyRightPrec(*P_,*Y);
513 lp_->applyOp(*Y,*V_);
517 for(i=0; i<numRHS_; i++) {
518 if (SCT::magnitude(rhatV[i]) < MT::sfmin())
524 alpha_[i] = rho_new[i] / rhatV[i];
528 axpy(one, *R_, alpha_, *V_, *S,
true);
531 if(lp_->isLeftPrec()) {
532 if(lp_->isRightPrec()) {
533 if(leftPrecVec == Teuchos::null) {
536 lp_->applyLeftPrec(*S,*leftPrecVec);
537 lp_->applyRightPrec(*leftPrecVec,*Z);
540 lp_->applyLeftPrec(*S,*Z);
543 else if(lp_->isRightPrec()) {
544 lp_->applyRightPrec(*S,*Z);
551 if(lp_->isLeftPrec()) {
552 if(leftPrecVec == Teuchos::null) {
555 if(leftPrecVec2 == Teuchos::null) {
558 lp_->applyLeftPrec(*T,*leftPrecVec2);
566 for(i=0; i<numRHS_; i++) {
567 if (SCT::magnitude(tT[i]) < MT::sfmin())
569 omega_[i] = SCT::zero();
573 omega_[i] = tS[i] / tT[i];
577 axpy(one, *X, alpha_, *Y, *X);
578 axpy(one, *X, omega_, *Z, *X);
581 axpy(one, *S, omega_, *T, *R_,
true);
591 template <
class ScalarType,
class MV,
class OP>
592 void BiCGStabIter<ScalarType,MV,OP>::axpy(
const ScalarType alpha,
const MV & A,
593 const std::vector<ScalarType> beta,
const MV& B, MV& mv,
bool minus)
595 Teuchos::RCP<const MV> A1, B1;
596 Teuchos::RCP<MV> mv1;
597 std::vector<int> index(1);
599 for(
int i=0; i<numRHS_; i++) {
601 A1 = MVT::CloneView(A,index);
602 B1 = MVT::CloneView(B,index);
603 mv1 = MVT::CloneViewNonConst(mv,index);
605 MVT::MvAddMv(alpha,*A1,-beta[i],*B1,*mv1);
608 MVT::MvAddMv(alpha,*A1,beta[i],*B1,*mv1);
Pure virtual base class which augments the basic interface for a conjugate gradient linear solver ite...
Belos header file which uses auto-configuration information to include necessary C++ headers.
Class which describes the linear problem to be solved by the iterative solver.
Templated virtual class for providing orthogonalization/orthonormalization methods with matrix-based ...
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.
void setBlockSize(int blockSize)
Set the blocksize.
Teuchos::ScalarTraits< ScalarType > SCT
bool breakdownDetected()
Has breakdown been detected in any linear system.
const LinearProblem< ScalarType, MV, OP > & getProblem() const
Get a constant reference to the linear problem.
Teuchos::RCP< MV > getCurrentUpdate() const
Get the current update to the linear system.
void resetNumIters(int iter=0)
Reset the iteration count.
BiCGStabIter(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)
BiCGStabIter constructor with linear problem, solver utilities, and parameter list of solver options.
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.
MultiVecTraits< ScalarType, MV > MVT
virtual ~BiCGStabIter()
Destructor.
void iterate()
This method performs BiCGStab iterations on each linear system until the status test indicates the ne...
int getNumIters() const
Get the current iteration count.
void initializeBiCGStab(BiCGStabIterationState< ScalarType, MV > &newstate)
Initialize the solver to an iterate, providing a complete state.
Teuchos::RCP< const MV > getNativeResiduals(std::vector< MagnitudeType > *) const
BiCGStabIterationState< ScalarType, MV > getState() const
Get the current state of the linear solver.
Teuchos::ScalarTraits< MagnitudeType > MT
SCT::magnitudeType MagnitudeType
void initialize()
Initialize the solver with the initial vectors from the linear problem or random data.
OperatorTraits< ScalarType, MV, OP > OPT
Iteration()
Default Constructor.
Traits class which defines basic operations on multivectors.
static ptrdiff_t GetGlobalLength(const MV &mv)
Return the number of rows in the given multivector mv.
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 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.
static void Assign(const MV &A, MV &mv)
mv := A
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...
A pure virtual class for defining the status tests for the Belos iterative solvers.
Structure to contain pointers to BiCGStabIteration state variables.
std::vector< ScalarType > omega
Teuchos::RCP< const MV > R
The current residual.
Teuchos::RCP< const MV > Rhat
The initial residual.
std::vector< ScalarType > rho_old
Teuchos::RCP< const MV > P
The first decent direction vector.
std::vector< ScalarType > alpha
Teuchos::RCP< const MV > V
A * M * the first decent direction vector.