Click here to read what the information on this page means.
Abbreviated names | MTW |
Degrees | \(k=3\) |
Polynomial subdegree | \(1\) |
Polynomial superdegree | \(k+d-2\) |
Reference elements | triangle, tetrahedron |
Polynomial set | \(\mathcal{Z}^{(17)}_{k}\) (triangle)
\(\mathcal{P}_{k}^d \oplus \mathcal{Z}^{(18)}_{k}\) (tetrahedron)
↓ Show polynomial set definitions ↓↑ Hide polynomial set definitions ↑\(\mathcal{Z}^{(17)}_k=\left\{\boldsymbol{p}\in\mathcal{P}_{k}^d\middle|\operatorname{div}\boldsymbol{p}\text{ is constant, and }\boldsymbol{p}\cdot\boldsymbol{n}_{e_i}\text{ is linear on each edge }e_i\right\}\)
\(\mathcal{P}_k=\operatorname{span}\left\{\prod_{i=1}^dx_i^{p_i}\middle|\sum_{i=1}^dp_i\leqslant k\right\}\)
\(\mathcal{Z}^{(18)}_k=\left\{\nabla\times(xyz(1-x-y-z)\boldsymbol{p})\middle|\boldsymbol{p}\in\mathcal{P}_{k}^d\right\}\) |
DOFs | On each facet (triangle): normal integral moments with an degree \(1\) Lagrange space, and tangent integral moments with an degree \(0\) Lagrange space
On each facet (tetrahedron): normal integral moments with an degree \(1\) Lagrange space, and integral moments with an degree \(1\) Nédélec (first kind) space |
Number of DOFs | triangle: \(9\) tetrahedron: \(24\) |
Mapping | contravariant Piola |
continuity | Components normal to facets are continuous |
Categories | Vector-valued elements, H(div) conforming elements |
This element is implemented in
FIAT ,
Symfem , and
(legacy) UFL.
↓ Show implementation detail ↓↑ Hide implementation detail ↑FIAT | FIAT.MardalTaiWinther ↓ 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 Mardal-Tai-Winther degree 3 element = FIAT.MardalTaiWinther(FIAT.ufc_cell("triangle"), 3) This implementation is correct for all the examples below that it supports. Correct: triangle,3 Not implemented: tetrahedron,3 Note: This implementation includes additional DOFs that are used then filtered out when mapping the element, as described in Kirby (2018). |
Symfem | "MTW" ↓ 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 Mardal-Tai-Winther degree 3 on a triangle element = symfem.create_element("triangle", "MTW", 3)
# Create Mardal-Tai-Winther degree 3 on a tetrahedron element = symfem.create_element("tetrahedron", "MTW", 3) This implementation is used to compute the examples below and verify other implementations. |
(legacy) UFL | "MTW" ↓ 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 Mardal-Tai-Winther degree 3 on a triangle element = ufl_legacy.FiniteElement("MTW", "triangle", 3) |
- Mardal, Kent-Andre, Tai, Xue-Cheng, and Winther, Ragner. A robust finite element method for Darcy–Stokes flow, SIAM Journal on Numerical Analysis 40(5), 1605–1631, 2002. [DOI: 10.1137/S0036142901383910] [BibTeX]
- Tai, Xue-Cheng and Winther, Ragner. A discrete de Rham complex with enhanced smoothness, Calcolo 43, 287–306, 2006. [DOI: 10.1007/s10092-006-0124-6] [BibTeX]
- Kirby, Robert C. A general approach to transforming finite elements, SMAI Journal of Computational Mathematics 4, 197–224, 2018. [DOI: 10.5802/smai-jcm.33] [BibTeX]
Element added | 09 January 2021 |
Element last updated | 27 September 2024 |