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

Convergence criterion which looks at the weighted absolute value of the residual. More...

#include <weightedresidreductioncriterion.hh>

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

Public Member Functions

 WeightedResidualReductionCriterion (const CollectiveCommunication &comm)
 
 WeightedResidualReductionCriterion (const CollectiveCommunication &comm, const Vector &residWeights, Scalar residualReductionTolerance, Scalar fixPointTolerance, Scalar absResidualTolerance=0.0, Scalar maxError=0.0)
 
void setResidualWeight (const Vector &residWeightVec)
 Sets the relative weight of each row of the residual. More...
 
Scalar residualWeight (size_t outerIdx, unsigned innerIdx) const
 Return the relative weight of a row of the residual. More...
 
void setResidualReductionTolerance (Scalar tol)
 Sets the residual reduction tolerance. More...
 
Scalar residualReductionTolerance () const
 Returns the tolerance of the residual reduction of the solution. More...
 
void setResidualTolerance (Scalar tol)
 Sets the maximum absolute tolerated residual. More...
 
Scalar absResidualTolerance () const
 Returns the maximum absolute tolerated residual. More...
 
Scalar residualAccuracy () const
 Returns the reduction of the weighted maximum of the residual compared to the initial solution. More...
 
void setFixPointTolerance (Scalar tol)
 Sets the fix-point tolerance. More...
 
Scalar fixPointTolerance () const
 Returns the maximum tolerated difference between two iterations to be met before a solution is considered to be converged. More...
 
Scalar fixPointAccuracy () const
 Returns the weighted maximum of the difference between the last two iterative solutions. More...
 
void setInitial (const Vector &curSol, const Vector &curResid)
 Set the initial solution of the linear system of equations. More...
 
void update (const Vector &curSol, const Vector &, const Vector &curResid)
 
bool converged () const
 Returns true if and only if the convergence criterion is met. More...
 
bool failed () const
 Returns true if the convergence criterion cannot be met anymore because the solver has broken down. More...
 
Scalar accuracy () const
 Returns the accuracy of the solution at the last update. More...
 
void printInitial (std::ostream &os=std::cout) const
 Prints the initial information about the convergence behaviour. More...
 
void print (Scalar iter, std::ostream &os=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::WeightedResidualReductionCriterion< Vector, CollectiveCommunication >

Convergence criterion which looks at the weighted absolute value of the residual.

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

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

where $r^k = \mathbf{A} x^k - b $ is the residual for the k-th iterative solution vector $x^k$ and $w_i$ is the weight of the $i$-th linear equation.

In addition, to the reduction of the maximum defect, the linear solver is also considered to be converged, if the defect goes below a given absolute limit.

Constructor & Destructor Documentation

◆ WeightedResidualReductionCriterion() [1/2]

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

◆ WeightedResidualReductionCriterion() [2/2]

template<class Vector , class CollectiveCommunication >
Opm::Linear::WeightedResidualReductionCriterion< Vector, CollectiveCommunication >::WeightedResidualReductionCriterion ( const CollectiveCommunication &  comm,
const Vector &  residWeights,
Scalar  residualReductionTolerance,
Scalar  fixPointTolerance,
Scalar  absResidualTolerance = 0.0,
Scalar  maxError = 0.0 
)
inline

Member Function Documentation

◆ absResidualTolerance()

template<class Vector , class CollectiveCommunication >
Scalar Opm::Linear::WeightedResidualReductionCriterion< Vector, CollectiveCommunication >::absResidualTolerance ( ) const
inline

Returns the maximum absolute tolerated residual.

◆ accuracy()

template<class Vector , class CollectiveCommunication >
Scalar Opm::Linear::WeightedResidualReductionCriterion< 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.

For the accuracy we only take the residual into account,

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

◆ converged()

◆ failed()

template<class Vector , class CollectiveCommunication >
bool Opm::Linear::WeightedResidualReductionCriterion< Vector, CollectiveCommunication >::failed ( ) const
inlinevirtual

Returns true if the convergence criterion cannot be met anymore because the solver has broken down.

Reimplemented from Opm::Linear::ConvergenceCriterion< Vector >.

References Opm::Linear::WeightedResidualReductionCriterion< Vector, CollectiveCommunication >::converged().

◆ fixPointAccuracy()

template<class Vector , class CollectiveCommunication >
Scalar Opm::Linear::WeightedResidualReductionCriterion< Vector, CollectiveCommunication >::fixPointAccuracy ( ) const
inline

Returns the weighted maximum of the difference between the last two iterative solutions.

Referenced by Opm::Linear::WeightedResidualReductionCriterion< Vector, CollectiveCommunication >::print().

◆ fixPointTolerance()

template<class Vector , class CollectiveCommunication >
Scalar Opm::Linear::WeightedResidualReductionCriterion< Vector, CollectiveCommunication >::fixPointTolerance ( ) const
inline

Returns the maximum tolerated difference between two iterations to be met before a solution is considered to be converged.

◆ print()

template<class Vector , class CollectiveCommunication >
void Opm::Linear::WeightedResidualReductionCriterion< Vector, CollectiveCommunication >::print ( Scalar  iter,
std::ostream &  os = std::cout 
) const
inlinevirtual

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 from Opm::Linear::ConvergenceCriterion< Vector >.

References Opm::Linear::WeightedResidualReductionCriterion< Vector, CollectiveCommunication >::fixPointAccuracy(), and Opm::Linear::WeightedResidualReductionCriterion< Vector, CollectiveCommunication >::residualAccuracy().

◆ printInitial()

template<class Vector , class CollectiveCommunication >
void Opm::Linear::WeightedResidualReductionCriterion< Vector, CollectiveCommunication >::printInitial ( std::ostream &  os = std::cout) const
inlinevirtual

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 from Opm::Linear::ConvergenceCriterion< Vector >.

◆ residualAccuracy()

template<class Vector , class CollectiveCommunication >
Scalar Opm::Linear::WeightedResidualReductionCriterion< Vector, CollectiveCommunication >::residualAccuracy ( ) const
inline

◆ residualReductionTolerance()

template<class Vector , class CollectiveCommunication >
Scalar Opm::Linear::WeightedResidualReductionCriterion< Vector, CollectiveCommunication >::residualReductionTolerance ( ) const
inline

Returns the tolerance of the residual reduction of the solution.

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

◆ residualWeight()

template<class Vector , class CollectiveCommunication >
Scalar Opm::Linear::WeightedResidualReductionCriterion< Vector, CollectiveCommunication >::residualWeight ( size_t  outerIdx,
unsigned  innerIdx 
) const
inline

Return the relative weight of a row of the residual.

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

◆ setFixPointTolerance()

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

Sets the fix-point tolerance.

◆ setInitial()

template<class Vector , class CollectiveCommunication >
void Opm::Linear::WeightedResidualReductionCriterion< Vector, CollectiveCommunication >::setInitial ( const Vector &  curSol,
const Vector &  curResid 
)
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 >.

◆ setResidualReductionTolerance()

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

Sets the residual reduction tolerance.

◆ setResidualTolerance()

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

Sets the maximum absolute tolerated residual.

◆ setResidualWeight()

template<class Vector , class CollectiveCommunication >
void Opm::Linear::WeightedResidualReductionCriterion< Vector, CollectiveCommunication >::setResidualWeight ( const Vector &  residWeightVec)
inline

Sets the relative weight of each row of the residual.

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

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

where $r^k = \mathbf{A} x^k - b $ is the residual for the k-th iterative solution vector $x^k$ and $w_i$ is the weight of the $i$-th linear equation.

This method is not part of the generic ConvergenceCriteria interface.

Parameters
residWeightVecA Dune::BlockVector<Dune::FieldVector<Scalar, n> > with the relative weights of the linear equations

◆ update()

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

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