This page was generated from unit-2.1.4-bddc/bddc.ipynb.
2.1.4 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 preconditioner 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:
Here, \(R\) is the averaging operator for the discontinuous 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_DOF
s 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_DOF
s 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:
[1]:
from ngsolve import *
from ngsolve.webgui import Draw
from ngsolve.la import EigenValues_Preconditioner
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.00052
1.03338
1.09991
...
4.55515
4.60067
4.79194
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.00023
1.02753
1.08363
...
4.06546
4.12581
4.22957
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 markedWIREBASKET_DOF
s). 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 markedWIREBASKET_DOF
and the remaining edge-dofs are markedINTERFACE_DOF
s).
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.00106
1.07035
1.22031
...
25.4058
25.406
26.9095
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]:
9.989385160527883
This is a slight improvement from the default.
[8]:
lams = TestPreconditioner (10, condense=True)
max(lams)/min(lams)
[8]:
13.974958344241331
Combine BDDC with AMG for large problems¶
coarsetype=h1amg
flag, we can ask BDDC to replace \({\,\widetilde{A}\,}^{-1}\) by an Algebraic MultiGrid (AMG) preconditioner. Since NGSolve’s h1amg
is well-suitedwb_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 = 747331
CG iteration 1, residual = 0.728413955096221
CG iteration 2, residual = 0.303806463397543
CG iteration 3, residual = 0.273157326519862
CG iteration 4, residual = 0.3216444680361296
CG iteration 5, residual = 0.33621375377763846
CG iteration 6, residual = 0.23091380686774454
CG iteration 7, residual = 0.1775817134464141
CG iteration 8, residual = 0.1457558490355257
CG iteration 9, residual = 0.10943345664815456
CG iteration 10, residual = 0.0798728796077207
CG iteration 11, residual = 0.05836078710530294
CG iteration 12, residual = 0.04630036174594882
CG iteration 13, residual = 0.037001482312351994
CG iteration 14, residual = 0.028875962787518954
CG iteration 15, residual = 0.021947291835534472
CG iteration 16, residual = 0.01688036939339385
CG iteration 17, residual = 0.013424016729028286
CG iteration 18, residual = 0.010407276786618043
CG iteration 19, residual = 0.007831119292639927
CG iteration 20, residual = 0.006068982634319486
CG iteration 21, residual = 0.004699791356400368
CG iteration 22, residual = 0.0035013081298882216
CG iteration 23, residual = 0.002688199833842728
CG iteration 24, residual = 0.0020913709940973877
CG iteration 25, residual = 0.0016265345742916685
CG iteration 26, residual = 0.0012407145818760833
CG iteration 27, residual = 0.0009742173263010199
CG iteration 28, residual = 0.0007649510758382971
CG iteration 29, residual = 0.0005840071544650607
CG iteration 30, residual = 0.00044828280645931625
CG iteration 31, residual = 0.00034441125654053115
CG iteration 32, residual = 0.0002614626901861619
CG iteration 33, residual = 0.00020256017312321543
CG iteration 34, residual = 0.00016069689398918858
CG iteration 35, residual = 0.0001246433042074301
CG iteration 36, residual = 9.349136831067515e-05
CG iteration 37, residual = 7.34444343353459e-05
CG iteration 38, residual = 6.474983270181665e-05
CG iteration 39, residual = 5.5173867890579765e-05
CG iteration 40, residual = 3.845130216038328e-05
CG iteration 41, residual = 2.9269583954166246e-05
CG iteration 42, residual = 2.285137661493088e-05
CG iteration 43, residual = 1.7602218994251856e-05
CG iteration 44, residual = 1.3885671650305907e-05
CG iteration 45, residual = 1.0813388493948484e-05
CG iteration 46, residual = 8.097364754102716e-06
CG iteration 47, residual = 6.141265221626117e-06
CG iteration 48, residual = 4.793277042950815e-06
CG iteration 49, residual = 3.6467780868827916e-06
CG iteration 50, residual = 2.830568735279166e-06
CG iteration 51, residual = 2.141470433460834e-06
CG iteration 52, residual = 1.6371318610893338e-06
CG iteration 53, residual = 1.2685213384086537e-06
CG iteration 54, residual = 9.687252649034577e-07
CG iteration 55, residual = 7.429325199359565e-07
CG iteration 56, residual = 5.718110256346721e-07
CG iteration 57, residual = 4.438687295366517e-07
CG iteration 58, residual = 3.4280742534723985e-07
CG iteration 59, residual = 2.691845408060228e-07
CG iteration 60, residual = 2.0761481848411183e-07
CG iteration 61, residual = 1.6004452368415978e-07
CG iteration 62, residual = 1.2386160608117877e-07
CG iteration 63, residual = 9.358105223660559e-08
CG iteration 64, residual = 7.218204286049857e-08
CG iteration 65, residual = 5.591433041798756e-08
CG iteration 66, residual = 4.243144589652259e-08
CG iteration 67, residual = 3.29804930095328e-08
CG iteration 68, residual = 2.7207879895608794e-08
CG iteration 69, residual = 2.546655816213624e-08
CG iteration 70, residual = 1.8332485280425107e-08
CG iteration 71, residual = 1.3331580858519546e-08
CG iteration 72, residual = 1.0162601199143588e-08
CG iteration 73, residual = 7.80326329596324e-09
CG iteration 74, residual = 6.026596381663069e-09
CG iteration 75, residual = 4.617040515499358e-09
CG iteration 76, residual = 3.5260073743319817e-09
CG iteration 77, residual = 2.705147668335879e-09
CG iteration 78, residual = 2.0910644448792158e-09
CG iteration 79, residual = 1.6587782476564174e-09
CG iteration 80, residual = 1.339031353371184e-09
CG iteration 81, residual = 1.073756443152403e-09
CG iteration 82, residual = 8.136011516325387e-10
CG iteration 83, residual = 6.186864023056255e-10
CG iteration 84, residual = 4.770264754423459e-10
CG iteration 85, residual = 3.690925048835036e-10
CG iteration 86, residual = 2.7715724122173027e-10
CG iteration 87, residual = 2.143529056200954e-10
CG iteration 88, residual = 1.684371707889656e-10
CG iteration 89, residual = 1.2970058530521201e-10
CG iteration 90, residual = 9.990284210179801e-11
CG iteration 91, residual = 7.72834462255674e-11
CG iteration 92, residual = 6.065228315588102e-11
CG iteration 93, residual = 4.965542424291387e-11
CG iteration 94, residual = 3.9095066734186765e-11
CG iteration 95, residual = 2.8790940511819408e-11
CG iteration 96, residual = 2.1915428443164043e-11
CG iteration 97, residual = 1.694521059391716e-11
CG iteration 98, residual = 1.33070313512233e-11
CG iteration 99, residual = 1.1325329132456237e-11
CG iteration 100, residual = 1.0242066152615846e-11
CG iteration 101, residual = 7.1371670769487215e-12
CG iteration 102, residual = 5.230924091303894e-12
CG iteration 103, residual = 3.946742077029183e-12
CG iteration 104, residual = 3.058083371122738e-12
CG iteration 105, residual = 2.3735675144565254e-12
CG iteration 106, residual = 1.828303943821735e-12
CG iteration 107, residual = 1.4502644382722582e-12
CG iteration 108, residual = 1.1301024967421625e-12
CG iteration 109, residual = 9.000558472790931e-13
CG iteration 110, residual = 6.798815243542712e-13
[9]:
BaseWebGuiScene
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 = 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],deformation=True, settings={"camera": {"transformations": [{"type": "rotateX", "angle": -45}]}})
Newton iteration 0
err = 0.3934614166999501
Newton iteration 1
err = 5.90497752440774e-15
[10]:
BaseWebGuiScene