This page was generated from unit-2.1.3-bddc/bddc.ipynb.

2.1.3 Element-wise BDDC Preconditioner

The element-wise BDDC (Balancing Domain Decomposition preconditioner with Constraints) preconditioner in NGSolve is a good general purpose preconditioner that works well both in the shared memory parallel mode as well as in distributed memory mode. In this tutorial, we discuss this preconditioner, related built-in options, and customization from python.

Let us start with a simple description of the element-wise BDDC in the context of Lagrange finite element space \(V\). Here the BDDC preconditoner is constructed on an auxiliary space \(\widetilde{V}\) obtained by connecting only element vertices (leaving edge and face shape functions discontinuous). Although larger, the auxiliary space allows local elimination of edge and face variables. Hence an analogue of the original matrix \(A\) on this space, named \(\widetilde A\), is less expensive to invert. This inverse is used to construct a preconditioner for \(A\) as follows:

\[C_{BDDC}^{-1} = R {\,\widetilde{A}\,}^{-1}\, R^t\]

Here, \(R\) is the averaging operator for the discontinous edge and face variables.

To construct a general purpose BDDC preconditioner, NGSolve generalizes this idea to all its finite element spaces by a classification of degrees of freedom. Dofs are classified into (condensable) LOCAL_DOFs that we saw in 1.4 and a remainder that includes these types:

WIREBASKET_DOF
INTERFACE_DOF

The original finite element space \(V\) is obtained by requiring conformity of both the above types of dofs, while the auxiliary space \(\widetilde{V}\) is obtained by requiring conformity of WIREBASKET_DOFs only.

Here is a figure of a typical function in the default \(\widetilde{V}\) (and the code to generate this is at the end of this tutorial) when \(V\) is the Lagrange space:

title

[1]:
import netgen.gui
from ngsolve import *
from ngsolve.la import EigenValues_Preconditioner
from netgen.csg import unit_cube
from netgen.geom2d import unit_square
SetHeapSize(100*1000*1000)
[2]:
mesh = Mesh(unit_cube.GenerateMesh(maxh=0.5))
# mesh = Mesh(unit_square.GenerateMesh(maxh=0.5))

Built-in options

Let us define a simple function to study how the spectrum of the preconditioned matrix changes with various options.

Effect of condensation

[3]:
def TestPreconditioner (p, condense=False, **args):
    fes = H1(mesh, order=p, **args)
    u,v = fes.TnT()
    a = BilinearForm(fes, eliminate_internal=condense)
    a += grad(u)*grad(v)*dx + u*v*dx
    c = Preconditioner(a, "bddc")
    a.Assemble()
    return EigenValues_Preconditioner(a.mat, c.mat)
[4]:
lams = TestPreconditioner(5)
print (lams[0:3], "...\n", lams[-3:])
 1.00409
 1.11188
 1.34814
 ...
  33.7436
 36.3195
 38.0012

Here is the effect of static condensation on the BDDC preconditioner.

[5]:
lams = TestPreconditioner(5, condense=True)
print (lams[0:3], "...\n", lams[-3:])
 1.00139
 1.04212
 1.12173
 ...
  6.95459
 7.03965
  7.3257

Tuning the auxiliary space

Next, let us study the effect of a few built-in flags for finite element spaces that are useful for tweaking the behavior of the BDDC preconditioner. The effect of these flags varies in two (2D) and three dimensions (3D), e.g.,

  • wb_fulledges=True: This option keeps all edge-dofs conforming (i.e. they are marked WIREBASKET_DOFs). This option is only meaningful in 3D. If used in 2D, the preconditioner becomes a direct solver.

  • wb_withedges=True: This option keeps only the first edge-dof conforming (i.e., the first edge-dof is marked WIREBASKET_DOF and the remaining edge-dofs are marked INTERFACE_DOFs).

The complete situation is a bit more complex due to the fact these options can take the three values True, False, or Undefined, the two options can be combined, and the space dimension can be 2 or 3. The default value of these flags in NGSolve is Undefined. Here is a table with the summary of the effect of these options:

wb_fulledges

wb_withedges

2D

3D

True

any value

all

all

False/Undefined

Undefined

none

first

False/Undefined

False

none

none

False/Undefined

True

first

first

An entry \(X \in\) {all, none, first} of the last two columns is to be read as follows: \(X\) of the edge-dofs is(are) WIREBASKET_DOF(s).

Here is an example of the effect of one of these flag values.

[6]:
lams = TestPreconditioner(5, condense=True,
                          wb_withedges=False)
print (lams[0:3], "...\n", lams[-3:])
 1.00468
 1.09359
 1.32949
 ...
  31.1877
 33.3524
 34.3632

Clearly, when moving from the default case (where the first of the edge dofs are wire basket dofs) to the case (where none of the edge dofs are wire basket dofs), the conditioning became less favorable.

Customize

From within python, we can change the types of degrees of freedom of finite element spaces, thus affecting the behavior of the BDDC preconditioner.

To depart from the default and commit the first two edge dofs to wire basket, we perform the next steps:

[7]:
fes = H1(mesh, order=10)
u,v = fes.TnT()

for ed in mesh.edges:
    dofs = fes.GetDofNrs(ed)
    for d in dofs:
        fes.SetCouplingType(d, COUPLING_TYPE.INTERFACE_DOF)

    # Set the first two edge dofs to be conforming
    fes.SetCouplingType(dofs[0], COUPLING_TYPE.WIREBASKET_DOF)
    fes.SetCouplingType(dofs[1], COUPLING_TYPE.WIREBASKET_DOF)

a = BilinearForm(fes, eliminate_internal=True)
a += grad(u)*grad(v)*dx + u*v*dx
c = Preconditioner(a, "bddc")
a.Assemble()

lams=EigenValues_Preconditioner(a.mat, c.mat)
max(lams)/min(lams)
[7]:
18.310200831463245

This is a slight improvement from the default.

[8]:
lams = TestPreconditioner (10, condense=True)
max(lams)/min(lams)
[8]:
27.55054816091915

Combine BDDC with AMG for large problems

The coarse inverse \({\,\widetilde{A}\,}^{-1}\) of BDDC is expensive on fine meshes. Using the option coarsetype=h1amg flag, we can ask BDDC to replace \({\,\widetilde{A}\,}^{-1}\) by an Algebraic MultiGrid (AMG) preconditioner. Since NGSolve’s h1amg is well-suited
for the lowest order Lagrange space, we use the option wb_withedges=False to ensure that \(\widetilde{A}\) is made solely with vertex unknowns.
[9]:
p = 5
mesh = Mesh(unit_cube.GenerateMesh(maxh=0.05))
fes = H1(mesh, order=p, dirichlet="left|bottom|back",
         wb_withedges=False)
print("NDOF = ", fes.ndof)
u,v = fes.TnT()
a = BilinearForm(fes)
a += grad(u)*grad(v)*dx
f = LinearForm(fes)
f += v*dx

with TaskManager():
    pre = Preconditioner(a, "bddc", coarsetype="h1amg")
    a.Assemble()
    f.Assemble()

    gfu = GridFunction(fes)
    solvers.CG(mat=a.mat, rhs=f.vec, sol=gfu.vec,
               pre=pre, maxsteps=500)
Draw(gfu)
NDOF =  1165206
iteration 0 error = 0.7051536867626971
iteration 1 error = 0.29102936294678683
iteration 2 error = 0.28528781035641354
iteration 3 error = 0.2917142728795673
iteration 4 error = 0.31907585981346387
iteration 5 error = 0.25435673384112695
iteration 6 error = 0.19165042832642773
iteration 7 error = 0.15435152936369143
iteration 8 error = 0.12622007530078627
iteration 9 error = 0.1013959211853625
iteration 10 error = 0.07418677746359831
iteration 11 error = 0.056486511340601925
iteration 12 error = 0.04498125252457271
iteration 13 error = 0.038014188902159025
iteration 14 error = 0.031577434682634783
iteration 15 error = 0.025217475114188027
iteration 16 error = 0.0193052360974247
iteration 17 error = 0.015541452232273545
iteration 18 error = 0.012412281391059155
iteration 19 error = 0.009974955403254984
iteration 20 error = 0.007771809754825912
iteration 21 error = 0.006076247217503247
iteration 22 error = 0.004842451370760304
iteration 23 error = 0.003893638861821674
iteration 24 error = 0.003120531484299098
iteration 25 error = 0.00251878779375597
iteration 26 error = 0.002018525693803704
iteration 27 error = 0.0016143532690698497
iteration 28 error = 0.0012916258860525749
iteration 29 error = 0.0010439900323139394
iteration 30 error = 0.0008376206442348002
iteration 31 error = 0.0006735089153655329
iteration 32 error = 0.0005413105814165692
iteration 33 error = 0.00042718920144069086
iteration 34 error = 0.00034294755957726817
iteration 35 error = 0.00027999675705078786
iteration 36 error = 0.0002235285405747343
iteration 37 error = 0.0001775030560273239
iteration 38 error = 0.00014043165560160363
iteration 39 error = 0.00011110211712027574
iteration 40 error = 8.931810657608277e-05
iteration 41 error = 7.114763395853235e-05
iteration 42 error = 5.677117370412992e-05
iteration 43 error = 4.5095377291277824e-05
iteration 44 error = 3.638908191289604e-05
iteration 45 error = 2.9053521919412603e-05
iteration 46 error = 2.3038045130388117e-05
iteration 47 error = 1.86262685750109e-05
iteration 48 error = 1.5028649777307698e-05
iteration 49 error = 1.2059390216946096e-05
iteration 50 error = 9.738150803348838e-06
iteration 51 error = 7.792649789413204e-06
iteration 52 error = 6.178945048384152e-06
iteration 53 error = 4.918593203198965e-06
iteration 54 error = 3.899572468267664e-06
iteration 55 error = 3.1251758611316553e-06
iteration 56 error = 2.5245280936060532e-06
iteration 57 error = 2.0377143394100593e-06
iteration 58 error = 1.6251994998389534e-06
iteration 59 error = 1.2786907526551923e-06
iteration 60 error = 1.037632368659669e-06
iteration 61 error = 8.320586378966094e-07
iteration 62 error = 6.63933575444813e-07
iteration 63 error = 5.30542117674907e-07
iteration 64 error = 4.2477996782203484e-07
iteration 65 error = 3.385813462325085e-07
iteration 66 error = 2.6912441809161477e-07
iteration 67 error = 2.1673091538655465e-07
iteration 68 error = 1.7316671895337035e-07
iteration 69 error = 1.394422416338081e-07
iteration 70 error = 1.1137097021044991e-07
iteration 71 error = 8.876246106641196e-08
iteration 72 error = 7.09780458563821e-08
iteration 73 error = 5.641008454984268e-08
iteration 74 error = 4.57388283505036e-08
iteration 75 error = 3.645831253793697e-08
iteration 76 error = 2.9577962759448863e-08
iteration 77 error = 2.362233991710163e-08
iteration 78 error = 1.871169684323783e-08
iteration 79 error = 1.5256534882983334e-08
iteration 80 error = 1.2259271368318211e-08
iteration 81 error = 9.737800285853582e-09
iteration 82 error = 7.819841489569124e-09
iteration 83 error = 6.256665763907396e-09
iteration 84 error = 4.963009754394482e-09
iteration 85 error = 3.982911697425353e-09
iteration 86 error = 3.1541784094214878e-09
iteration 87 error = 2.5428426183637242e-09
iteration 88 error = 2.0713012262709986e-09
iteration 89 error = 1.6953420504644155e-09
iteration 90 error = 1.3920527361051839e-09
iteration 91 error = 1.128318909918821e-09
iteration 92 error = 8.947990860511049e-10
iteration 93 error = 7.194970597370553e-10
iteration 94 error = 5.769398329678581e-10
iteration 95 error = 4.5414287437165813e-10
iteration 96 error = 3.646145049515405e-10
iteration 97 error = 2.927891800935106e-10
iteration 98 error = 2.3316516368914693e-10
iteration 99 error = 1.8754567519706953e-10
iteration 100 error = 1.4988320658187106e-10
iteration 101 error = 1.1970921563271457e-10
iteration 102 error = 9.57399443829608e-11
iteration 103 error = 7.715356126780888e-11
iteration 104 error = 6.219996922068764e-11
iteration 105 error = 5.038159770574872e-11
iteration 106 error = 4.116890777863329e-11
iteration 107 error = 3.3581810340181665e-11
iteration 108 error = 2.6592050553271824e-11
iteration 109 error = 2.1513120273965233e-11
iteration 110 error = 1.7403882946364987e-11
iteration 111 error = 1.39013314895141e-11
iteration 112 error = 1.1104965885032501e-11
iteration 113 error = 8.826115013208097e-12
iteration 114 error = 6.991651431969361e-12
iteration 115 error = 5.569051958762189e-12
iteration 116 error = 4.4518732541648345e-12
iteration 117 error = 3.5419066040742178e-12
iteration 118 error = 2.8662879978211275e-12
iteration 119 error = 2.2995099546138864e-12
iteration 120 error = 1.853344415170077e-12
iteration 121 error = 1.5443247550707785e-12
iteration 122 error = 1.2514318112614145e-12
iteration 123 error = 9.803718010382852e-13
iteration 124 error = 7.786654485062723e-13
iteration 125 error = 6.160483827643422e-13

Postscript

By popular demand, here is the code to draw the figure found at the beginning of this tutorial:

[10]:
from netgen.geom2d import unit_square
mesh = Mesh(unit_square.GenerateMesh(maxh=0.1))
fes_ho = Discontinuous(H1(mesh, order=10))
fes_lo = H1(mesh, order=1, dirichlet=".*")
fes_lam = Discontinuous(H1(mesh, order=1))
fes = FESpace([fes_ho, fes_lo, fes_lam])
uho, ulo, lam = fes.TrialFunction()

a = BilinearForm(fes)
a += Variation(0.5 * grad(uho)*grad(uho)*dx
               - 1*uho*dx
               + (uho-ulo)*lam*dx(element_vb=BBND))
gfu = GridFunction(fes)
solvers.Newton(a=a, u=gfu)
Draw(gfu.components[0])
Newton iteration  0
err =  0.39392955469293095
Newton iteration  1
err =  3.99279175725607e-15