KASKADE 7 development version
|
Preconditioned chebyshev semi iteration. More...
#include <chebyshev.hh>
Preconditioned chebyshev semi iteration.
Standard implementation based on a three-term recurrence and explicit computation of the residuals to avoid accumulation of round of errors in the computation of the residuals. In contrast to Krylov methods this does not slow down convergence (see Gutknecht,Röllin: The Chebyshev iteration revisited).
The Chebyshev semi-iteration requires bounds on the spectrum! Currently we only provide these bounds for tetrahedral elements (3D!) with piecewise linear ansatz functions. For more bounds see a.o. "Wathen: Realistic eigenvalue bounds for the Galerkin mass matrix the spectrum of the preconditioned mass matrix". For reason of consistency with the published bounds the implementation is restricted to a one-step Jacobi preconditioner.
Definition at line 53 of file chebyshev.hh.
Public Types | |
typedef Operator::Domain | Domain |
typedef Operator::Range | Range |
Public Member Functions | |
ChebyshevSemiIteration (Operator const &A_, double tol, size_t steps, size_t verbose_=0) | |
virtual | ~ChebyshevSemiIteration () |
virtual void | pre (Domain &, Range &) |
virtual void | post (Domain &) |
virtual void | apply (Domain &x, Range const &y) |
void | initForMassMatrix_TetrahedralQ1Elements () |
Init spectral bounds for the mass matrix arising from tetrahedral discretization of the domain and linear elements. More... | |
void | setSteps (size_t steps) |
typedef Operator::Domain Kaskade::ChebyshevSemiIteration< Operator, isPreconditioner >::Domain |
Definition at line 56 of file chebyshev.hh.
typedef Operator::Range Kaskade::ChebyshevSemiIteration< Operator, isPreconditioner >::Range |
Definition at line 57 of file chebyshev.hh.
|
inline |
Definition at line 59 of file chebyshev.hh.
|
inlinevirtual |
Definition at line 63 of file chebyshev.hh.
|
inlinevirtual |
Definition at line 68 of file chebyshev.hh.
Referenced by Kaskade::NormalStepPreconditioner< Functional, Assembler >::apply(), Kaskade::TangentSpacePreconditioner< Functional, Assembler, components >::apply(), Kaskade::TangentSpacePreconditioner2< Functional, Assembler, components >::apply(), Kaskade::NormalStepPreconditioner3< Functional, Assembler, components >::apply(), Kaskade::InexactTangentSpacePreconditioner< Functional, Assembler, components, exactConstraint >::apply(), Kaskade::InexactTangentSpacePreconditionerILU< Functional, Assembler, components, exactConstraint >::apply(), Kaskade::YetAnotherHBErrorEstimator< Functional, VariableSetDescription, ExtensionVariableSetDescription, ExtensionSpace, NormFunctional, RefinementStrategy, lump, components, ReferenceSolution, ReferenceOperator >::operator()(), and Kaskade::YetAnotherHBErrorEstimator_Elasticity< Functional, VariableSetDescription, ExtensionVariableSetDescription, ExtensionSpace, NormFunctional, RefinementStrategy, lump, components >::operator()().
|
inline |
Init spectral bounds for the mass matrix arising from tetrahedral discretization of the domain and linear elements.
According to "Wathen: Realistic eigenvalue bounds for the Galerkin mass matrix" the spectrum of the preconditioned mass matrix is contained in \([\frac{1}{2},\frac{5}{2}]\). For chebyshev semi-iteration the bounds on the spectrum are in general given in the form \( [a-b,a+b] \) yielding the coefficients \(a=1.5,\ b=1\). See Wathen,Rees: "Chebyshev semi-iteration in preconditioning for problems including the mass matrix".
Definition at line 133 of file chebyshev.hh.
Referenced by Kaskade::InexactTangentSpacePreconditioner< Functional, Assembler, components, exactConstraint >::InexactTangentSpacePreconditioner(), Kaskade::InexactTangentSpacePreconditionerILU< Functional, Assembler, components, exactConstraint >::InexactTangentSpacePreconditionerILU(), Kaskade::NormalStepPreconditioner< Functional, Assembler >::NormalStepPreconditioner(), Kaskade::NormalStepPreconditioner3< Functional, Assembler, components >::NormalStepPreconditioner3(), Kaskade::YetAnotherHBErrorEstimator< Functional, VariableSetDescription, ExtensionVariableSetDescription, ExtensionSpace, NormFunctional, RefinementStrategy, lump, components, ReferenceSolution, ReferenceOperator >::operator()(), Kaskade::YetAnotherHBErrorEstimator_Elasticity< Functional, VariableSetDescription, ExtensionVariableSetDescription, ExtensionSpace, NormFunctional, RefinementStrategy, lump, components >::operator()(), Kaskade::TangentSpacePreconditioner< Functional, Assembler, components >::TangentSpacePreconditioner(), and Kaskade::TangentSpacePreconditioner2< Functional, Assembler, components >::TangentSpacePreconditioner2().
|
inlinevirtual |
Definition at line 66 of file chebyshev.hh.
|
inlinevirtual |
Definition at line 65 of file chebyshev.hh.
|
inline |
Definition at line 143 of file chebyshev.hh.
Referenced by Kaskade::NormalStepPreconditioner3< Functional, Assembler, components >::setParameter().