Click here to read what the information on this page means.
Degrees | \(1\leqslant k\) where \(k\) is the Polynomial superdegree |
Polynomial subdegree | \(k\) |
Polynomial superdegree | \(k\) |
Reference elements | triangle, tetrahedron |
Polynomial set | \(\mathcal{Z}^{(2)}_{k}\) ↓ Show polynomial set definitions ↓↑ Hide polynomial set definitions ↑\(\mathcal{Z}^{(2)}_k=\left\{\mathbf{M}\in\mathcal{P}_{k}^{d\times d}\middle|\mathbf{M}^t=\mathbf{M}\right\}\)
\(\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 edge: point evaluations of inner products with direction of edge
On each face: point evaluations of inner products with direction of edges
On each volume: point evaluations of inner products with direction of edges |
Number of DOFs | triangle: \(3(k+1)(k+2)/2\) (A045943) tetrahedron: \((k+1)(k+2)(k+3)\) (A007531) |
Mapping | double covariant Piola |
continuity | Inner products with tangents to facets are continuous |
Categories | Matrix-valued elements |
This element is implemented in
Basix ,
Basix.UFL ,
FIAT ,
Symfem , and
(legacy) UFL.
↓ Show implementation detail ↓↑ Hide implementation detail ↑Basix | basix.ElementFamily.Regge ↓ 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 Regge degree 1 on a triangle element = basix.create_element(basix.ElementFamily.Regge, basix.CellType.triangle, 1)
# Create Regge degree 2 on a triangle element = basix.create_element(basix.ElementFamily.Regge, basix.CellType.triangle, 2) This implementation is correct for all the examples below. |
Basix.UFL | basix.ElementFamily.Regge ↓ Show Basix.UFL examples ↓↑ Hide Basix.UFL examples ↑Before running this example, you must install Basix.UFL: pip3 install fenics-basix fenics-ufl This element can then be created with the following lines of Python: import basix import basix.ufl
# Create Regge degree 1 on a triangle element = basix.ufl.element(basix.ElementFamily.Regge, basix.CellType.triangle, 1)
# Create Regge degree 2 on a triangle element = basix.ufl.element(basix.ElementFamily.Regge, basix.CellType.triangle, 2) This implementation is correct for all the examples below. |
FIAT | FIAT.Regge ↓ 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 Regge degree 1 element = FIAT.Regge(FIAT.ufc_cell("triangle"), 1)
# Create Regge degree 2 element = FIAT.Regge(FIAT.ufc_cell("triangle"), 2) This implementation is correct for all the examples below. |
Symfem | "Regge" ↓ 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 Regge degree 1 on a triangle element = symfem.create_element("triangle", "Regge", 1)
# Create Regge degree 2 on a triangle element = symfem.create_element("triangle", "Regge", 2) This implementation is used to compute the examples below and verify other implementations. |
(legacy) UFL | "Regge" ↓ 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 Regge degree 1 on a triangle element = ufl_legacy.FiniteElement("Regge", "triangle", 1)
# Create Regge degree 2 on a triangle element = ufl_legacy.FiniteElement("Regge", "triangle", 2) |
- Regge, Tullio. General relativity without coordinates, Il Nuovo Cimento 19(3), 558–571, 1961. [DOI: 10.1007/BF02733251] [BibTeX]
- Christiansen, Snorre H. On the linearization of Regge calculus, Numerische Mathematik 119(4), 613–640, 2011. [DOI: 10.1007/s00211-011-0394-z] [BibTeX]
Element added | 01 January 2021 |
Element last updated | 27 September 2024 |