Click here to read what the information on this page means.
Abbreviated names | BDFM |
Degrees | \(1\leqslant k\) where \(k\) is the Lagrange superdegree |
Polynomial subdegree | \(k-1\) |
Polynomial superdegree | \(k\) |
Lagrange subdegree | triangle: \(k-1\) tetrahedron: \(k-1\) quadrilateral: \(\operatorname{floor}(k/d)\) hexahedron: \(\operatorname{floor}(k/d)\) |
Lagrange superdegree | \(k\) |
Reference elements | triangle, quadrilateral, tetrahedron, hexahedron |
Polynomial set | \(\mathcal{Z}^{(10)}_{k}\) ↓ Show polynomial set definitions ↓↑ Hide polynomial set definitions ↑\(\mathcal{Z}^{(10)}_k=\left\{\boldsymbol{p}\in\mathcal{P}_{k}^d\middle|\boldsymbol{p}\cdot\boldsymbol{n}\in\mathcal{P}_{k-1}\text{ on each facet}\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 facet (triangle): normal integral moments with an degree \(k-1\) Lagrange space
On each facet (tetrahedron): normal integral moments with an degree \(k-1\) Lagrange space
On each facet (quadrilateral): normal integral moments with an degree \(k-1\) dPc space
On each facet (hexahedron): normal integral moments with an degree \(k-1\) dPc space
On the interior of the reference element (triangle): integral moments with an degree \(k-1\) Nédélec (first kind) space
On the interior of the reference element (tetrahedron): integral moments with an degree \(k-1\) Nédélec (first kind) space
On the interior of the reference element (quadrilateral): integral moments with an degree \(k-2\) vector dPc space
On the interior of the reference element (hexahedron): integral moments with an degree \(k-2\) vector dPc space |
Number of DOFs | triangle: \(k^2+3k-1\) quadrilateral: \(k(k+3)\) (A028552) tetrahedron: \((k+1)(k^2+5k-2)/2\) hexahedron: \(k(k+1)(k+5)/2\) |
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.BrezziDouglasFortinMarini ↓ 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 Brezzi-Douglas-Fortin-Marini degree 2 element = FIAT.BrezziDouglasFortinMarini(FIAT.ufc_cell("triangle"), 2) This implementation is correct for all the examples below that it supports. Correct: triangle,2 Not implemented: triangle,1; quadrilateral,1; quadrilateral,2; tetrahedron,2; hexahedron,2 |
Symfem | "BDFM" ↓ 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 Brezzi-Douglas-Fortin-Marini degree 1 on a triangle element = symfem.create_element("triangle", "BDFM", 1)
# Create Brezzi-Douglas-Fortin-Marini degree 2 on a triangle element = symfem.create_element("triangle", "BDFM", 2)
# Create Brezzi-Douglas-Fortin-Marini degree 1 on a quadrilateral element = symfem.create_element("quadrilateral", "BDFM", 1)
# Create Brezzi-Douglas-Fortin-Marini degree 2 on a quadrilateral element = symfem.create_element("quadrilateral", "BDFM", 2)
# Create Brezzi-Douglas-Fortin-Marini degree 2 on a tetrahedron element = symfem.create_element("tetrahedron", "BDFM", 2)
# Create Brezzi-Douglas-Fortin-Marini degree 2 on a hexahedron element = symfem.create_element("hexahedron", "BDFM", 2) This implementation is used to compute the examples below and verify other implementations. |
(legacy) UFL | "BDFM" ↓ 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 Brezzi-Douglas-Fortin-Marini degree 1 on a triangle element = ufl_legacy.FiniteElement("BDFM", "triangle", 1)
# Create Brezzi-Douglas-Fortin-Marini degree 2 on a triangle element = ufl_legacy.FiniteElement("BDFM", "triangle", 2)
# Create Brezzi-Douglas-Fortin-Marini degree 2 on a tetrahedron element = ufl_legacy.FiniteElement("BDFM", "tetrahedron", 2) |
- Brezzi, Franco, Douglas, Jim, Fortin, Michel, and Marini, L. Donatella. Efficient rectangular mixed finite elements in two and three space variables, ESAIM: Mathematical Modelling and Numerical Analysis 21, 581–604, 1987. [DOI: 10.1051/m2an/1987210405811] [BibTeX]
- Brezzi, Franco and Fortin, Michel. Function spaces and finite element approximations, in Mixed and hybrid finite element methods (eds: Brezzi, Franco and Fortin, Michel), 1991. [DOI: 10.1007/978-1-4612-3172-1_3] [BibTeX]
Element added | 30 January 2021 |
Element last updated | 16 October 2024 |