42 #ifndef BELOS_FIXEDPOINT_ITER_HPP 
   43 #define BELOS_FIXEDPOINT_ITER_HPP 
   76 template<
class ScalarType, 
class MV, 
class OP>
 
  233   bool stateStorageInitialized_;
 
  252   template<
class ScalarType, 
class MV, 
class OP>
 
  262     stateStorageInitialized_(false),
 
  270   template <
class ScalarType, 
class MV, 
class OP>
 
  273     if (!stateStorageInitialized_) {
 
  277       if (lhsMV == Teuchos::null && rhsMV == Teuchos::null) {
 
  278         stateStorageInitialized_ = 
false;
 
  285         if (R_ == Teuchos::null) {
 
  289                              "Belos::FixedPointIter::setStateSize(): linear problem does not specify multivectors to clone from.");
 
  290           R_ = MVT::Clone( *tmp, numRHS_ );
 
  291           Z_ = MVT::Clone( *tmp, numRHS_ );
 
  295         stateStorageInitialized_ = 
true;
 
  302   template <
class ScalarType, 
class MV, 
class OP>
 
  308     TEUCHOS_TEST_FOR_EXCEPTION(blockSize != MVT::GetNumberVecs(*lp_->getCurrRHSVec()), std::invalid_argument, 
"Belos::FixedPointIter::setBlockSize size must match # RHS.");
 
  310     TEUCHOS_TEST_FOR_EXCEPTION(blockSize <= 0, std::invalid_argument, 
"Belos::FixedPointIter::setBlockSize was passed a non-positive argument.");
 
  311     if (blockSize == numRHS_) {
 
  316     if (blockSize!=numRHS_)
 
  317       stateStorageInitialized_ = 
false;
 
  321     initialized_ = 
false;
 
  330   template <
class ScalarType, 
class MV, 
class OP>
 
  334     if (!stateStorageInitialized_)
 
  338                        "Belos::FixedPointIter::initialize(): Cannot initialize state storage!");
 
  342     std::string errstr(
"Belos::FixedPointIter::initialize(): Specified multivectors must have a consistent length and width.");
 
  348     if (newstate.
R != Teuchos::null) {
 
  350                                   std::invalid_argument, errstr );
 
  353                           std::invalid_argument, errstr );
 
  355                           std::invalid_argument, errstr );
 
  358       if (newstate.
R != R_) {
 
  360         MVT::MvAddMv( one, *newstate.
R, zero, *newstate.
R, *R_ );
 
  366                          "Belos::FixedPointIter::initialize(): FixedPointIterationState does not have initial residual.");
 
  376   template <
class ScalarType, 
class MV, 
class OP>
 
  382     if (initialized_ == 
false) {
 
  396     if (lp_->getRightPrec() != Teuchos::null) {
 
  401       MVT::MvInit( *Z_, zero );
 
  406       while (stest_->checkStatus(
this) != 
Passed) {
 
  412         lp_->applyRightPrec( *R_, *tmp );
 
  415         MVT::MvAddMv( one, *cur_soln_vec, one, *tmp, *cur_soln_vec );
 
  416         lp_->updateSolution();
 
  419         MVT::MvAddMv( one, *Z_, one, *tmp, *Z_ );
 
  422         lp_->applyOp (*Z_, *tmp );
 
  423         MVT::MvAddMv( one, *rhs, -one, *tmp, *R_ );
 
  433       while (stest_->checkStatus(
this) != 
Passed) {
 
  439         if ( lp_->getLeftPrec() != Teuchos::null ) {
 
  440           lp_->applyLeftPrec( *R_, *Z_ );
 
  447         MVT::MvAddMv(one,*cur_soln_vec,one,*Z_,*cur_soln_vec);
 
  448         lp_->updateSolution();
 
  451         lp_->applyOp(*cur_soln_vec,*tmp);
 
  452         MVT::MvAddMv(one,*rhs,-one,*tmp,*R_);
 
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...
Class which manages the output and verbosity of the Belos solvers. 
MultiVecTraits< ScalarType, MV > MVT
T & get(const std::string &name, T def_value)
void setBlockSize(int blockSize)
Set the blocksize to be used by the iterative solver in solving this linear problem. 
Pure virtual base class for defining the status testing capabilities of Belos. 
virtual ~FixedPointIter()
Destructor. 
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Teuchos::ScalarTraits< ScalarType > SCT
void initializeFixedPoint(FixedPointIterationState< ScalarType, MV > &newstate)
Initialize the solver to an iterate, providing a complete state. 
Declaration of basic traits for the multivector type. 
int getBlockSize() const 
Get the blocksize to be used by the iterative solver in solving this linear problem. 
void iterate()
This method performs Fixed Point iterations until the status test indicates the need to stop or an er...
A pure virtual class for defining the status tests for the Belos iterative solvers. 
Class which defines basic traits for the operator type. 
static int GetNumberVecs(const MV &mv)
Obtain the number of vectors in mv. 
int getNumIters() const 
Get the current iteration count. 
Traits class which defines basic operations on multivectors. 
FixedPointIterationState< ScalarType, MV > getState() const 
Get the current state of the linear solver. 
void resetNumIters(int iter=0)
Reset the iteration count. 
bool isInitialized()
States whether the solver has been initialized or not. 
A linear system to solve, and its associated information. 
Class which describes the linear problem to be solved by the iterative solver. 
Teuchos::RCP< const MV > getNativeResiduals(std::vector< MagnitudeType > *) const 
Get the norms of the residuals native to the solver. 
FixedPointIter(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)
FixedPointIter constructor with linear problem, solver utilities, and parameter list of solver option...
Teuchos::RCP< const MV > R
The current residual. 
const LinearProblem< ScalarType, MV, OP > & getProblem() const 
Get a constant reference to the linear problem. 
SCT::magnitudeType MagnitudeType
Structure to contain pointers to FixedPointIteration state variables. 
Teuchos::RCP< MV > getCurrentUpdate() const 
Get the current update to the linear system. 
void initialize()
Initialize the solver with the initial vectors from the linear problem or random data. 
Class which defines basic traits for the operator type. 
Teuchos::RCP< const MV > Z
The current preconditioned residual. 
Pure virtual base class which augments the basic interface for a fixed point linear solver iteration...
Belos header file which uses auto-configuration information to include necessary C++ headers...
This class implements the preconditioned fixed point iteration. 
OperatorTraits< ScalarType, MV, OP > OPT