Click here to read what the information on this page means.
De Rham complex families | \(\left[S_{2,k}^\square\right]_{d-1}\) or \(\mathcal{S}^-_{k}\Lambda^{d-1}(\square_d)\) |
Degrees | \(1\leqslant k\) where \(k\) is the Lagrange superdegree |
Polynomial subdegree | \(k-1\) |
Polynomial superdegree | \(k\) |
Lagrange subdegree | \(\operatorname{floor}((k+1)/(d+1))\) |
Lagrange superdegree | \(k\) |
Reference elements | quadrilateral, hexahedron |
Polynomial set | \(\mathcal{P}_{k-1}^d \oplus \mathcal{Z}^{(47)}_{k} \oplus \mathcal{Z}^{(48)}_{k}\) (quadrilateral)
\(\mathcal{P}_{k-1}^d \oplus \mathcal{Z}^{(47)}_{k} \oplus \mathcal{Z}^{(49)}_{k} \oplus \mathcal{Z}^{(50)}_{k} \oplus \mathcal{Z}^{(51)}_{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}^{(47)}_k=\operatorname{span}\left\{\boldsymbol{x}p\middle|p\in\hat{\mathcal{P}}_{k-1}\right\}\)
\(\hat{\mathcal{P}}_k=\operatorname{span}\left\{\prod_{i=1}^dx_i^{p_i}\middle|\sum_{i=1}^dp_i=k\right\}=\mathcal{P}_k\setminus\mathcal{P}_{k-1}\)
\(\mathcal{Z}^{(48)}_k=\operatorname{span}\left\{\left(\begin{array}{c}kxy^{k-1}\\-y^k\end{array}\right),\left(\begin{array}{c}-x^k\\kx^{k-1}y\end{array}\right)\right\}\)
\(\mathcal{Z}^{(49)}_k=\operatorname{span}\left\{\operatorname{curl}\left(\left(\begin{array}{c}0\\-z\\y\end{array}\right)xy^az^{k-1-a}\right)\middle|a=0,1,...,a-1\right\}\)
\(\mathcal{Z}^{(50)}_k=\operatorname{span}\left\{\operatorname{curl}\left(\left(\begin{array}{c}-z\\0\\x\end{array}\right)yx^az^{k-1-a}\right)\middle|a=0,1,...,a-1\right\}\)
\(\mathcal{Z}^{(51)}_k=\operatorname{span}\left\{\operatorname{curl}\left(\left(\begin{array}{c}-y\\x\\0\end{array}\right)zx^ay^{k-1-a}\right)\middle|a=0,1,...,a-1\right\}\) |
DOFs | On each facet: normal integral moments with an degree \(k-1\) dPc space
On the interior of the reference element: integral moments with an degree \(k-3\) vector dPc space, and integral moments with \(\left\{\nabla(p)\middle|p\text{ is a degree \(k-1\) monomial}\right\}\) |
Mapping | contravariant Piola |
continuity | Components normal to facets are continuous |
Categories | Vector-valued elements, H(div) conforming elements |
This element is implemented in
FIAT and
Symfem .
↓ Show implementation detail ↓↑ Hide implementation detail ↑FIAT | FIAT.TrimmedSerendipityDiv ↓ 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 trimmed serendipity H(div) degree 1 element = FIAT.TrimmedSerendipityDiv(FIAT.reference_element.UFCQuadrilateral(), 1)
# Create trimmed serendipity H(div) degree 2 element = FIAT.TrimmedSerendipityDiv(FIAT.reference_element.UFCQuadrilateral(), 2)
# Create trimmed serendipity H(div) degree 3 element = FIAT.TrimmedSerendipityDiv(FIAT.reference_element.UFCQuadrilateral(), 3)
# Create trimmed serendipity H(div) degree 1 element = FIAT.TrimmedSerendipityDiv(FIAT.reference_element.UFCHexahedron(), 1)
# Create trimmed serendipity H(div) degree 2 element = FIAT.TrimmedSerendipityDiv(FIAT.reference_element.UFCHexahedron(), 2)
# Create trimmed serendipity H(div) degree 3 element = FIAT.TrimmedSerendipityDiv(FIAT.reference_element.UFCHexahedron(), 3) This implementation is correct for all the examples below. |
Symfem | "TSdiv" ↓ 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 trimmed serendipity H(div) degree 1 on a quadrilateral element = symfem.create_element("quadrilateral", "TSdiv", 1)
# Create trimmed serendipity H(div) degree 2 on a quadrilateral element = symfem.create_element("quadrilateral", "TSdiv", 2)
# Create trimmed serendipity H(div) degree 3 on a quadrilateral element = symfem.create_element("quadrilateral", "TSdiv", 3)
# Create trimmed serendipity H(div) degree 1 on a hexahedron element = symfem.create_element("hexahedron", "TSdiv", 1)
# Create trimmed serendipity H(div) degree 2 on a hexahedron element = symfem.create_element("hexahedron", "TSdiv", 2)
# Create trimmed serendipity H(div) degree 3 on a hexahedron element = symfem.create_element("hexahedron", "TSdiv", 3) This implementation is used to compute the examples below and verify other implementations. |
- Cockburn, Bernardo and Fu, Guosheng. A systematic construction of finite element commuting exact sequences, SIAM Journal of Numerical Analysis 55(4), 1650–1688, 2017. [DOI: 10.1137/16M1073352] [BibTeX]
- Gillette, Andrew and Kloefkorn, Tyler. Trimmed serendipity finite element differential forms, Mathematics of Computation 88, 583–606, 2019. [DOI: 10.1090/mcom/3354] [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 | 07 October 2021 |
Element last updated | 16 October 2024 |