Scalar and vectorial elements in NGSolve:
Standard continuous H1 elements:
Nedelec's tangentially-continuous H(curl)-conforming edge elements:
Raviart-Thomas normally-continuous H(div)-conforming face elements:
Discontinuous L2 elements:
These vector-valued spaces allow to represent physical quantities which are either normally or tangentially continuous.
The finite element spaces are related by the de Rham complex:
H1grad⟶H(curl)curl⟶H(div)div⟶L2⋃⋃⋃⋃Whgrad⟶Vhcurl⟶Qhdiv⟶ShNGSolve supports these elements of arbitrary order, on all common element shapes (trigs, quads, tets, prisms, pyramids, hexes). Elements may be curved.
import netgen.gui
%gui tk
from ngsolve import *
from netgen.geom2d import unit_square
from netgen.csg import unit_cube
mesh = Mesh(unit_square.GenerateMesh(maxh=0.3))
Generate a higher order H1-space. We first explore its different types of basis functions.
order=3
fes = H1(mesh, order=order)
gfu = GridFunction(fes)
Draw(gfu)
The first #vertices basis functions are hat-functions. By setting the solution vector to a unit-vector, we may look at the individual basis functions:
SetVisualization(min=0, max=1)
gfu.vec[:] = 0
# vertex nr:
gfu.vec[17] = 1
Redraw()
The next are edge-bubbles, where we have (order−1) basis functions per edge. A NodeId
object refers to a particular vertex, edge, face or cell node in the mesh. We can ask for the degrees of freedom on a node:
# basis functions on edge nr:
edge_dofs = fes.GetDofNrs(NodeId(EDGE,10))
print("edge_dofs =", edge_dofs)
SetVisualization(min=-0.05, max=0.05)
gfu.vec[:] = 0
gfu.vec[edge_dofs[1]] = 1
Redraw()
Finally, we have (p−1)(p−2)/2 inner basis functions on every triangle:
trig_dofs = fes.GetDofNrs(NodeId(FACE,0))
print("trig_dofs = ", trig_dofs)
SetVisualization(min=0, max=0.03)
gfu.vec[:] = 0
gfu.vec[trig_dofs[0]] = 1
Redraw()
The FESpace
also maintains information about local dofs, interface dofs and wire-basket dofs for the BDDC preconditioner:
for i in range(fes.ndof):
print (i,":", fes.CouplingType(i))
In NGSolve we use hierarchical high order finite element basis functions with node-wise exact sequences. The lowest order space Wl.o is the edge-element space:
Whp=Wp=1+∑EWE+∑FWF+∑CWC⊂H1Vhp=Wl.o+∑EVE+∑FVF+∑CVC⊂H(curl)where the edge, face and cell blocks are compatible in the sense that
∇WE=VE,∇WF⊂VF,∇WC⊂VCWe obtain this by using gradients of H1 basis functions as H(curl) basis functions, and some more (see thesis Sabine Zaglmayr):
VE=span{∇φH1E,i}VF=span{∇φH1F,i∪˜φH(curl)F,i}VC=span{∇φH1C,i∪˜φH(curl)C,i}fes = HCurl(mesh, order=2)
uc = GridFunction(fes, name="uc")
Draw (uc)
Draw (curl(uc), mesh, "curl")
edge_dofs = fes.GetDofNrs(NodeId(EDGE,10))
print ("edgedofs: ", edge_dofs)
uc.vec[:] = 0
uc.vec[edge_dofs[1]] = 1
Redraw()
look at them by activating Draw Surface Vectors.
face_dofs = fes.GetDofNrs(NodeId(FACE,10))
print ("facedofs: ", face_dofs)
uc.vec[:] = 0
uc.vec[face_dofs[1]] = 1
Redraw()