Click here to read what the information on this page means.
Alternative names | Brezzi–Douglas–Marini cubical H(div) (quadrilateral), Arnold–Awanou H(div) (hexahedron) |
De Rham complex families | \(\left[S_{1,k}^\square\right]_{d-1}\) or \(\mathcal{S}_{k}\Lambda^{d-1}(\square_d)\) |
Abbreviated names | BDMcf (quadrilateral), AAf (hexahedron) |
Degrees | \(1\leqslant k\) where \(k\) is the Polynomial subdegree |
Polynomial subdegree | \(k\) |
Polynomial superdegree | \(k+1\) |
Lagrange subdegree | \(\operatorname{floor}(k/d)\) |
Lagrange superdegree | \(k+1\) |
Reference elements | quadrilateral, hexahedron |
Polynomial set | \(\mathcal{P}_{k}^d \oplus \mathcal{Z}^{(28)}_{k}\) (quadrilateral)
\(\mathcal{P}_{k}^d \oplus \mathcal{Z}^{(29)}_{k}\) (hexahedron)
↓ 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\}\)
\(\mathcal{Z}^{(28)}_k=\operatorname{span}\left\{\left(\begin{array}{c}x^{k+1}\\(k+1)x^ky\end{array}\right),\left(\begin{array}{c}(k+1)xy^k\\y^{k+1}\end{array}\right)\right\}\)
\(\mathcal{Z}^{(29)}_k=\left\{\nabla\times\boldsymbol{p}\middle|\boldsymbol{p}\in\mathcal{A}_{k}\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\}\) |
DOFs | On each facet: normal integral moments with an degree \(k\) dPc space
On the interior of the reference element: integral moments with an degree \(k-2\) vector dPc space |
Number of DOFs | quadrilateral: \(k^2+3k+4\) (A014206) hexahedron: \((k+1)(k^2+5k+12)/2\) |
Mapping | contravariant Piola |
continuity | Components normal to facets are continuous |
Categories | Vector-valued elements, H(div) conforming elements |
This element is implemented in
Symfem and
(legacy) UFL.
↓ Show implementation detail ↓↑ Hide implementation detail ↑Symfem | "Sdiv" ↓ 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 serendipity H(div) degree 1 on a quadrilateral element = symfem.create_element("quadrilateral", "Sdiv", 1)
# Create serendipity H(div) degree 2 on a quadrilateral element = symfem.create_element("quadrilateral", "Sdiv", 2)
# Create serendipity H(div) degree 1 on a hexahedron element = symfem.create_element("hexahedron", "Sdiv", 1)
# Create serendipity H(div) degree 2 on a hexahedron element = symfem.create_element("hexahedron", "Sdiv", 2) This implementation is used to compute the examples below and verify other implementations. |
(legacy) UFL | "BDMCF" (quadrilateral)
"AAF" (hexahedron) ↓ 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 serendipity H(div) degree 1 on a quadrilateral element = ufl_legacy.FiniteElement("BDMCF", "quadrilateral", 1)
# Create serendipity H(div) degree 2 on a quadrilateral element = ufl_legacy.FiniteElement("BDMCF", "quadrilateral", 2)
# Create serendipity H(div) degree 1 on a hexahedron element = ufl_legacy.FiniteElement("AAF", "hexahedron", 1)
# Create serendipity H(div) degree 2 on a hexahedron element = ufl_legacy.FiniteElement("AAF", "hexahedron", 2) |
- Arnold, Douglas N. and Awanou, Gerard. Finite element differential forms on cubical meshes, Mathematics of computation 83, 1551–5170, 2014. [BibTeX]
- Brezzi, Franco, Douglas, Jim, and Marini, L. Donatella. Two families of mixed finite elements for second order elliptic problems, Numerische Mathematik 47, 217–235, 1985. [DOI: 10.1007/BF01389710] [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 | 30 December 2020 |
Element last updated | 27 September 2024 |