Click here to read what the information on this page means.
| Alternative names | Brezzi–Douglas–Marini cubical H(curl) (quadrilateral), Arnold–Awanou H(curl) (hexahedron) | 
| De Rham complex families | \(\left[S_{1,k}^\square\right]_{1}\) or \(\mathcal{S}_{k}\Lambda^{1}(\square_d)\) | 
| Abbreviated names | BDMce (quadrilateral), AAe (hexahedron) | 
| Degrees | \(1\leqslant k\) where \(k\) is the polynomial subdegree
 | 
| Polynomial subdegree | \(k\) | 
| Polynomial superdegree | \(k + d - 1\) | 
| Lagrange subdegree | quadrilateral: \(\operatorname{floor}(k/d)\) hexahedron: \(1\) (degree=2), \(\operatorname{floor}(k/d)\) (otherwise)
 | 
| Lagrange superdegree | \(k+1\) | 
| Reference cells | quadrilateral, hexahedron | 
| Finite dimensional space | \(\mathcal{P}_{k}^d \oplus \mathcal{Z}^{(26)}_{k}\) (quadrilateral) \(\mathcal{P}_{k}^d \oplus \mathcal{A}_{k-1} \oplus \mathcal{Z}^{(27)}_{k}\) (hexahedron)
 ↓ Show set definitions ↓↑ Hide set definitions ↑
 \(\mathcal{P}_k=\operatorname{span}\left\{\prod_{i=1}^dx_i^{p_i}\middle|\sum_{i=1}^dp_i\leqslant k\right\}\)
 \(\mathcal{Z}^{(26)}_k=\operatorname{span}\left\{\left(\begin{array}{c}(k+1)x^ky\\-x^{k+1}\end{array}\right),\left(\begin{array}{c}y^{k+1}\\-(k+1)xy^k\end{array}\right)\right\}\)
 
 \(\mathcal{A}_k=\operatorname{span}\bigcup_{p\in\mathcal{P}_k(\mathbb{R}^{d-1})}\left\{\left(\begin{array}{c}0\\xzp(y,z)\\-xyp(y,z)\end{array}\right),\left(\begin{array}{c}yzp(x,z)\\0\\-xyp(x,z)\end{array}\right),\left(\begin{array}{c}yzp(x,y)\\-xzp(x,y)\\0\end{array}\right)\right\}\)
 
 \(\mathcal{Z}^{(27)}_k=\left\{\nabla p\middle|p\in\mathcal{X}_{k+1}\right\}\)
 
 \(\mathcal{X}_k=\operatorname{span}\left\{\prod_{i=1}^dx_i^{p_i}\middle|k<\sum_{i=1}^da_i\leqslant k+\#\left\{i\in\{1,\dots,d\}\middle| a_i=1\right\}\right\}\)
 | 
| DOFs | On each edge: tangent integral moments with a degree \(k\) dPc space On each face: integral moments with a degree \(k-2\) vector dPc space
 On each volume: integral moments with a degree \(k-4\) vector dPc space
 | 
| Number of DOFs | quadrilateral: \(k^2+3k+4\) (A014206) hexahedron: \(\begin{cases}6(k^2+k+2)&k=1,2,3\\k(k+1)(k-1)/2 + 3k^2 + 12k + 9&k > 3\end{cases}\)
 | 
| Mapping | covariant Piola | 
| continuity | Components tangential to facets are continuous | 
| Categories | Vector-valued elements, H(curl) conforming elements | 
This element is implemented in 
Symfem  and 
(legacy) UFL.
↓ Show implementation detail ↓↑ Hide implementation detail ↑| Symfem | "Scurl"↓ Show Symfem examples ↓↑ Hide Symfem examples ↑
 This implementation is used to compute the examples below and verify other implementations.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 serendipity H(curl) degree 1 on a quadrilateral
 element = symfem.create_element("quadrilateral", "Scurl", 1)
 
 # Create serendipity H(curl) degree 2 on a quadrilateral
 element = symfem.create_element("quadrilateral", "Scurl", 2)
 
 # Create serendipity H(curl) degree 1 on a hexahedron
 element = symfem.create_element("hexahedron", "Scurl", 1)
 
 # Create serendipity H(curl) degree 2 on a hexahedron
 element = symfem.create_element("hexahedron", "Scurl", 2)
 | 
| (legacy) UFL | "BDMCE"(quadrilateral)
 "AAE"(hexahedron)↓ Show (legacy) UFL examples ↓↑ Hide (legacy) UFL examples ↑
 Before running this example, you must install (legacy) UFL : pip3 install setuptoolspip3 install fenics-ufl-legacy
 This element can then be created with the following lines of Python: import ufl_legacy
 # Create serendipity H(curl) degree 1 on a quadrilateral
 element = ufl_legacy.FiniteElement("BDMCE", "quadrilateral", 1)
 
 # Create serendipity H(curl) degree 2 on a quadrilateral
 element = ufl_legacy.FiniteElement("BDMCE", "quadrilateral", 2)
 
 # Create serendipity H(curl) degree 1 on a hexahedron
 element = ufl_legacy.FiniteElement("AAE", "hexahedron", 1)
 
 # Create serendipity H(curl) degree 2 on a hexahedron
 element = ufl_legacy.FiniteElement("AAE", "hexahedron", 2)
 | 
 
- Arnold, Douglas N. and Awanou, Gerard. Finite element differential forms on cubical meshes, Mathematics of computation 83, 1551–5170, 2014. [DOI: 10.1090/S0025-5718-2013-02783-4] [BibTeX]
- Arnold, Douglas N. and Logg, Anders. Periodic table of the finite elements, SIAM News 47, 2014. [www.siam.org/publications/siam-news/issues/volume-47-number-09-november-2014/] [BibTeX]
- Cockburn, Bernardo and Fu, Guosheng. A systematic construction of finite element commuting exact sequences, SIAM journal on numerical analysis 55, 1650–1688, 2017. [DOI: 10.1137/16M1073352] [BibTeX]
| Element added | 31 December 2020 | 
| Element last updated | 04 June 2025 |