Click here to read what the information on this page means.
Degrees | triangle: \(k=1\)
tetrahedron: \(1\leqslant k\leqslant 2\) |
Reference elements | triangle, tetrahedron |
DOFs | On each vertex: point evaluations in coordinate directions
On each edge: (if \(k>1\)) point evaluations in coordinate directions at midpoints
On each facet: normal integral moments with an degree \(0\) Lagrange space
On the interior of the reference element: point evaluations in coordinate directions at the midpoint, and (if \(k>1\)) point evaluations in coordinate directions at midpoints of internal edges |
Number of DOFs | triangle: \(11\) tetrahedron: \(\begin{cases}19&k=1\\49&k=2\end{cases}\) |
Mapping | contravariant Piola |
continuity | Function values are continuous. |
Categories | Vector-valued elements, Macro elements |
This element is implemented in
FIAT and
Symfem .
↓ Show implementation detail ↓↑ Hide implementation detail ↑FIAT | FIAT.GuzmanNeilanSecondKindH1 ↓ Show FIAT examples ↓↑ Hide FIAT examples ↑Before running this example, you must install FIAT: pip3 install git+https://github.com/firedrakeproject/fiat.git This element can then be created with the following lines of Python: import FIAT
# Create Guzman-Neilan (second kind) degree 1 element = FIAT.GuzmanNeilanSecondKindH1(FIAT.ufc_cell("triangle"), 1)
# Create Guzman-Neilan (second kind) degree 1 element = FIAT.GuzmanNeilanSecondKindH1(FIAT.ufc_cell("tetrahedron"), 1)
# Create Guzman-Neilan (second kind) degree 2 element = FIAT.GuzmanNeilanSecondKindH1(FIAT.ufc_cell("tetrahedron"), 2) This implementation is correct for some of the examples below. Correct: triangle,1 Incorrect: tetrahedron,1; tetrahedron,2 |
Symfem | "Guzman-Neilan second kind" ↓ Show Symfem examples ↓↑ Hide Symfem examples ↑Before running this example, you must install Symfem: pip3 install symfem This element can then be created with the following lines of Python: import symfem
# Create Guzman-Neilan (second kind) degree 1 on a triangle element = symfem.create_element("triangle", "Guzman-Neilan second kind", 1)
# Create Guzman-Neilan (second kind) degree 1 on a tetrahedron element = symfem.create_element("tetrahedron", "Guzman-Neilan second kind", 1)
# Create Guzman-Neilan (second kind) degree 2 on a tetrahedron element = symfem.create_element("tetrahedron", "Guzman-Neilan second kind", 2) This implementation is used to compute the examples below and verify other implementations. |
- Guzmán, Johnny and Neilan, Michael. Inf-sup stable finite elements on barycentric refinements producing divergence-free approximations in arbitrary dimensions, SIAM Journal on Numerical Analysis 56, 2826–2844, 2018. [DOI: 10.1137/17M1153467] [BibTeX]
Element added | 16 October 2024 |
Element last updated | 16 October 2024 |