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:
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_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
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.00229
1.04141
1.13674
...
7.35531
7.49161
7.66495
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.00137
1.03997
1.11641
...
6.92606
7.01828
7.29664
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.00412
1.09465
1.334
...
31.0839
33.0367
33.6481
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.17035809144672
This is a slight improvement from the default.
[8]:
lams = TestPreconditioner (10, condense=True)
max(lams)/min(lams)
[8]:
27.352602492150286
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 = 1175801
WARNING: kwarg 'coarsetype' is an undocumented flags option for class <class 'ngsolve.comp.Preconditioner'>, maybe there is a typo?
H1AMG: level = 0, num_edges = 58879, nv = 1175801
H1AMG: level = 1, num_edges = 31014, nv = 4852
H1AMG: level = 2, num_edges = 16548, nv = 2536
H1AMG: level = 3, num_edges = 8627, nv = 1328
H1AMG: level = 4, num_edges = 4447, nv = 694
H1AMG: level = 5, num_edges = 2260, nv = 364
H1AMG: level = 6, num_edges = 1137, nv = 191
H1AMG: level = 7, num_edges = 547, nv = 101
H1AMG: level = 8, num_edges = 254, nv = 53
H1AMG: level = 9, num_edges = 118, nv = 27
H1AMG: level = 10, num_edges = 50, nv = 14
WARNING: maxsteps is deprecated, use maxiter instead!
CG iteration 1, residual = 0.6993046332934342
CG iteration 2, residual = 0.2816405683923466
CG iteration 3, residual = 0.2534270702488981
CG iteration 4, residual = 0.29374588782750743
CG iteration 5, residual = 0.24979838131887253
CG iteration 6, residual = 0.19307143208361677
CG iteration 7, residual = 0.14507027415991997
CG iteration 8, residual = 0.11425663721430869
CG iteration 9, residual = 0.08437842455990938
CG iteration 10, residual = 0.06181892970837919
CG iteration 11, residual = 0.04329458652656803
CG iteration 12, residual = 0.035884357865296646
CG iteration 13, residual = 0.030052509278868603
CG iteration 14, residual = 0.021520138691667457
CG iteration 15, residual = 0.01638141038753024
CG iteration 16, residual = 0.01284355058668398
CG iteration 17, residual = 0.009972266196637125
CG iteration 18, residual = 0.007817539069370404
CG iteration 19, residual = 0.005586155133694766
CG iteration 20, residual = 0.004194069760791714
CG iteration 21, residual = 0.003254149961144854
CG iteration 22, residual = 0.0025160733693793563
CG iteration 23, residual = 0.001938975232682286
CG iteration 24, residual = 0.0015342219235301557
CG iteration 25, residual = 0.0011592922349818555
CG iteration 26, residual = 0.0008947546311729194
CG iteration 27, residual = 0.000694818282673506
CG iteration 28, residual = 0.0005250840092673742
CG iteration 29, residual = 0.0004086187156543107
CG iteration 30, residual = 0.00031345299930079167
CG iteration 31, residual = 0.0002388111691206964
CG iteration 32, residual = 0.00018847457263429927
CG iteration 33, residual = 0.00014361299332016895
CG iteration 34, residual = 0.00011134704615984277
CG iteration 35, residual = 8.751774966168572e-05
CG iteration 36, residual = 6.515155202704007e-05
CG iteration 37, residual = 4.938447762218151e-05
CG iteration 38, residual = 3.874726887941691e-05
CG iteration 39, residual = 2.9928231414717188e-05
CG iteration 40, residual = 2.2859641924446893e-05
CG iteration 41, residual = 1.7333872753610214e-05
CG iteration 42, residual = 1.3547429120723057e-05
CG iteration 43, residual = 1.0459904675199333e-05
CG iteration 44, residual = 7.988330244467361e-06
CG iteration 45, residual = 6.090382630526065e-06
CG iteration 46, residual = 4.718243021902477e-06
CG iteration 47, residual = 3.620566517598434e-06
CG iteration 48, residual = 2.795568140007247e-06
CG iteration 49, residual = 2.1470877720153667e-06
CG iteration 50, residual = 1.6186115780272492e-06
CG iteration 51, residual = 1.2272942165706068e-06
CG iteration 52, residual = 9.308515571295224e-07
CG iteration 53, residual = 7.125134709266638e-07
CG iteration 54, residual = 5.353433795307413e-07
CG iteration 55, residual = 4.101157329993948e-07
CG iteration 56, residual = 3.149076991308879e-07
CG iteration 57, residual = 2.379943913561529e-07
CG iteration 58, residual = 1.8386367447047412e-07
CG iteration 59, residual = 1.4364864692883897e-07
CG iteration 60, residual = 1.0886099401070698e-07
CG iteration 61, residual = 8.334836077168122e-08
CG iteration 62, residual = 6.450942015892669e-08
CG iteration 63, residual = 4.90650252979628e-08
CG iteration 64, residual = 3.718256250710131e-08
CG iteration 65, residual = 2.852501802801346e-08
CG iteration 66, residual = 2.1912103226930428e-08
CG iteration 67, residual = 1.6567078834605155e-08
CG iteration 68, residual = 1.2671328071097729e-08
CG iteration 69, residual = 9.691497086979246e-09
CG iteration 70, residual = 7.309878775823902e-09
CG iteration 71, residual = 5.735969674274909e-09
CG iteration 72, residual = 4.350285885579361e-09
CG iteration 73, residual = 3.292240769686085e-09
CG iteration 74, residual = 2.5414378881654083e-09
CG iteration 75, residual = 1.9644083701698676e-09
CG iteration 76, residual = 1.5268005760274554e-09
CG iteration 77, residual = 1.229642820947494e-09
CG iteration 78, residual = 9.747784638925182e-10
CG iteration 79, residual = 7.391644957793579e-10
CG iteration 80, residual = 5.675036029900909e-10
CG iteration 81, residual = 4.25442788954528e-10
CG iteration 82, residual = 3.2470566810389147e-10
CG iteration 83, residual = 2.476138899317214e-10
CG iteration 84, residual = 1.9010848406383188e-10
CG iteration 85, residual = 1.5005465239941174e-10
CG iteration 86, residual = 1.2084418621256235e-10
CG iteration 87, residual = 9.489772793025294e-11
CG iteration 88, residual = 7.030542848249846e-11
CG iteration 89, residual = 5.381494944947093e-11
CG iteration 90, residual = 4.0508675930154565e-11
CG iteration 91, residual = 3.033960261991617e-11
CG iteration 92, residual = 2.3311662414042465e-11
CG iteration 93, residual = 1.8177836733798322e-11
CG iteration 94, residual = 1.4417457995704023e-11
CG iteration 95, residual = 1.1470230125188796e-11
CG iteration 96, residual = 8.76851219779089e-12
CG iteration 97, residual = 6.538714241188641e-12
CG iteration 98, residual = 4.9327891110525316e-12
CG iteration 99, residual = 3.7355349937260855e-12
CG iteration 100, residual = 2.8416145418870637e-12
CG iteration 101, residual = 2.157092846078873e-12
CG iteration 102, residual = 1.6311977816195852e-12
CG iteration 103, residual = 1.2530499385800884e-12
CG iteration 104, residual = 9.613502496971098e-13
CG iteration 105, residual = 7.337710180528697e-13
CG iteration 106, residual = 5.746015905728054e-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)
Newton iteration 0
err = 0.39346141669994983
Newton iteration 1
err = 1.8488788182878563e-15
[10]:
BaseWebGuiScene
[ ]: