KASKADE 7 development version
Public Types | Public Member Functions | List of all members
Kaskade::TCG< X, Xstar, Regularization > Class Template Reference

preconditioned conjugate gradient method More...

#include <tcg.hh>

Detailed Description

template<class X, class Xstar, class Regularization = NoRegularization>
class Kaskade::TCG< X, Xstar, Regularization >

preconditioned conjugate gradient method

This implements a preconditioned IterateType::CG iteration for an operator \( A: X\to x^* \), preconditioned by a preconditioner \( B^{-1}: X^* \to X \). The termination is based on an estimate of the absolute energy error.

The implementation follows Deuflhard/Weiser, Section 5.3.3.

Definition at line 79 of file tcg.hh.

Inheritance diagram for Kaskade::TCG< X, Xstar, Regularization >:

Public Types

enum class  Result { Converged , Failed , EncounteredNonConvexity }
 
typedef ScalarTraits< typenameGetScalar< X >::type >::Real Real
 the real field type corresponding to X::field_type More...
 

Public Member Functions

template<class Int = int, class enable = typename std::enable_if<!std::is_same<Regularization,NoRegularization>::value && std::is_same<Int,int>::value>::type>
 TCG (Dune::LinearOperator< X, Xstar > &op_, Dune::Preconditioner< X, Xstar > &prec_, DualPairing< X, Xstar > const &dp_, Regularization &regularization_, double relTol_=1e-3, size_t maxSteps_=100, Int verbose_=0)
 Set up conjugate gradient solver with absolute energy error termination criterion. More...
 
template<class Int = int, class enable = typename std::enable_if<std::is_same<Regularization,NoRegularization>::value && std::is_same<Int,int>::value>::type>
 TCG (Dune::LinearOperator< X, Xstar > &op_, Dune::Preconditioner< X, Xstar > &prec_, DualPairing< X, Xstar > const &dp_, double relTol_=1e-3, size_t maxSteps_=100, Int verbose_=0)
 Set up conjugate gradient solver with absolute energy error termination criterion. More...
 
virtual int apply (X &u, X &q, Xstar &b, Dune::InverseOperatorResult &res)
 Apply inverse operator. More...
 
virtual void apply (X &u, Xstar &b, Dune::InverseOperatorResult &res)
 Apply tcg and possibly forget second descent direction. More...
 
void apply (X &u, Xstar &b)
 
virtual void apply (X &x, X &b, double relTol, Dune::InverseOperatorResult &res)
 Apply inverse operator with given absolute tolerance. More...
 
void setRelativeAccuracy (double relTol_)
 
void setMaxSteps (size_t maxSteps_)
 
bool localConvergenceLikely () const
 
bool encounteredNonConvexity () const
 

Member Typedef Documentation

◆ Real

template<class X , class Xstar , class Regularization = NoRegularization>
typedef ScalarTraits<typenameGetScalar<X>::type>::Real Kaskade::TCG< X, Xstar, Regularization >::Real

the real field type corresponding to X::field_type

Definition at line 85 of file tcg.hh.

Member Enumeration Documentation

◆ Result

template<class X , class Xstar , class Regularization = NoRegularization>
enum class Kaskade::TCG::Result
strong
Enumerator
Converged 
Failed 
EncounteredNonConvexity 

Definition at line 81 of file tcg.hh.

Constructor & Destructor Documentation

◆ TCG() [1/2]

template<class X , class Xstar , class Regularization = NoRegularization>
template<class Int = int, class enable = typename std::enable_if<!std::is_same<Regularization,NoRegularization>::value && std::is_same<Int,int>::value>::type>
Kaskade::TCG< X, Xstar, Regularization >::TCG ( Dune::LinearOperator< X, Xstar > &  op_,
Dune::Preconditioner< X, Xstar > &  prec_,
DualPairing< X, Xstar > const &  dp_,
Regularization &  regularization_,
double  relTol_ = 1e-3,
size_t  maxSteps_ = 100,
Int  verbose_ = 0 
)
inline

Set up conjugate gradient solver with absolute energy error termination criterion.

Parameters
verbose

Definition at line 93 of file tcg.hh.

◆ TCG() [2/2]

template<class X , class Xstar , class Regularization = NoRegularization>
template<class Int = int, class enable = typename std::enable_if<std::is_same<Regularization,NoRegularization>::value && std::is_same<Int,int>::value>::type>
Kaskade::TCG< X, Xstar, Regularization >::TCG ( Dune::LinearOperator< X, Xstar > &  op_,
Dune::Preconditioner< X, Xstar > &  prec_,
DualPairing< X, Xstar > const &  dp_,
double  relTol_ = 1e-3,
size_t  maxSteps_ = 100,
Int  verbose_ = 0 
)
inline

Set up conjugate gradient solver with absolute energy error termination criterion.

Parameters
verbose

Definition at line 110 of file tcg.hh.

Member Function Documentation

◆ apply() [1/4]

template<class X , class Xstar , class Regularization = NoRegularization>
virtual int Kaskade::TCG< X, Xstar, Regularization >::apply ( X &  u,
X &  q,
Xstar &  b,
Dune::InverseOperatorResult &  res 
)
inlinevirtual

Apply inverse operator.

Parameters
uthe initial value (starting iterate)
bthe right hand side

Definition at line 127 of file tcg.hh.

Referenced by Kaskade::TCG< X, Xstar, Regularization >::apply().

◆ apply() [2/4]

template<class X , class Xstar , class Regularization = NoRegularization>
void Kaskade::TCG< X, Xstar, Regularization >::apply ( X &  u,
Xstar &  b 
)
inline

Definition at line 257 of file tcg.hh.

◆ apply() [3/4]

template<class X , class Xstar , class Regularization = NoRegularization>
virtual void Kaskade::TCG< X, Xstar, Regularization >::apply ( X &  u,
Xstar &  b,
Dune::InverseOperatorResult &  res 
)
inlinevirtual

Apply tcg and possibly forget second descent direction.

Definition at line 252 of file tcg.hh.

◆ apply() [4/4]

template<class X , class Xstar , class Regularization = NoRegularization>
virtual void Kaskade::TCG< X, Xstar, Regularization >::apply ( X &  x,
X &  b,
double  relTol,
Dune::InverseOperatorResult &  res 
)
inlinevirtual

Apply inverse operator with given absolute tolerance.

Definition at line 266 of file tcg.hh.

◆ encounteredNonConvexity()

template<class X , class Xstar , class Regularization = NoRegularization>
bool Kaskade::TCG< X, Xstar, Regularization >::encounteredNonConvexity ( ) const
inline

Definition at line 277 of file tcg.hh.

◆ localConvergenceLikely()

template<class X , class Xstar , class Regularization = NoRegularization>
bool Kaskade::TCG< X, Xstar, Regularization >::localConvergenceLikely ( ) const
inline

Definition at line 275 of file tcg.hh.

◆ setMaxSteps()

template<class X , class Xstar , class Regularization = NoRegularization>
void Kaskade::TCG< X, Xstar, Regularization >::setMaxSteps ( size_t  maxSteps_)
inline

Definition at line 273 of file tcg.hh.

◆ setRelativeAccuracy()

template<class X , class Xstar , class Regularization = NoRegularization>
void Kaskade::TCG< X, Xstar, Regularization >::setRelativeAccuracy ( double  relTol_)
inline

Definition at line 271 of file tcg.hh.


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