Click here to read what the information on this page means.
| Degrees | \(k=1\) where \(k\) is the polynomial subdegree
 | 
| Polynomial subdegree | \(k\) | 
| Polynomial superdegree | \(k+1\) | 
| Reference cells | triangle | 
| Finite dimensional space | \(\mathcal{Z}^{(23)}_{k+1}\) ↓ Show set definitions ↓↑ Hide set definitions ↑
 \(\mathcal{Z}^{(23)}_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 a degree \(1\) Lagrange space On each face: integral moments of three components with a 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 ↑
 This implementation is correct for all the examples below.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 1
 element = FIAT.ArnoldWintherNC(FIAT.ufc_cell("triangle"), 2)
Note: This implementation includes additional DOFs that are used then filtered out when mapping the element, as described in Kirby (2018). Note: This element uses the Lagrange superdegree as the canonical degree of this element | 
| Symfem | "nonconforming AW"↓ 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 nonconforming Arnold-Winther degree 1 on a triangle
 element = symfem.create_element("triangle", "nonconforming AW", 1)
 | 
| (legacy) UFL | "AWnc"↓ 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 nonconforming Arnold-Winther degree 1 on a triangle
 element = ufl_legacy.FiniteElement("AWnc", "triangle", 2)
Note: This element uses the Lagrange superdegree as the canonical degree of this element | 
 
- 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 | 04 June 2025 |