Click here to read what the information on this page means.
Alternative names | non-conforming Crouzeix–Raviart |
Abbreviated names | CR |
Degrees | \(k=1\) where \(k\) is the Polynomial superdegree |
Polynomial subdegree | \(k\) |
Polynomial superdegree | \(k\) |
Reference elements | triangle, tetrahedron |
Polynomial set | \(\mathcal{P}_{k}\) ↓ 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\}\) |
DOFs | On each facet: point evaluation at midpoint |
Number of DOFs | triangle: \(3\) tetrahedron: \(4\) |
Mapping | identity |
continuity | Discontinuous. |
Categories | Scalar-valued elements |
This element is implemented in
Basix ,
Basix.UFL ,
FIAT ,
Symfem , and
(legacy) UFL.
↓ Show implementation detail ↓↑ Hide implementation detail ↑Basix | basix.ElementFamily.CR ↓ Show Basix examples ↓↑ Hide Basix examples ↑Before running this example, you must install Basix: pip3 install fenics-basix This element can then be created with the following lines of Python: import basix
# Create Crouzeix-Raviart degree 1 on a triangle element = basix.create_element(basix.ElementFamily.CR, basix.CellType.triangle, 1) This implementation is correct for all the examples below. |
Basix.UFL | basix.ElementFamily.CR ↓ Show Basix.UFL examples ↓↑ Hide Basix.UFL examples ↑Before running this example, you must install Basix.UFL: pip3 install fenics-basix fenics-ufl This element can then be created with the following lines of Python: import basix import basix.ufl
# Create Crouzeix-Raviart degree 1 on a triangle element = basix.ufl.element(basix.ElementFamily.CR, basix.CellType.triangle, 1) This implementation is correct for all the examples below. |
FIAT | FIAT.CrouzeixRaviart ↓ 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 Crouzeix-Raviart degree 1 element = FIAT.CrouzeixRaviart(FIAT.ufc_cell("triangle"), 1) This implementation is correct for all the examples below. |
Symfem | "Crouzeix-Raviart" ↓ 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 Crouzeix-Raviart degree 1 on a triangle element = symfem.create_element("triangle", "Crouzeix-Raviart", 1) This implementation is used to compute the examples below and verify other implementations. |
(legacy) UFL | "Crouzeix-Raviart" ↓ 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 Crouzeix-Raviart degree 1 on a triangle element = ufl_legacy.FiniteElement("Crouzeix-Raviart", "triangle", 1) |
- Crouzeix, Michel and Raviart, Pierre-Arnaud. Conforming and nonconforming finite element methods for solving the stationary Stokes equations, Revue Française d'Automatique, Informatique et Recherche Opérationnelle 3, 33–75, 1973. [DOI: 10.1051/m2an/197307R300331] [BibTeX]
Element added | 01 January 2021 |
Element last updated | 27 September 2024 |