87. Algebraic Multigrid Methods#

Algebraic multigrid methods (AMG) build a multigrid hierarchy from the given matrix. In contrast to geometric multigrid methods, they do not need a mesh hierarchy. Just one finite element mesh is enough. This makes them useful for a wide range of of applications.

AMG takes profit from providing the type of problem (Poisson equation, elasticity, Maxwell, …).

NGSolve comes with builtin AMG solvers for scalar equations, and for Maxwell equations. It provides also interfaces to external, parallel AMG solvers (hypre, gamg, …)

87.1. Agglomeration based scalar AMG#

The coarsening of degrees of freedom is steered by the strength of connections between dofs, one may think of a network of resistors. For this, one finds edge-based weights \(w_E\) such that the energy norm is equivalent to the weighted sum of squared differences:

\[ u^T A u \approx \sum_{{\text edges} E} w_E \, (u_{E_1} - u_{E_2})^2, \]

where \(w_E\) is the edge-weight (i.e. the conductivity of each resistor), and \(E_1\) and \(E_2\) are the vertex numbers of the end-points of the edge \(E\). The right hand side is a norm represented by a surrogate matrix \(\tilde A\).

Coarse grid vertices are defined by the index mapping

\[ Ind : \text{Vertex} \rightarrow \text{Cluster} \]

All vertices having the same index are combined to a cluster. Each cluster is a degree of freedom on the coarse grid.

Alternative text

How to select the agglomerates ? The basic principle is to combine vertices to clusters which are connected by a strong conductivity. For each edge we compute the coarsening weight

\[ \frac { w_{E_{ij}} ^2 } { \sum_{k \neq j} w_{E_{ik}} \sum_{k \neq i} w_{E_{ik}} } \]

which is the ratio of the edge-weight, to the geometric mean of the sum of all edge-weights connected to vertices i, and j. One sorts the edges by the coarsening weights, and selects all edges for collapsing whose vertices don’t belong to edges selected before.

87.2. The builtin h1amg preconditioner#

The h1amg preconditioner works for symmetric, scalar problems with nodal degrees of freedom. It uses unsmoothed agglomeration for the generation of coarse spaces.

The first task is to determine the edge-weights. If one has access to element-matrices (instead of the assembled matrix), one has better possibilities. One may compute Schur complements with respect to all edges of each element, which gives a surrogate matrix for each element. Then sum up the weights (conductivities) of all elements sharing the edge.

To have access to element matrices, the setup of the surrogate matrix is included into the assembling loop. Thus, the workflow is to

  • define the biliear-form

  • define the h1amg preconditioner, which registers at the bilinear-form

  • finally assemble the bilinear-form, which also runs the setup of the preconditioner

from ngsolve import *
from ngsolve.la import EigenValues_Preconditioner

with TaskManager():
    mesh = Mesh(unit_cube.GenerateMesh(maxh=0.1))
    for l in range(3): mesh.Refine()    
copy mesh complete
# ngsglobals.msg_level=3
fes = FESpace("nodal", mesh, order=1)
print ("ndof=", fes.ndof)
u,v = fes.TnT()
a = BilinearForm(grad(u)*grad(v)*dx + 1e-3*u*v*dx)
pre = Preconditioner(a, "h1amg")
with TaskManager():
    a.Assemble();
    lam = EigenValues_Preconditioner(a.mat, pre.mat)
    print (list(lam[0:3]), '...', list(lam[-3:-1]))
ndof= 442577
[0.18071473400353516, 0.20112839471453497, 0.24113984238174407] ... [0.9648042230973743, 0.9926905688031753]

87.3. low order AMG, high order smoothing#

from ngsolve import *
from ngsolve.la import EigenValues_Preconditioner
ngsglobals.msg_level=0

with TaskManager():
    mesh = Mesh(unit_cube.GenerateMesh(maxh=0.05))
# ngsglobals.msg_level=3

fes = H1(mesh, order=3) 
print ("ndof=", fes.ndof)
u,v = fes.TnT()
a = BilinearForm(grad(u)*grad(v)*dx + 1e-3*u*v*dx)
pre = Preconditioner(a, "multigrid", coarsetype="h1amg")
with TaskManager():
    a.Assemble()
    lam = EigenValues_Preconditioner(a.mat, pre.mat)
    print (list(lam[0:3]), '...', list(lam[-3:-1]))
ndof= 166262
[0.23482446150988417, 0.24293911496235493, 0.2586032719217159] ... [0.9700774892959718, 0.9906464426157355]

Exercises:

  • test the h1amg preconditioner with more complicated domains, varying coefficients across sub-domains, and additional \(L_2(\Omega)\) and \(L_2(\Gamma)\)-terms