10 #ifndef BELOS_MINRES_ITER_HPP
11 #define BELOS_MINRES_ITER_HPP
61 template<
class ScalarType,
class MV,
class OP>
150 throw std::logic_error(
"getState() cannot be called unless "
151 "the state has been initialized");
181 std::vector<MagnitudeType>& theNorms = *norms;
182 if (theNorms.size() < 1)
184 theNorms[0] = phibar_;
186 return Teuchos::null;
195 void symOrtho( ScalarType a, ScalarType b, ScalarType *c, ScalarType *s, ScalarType *r );
211 "Belos::MinresIter::setBlockSize(): Cannot use a block size that is not one.");
251 bool stateStorageInitialized_;
292 template<
class ScalarType,
class MV,
class OP>
301 stateStorageInitialized_(false),
309 template <
class ScalarType,
class MV,
class OP>
312 if (!stateStorageInitialized_) {
317 if (lhsMV == Teuchos::null && rhsMV == Teuchos::null) {
318 stateStorageInitialized_ =
false;
325 if (Y_ == Teuchos::null) {
329 std::invalid_argument,
330 "Belos::MinresIter::setStateSize(): linear problem does not specify multivectors to clone from.");
331 Y_ = MVT::Clone( *tmp, 1 );
332 R1_ = MVT::Clone( *tmp, 1 );
333 R2_ = MVT::Clone( *tmp, 1 );
334 W_ = MVT::Clone( *tmp, 1 );
335 W1_ = MVT::Clone( *tmp, 1 );
336 W2_ = MVT::Clone( *tmp, 1 );
339 stateStorageInitialized_ =
true;
347 template <
class ScalarType,
class MV,
class OP>
351 if (!stateStorageInitialized_)
355 std::invalid_argument,
356 "Belos::MinresIter::initialize(): Cannot initialize state storage!" );
359 std::invalid_argument,
360 "Belos::MinresIter::initialize(): MinresIterationState does not have initial residual.");
362 std::string errstr(
"Belos::MinresIter::initialize(): Specified multivectors must have a consistent length and width.");
364 std::invalid_argument,
367 std::invalid_argument,
371 const ScalarType one = SCT::one();
377 MVT::Assign( *newstate.
Y, *R2_ );
378 MVT::Assign( *newstate.
Y, *R1_ );
382 MVT::MvInit ( *W2_ );
384 if ( lp_->getLeftPrec() != Teuchos::null ) {
385 lp_->applyLeftPrec( *newstate.
Y, *Y_ );
386 if ( lp_->getRightPrec() != Teuchos::null ) {
388 lp_->applyRightPrec( *tmp, *Y_ );
391 else if ( lp_->getRightPrec() != Teuchos::null ) {
392 lp_->applyRightPrec( *newstate.
Y, *Y_ );
395 if (newstate.
Y != Y_) {
397 MVT::Assign( *newstate.
Y, *Y_ );
403 MVT::MvTransMv( one, *newstate.
Y, *Y_, beta1_ );
406 std::invalid_argument,
407 "The preconditioner is not positive definite." );
409 if( SCT::magnitude(beta1_(0,0)) == m_zero )
413 MVT::MvInit( *cur_soln_vec );
416 beta1_(0,0) = SCT::squareroot( beta1_(0,0) );
425 template <
class ScalarType,
class MV,
class OP>
431 if (initialized_ ==
false) {
436 const ScalarType one = SCT::one();
437 const ScalarType zero = SCT::zero();
446 ScalarType oldBeta = zero;
447 ScalarType epsln = zero;
448 ScalarType cs = -one;
449 ScalarType sn = zero;
450 ScalarType dbar = zero;
469 "Belos::MinresIter::iterate(): current linear system has more than one vector!" );
474 while (stest_->checkStatus(
this) !=
Passed) {
481 MVT::MvAddMv (one / beta(0,0), *Y_, zero, *Y_, *V);
484 lp_->applyOp (*V, *Y_);
487 MVT::MvAddMv (one, *Y_, -beta(0,0)/oldBeta, *R1_, *Y_);
490 MVT::MvTransMv (one, *V, *Y_, alpha);
493 MVT::MvAddMv (one, *Y_, -alpha(0,0)/beta(0,0), *R2_, *Y_);
503 if ( lp_->getLeftPrec() != Teuchos::null ) {
504 lp_->applyLeftPrec( *R2_, *Y_ );
505 if ( lp_->getRightPrec() != Teuchos::null ) {
507 lp_->applyRightPrec( *tmp, *Y_ );
510 else if ( lp_->getRightPrec() != Teuchos::null ) {
511 lp_->applyRightPrec( *R2_, *Y_ );
514 MVT::Assign( *R2_, *Y_ );
519 MVT::MvTransMv( one, *R2_, *Y_, beta );
533 "Belos::MinresIter::iterate(): Encountered negative "
534 "value " << beta(0,0) <<
" for r2^H*M*r2 at itera"
535 "tion " << iter_ <<
": MINRES cannot continue." );
536 beta(0,0) = SCT::squareroot( beta(0,0) );
544 delta = cs*dbar + sn*alpha(0,0);
545 gbar = sn*dbar - cs*alpha(0,0);
546 epsln = sn*beta(0,0);
547 dbar = - cs*beta(0,0);
550 this->symOrtho(gbar, beta(0,0), &cs, &sn, &gamma);
557 MVT::Assign( *W_, *W1_ );
564 MVT::MvAddMv( one, *V, -oldeps, *W1_, *W_ );
565 MVT::MvAddMv( one, *W_, -delta, *W2_, *W_ );
566 MVT::MvScale( *W_, one / gamma );
570 MVT::MvAddMv( one, *cur_soln_vec, phi, *W_, *cur_soln_vec );
571 lp_->updateSolution();
581 template <
class ScalarType,
class MV,
class OP>
583 ScalarType *c, ScalarType *s, ScalarType *r
586 const ScalarType one = SCT::one();
587 const ScalarType zero = SCT::zero();
591 if ( absB == m_zero ) {
594 if ( absA == m_zero )
598 }
else if ( absA == m_zero ) {
602 }
else if ( absB >= absA ) {
603 ScalarType tau = a / b;
605 *s = -one / SCT::squareroot( one+tau*tau );
607 *s = one / SCT::squareroot( one+tau*tau );
611 ScalarType tau = b / a;
613 *c = -one / SCT::squareroot( one+tau*tau );
615 *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.