DUNE-FEM (unstable)

Dune::Fem::CGInverseOperator< DiscreteFunction, Op > Class Template Reference

Inverse operator base on CG method. Uses a runtime parameter fem.preconditioning which enables diagonal preconditioning if diagonal matrix entries are available, i.e., Op :: assembled is true. More...

#include <dune/fem/solver/cginverseoperator.hh>

Public Types

typedef Op OperatorType
 type of operator
 
typedef DomainFunction::RangeFieldType DomainFieldType
 field type of the operator's domain
 

Public Member Functions

 CGInverseOperator (RealType redEps, RealType absLimit, unsigned int maxIter, bool verbose, const ParameterReader &parameter=Parameter::container())
 constructor of CGInverseOperator More...
 
 CGInverseOperator (RealType redEps, RealType absLimit, unsigned int maxIter, const ParameterReader &parameter=Parameter::container())
 constructor of CGInverseOperator More...
 
template<class LinearOperator , std::enable_if_t< std::is_base_of< OperatorType, LinearOperator >::value, int > = 0>
 CGInverseOperator (const LinearOperator &op, RealType redEps, RealType absLimit, unsigned int maxIter, bool verbose, const ParameterReader &parameter=Parameter::container())
 constructor of CGInverseOperator More...
 
template<class LinearOperator , std::enable_if_t< std::is_base_of< OperatorType, LinearOperator >::value, int > = 0>
 CGInverseOperator (const LinearOperator &op, RealType redEps, RealType absLimit, unsigned int maxIter, const ParameterReader &parameter=Parameter::container())
 constructor of CGInverseOperator More...
 
 CGInverseOperator (const OperatorType &op, const PreconditioningType &precond, RealType redEps, RealType absLimit, unsigned int maxIter, bool verbose, const ParameterReader &parameter=Parameter::container())
 constructor of CGInverseOperator More...
 
 CGInverseOperator (const OperatorType &op, const PreconditioningType &precond, RealType redEps, RealType absLimit, unsigned int maxIter, const ParameterReader &parameter=Parameter::container())
 constructor of CGInverseOperator More...
 
virtual void operator() (const DomainFunctionType &arg, RangeFunctionType &dest) const
 application operator More...
 
virtual void apply (const DomainFunctionType &arg, RangeFunctionType &dest) const
 application operator More...
 
unsigned int iterations () const
 number of iterations needed for last solve
 
double averageCommTime () const
 return average communication time during last solve
 
virtual void finalize ()
 finalization of operator More...
 
virtual bool nonlinear () const
 

Detailed Description

template<class DiscreteFunction, class Op = Fem::Operator< DiscreteFunction, DiscreteFunction >>
class Dune::Fem::CGInverseOperator< DiscreteFunction, Op >

Inverse operator base on CG method. Uses a runtime parameter fem.preconditioning which enables diagonal preconditioning if diagonal matrix entries are available, i.e., Op :: assembled is true.

Constructor & Destructor Documentation

◆ CGInverseOperator() [1/6]

template<class DiscreteFunction , class Op = Fem::Operator< DiscreteFunction, DiscreteFunction >>
Dune::Fem::CGInverseOperator< DiscreteFunction, Op >::CGInverseOperator ( RealType  redEps,
RealType  absLimit,
unsigned int  maxIter,
bool  verbose,
const ParameterReader &  parameter = Parameter::container() 
)
inline

constructor of CGInverseOperator

Parameters
[in]redEpsreduction epsilon
[in]absLimitabsolut limit of residual
[in]maxItermaximum number of iteration steps
[in]verboseverbosity

◆ CGInverseOperator() [2/6]

template<class DiscreteFunction , class Op = Fem::Operator< DiscreteFunction, DiscreteFunction >>
Dune::Fem::CGInverseOperator< DiscreteFunction, Op >::CGInverseOperator ( RealType  redEps,
RealType  absLimit,
unsigned int  maxIter,
const ParameterReader &  parameter = Parameter::container() 
)
inline

constructor of CGInverseOperator

Parameters
[in]redEpsreduction epsilon
[in]absLimitabsolut limit of residual
[in]maxItermaximum number of iteration steps

◆ CGInverseOperator() [3/6]

template<class DiscreteFunction , class Op = Fem::Operator< DiscreteFunction, DiscreteFunction >>
template<class LinearOperator , std::enable_if_t< std::is_base_of< OperatorType, LinearOperator >::value, int > = 0>
Dune::Fem::CGInverseOperator< DiscreteFunction, Op >::CGInverseOperator ( const LinearOperator op,
RealType  redEps,
RealType  absLimit,
unsigned int  maxIter,
bool  verbose,
const ParameterReader &  parameter = Parameter::container() 
)
inline

constructor of CGInverseOperator

Parameters
[in]opoperator to invert
[in]redEpsreduction epsilon
[in]absLimitabsolut limit of residual
[in]maxItermaximum number of iteration steps
[in]verboseverbosity

◆ CGInverseOperator() [4/6]

template<class DiscreteFunction , class Op = Fem::Operator< DiscreteFunction, DiscreteFunction >>
template<class LinearOperator , std::enable_if_t< std::is_base_of< OperatorType, LinearOperator >::value, int > = 0>
Dune::Fem::CGInverseOperator< DiscreteFunction, Op >::CGInverseOperator ( const LinearOperator op,
RealType  redEps,
RealType  absLimit,
unsigned int  maxIter,
const ParameterReader &  parameter = Parameter::container() 
)
inline

constructor of CGInverseOperator

Parameters
[in]opoperator to invert
[in]redEpsreduction epsilon
[in]absLimitabsolut limit of residual
[in]maxItermaximum number of iteration steps

◆ CGInverseOperator() [5/6]

template<class DiscreteFunction , class Op = Fem::Operator< DiscreteFunction, DiscreteFunction >>
Dune::Fem::CGInverseOperator< DiscreteFunction, Op >::CGInverseOperator ( const OperatorType op,
const PreconditioningType precond,
RealType  redEps,
RealType  absLimit,
unsigned int  maxIter,
bool  verbose,
const ParameterReader &  parameter = Parameter::container() 
)
inline

constructor of CGInverseOperator

Parameters
[in]opoperator to invert
[in]precondprecondition operator
[in]redEpsreduction epsilon
[in]absLimitabsolut limit of residual
[in]maxItermaximum number of iteration steps

◆ CGInverseOperator() [6/6]

template<class DiscreteFunction , class Op = Fem::Operator< DiscreteFunction, DiscreteFunction >>
Dune::Fem::CGInverseOperator< DiscreteFunction, Op >::CGInverseOperator ( const OperatorType op,
const PreconditioningType precond,
RealType  redEps,
RealType  absLimit,
unsigned int  maxIter,
const ParameterReader &  parameter = Parameter::container() 
)
inline

constructor of CGInverseOperator

Parameters
[in]opoperator to invert
[in]precondprecondition operator
[in]redEpsreduction epsilon
[in]absLimitabsolut limit of residual
[in]maxItermaximum number of iteration steps

Member Function Documentation

◆ apply()

template<class DiscreteFunction >
virtual void Dune::Fem::Solver::CGInverseOperator< DiscreteFunction >::apply ( const DomainFunctionType &  arg,
RangeFunctionType &  dest 
) const
inlinevirtualinherited

application operator

The application operator actually solves the linear system \(op(dest) = arg\) using the CG method.

Parameters
[in]argargument discrete function
[out]destdestination discrete function

References Dune::Fem::ConjugateGradientSolver< Operator >::solve().

Referenced by Dune::Fem::Solver::CGInverseOperator< DiscreteFunction >::operator()().

◆ finalize()

virtual void Dune::Fem::Operator< DiscreteFunction , DiscreteFunction >::finalize ( )
inlinevirtualinherited

finalization of operator

Note
The default implementation is empty.

◆ nonlinear()

virtual bool Dune::Fem::Operator< DiscreteFunction , DiscreteFunction >::nonlinear ( ) const
inlinevirtualinherited

Return true if the Operator is nonlinear and false otherwise (default is true).

◆ operator()()

template<class DiscreteFunction >
virtual void Dune::Fem::Solver::CGInverseOperator< DiscreteFunction >::operator() ( const DomainFunctionType &  arg,
RangeFunctionType &  dest 
) const
inlinevirtualinherited

application operator

The application operator actually solves the linear system \(op(dest) = arg\) using the CG method.

Parameters
[in]argargument discrete function
[out]destdestination discrete function

Implements Dune::Fem::Operator< DiscreteFunction, DiscreteFunction >.

References Dune::Fem::Solver::CGInverseOperator< DiscreteFunction >::apply(), and Dune::Fem::Operator< DomainFunction, RangeFunction >::finalize().


The documentation for this class was generated from the following file:
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 13, 23:29, 2024)