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

preconditioned conjugate gradient method More...

#include <apcg.hh>

Detailed Description

template<class X, class Xstar>
class Kaskade::Pcg< X, Xstar >

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 policy-based.

The implementation follows Deuflhard/Weiser, Section 5.3.3.

Template Parameters
Xthe type of vectors from the primal space
Xstarthe type of vectors from the dual space
See also
PCGEnergyErrorTerminationCriterion

Definition at line 280 of file linalg/apcg.hh.

Inheritance diagram for Kaskade::Pcg< X, Xstar >:

Public Types

typedef X domain_type
 The domain type of the operator to be inverted. More...
 
typedef Xstar range_type
 The range type of the operator to be inverted. More...
 
typedef X::field_type field_type
 The field type of the operator to be inverted. More...
 
typedef ScalarTraits< field_type >::Real Real
 the real field type corresponding to field_type More...
 

Public Member Functions

 Pcg (SymmetricLinearOperator< X, Xstar > &op_, SymmetricPreconditioner< X, Xstar > &prec_, PCGTerminationCriterion< Real > &terminate_, int verbose_=0)
 Set up conjugate gradient solver. More...
 
 Pcg (Dune::LinearOperator< X, Xstar > &op_, Dune::Preconditioner< X, Xstar > &prec_, DualPairing< X, Xstar > const &dp_, PCGTerminationCriterion< Real > &terminate_, int verbose_=0)
 Set up conjugate gradient solver. More...
 
virtual void apply (X &u, Xstar &b, Dune::InverseOperatorResult &res)
 Apply inverse operator by performing a number of IterateType::PCG iterations. More...
 
virtual void apply (X &x, Xstar &b, double tol, Dune::InverseOperatorResult &res)
 Apply inverse operator with given tolerance. More...
 
virtual Dune::SolverCategory::Category category () const override
 returns the category of the operator More...
 

Member Typedef Documentation

◆ domain_type

template<class X , class Xstar >
typedef X Kaskade::Pcg< X, Xstar >::domain_type

The domain type of the operator to be inverted.

Definition at line 284 of file linalg/apcg.hh.

◆ field_type

template<class X , class Xstar >
typedef X::field_type Kaskade::Pcg< X, Xstar >::field_type

The field type of the operator to be inverted.

Definition at line 288 of file linalg/apcg.hh.

◆ range_type

template<class X , class Xstar >
typedef Xstar Kaskade::Pcg< X, Xstar >::range_type

The range type of the operator to be inverted.

Definition at line 286 of file linalg/apcg.hh.

◆ Real

template<class X , class Xstar >
typedef ScalarTraits<field_type>::Real Kaskade::Pcg< X, Xstar >::Real

the real field type corresponding to field_type

Definition at line 293 of file linalg/apcg.hh.

Constructor & Destructor Documentation

◆ Pcg() [1/2]

template<class X , class Xstar >
Kaskade::Pcg< X, Xstar >::Pcg ( SymmetricLinearOperator< X, Xstar > &  op_,
SymmetricPreconditioner< X, Xstar > &  prec_,
PCGTerminationCriterion< Real > &  terminate_,
int  verbose_ = 0 
)
inline

Set up conjugate gradient solver.

Parameters
opthe operator
precthe preconditioner
terminatethe termination criterion. The object has to exist during the lifetime of the pcg object as it is referenced.
verbosecontrols the verbosity of logging to std::cout. 0 means no output at all (default). Values in 0,1,2 are valid.

Definition at line 303 of file linalg/apcg.hh.

◆ Pcg() [2/2]

template<class X , class Xstar >
Kaskade::Pcg< X, Xstar >::Pcg ( Dune::LinearOperator< X, Xstar > &  op_,
Dune::Preconditioner< X, Xstar > &  prec_,
DualPairing< X, Xstar > const &  dp_,
PCGTerminationCriterion< Real > &  terminate_,
int  verbose_ = 0 
)
inline

Set up conjugate gradient solver.

Parameters
opthe operator
precthe preconditioner
dpthe dual pairing with respect to which the operator is symmetric
terminatethe termination criterion. The object has to exist during the lifetime of the pcg object as it is referenced.
verbosecontrols the verbosity of logging to std::cout. 0 means no output at all (default). Values in 0,1,2 are valid.

Definition at line 321 of file linalg/apcg.hh.

Member Function Documentation

◆ apply() [1/2]

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

Apply inverse operator by performing a number of IterateType::PCG iterations.

Parameters
uthe initial value (starting iterate)
bthe right hand side (which will not be modified)

Definition at line 334 of file linalg/apcg.hh.

◆ apply() [2/2]

template<class X , class Xstar >
virtual void Kaskade::Pcg< X, Xstar >::apply ( X &  x,
Xstar &  b,
double  tol,
Dune::InverseOperatorResult &  res 
)
inlinevirtual

Apply inverse operator with given tolerance.

This method is equivalent to setting first the tolerance in the termination criterion, then calling apply().

Definition at line 481 of file linalg/apcg.hh.

◆ category()

template<class X , class Xstar >
virtual Dune::SolverCategory::Category Kaskade::Pcg< X, Xstar >::category ( ) const
inlineoverridevirtual

returns the category of the operator

From the Dune doxygen documentation it is unclear what this is supposed to mean. We return a dummy here.

Definition at line 493 of file linalg/apcg.hh.


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