Dune Core Modules (unstable)
A fast (sequential) algebraic multigrid based on agglomeration that saves memory bandwidth. More...
#include <dune/istl/paamg/fastamg.hh>
Public Types | |
| typedef M | Operator |
| The matrix operator type. | |
| typedef PI | ParallelInformation |
| The type of the parallel information. Either OwnerOverlapCommunication or another type describing the parallel data distribution and providing communication methods. | |
| typedef MatrixHierarchy< M, ParallelInformation, A > | OperatorHierarchy |
| The operator hierarchy type. | |
| typedef OperatorHierarchy::ParallelInformationHierarchy | ParallelInformationHierarchy |
| The parallal data distribution hierarchy type. | |
| typedef X | Domain |
| The domain type. | |
| typedef X | Range |
| The range type. | |
| typedef InverseOperator< X, X > | CoarseSolver |
| the type of the coarse solver. | |
| typedef X | domain_type |
| The domain type of the preconditioner. | |
| typedef X | range_type |
| The range type of the preconditioner. | |
| typedef X::field_type | field_type |
| The field type of the preconditioner. | |
Public Member Functions | |
| FastAMG (OperatorHierarchy &matrices, CoarseSolver &coarseSolver, const Parameters &parms, bool symmetric=true) | |
| Construct a new amg with a specific coarse solver. More... | |
| template<class C > | |
| FastAMG (std::shared_ptr< const Operator > fineOperator, const C &criterion, const Parameters &parms=Parameters(), bool symmetric=true, const ParallelInformation &pinfo=ParallelInformation()) | |
| Construct an AMG with an inexact coarse solver based on the smoother. More... | |
| template<class C > | |
| FastAMG (const Operator &fineOperator, const C &criterion, const Parameters &parms=Parameters(), bool symmetric=true, const ParallelInformation &pinfo=ParallelInformation()) | |
| Construct an AMG with an inexact coarse solver based on the smoother. More... | |
| FastAMG (const FastAMG &amg) | |
| Copy constructor. | |
| void | pre (Domain &x, Range &b) |
| Prepare the preconditioner. More... | |
| void | apply (Domain &v, const Range &d) |
| Apply one step of the preconditioner to the system A(v)=d. More... | |
| virtual SolverCategory::Category | category () const |
| Category of the preconditioner (see SolverCategory::Category) | |
| void | post (Domain &x) |
| Clean up. More... | |
| template<class A1 > | |
| void | getCoarsestAggregateNumbers (std::vector< std::size_t, A1 > &cont) |
| Get the aggregate number of each unknown on the coarsest level. More... | |
| void | recalculateHierarchy () |
| Recalculate the matrix hierarchy. More... | |
| bool | usesDirectCoarseLevelSolver () const |
| Check whether the coarse solver used is a direct solver. More... | |
| virtual void | pre (X &x, X &b)=0 |
| Prepare the preconditioner. More... | |
| virtual void | apply (X &v, const X &d)=0 |
| Apply one step of the preconditioner to the system A(v)=d. More... | |
Detailed Description
class Dune::Amg::FastAMG< M, X, PI, A >
A fast (sequential) algebraic multigrid based on agglomeration that saves memory bandwidth.
It combines one Gauss-Seidel smoothing sweep with the defect calculation to reduce memory transfers.
- Template Parameters
-
M The matrix-operator type, e.g. MatrixAdapter or AssembledLinearOperator X The vector type PI Currently ignored. A An allocator for X
Member Function Documentation
◆ apply()
|
pure virtualinherited |
Apply one step of the preconditioner to the system A(v)=d.
On entry v=0 and d=b-A(x) (although this might not be computed in that way. On exit v contains the update, i.e one step computes \( v = M^{-1} d \) where \( M \) is the approximate inverse of the operator \( A \) characterizing the preconditioner.
- Parameters
-
[out] v The update to be computed d The current defect.
Implemented in Dune::SeqOverlappingSchwarz< M, X, TM, TD, TA >, Dune::SeqOverlappingSchwarz< M, X, TM, TS, TA >, Dune::SeqOverlappingSchwarz< M, X, TM, TD, TA >, and Dune::SeqOverlappingSchwarz< M, X, TM, TS, TA >.
◆ pre()
|
pure virtualinherited |
Prepare the preconditioner.
A solver solves a linear operator equation A(x)=b by applying one or several steps of the preconditioner. The method pre() is called before the first apply operation. b and x are right hand side and solution vector of the linear system respectively. It may. e.g., scale the system, allocate memory or compute a (I)LU decomposition. Note: The ILU decomposition could also be computed in the constructor or with a separate method of the derived method if several linear systems with the same matrix are to be solved.
- Note
- if a preconditioner is copied (e.g. for a second thread) again the pre() method has to be called to ensure proper memory management.
- Parameters
-
x The left hand side of the equation. b The right hand side of the equation.
Implemented in Dune::SeqOverlappingSchwarz< M, X, TM, TD, TA >, and Dune::SeqOverlappingSchwarz< M, X, TM, TS, TA >.
The documentation for this class was generated from the following file:
- dune/istl/paamg/fastamg.hh
|
Legal Statements / Impressum |
Hosted by TU Dresden & Uni Heidelberg |
generated with Hugo v0.111.3
(Nov 2, 23:43, 2025)