Click here to read what the information on this page means.
Degrees | \(k=2\) where \(k\) is the Polynomial superdegree |
Polynomial subdegree | \(k-1\) |
Polynomial superdegree | \(k\) |
Reference elements | triangle |
Polynomial set | \(\mathcal{Z}^{(21)}_{k}\) ↓ Show polynomial set definitions ↓↑ Hide polynomial set definitions ↑\(\mathcal{Z}^{(21)}_k=\left\{\mathbf{M}\in\mathcal{P}_{k}^{d\times d}\middle|\mathbf{M}^t=\mathbf{M}\text{ and for every edge }\hat{\mathbf{n}}_i^t\mathbf{M}\hat{\mathbf{n}}_i\text{ is linear }\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: integral moments of normal-normal and normal-tangent inner products with an degree \(1\) Lagrange space
On each face: integral moments of three components with an degree \(0\) Lagrange space |
Number of DOFs | triangle: \(15\) |
Mapping | double contravariant Piola |
continuity | Inner products with normals to facets are continuous |
Categories | Matrix-valued elements |
This element is implemented in
FIAT ,
Symfem , and
(legacy) UFL.
↓ Show implementation detail ↓↑ Hide implementation detail ↑FIAT | FIAT.ArnoldWintherNC ↓ 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 nonconforming Arnold-Winther degree 2 element = FIAT.ArnoldWintherNC(FIAT.ufc_cell("triangle"), 2) This implementation is correct for all the examples below.Note: This implementation includes additional DOFs that are used then filtered out when mapping the element, as described in Kirby (2018). |
Symfem | "nonconforming AW" ↓ 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 nonconforming Arnold-Winther degree 2 on a triangle element = symfem.create_element("triangle", "nonconforming AW", 2) This implementation is used to compute the examples below and verify other implementations. |
(legacy) UFL | "AWnc" ↓ 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 nonconforming Arnold-Winther degree 2 on a triangle element = ufl_legacy.FiniteElement("AWnc", "triangle", 2) |
- Arnold, Douglas N. and Winther, Ragnar. Nonconforming mixed elements for elasticity, Numerische Mathematik 13(3), 295–307, 2003. [DOI: 10.1142/S0218202503002507] [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 | 02 August 2021 |
Element last updated | 16 October 2024 |