This page was generated from unit-2.1.6-highorder-mg/highorder-multigrid.ipynb.
2.1.6. Multigrid for high order finite element spaces¶
Officially released with NGSolve-2504, most high order spaces provide high order polynomial preserving prolongation operators. Current restrictions:
require simplicial meshes
refinement by bisection
require uniform polynomial order
[1]:
from ngsolve import *
from ngsolve.webgui import Draw
from ngsolve.la import EigenValues_Preconditioner
Many FESpaces provide now a high-order accurate prolongation operator. It has to be enabled by the flag hoprolongation=True, maybe becoming default in future.
[2]:
mesh = Mesh(unit_square.GenerateMesh(maxh=0.2))
gfu = GridFunction(H1(mesh, order=10, hoprolongation=True))
gfu.Set (sin(50*x*y))
Draw (gfu, order=10);
[3]:
for l in range(2):
mesh.Refine()
Draw (gfu, order=5);
Multigrid preconditioners for high order spaces¶
If the high order prolongation is enabled, the multigrid preconditioner uses the high order discretization on the mesh hierarchy. If not, the coarse grid spaces use the lowest order spaces in the background.
[4]:
mesh = Mesh(unit_square.GenerateMesh(maxh=0.2))
fes = H1(mesh, order=5, hoprolongation=True, dirichlet=".*")
u,v = fes.TnT()
a = BilinearForm(grad(u)*grad(v)*dx).Assemble()
pre = preconditioners.MultiGrid(a, blocktype=["vertexpatch"])
# pre = preconditioners.MultiGrid(a, smoother="block", blocktype=["vertexedge","face"])
# pre = preconditioners.MultiGrid(a, smoother="block", blocktype=["vertpatch"])
for l in range(4):
mesh.Refine()
a.Assemble()
lam = EigenValues_Preconditioner(a.mat, pre)
print (mesh.levels, fes.ndof, lam[0], lam[-1])
2 2801 0.9737767047975869 0.9999650900439694
3 11001 0.973189334176511 0.9999569893695426
4 43601 0.9738508901099723 0.9999593243775629
5 173601 0.9740047338416953 0.9999581308397345
[5]:
f = LinearForm(x*v*dx).Assemble()
gfu = GridFunction(fes)
Solve(a * gfu == f, dirichlet=x*(1-x), lin_solver=solvers.CGSolver, pre=pre, printrates=True)
CG iteration 1, residual = 2.5006659907705164
CG iteration 2, residual = 0.02438398776423324
CG iteration 3, residual = 0.00018799791819351172
CG iteration 4, residual = 1.1484794277917987e-06
CG iteration 5, residual = 9.252209912922827e-09
CG iteration 6, residual = 6.018344774975674e-11
CG iteration 7, residual = 4.580225539555759e-13
[6]:
ea = { "euler_angles" : (-70, 0,-55) }
Draw (gfu, deformation=True, **ea);
High order with static condensation¶
For high order methods, static condensation may save a lot of computation. However, canonical prolongation for the skeleton variables only does not preserve high order polynomials. Here the HarmonicProlongation comes into play: It prolongates functions on the boundaries of the coarse grid elements, and then solves local Dirichlet problems for the dofs on the fine-grid skeleton inside the coarse grid elements.
The Dirichlet problem is solved for the problem-specific bilinear form, which has to be provided when enabling the HarmonicProlongation:
[7]:
mesh = Mesh(unit_square.GenerateMesh(maxh=0.2))
fes = H1(mesh, order=5, hoprolongation=True, dirichlet=".*")
u,v = fes.TnT()
a = BilinearForm(grad(u)*grad(v)*dx, condense=True).Assemble()
fes.SetHarmonicProlongation(a)
pre = preconditioners.MultiGrid(a, blocktype="vertexpatch")
for l in range(5):
mesh.Refine()
a.Assemble()
lam = EigenValues_Preconditioner(a.mat, pre)
print (mesh.levels, fes.ndof, lam[0], lam[-1])
2 2801 0.8873455539760273 0.9998975910633989
3 11001 0.8957401274225063 0.9998894014559885
4 43601 0.8945786078982367 0.999886328500714
5 173601 0.8943861692924868 0.9998960435342539
6 692801 0.8953852278217913 0.9998858093863499
[8]:
f = LinearForm(x*v*dx).Assemble()
gfu = GridFunction(fes)
Solve(a * gfu == f, dirichlet=x*(1-x), lin_solver=solvers.CGSolver, pre=pre, printrates=True)
CG iteration 1, residual = 3.488608797786568
CG iteration 2, residual = 0.12345444219082123
CG iteration 3, residual = 0.003639088091883729
CG iteration 4, residual = 9.080367025015463e-05
CG iteration 5, residual = 2.783128212652323e-06
CG iteration 6, residual = 7.69630828255892e-08
CG iteration 7, residual = 2.1727680277905276e-09
CG iteration 8, residual = 6.173415107276046e-11
CG iteration 9, residual = 1.8887182888427874e-12
[9]:
Draw (gfu, deformation=True, **ea);
This example shows the result of a HarmonicProlongation:
[10]:
mesh = Mesh(unit_square.GenerateMesh(maxh=2))
fes = H1(mesh, order=8, hoprolongation=True)
u,v = fes.TnT()
a = BilinearForm(grad(u)*grad(v)*dx, condense=True).Assemble()
fes.SetHarmonicProlongation(a)
pre = preconditioners.MultiGrid(a, smoother="block", blocktype="vertexpatch")
[11]:
gfu = GridFunction(fes)
gfu.Set(sin(10*x))
Draw(gfu, order=5)
mesh.Refine()
Draw(gfu, order=5)
gfu.vec.data += a.harmonic_extension * gfu.vec
Draw (gfu, order=5);
Nearly incompressible elasticity and Stokes¶
The Scott-Vogelius element:
[12]:
mesh = Mesh(unit_square.GenerateMesh(maxh=0.2))
fes = VectorH1(mesh, order=4, hoprolongation=True, dirichlet=".*")
u,v = fes.TnT()
eq = InnerProduct(Grad(u),Grad(v)) * dx + 1e4 * div(u)*div(v)* dx
a = BilinearForm(eq).Assemble()
# a = BilinearForm( (Grad(u)|Grad(v)) * dx).Assemble()
pre = preconditioners.MultiGrid(a, blocktype=["vertexpatch"])
with TaskManager():
for l in range(5):
mesh.Refine()
a.Assemble()
lam = EigenValues_Preconditioner(a.mat, pre)
print (mesh.levels, fes.ndof, lam[0], lam[-1])
2 3618 0.3132568089232025 0.999900784186339
3 14146 0.19750181032346975 0.999910050962979
4 55938 0.1954916126779867 0.9999047978307694
5 222466 0.19549539298190566 0.9998999587479178
6 887298 0.19556232072451185 0.9998955800582734
A robust method with reduced integration¶
We define the bilinear form as
where \(P_{L_2}^Q\) is an element-wise \(L_2\)-projector into a lower order space.
[13]:
mesh = Mesh(unit_square.GenerateMesh(maxh=0.2))
fesp = L2(mesh, order=0)
fes = VectorH1(mesh, order=2, hoprolongation=True, dirichlet=".*")
u,v = fes.TnT()
eq = InnerProduct(Grad(u),Grad(v)) * dx + 1e4 * Interpolate(div(u),fesp)*div(v)* dx
a = BilinearForm(eq).Assemble()
# a = BilinearForm( (Grad(u)|Grad(v)) * dx).Assemble()
fes.SetHarmonicProlongation(a)
pre = preconditioners.MultiGrid(a, blocktype=["vertexpatch"])
for l in range(4):
mesh.Refine()
a.Assemble()
lam = EigenValues_Preconditioner(a.mat, pre)
print (mesh.levels, fes.ndof, lam[0], lam[-1])
2 946 0.5451163631462156 1.2567450360094283
3 3618 0.3588082330609003 1.4293139465312101
4 14146 0.3055401019237036 1.5320915176935577
5 55938 0.2798589825176979 1.6946834108518791
[14]:
f = LinearForm((x-0.5)*v[1]*dx).Assemble()
gfu = GridFunction(fes)
Solve (a*gfu==f, pre, solvers.CGSolver, printrates=True)
CG iteration 1, residual = 0.020553570901516955
CG iteration 2, residual = 0.0024505017134841683
CG iteration 3, residual = 0.0006507984597058582
CG iteration 4, residual = 0.00020038866210218308
CG iteration 5, residual = 7.105863986260101e-05
CG iteration 6, residual = 3.028758180103971e-05
CG iteration 7, residual = 1.3957421251326409e-05
CG iteration 8, residual = 6.026536683793543e-06
CG iteration 9, residual = 2.743921037373388e-06
CG iteration 10, residual = 1.218630532456111e-06
CG iteration 11, residual = 4.964351575174683e-07
CG iteration 12, residual = 1.9398388083208307e-07
CG iteration 13, residual = 7.686678702176261e-08
CG iteration 14, residual = 3.304683227185746e-08
CG iteration 15, residual = 1.4352549305264776e-08
CG iteration 16, residual = 6.021165087793069e-09
CG iteration 17, residual = 2.5562823884468053e-09
CG iteration 18, residual = 1.0915843646613853e-09
CG iteration 19, residual = 4.5195816561668957e-10
CG iteration 20, residual = 1.7588626620383023e-10
CG iteration 21, residual = 6.957906924142435e-11
CG iteration 22, residual = 2.8502107548652563e-11
CG iteration 23, residual = 1.1345979806385866e-11
CG iteration 24, residual = 4.16916211928785e-12
CG iteration 25, residual = 1.50906111268178e-12
CG iteration 26, residual = 5.394782139787027e-13
CG iteration 27, residual = 1.9871234166889255e-13
CG iteration 28, residual = 7.249987931780704e-14
CG iteration 29, residual = 2.7920178284556712e-14
CG iteration 30, residual = 1.107782882396837e-14
[15]:
Draw (gfu);
[16]:
from netgen.occ import *
shape = MoveTo(0,0).LineTo(1,0,"in").LineTo(1,1).LineTo(2,1).LineTo(3,0).LineTo(4,1).LineTo(5,1) \
.LineTo(5,2,"out").LineTo(4,2).LineTo(3,1).LineTo(2,2).LineTo(1,2).Rotate(180).Arc(1,90).Close().Face()
mesh = shape.GenerateMesh(dim=2, maxh=0.25).Curve(3)
Draw (mesh)
print (mesh.GetBoundaries())
fesp = L2(mesh, order=0)
fes = VectorH1(mesh, order=2, hoprolongation=True, dirichlet="in|default")
u,v = fes.TnT()
eq = InnerProduct(Grad(u),Grad(v)) * dx + 1e4 * Interpolate(div(u),fesp)*div(v)* dx
a = BilinearForm(eq).Assemble()
# a = BilinearForm( (Grad(u)|Grad(v)) * dx).Assemble()
fes.SetHarmonicProlongation(a)
pre = preconditioners.MultiGrid(a, blocktype=["vertexpatch"])
for l in range(4):
mesh.Refine()
a.Assemble()
print (mesh.levels, fes.ndof)
('in', 'default', 'default', 'default', 'default', 'default', 'out', 'default', 'default', 'default', 'default', 'default', 'default')
2 3482
3 13426
4 52706
5 208834
from Dissertation J. Schöberl, page 116.
[17]:
f = LinearForm(fes)
gfu = GridFunction(fes)
Solve (a*gfu==f, pre, solvers.CGSolver, dirichlet=CF((0,x*(1-x))) | mesh.Boundaries("in"), printrates=True)
CG iteration 1, residual = 156.02169771864763
CG iteration 2, residual = 4.69699351846272
CG iteration 3, residual = 0.44126979160195273
CG iteration 4, residual = 0.06908617512247596
CG iteration 5, residual = 0.02943344479984058
CG iteration 6, residual = 0.014480167237763962
CG iteration 7, residual = 0.008734149812355993
CG iteration 8, residual = 0.005471780853094519
CG iteration 9, residual = 0.0020406404664476004
CG iteration 10, residual = 0.000867522773548497
CG iteration 11, residual = 0.00030334263942497283
CG iteration 12, residual = 0.0001495584892510743
CG iteration 13, residual = 5.951615094267259e-05
CG iteration 14, residual = 2.7117964901232846e-05
CG iteration 15, residual = 1.0442357021417666e-05
CG iteration 16, residual = 5.095531021683239e-06
CG iteration 17, residual = 2.342587632655926e-06
CG iteration 18, residual = 9.935176562377645e-07
CG iteration 19, residual = 3.816539545954556e-07
CG iteration 20, residual = 1.5921102516389785e-07
CG iteration 21, residual = 7.002462032025954e-08
CG iteration 22, residual = 3.0241428333010755e-08
CG iteration 23, residual = 1.4336012131891562e-08
CG iteration 24, residual = 6.296140098357637e-09
CG iteration 25, residual = 2.5090395542746937e-09
CG iteration 26, residual = 1.0887999974377092e-09
CG iteration 27, residual = 4.76529711649176e-10
CG iteration 28, residual = 2.1085634875243078e-10
CG iteration 29, residual = 1.1199062925516361e-10
[18]:
Draw (gfu);
[ ]: