DUNE-FEM (unstable)
Dune::Fem::ConjugateGradientSolver< Operator > Class Template Reference
linear solver using the CG algorithm More...
#include <dune/fem/solver/cginverseoperator.hh>
Public Types | |
| typedef Operator | OperatorType | 
| type of the operators to invert  | |
| typedef OperatorType::DomainFieldType | DomainFieldType | 
| field type of the operator's domain vectors  | |
| typedef OperatorType::RangeFieldType | RangeFieldType | 
| field type of the operator's range vectors  | |
| typedef OperatorType::DomainFunctionType | DomainFunctionType | 
| type of the operator's domain vectors  | |
| typedef OperatorType::RangeFunctionType | RangeFunctionType | 
| type of the operator's range vectors  | |
| typedef Fem::Operator< RangeFunctionType, DomainFunctionType > | PreconditionerType | 
| type of the preconditioner, maps from the range of the operator (the dual space) in it's domain  | |
Public Member Functions | |
| ConjugateGradientSolver (const RealType &epsilon, unsigned int maxIterations, int errorMeasure, bool verbose) | |
| constructor  More... | |
| ConjugateGradientSolver (RealType epsilon, unsigned int maxIterations, const ParameterReader ¶meter=Parameter::container()) | |
| constructor  More... | |
| void | solve (const OperatorType &op, const RangeFunctionType &b, DomainFunctionType &x) const | 
| solve \(op( x ) = b\)  More... | |
| void | solve (const OperatorType &op, const PreconditionerType &p, const RangeFunctionType &b, DomainFunctionType &x) const | 
| solve \(op( x ) = b\)  More... | |
| unsigned int | iterations () const | 
| number of iterations needed for last solve  | |
| double | averageCommTime () const | 
| return average communication time during last solve  | |
Detailed Description
template<class Operator>
class Dune::Fem::ConjugateGradientSolver< Operator >
class Dune::Fem::ConjugateGradientSolver< Operator >
linear solver using the CG algorithm
- Parameters
 - 
  
Operator type of the operator to invert  
Constructor & Destructor Documentation
◆ ConjugateGradientSolver() [1/2]
template<class Operator > 
      
  | 
  inline | 
constructor
- Parameters
 - 
  
[in] epsilon tolerance [in] maxIterations maximum number of CG iterations [in] errorMeasure use absolute (0) or relative (1) error [in] verbose verbose output  
◆ ConjugateGradientSolver() [2/2]
template<class Operator > 
      
  | 
  inline | 
constructor
- Parameters
 - 
  
[in] epsilon tolerance [in] maxIterations maximum number of CG iterations  
Member Function Documentation
◆ solve() [1/2]
template<class Operator > 
      
  | 
  inline | 
solve \(op( x ) = b\)
- Note
 - The CG algorithm also works for positive semidefinite operators. In this case, \(x \cdot v = b \cdot v\) for all \(v\) in the operator's kernel.
 
- Parameters
 - 
  
[in] op linear operator to invert (must be symmetic and positive definite) [in] p (lef) preconditioning operator [in] b right hand side x solution (must be initialized to a start value)  
◆ solve() [2/2]
template<class Operator > 
      
  | 
  inline | 
solve \(op( x ) = b\)
- Note
 - The CG algorithm also works for positive semidefinite operators. In this case, \(x \cdot v = b \cdot v\) for all \(v\) in the operator's kernel.
 
- Parameters
 - 
  
[in] op linear operator to invert (must be symmetic and positive definite) [in] b right hand side x solution (must be initialized to a start value)  
Referenced by Dune::Fem::Solver::CGInverseOperator< DiscreteFunction >::apply().
The documentation for this class was generated from the following file:
- dune/fem/solver/cginverseoperator.hh
 
   | 
                                Legal Statements / Impressum  | 
                                Hosted by  TU Dresden & Uni Heidelberg  | 
				  generated with Hugo v0.111.3
								(Nov 3, 23:36, 2025)