Opm::Linear::FixPointCriterion< Vector, CollectiveCommunication > Class Template Reference

Provides a convergence criterion for the linear solvers which looks at the weighted maximum of the difference between two iterations. More...

#include <fixpointcriterion.hh>

Inheritance diagram for Opm::Linear::FixPointCriterion< Vector, CollectiveCommunication >:
Inheritance graph

Public Member Functions

 FixPointCriterion (const CollectiveCommunication &comm)
 
 FixPointCriterion (const CollectiveCommunication &comm, const Vector &weightVec, Scalar reduction)
 
void setWeight (const Vector &weightVec)
 Sets the relative weight of a primary variable. More...
 
Scalar weight (int outerIdx, int innerIdx) const
 Return the relative weight of a primary variable. More...
 
void setTolerance (Scalar tol)
 Set the maximum allowed weighted maximum difference between two iterations. More...
 
Scalar tolerance () const
 Return the maximum allowed weighted difference between two iterations for the solution considered to be converged. More...
 
void setInitial (const Vector &curSol, const Vector &)
 Set the initial solution of the linear system of equations. More...
 
void update (const Vector &curSol, const Vector &, const Vector &)
 Update the internal members of the convergence criterion with the current solution. More...
 
bool converged () const
 Returns true if and only if the convergence criterion is met. More...
 
Scalar accuracy () const
 Returns the accuracy of the solution at the last update. More...
 
virtual bool failed () const
 Returns true if the convergence criterion cannot be met anymore because the solver has broken down. More...
 
virtual void printInitial (std::ostream &=std::cout) const
 Prints the initial information about the convergence behaviour. More...
 
virtual void print (Scalar, std::ostream &=std::cout) const
 Prints the information about the convergence behaviour for the current iteration. More...
 

Detailed Description

template<class Vector, class CollectiveCommunication>
class Opm::Linear::FixPointCriterion< Vector, CollectiveCommunication >

Provides a convergence criterion for the linear solvers which looks at the weighted maximum of the difference between two iterations.

For the FixPointCriterion, the error of the solution is defined as

\[ e^k = \max_i\{ \left| w_i \delta^k_i \right| \}\;, \]

where $\delta = x^k - x^{k + 1} $ is the difference between two consequtive iterative solution vectors $x^k$ and $x^{k + 1}$ and $w_i$ is the weight of the $i$-th degree of freedom.

This criterion requires that the block type of the vector is a Dune::FieldVector

Constructor & Destructor Documentation

◆ FixPointCriterion() [1/2]

template<class Vector , class CollectiveCommunication >
Opm::Linear::FixPointCriterion< Vector, CollectiveCommunication >::FixPointCriterion ( const CollectiveCommunication &  comm)
inline

◆ FixPointCriterion() [2/2]

template<class Vector , class CollectiveCommunication >
Opm::Linear::FixPointCriterion< Vector, CollectiveCommunication >::FixPointCriterion ( const CollectiveCommunication &  comm,
const Vector &  weightVec,
Scalar  reduction 
)
inline

Member Function Documentation

◆ accuracy()

template<class Vector , class CollectiveCommunication >
Scalar Opm::Linear::FixPointCriterion< Vector, CollectiveCommunication >::accuracy ( ) const
inlinevirtual

Returns the accuracy of the solution at the last update.

A value of zero means that the solution was exact.

Implements Opm::Linear::ConvergenceCriterion< Vector >.

Referenced by Opm::Linear::FixPointCriterion< Vector, CollectiveCommunication >::converged().

◆ converged()

template<class Vector , class CollectiveCommunication >
bool Opm::Linear::FixPointCriterion< Vector, CollectiveCommunication >::converged ( ) const
inlinevirtual

◆ failed()

template<class Vector >
virtual bool Opm::Linear::ConvergenceCriterion< Vector >::failed ( ) const
inlinevirtualinherited

◆ print()

template<class Vector >
virtual void Opm::Linear::ConvergenceCriterion< Vector >::print ( Scalar  ,
std::ostream &  = std::cout 
) const
inlinevirtualinherited

Prints the information about the convergence behaviour for the current iteration.

Parameters
iterThe iteration number. The semantics of this parameter are chosen by the linear solver.
osThe output stream to which the message gets written.

Reimplemented in Opm::Linear::ResidReductionCriterion< Vector >, and Opm::Linear::WeightedResidualReductionCriterion< Vector, CollectiveCommunication >.

Referenced by Opm::Linear::BiCGStabSolver< LinearOperator, Vector, Preconditioner >::apply().

◆ printInitial()

template<class Vector >
virtual void Opm::Linear::ConvergenceCriterion< Vector >::printInitial ( std::ostream &  = std::cout) const
inlinevirtualinherited

Prints the initial information about the convergence behaviour.

This method is called after setInitial() if the solver thinks it's a good idea to be verbose. In practice, "printing the initial information" means printing column headers and the initial state.

Parameters
osThe output stream to which the message gets written.

Reimplemented in Opm::Linear::ResidReductionCriterion< Vector >, Opm::Linear::WeightedResidualReductionCriterion< Vector, CollectiveCommunication >, and Opm::Linear::CombinedCriterion< Vector, CollectiveCommunication >.

Referenced by Opm::Linear::BiCGStabSolver< LinearOperator, Vector, Preconditioner >::apply().

◆ setInitial()

template<class Vector , class CollectiveCommunication >
void Opm::Linear::FixPointCriterion< Vector, CollectiveCommunication >::setInitial ( const Vector &  curSol,
const Vector &   
)
inlinevirtual

Set the initial solution of the linear system of equations.

This version of the method does NOT take the two-norm of the residual as argument. If the two-norm of the defect is available for the linear solver, the version of the update() method with it should be called.

Parameters
curSolThe current iterative solution of the linear system of equations
curResidThe residual vector of the current iterative solution of the linear system of equations

Implements Opm::Linear::ConvergenceCriterion< Vector >.

◆ setTolerance()

template<class Vector , class CollectiveCommunication >
void Opm::Linear::FixPointCriterion< Vector, CollectiveCommunication >::setTolerance ( Scalar  tol)
inline

Set the maximum allowed weighted maximum difference between two iterations.

Set the maximum allowed maximum difference between two iterationsfor the solution considered to be converged.

◆ setWeight()

template<class Vector , class CollectiveCommunication >
void Opm::Linear::FixPointCriterion< Vector, CollectiveCommunication >::setWeight ( const Vector &  weightVec)
inline

Sets the relative weight of a primary variable.

For the FixPointCriterion, the error of the solution is defined as

\[ e^k = \max_i\{ \left| w_i \delta^k_i \right| \}\;, \]

where $\delta = x^k - x^{k + 1} $ is the difference between two consequtive iterative solution vectors $x^k$ and $x^{k + 1}$ and $w_i$ is the weight of the $i$-th degree of freedom.

This method is specific to the FixPointCriterion.

Parameters
weightVecA Dune::BlockVector<Dune::FieldVector<Scalar, n> > with the relative weights of the degrees of freedom

◆ tolerance()

template<class Vector , class CollectiveCommunication >
Scalar Opm::Linear::FixPointCriterion< Vector, CollectiveCommunication >::tolerance ( ) const
inline

Return the maximum allowed weighted difference between two iterations for the solution considered to be converged.

Referenced by Opm::Linear::FixPointCriterion< Vector, CollectiveCommunication >::converged().

◆ update()

template<class Vector , class CollectiveCommunication >
void Opm::Linear::FixPointCriterion< Vector, CollectiveCommunication >::update ( const Vector &  curSol,
const Vector &  ,
const Vector &   
)
inlinevirtual

Update the internal members of the convergence criterion with the current solution.

This version of the method does NOT take the two-norm of the residual as argument. If the two-norm of the defect is available for the linear solver, the version of the update() method with it should be called.

Parameters
curSolThe current iterative solution of the linear system of equations
changeIndicatorA vector where all non-zero values indicate that the solution has changed since the last iteration.
curResidThe residual vector of the current iterative solution of the linear system of equations

Implements Opm::Linear::ConvergenceCriterion< Vector >.

References Opm::Linear::FixPointCriterion< Vector, CollectiveCommunication >::weight().

◆ weight()

template<class Vector , class CollectiveCommunication >
Scalar Opm::Linear::FixPointCriterion< Vector, CollectiveCommunication >::weight ( int  outerIdx,
int  innerIdx 
) const
inline

Return the relative weight of a primary variable.

For the FixPointCriterion, the error of the solution is defined as

\[ e^k = \max_i\{ \left| w_i \delta^k_i \right| \}\;, \]

where $\delta = x^k - x^{k + 1} $ is the difference between two consequtive iterative solution vectors $x^k$ and $x^{k + 1}$ and $w_i$ is the weight of the $i$-th degree of freedom.

This method is specific to the FixPointCriterion.

Parameters
outerIdxThe index of the outer vector (i.e. Dune::BlockVector)
innerIdxThe index of the inner vector (i.e. Dune::FieldVector)

Referenced by Opm::Linear::FixPointCriterion< Vector, CollectiveCommunication >::update().


The documentation for this class was generated from the following file: