This page was generated from unit-2.7-hybrid/hybrid.ipynb.
2.7 Facet spaces and hybrid methods¶
Mixed methods for second order problems lead to saddle point problems, and indefinite matrices. By hybridization one obtains a positive definite system again. It’s structure is similar to the non-conforming \(P^1\) method, but hybridization works for any order. See text-book by Brezzi and Fortin.
One skips the normal-continuity of the \(H(div)\) variable, and reinforces it by a Lagrange parameter. This leads to the following discrete system:
Find \(\sigma, u, \widehat u \in \Sigma_h \times V_h \times F_h\):
where \(\Sigma_h\) is an discontinuous \(H(div)\) finite element space, \(V_h\) a sub-space of \(L_2\), and \(F_h\) consists of polynomials on every edge.
[1]:
from ngsolve import *
from ngsolve.webgui import Draw
mesh = Mesh(unit_square.GenerateMesh(maxh=0.2))
same example as in 'mixed':
[2]:
source = sin(pi*x)
ud = CF(5)
g = mesh.BoundaryCF( {"left" : y*(1-y)}, default=0)
lam = 10
define spaces:
The Discontinuous FESpace-wrapper generates an element-wise discontinuous space
FacetFESpace lives only on facets (i.e. faces in 3D, edges in 2D, points in 1D)
Boundary conditions are now posed for the facet-space
[3]:
order = 3
V = Discontinuous (HDiv(mesh, order=order))
Q = L2(mesh, order=order-1)
F = FacetFESpace(mesh, order=order, dirichlet="bottom")
X = V*Q*F
print ("sigmadofs:", X.Range(0))
print ("udofs: ", X.Range(1))
print ("uhatdofs: ", X.Range(2))
sigmadofs: [0,1120)
udofs: [1120,1456)
uhatdofs: [1456,1832)
Assemble forms. The jump-term is rewritten as
[4]:
sigma,u,uhat = X.TrialFunction()
tau,v,vhat = X.TestFunction()
a = BilinearForm(X, condense=True)
a += (1/lam * sigma*tau + div(sigma)*v + div(tau)*u) * dx
n = specialcf.normal(mesh.dim)
a += (-sigma*n*vhat-tau*n*uhat) * dx(element_boundary=True)
c = Preconditioner(a, "bddc")
f = LinearForm(X)
f += -source*v * dx - g*vhat.Trace() * ds
a.Assemble()
print ("A non-zero", a.mat.nze)
gfu = GridFunction(X)
A non-zero 6880
Solve system. Either we leave everything to the sparse direct solver, or use CG
[5]:
f.Assemble()
gfu.components[2].Set(ud, BND)
if a.condense:
fmod = (f.vec + a.harmonic_extension_trans * f.vec).Evaluate()
solvers.CG(mat=a.mat, pre=c.mat, rhs=fmod, sol=gfu.vec, initialize=False)
gfu.vec.data += a.harmonic_extension * gfu.vec
gfu.vec.data += a.inner_solve * f.vec
else:
r = f.vec - a.mat * gfu.vec
inv = a.mat.Inverse(freedofs=X.FreeDofs())
gfu.vec.data += inv * r
CG iteration 1, residual = 136.14101710171593
CG iteration 2, residual = 59.971998967057004
CG iteration 3, residual = 24.09708229776067
CG iteration 4, residual = 12.365259367332888
CG iteration 5, residual = 4.500227573740079
CG iteration 6, residual = 1.5436222690215968
CG iteration 7, residual = 0.6253606118634968
CG iteration 8, residual = 0.23890072286213285
CG iteration 9, residual = 0.09410401770377796
CG iteration 10, residual = 0.03045052624547956
CG iteration 11, residual = 0.014237311920756503
CG iteration 12, residual = 0.004932886745754303
CG iteration 13, residual = 0.0021466580665278276
CG iteration 14, residual = 0.0008801579888635684
CG iteration 15, residual = 0.00040538011416669064
CG iteration 16, residual = 0.0001073184650755901
CG iteration 17, residual = 4.241130888841876e-05
CG iteration 18, residual = 1.6149940767707053e-05
CG iteration 19, residual = 7.046696006269176e-06
CG iteration 20, residual = 2.5136423160207232e-06
CG iteration 21, residual = 6.758564102270889e-07
CG iteration 22, residual = 3.182144556291163e-07
CG iteration 23, residual = 1.105120930565134e-07
CG iteration 24, residual = 5.4002776862598917e-08
CG iteration 25, residual = 1.890113451436262e-08
CG iteration 26, residual = 7.3118615775742255e-09
CG iteration 27, residual = 2.69270615226059e-09
CG iteration 28, residual = 1.028485720121852e-09
CG iteration 29, residual = 4.148508907053663e-10
CG iteration 30, residual = 1.420499282388908e-10
CG iteration 31, residual = 5.2792155727692876e-11
[6]:
Draw (gfu.components[0], mesh, "sigma")
Draw (gfu.components[1], mesh, "u");
[ ]:
[ ]: