KASKADE 7 development version
|
#include <multiGridSolver.hh>
MultigridSolver realizes a classical v-cycle multigrid solver
truncation criterion is relative defect tolerance ||Ax-b||_2/||b||_2 < maxTol
Definition at line 377 of file multiGridSolver.hh.
Classes | |
struct | Parameter |
Public Member Functions | |
MultigridSolver (Dune::BCRSMatrix< MatrixBlock > const &A, Grid const &grid, Parameter params=Parameter()) | |
template<class Assembler , int row, int col, class BlockFilter > | |
MultigridSolver (AssembledGalerkinOperator< Assembler, row, row+1, col, col+1, BlockFilter > const &A, Grid const &grid, Parameter params=Parameter(), bool transposed=false) | |
void | apply (CoeffVector &solution, CoeffVector &rightHand) |
void | setParameter (Parameter p) |
|
inline |
A | const reference to the stiffness matrix on leaf level. |
Definition at line 390 of file multiGridSolver.hh.
|
inline |
Definition at line 397 of file multiGridSolver.hh.
|
inline |
Definition at line 406 of file multiGridSolver.hh.
Referenced by 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::TangentSpacePreconditioner2< Functional, Assembler, components >::applyAdjointPreconditioner(), Kaskade::NormalStepPreconditioner3< Functional, Assembler, components >::applyAdjointPreconditioner(), Kaskade::TangentSpacePreconditioner2< Functional, Assembler, components >::applyStatePreconditioner(), Kaskade::NormalStepPreconditioner3< Functional, Assembler, components >::applyStatePreconditioner(), Kaskade::InexactTangentSpacePreconditioner< Functional, Assembler, components, exactConstraint >::applyStatePreconditioner(), Kaskade::InexactTangentSpacePreconditionerILU< Functional, Assembler, components, exactConstraint >::applyStatePreconditioner(), 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 |
Definition at line 427 of file multiGridSolver.hh.
Referenced by Kaskade::NormalStepPreconditioner3< Functional, Assembler, components >::setParameter().