Skip to content

Commit 1ea0f4e

Browse files
committed
Merge branch 'main' of github.com:mscroggs/feast into main
2 parents bb0f2b9 + e0649d3 commit 1ea0f4e

File tree

3 files changed

+51
-0
lines changed

3 files changed

+51
-0
lines changed

symfem/create.py

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -143,6 +143,7 @@ def create_element(cell_type, element_type, order, variant="equispaced"):
143143
matrix discontinuous Lagrange,
144144
symmetric matrix discontinuous Lagrange,
145145
Crouzeix-Raviart, CR, Crouzeix-Falk, CF,
146+
conforming Crouzeix-Raviart, conforming CR,
146147
serendipity, S,
147148
serendipity Hcurl, Scurl, BDMCE, AAE,
148149
serendipity Hdiv, Sdiv, BDMCF, AAF,
Lines changed: 49 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,49 @@
1+
"""Conforming Crouzeix-Raviart elements on simplices.
2+
3+
This element's definition appears in https://doi.org/10.1051/m2an/197307R300331
4+
(Crouzeix, Raviart, 1973)
5+
"""
6+
7+
import sympy
8+
from ..finite_element import CiarletElement
9+
from ..polynomials import polynomial_set
10+
from ..functionals import PointEvaluation
11+
from ..symbolic import x
12+
13+
14+
class ConformingCrouzeixRaviart(CiarletElement):
15+
"""Conforming Crouzeix-Raviart finite element."""
16+
17+
def __init__(self, reference, order, variant):
18+
assert reference.name == "triangle"
19+
20+
poly = polynomial_set(reference.tdim, 1, order)
21+
22+
poly += [
23+
x[0] ** i * x[1] ** (order - i) * (x[0] + x[1])
24+
for i in range(1, order)
25+
]
26+
27+
dofs = []
28+
for i, v in enumerate(reference.vertices):
29+
dofs.append(PointEvaluation(v, entity=(0, i)))
30+
if order >= 2:
31+
for i, edge in enumerate(reference.edges):
32+
for p in range(1, order):
33+
v = tuple(sympy.Rational((order - p) * a + p * b, order) for a, b in zip(
34+
reference.vertices[edge[0]], reference.vertices[edge[1]]))
35+
dofs.append(PointEvaluation(v, entity=(1, i)))
36+
for i in range(1, order):
37+
for j in range(1, order + 1 - i):
38+
point = (
39+
sympy.Rational(3 * i - 1, 3 * order),
40+
sympy.Rational(3 * j - 1, 3 * order)
41+
)
42+
dofs.append(PointEvaluation(point, entity=(2, 0)))
43+
44+
super().__init__(reference, order, poly, dofs, reference.tdim, 1)
45+
46+
names = ["conforming Crouzeix-Raviart", "conforming CR"]
47+
references = ["triangle"]
48+
min_order = 1
49+
continuity = "L2"

test/utils.py

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -32,6 +32,7 @@
3232
"bubble enriched vector Lagrange":
3333
{"equispaced": range(1, 3), "lobatto": range(1, 3)},
3434
"CR": {"equispaced": [1, 3, 5]},
35+
"conforming CR": {"equispaced": range(1, 6)},
3536
"HHJ": {"equispaced": range(3)},
3637
"Bell": {"equispaced": [5]},
3738
"Morley": {None: [2]},

0 commit comments

Comments
 (0)