42 #ifndef BELOS_MINRES_ITER_HPP
43 #define BELOS_MINRES_ITER_HPP
94 template<
class ScalarType,
class MV,
class OP>
183 throw std::logic_error(
"getState() cannot be called unless "
184 "the state has been initialized");
214 std::vector<MagnitudeType>& theNorms = *norms;
215 if (theNorms.size() < 1)
217 theNorms[0] = phibar_;
219 return Teuchos::null;
228 void symOrtho( ScalarType a, ScalarType b, ScalarType *c, ScalarType *s, ScalarType *r );
244 "Belos::MinresIter::setBlockSize(): Cannot use a block size that is not one.");
284 bool stateStorageInitialized_;
325 template<
class ScalarType,
class MV,
class OP>
334 stateStorageInitialized_(false),
342 template <
class ScalarType,
class MV,
class OP>
345 if (!stateStorageInitialized_) {
350 if (lhsMV == Teuchos::null && rhsMV == Teuchos::null) {
351 stateStorageInitialized_ =
false;
358 if (Y_ == Teuchos::null) {
362 std::invalid_argument,
363 "Belos::MinresIter::setStateSize(): linear problem does not specify multivectors to clone from.");
364 Y_ = MVT::Clone( *tmp, 1 );
365 R1_ = MVT::Clone( *tmp, 1 );
366 R2_ = MVT::Clone( *tmp, 1 );
367 W_ = MVT::Clone( *tmp, 1 );
368 W1_ = MVT::Clone( *tmp, 1 );
369 W2_ = MVT::Clone( *tmp, 1 );
372 stateStorageInitialized_ =
true;
380 template <
class ScalarType,
class MV,
class OP>
384 if (!stateStorageInitialized_)
388 std::invalid_argument,
389 "Belos::MinresIter::initialize(): Cannot initialize state storage!" );
392 std::invalid_argument,
393 "Belos::MinresIter::initialize(): MinresIterationState does not have initial residual.");
395 std::string errstr(
"Belos::MinresIter::initialize(): Specified multivectors must have a consistent length and width.");
397 std::invalid_argument,
400 std::invalid_argument,
404 const ScalarType one = SCT::one();
410 MVT::MvAddMv( one, *newstate.
Y, zero, *newstate.
Y, *R2_ );
411 MVT::MvAddMv( one, *newstate.
Y, zero, *newstate.
Y, *R1_ );
415 MVT::MvInit ( *W2_ );
417 if ( lp_->getLeftPrec() != Teuchos::null ) {
418 lp_->applyLeftPrec( *newstate.
Y, *Y_ );
421 if (newstate.
Y != Y_) {
423 MVT::MvAddMv( one, *newstate.
Y, zero, *newstate.
Y, *Y_ );
429 MVT::MvTransMv( one, *newstate.
Y, *Y_, beta1_ );
432 std::invalid_argument,
433 "The preconditioner is not positive definite." );
435 if( SCT::magnitude(beta1_(0,0)) == zero )
439 MVT::MvInit( *cur_soln_vec );
442 beta1_(0,0) = SCT::squareroot( beta1_(0,0) );
451 template <
class ScalarType,
class MV,
class OP>
457 if (initialized_ ==
false) {
464 const ScalarType one = SCT::one();
471 ScalarType shift = zero;
474 ScalarType oldBeta = zero;
475 ScalarType epsln = zero;
476 ScalarType cs = -one;
477 ScalarType sn = zero;
478 ScalarType dbar = zero;
497 "Belos::MinresIter::iterate(): current linear system has more than one vector!" );
502 while (stest_->checkStatus(
this) !=
Passed) {
509 MVT::MvAddMv (one / beta(0,0), *Y_, zero, *Y_, *V);
512 lp_->applyOp (*V, *Y_);
516 MVT::MvAddMv (one, *Y_, -shift, *V, *Y_);
519 MVT::MvAddMv (one, *Y_, -beta(0,0)/oldBeta, *R1_, *Y_);
522 MVT::MvTransMv (one, *V, *Y_, alpha);
525 MVT::MvAddMv (one, *Y_, -alpha(0,0)/beta(0,0), *R2_, *Y_);
535 if ( lp_->getLeftPrec() != Teuchos::null ) {
536 lp_->applyLeftPrec( *R2_, *Y_ );
539 MVT::MvAddMv( one, *R2_, zero, *R2_, *Y_ );
544 MVT::MvTransMv( one, *R2_, *Y_, beta );
558 "Belos::MinresIter::iterate(): Encountered nonpositi"
559 "ve value " << beta(0,0) <<
" for r2^H*M*r2 at itera"
560 "tion " << iter_ <<
": MINRES cannot continue." );
561 beta(0,0) = SCT::squareroot( beta(0,0) );
569 delta = cs*dbar + sn*alpha(0,0);
570 gbar = sn*dbar - cs*alpha(0,0);
571 epsln = sn*beta(0,0);
572 dbar = - cs*beta(0,0);
575 this->symOrtho(gbar, beta(0,0), &cs, &sn, &gamma);
582 MVT::MvAddMv( one, *W_, zero, *W_, *W1_ );
589 MVT::MvAddMv( one, *V, -oldeps, *W1_, *W_ );
590 MVT::MvAddMv( one, *W_, -delta, *W2_, *W_ );
591 MVT::MvScale( *W_, one / gamma );
595 MVT::MvAddMv( one, *cur_soln_vec, phi, *W_, *cur_soln_vec );
596 lp_->updateSolution();
606 template <
class ScalarType,
class MV,
class OP>
608 ScalarType *c, ScalarType *s, ScalarType *r
611 const ScalarType one = SCT::one();
612 const ScalarType zero = SCT::zero();
616 if ( absB == m_zero ) {
619 if ( absA == m_zero )
623 }
else if ( absA == m_zero ) {
627 }
else if ( absB >= absA ) {
628 ScalarType tau = a / b;
630 *s = -one / SCT::squareroot( one+tau*tau );
632 *s = one / SCT::squareroot( one+tau*tau );
636 ScalarType tau = b / a;
638 *c = -one / SCT::squareroot( one+tau*tau );
640 *c = one / SCT::squareroot( one+tau*tau );
Teuchos::RCP< const MV > R1
Previous residual.
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.
Pure virtual base class for defining the status testing capabilities of Belos.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
void resetNumIters(int iter=0)
Reset the iteration count.
Declaration of basic traits for the multivector type.
Teuchos::RCP< const MV > getNativeResiduals(std::vector< MagnitudeType > *norms) const
Get the norms of the residuals native to the solver.
A pure virtual class for defining the status tests for the Belos iterative solvers.
SCT::magnitudeType MagnitudeType
Class which defines basic traits for the operator type.
virtual ~MinresIter()
Destructor.
MinresIter(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem, const Teuchos::RCP< OutputManager< ScalarType > > &printer, const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &tester, const Teuchos::ParameterList ¶ms)
Constructor.
Traits class which defines basic operations on multivectors.
bool isInitialized()
States whether the solver has been initialized or not.
void iterate()
Perform MINRES iterations until convergence or error.
Teuchos::RCP< const MV > R2
Previous residual.
MinresIterateFailure is thrown when the MinresIteration object is unable to compute the next iterate ...
void symOrtho(ScalarType a, ScalarType b, ScalarType *c, ScalarType *s, ScalarType *r)
A linear system to solve, and its associated information.
Class which describes the linear problem to be solved by the iterative solver.
Structure to contain pointers to MinresIteration state variables.
Teuchos::ScalarTraits< MagnitudeType > SMT
MultiVecTraits< ScalarType, MV > MVT
MinresIterationState< ScalarType, MV > getState() const
Get the current state of the linear solver.
Teuchos::RCP< const MV > W2
Previous direction vector.
static magnitudeType magnitude(T a)
Pure virtual base class which augments the basic interface for a minimal residual linear solver itera...
OperatorTraits< ScalarType, MV, OP > OPT
void initialize()
Initialize the solver.
Teuchos::RCP< const MV > Y
The current residual.
Teuchos::RCP< MV > getCurrentUpdate() const
Get the current update to the linear system.
void initializeMinres(const MinresIterationState< ScalarType, MV > &newstate)
Initialize the solver to an iterate, providing a complete state.
int getBlockSize() const
Get the blocksize to be used by the iterative solver in solving this linear problem.
int getNumIters() const
Get the current iteration count.
void setBlockSize(int blockSize)
Set the blocksize to be used by the iterative solver in solving this linear problem.
Class which defines basic traits for the operator type.
const LinearProblem< ScalarType, MV, OP > & getProblem() const
Get a constant reference to the linear problem.
Teuchos::ScalarTraits< ScalarType > SCT
bool isInitialized() const
States whether the solver has been initialized or not.
Belos header file which uses auto-configuration information to include necessary C++ headers...
Teuchos::RCP< const MV > W1
Previous direction vector.
Teuchos::RCP< const MV > W
The current direction vector.