Click here to read what the information on this page means.
Degrees | \(k=3\) where \(k\) is the Polynomial superdegree |
Polynomial subdegree | \(k\) |
Polynomial superdegree | \(k\) |
Reference elements | interval, triangle, tetrahedron |
Polynomial set | \(\mathcal{P}_{k}\) ↓ Show polynomial set definitions ↓↑ Hide polynomial set definitions ↑\(\mathcal{P}_k=\operatorname{span}\left\{\prod_{i=1}^dx_i^{p_i}\middle|\sum_{i=1}^dp_i\leqslant k\right\}\) |
DOFs | On each vertex: point evaluations, and point evaluations of derivatives in coordinate directions
On each face: point evaluations at midpoints |
Number of DOFs | interval: \(4\) triangle: \(10\) tetrahedron: \(20\) |
Mapping | identity |
continuity | Function values are continuous. |
Notes | The derivatives of the basis functions are continuous between cells at the vertices of the element |
Categories | Scalar-valued elements |
This element is implemented in
Basix ,
FIAT ,
Symfem , and
(legacy) UFL.
↓ Show implementation detail ↓↑ Hide implementation detail ↑Basix | basix.ElementFamily.Hermite ↓ Show Basix examples ↓↑ Hide Basix examples ↑Before running this example, you must install Basix: pip3 install fenics-basix This element can then be created with the following lines of Python: import basix
# Create Hermite degree 3 on a interval element = basix.create_element(basix.ElementFamily.Hermite, basix.CellType.interval, 3)
# Create Hermite degree 3 on a triangle element = basix.create_element(basix.ElementFamily.Hermite, basix.CellType.triangle, 3)
# Create Hermite degree 3 on a tetrahedron element = basix.create_element(basix.ElementFamily.Hermite, basix.CellType.tetrahedron, 3) This implementation is correct for all the examples below. |
FIAT | FIAT.CubicHermite ↓ 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 Hermite degree 3 element = FIAT.CubicHermite(FIAT.ufc_cell("interval"), 3)
# Create Hermite degree 3 element = FIAT.CubicHermite(FIAT.ufc_cell("triangle"), 3)
# Create Hermite degree 3 element = FIAT.CubicHermite(FIAT.ufc_cell("tetrahedron"), 3) This implementation is correct for all the examples below. |
Symfem | "Hermite" ↓ 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 Hermite degree 3 on a interval element = symfem.create_element("interval", "Hermite", 3)
# Create Hermite degree 3 on a triangle element = symfem.create_element("triangle", "Hermite", 3)
# Create Hermite degree 3 on a tetrahedron element = symfem.create_element("tetrahedron", "Hermite", 3) This implementation is used to compute the examples below and verify other implementations. |
(legacy) UFL | "Hermite" ↓ Show (legacy) UFL examples ↓↑ Hide (legacy) UFL examples ↑Before running this example, you must install (legacy) UFL: pip3 install setuptools pip3 install fenics-ufl-legacy This element can then be created with the following lines of Python: import ufl_legacy
# Create Hermite degree 3 on a interval element = ufl_legacy.FiniteElement("Hermite", "interval", 3)
# Create Hermite degree 3 on a triangle element = ufl_legacy.FiniteElement("Hermite", "triangle", 3)
# Create Hermite degree 3 on a tetrahedron element = ufl_legacy.FiniteElement("Hermite", "tetrahedron", 3) |
- Ciarlet, Philippe G. and Raviart, Pierre-Arnaud. Interpolation theory over curved elements, with applications to finite element methods, Computer Methods in Applied Mechanics and Engineering 1(2), 217–249, 1972. [DOI: 10.1016/0045-7825(72)90006-0] [BibTeX]
Element added | 09 January 2021 |
Element last updated | 27 September 2024 |