an algebraic multigrid method using a Krylov-cycle.
More...
#include <dune/istl/paamg/amg.hh>
|
virtual SolverCategory::Category | category () const |
| Category of the preconditioner (see SolverCategory::Category) More...
|
|
| KAMG (const OperatorHierarchy &matrices, CoarseSolver &coarseSolver, const SmootherArgs &smootherArgs, std::size_t gamma, std::size_t preSmoothingSteps=1, std::size_t postSmoothingSteps=1, std::size_t maxLevelKrylovSteps=3, double minDefectReduction=1e-1) |
| Construct a new amg with a specific coarse solver. More...
|
|
| KAMG (const OperatorHierarchy &matrices, CoarseSolver &coarseSolver, const SmootherArgs &smootherArgs, const Parameters &parms, std::size_t maxLevelKrylovSteps=3, double minDefectReduction=1e-1) |
| Construct a new amg with a specific coarse solver. More...
|
|
template<class C > |
| KAMG (const Operator &fineOperator, const C &criterion, const SmootherArgs &smootherArgs, std::size_t gamma, std::size_t preSmoothingSteps=1, std::size_t postSmoothingSteps=1, std::size_t maxLevelKrylovSteps=3, double minDefectReduction=1e-1, const ParallelInformation &pinfo=ParallelInformation()) |
| Construct an AMG with an inexact coarse solver based on the smoother. More...
|
|
template<class C > |
| KAMG (const Operator &fineOperator, const C &criterion, const SmootherArgs &smootherArgs=SmootherArgs(), std::size_t maxLevelKrylovSteps=3, double minDefectReduction=1e-1, const ParallelInformation &pinfo=ParallelInformation()) |
| Construct an AMG with an inexact coarse solver based on the smoother. More...
|
|
void | pre (Domain &x, Range &b) |
| Prepare the preconditioner. More...
|
|
void | post (Domain &x) |
| Clean up. More...
|
|
void | apply (Domain &v, const Range &d) |
| Apply one step of the preconditioner to the system A(v)=d. More...
|
|
std::size_t | maxlevels () |
|
template<class M, class X, class S, class PI = SequentialInformation, class K = GeneralizedPCGSolver<X>, class A = std::allocator<X>>
class Dune::Amg::KAMG< M, X, S, PI, K, A >
an algebraic multigrid method using a Krylov-cycle.
The implementation is based on the paper [Notay and Vassilevski, 2007]
- Template Parameters
-
M | The type of the linear operator. |
X | The type of the range and domain. |
PI | The parallel information object. Use SequentialInformation (default) for a sequential AMG, OwnerOverlapCopyCommunication for the parallel case. |
K | The type of the Krylov method to use for the cycle. |
A | The type of the allocator to use. |
◆ Amg
template<class M , class X , class S , class PI = SequentialInformation, class K = GeneralizedPCGSolver<X>, class A = std::allocator<X>>
The type of the underlying AMG.
◆ CoarseSolver
template<class M , class X , class S , class PI = SequentialInformation, class K = GeneralizedPCGSolver<X>, class A = std::allocator<X>>
The type of the coarse solver.
◆ Domain
template<class M , class X , class S , class PI = SequentialInformation, class K = GeneralizedPCGSolver<X>, class A = std::allocator<X>>
◆ domain_type
The domain type of the preconditioner.
◆ field_type
The field type of the preconditioner.
◆ KrylovSolver
template<class M , class X , class S , class PI = SequentialInformation, class K = GeneralizedPCGSolver<X>, class A = std::allocator<X>>
The type of the Krylov solver for the cycle.
◆ Operator
template<class M , class X , class S , class PI = SequentialInformation, class K = GeneralizedPCGSolver<X>, class A = std::allocator<X>>
the type of the lineatr operator.
◆ OperatorHierarchy
template<class M , class X , class S , class PI = SequentialInformation, class K = GeneralizedPCGSolver<X>, class A = std::allocator<X>>
The type of the hierarchy of operators.
◆ ParallelInformation
template<class M , class X , class S , class PI = SequentialInformation, class K = GeneralizedPCGSolver<X>, class A = std::allocator<X>>
the type of the parallelinformation to use.
◆ ParallelInformationHierarchy
template<class M , class X , class S , class PI = SequentialInformation, class K = GeneralizedPCGSolver<X>, class A = std::allocator<X>>
The type of the hierarchy of parallel information.
◆ Range
template<class M , class X , class S , class PI = SequentialInformation, class K = GeneralizedPCGSolver<X>, class A = std::allocator<X>>
◆ range_type
The range type of the preconditioner.
◆ ScalarProduct
template<class M , class X , class S , class PI = SequentialInformation, class K = GeneralizedPCGSolver<X>, class A = std::allocator<X>>
The type of the scalar product.
◆ SmootherArgs
template<class M , class X , class S , class PI = SequentialInformation, class K = GeneralizedPCGSolver<X>, class A = std::allocator<X>>
The type of the arguments for construction of the smoothers.
◆ category()
template<class M , class X , class S , class PI = SequentialInformation, class K = GeneralizedPCGSolver<X>, class A = std::allocator<X>>
The documentation for this class was generated from the following files: