Skip to content

Commit

Permalink
Tidy up BDMDualSet
Browse files Browse the repository at this point in the history
  • Loading branch information
pbrubeck committed Dec 30, 2023
1 parent 084f285 commit 861be97
Show file tree
Hide file tree
Showing 2 changed files with 25 additions and 30 deletions.
54 changes: 25 additions & 29 deletions FIAT/brezzi_douglas_marini.py
Original file line number Diff line number Diff line change
Expand Up @@ -9,6 +9,7 @@
polynomial_set, nedelec)
from FIAT.check_format_variant import check_format_variant
from FIAT.quadrature_schemes import create_quadrature
from FIAT.quadrature import FacetQuadratureRule


class BDMDualSet(dual_set.DualSet):
Expand All @@ -20,51 +21,46 @@ def __init__(self, ref_el, degree, variant, interpolant_deg):
nodes = []

sd = ref_el.get_spatial_dimension()
t = ref_el.get_topology()
top = ref_el.get_topology()

if variant == "integral":
facet = ref_el.get_facet_element()
# Facet nodes are \int_F v\cdot n p ds where p \in P_{q-1}
# degree is q - 1
Q = create_quadrature(facet, interpolant_deg + degree)
# Facet nodes are \int_F v\cdot n p ds where p \in P_{q}
# degree is q
Q_ref = create_quadrature(facet, interpolant_deg + degree)
Pq = polynomial_set.ONPolynomialSet(facet, degree)
Pq_at_qpts = Pq.tabulate(Q.get_points())[(0,)*(sd - 1)]
nodes.extend(functional.IntegralMomentOfScaledNormalEvaluation(ref_el, Q, phi, f)
for f in range(len(t[sd - 1]))
for phi in Pq_at_qpts)

# internal nodes
if degree > 1:
Q = create_quadrature(ref_el, interpolant_deg + degree - 1)
qpts = Q.get_points()
Nedel = nedelec.Nedelec(ref_el, degree - 1, variant)
Nedfs = Nedel.get_nodal_basis()
Ned_at_qpts = Nedfs.tabulate(qpts)[(0,) * sd]
Pq_at_qpts = Pq.tabulate(Q_ref.get_points())[(0,)*(sd - 1)]
for f in top[sd - 1]:
Q = FacetQuadratureRule(ref_el, sd - 1, f, Q_ref)
Jdet = Q.jacobian_determinant()
n = ref_el.compute_scaled_normal(f) / Jdet
phis = n[None, :, None] * Pq_at_qpts[:, None, :]
nodes.extend(functional.FrobeniusIntegralMoment(ref_el, Q, phi)
for phi in Ned_at_qpts)
for phi in phis)

elif variant == "point":
# Define each functional for the dual set
# codimension 1 facets
for i in range(len(t[sd - 1])):
for i in top[sd - 1]:
pts_cur = ref_el.make_points(sd - 1, i, sd + degree)
nodes.extend(functional.PointScaledNormalEvaluation(ref_el, i, pt)
for pt in pts_cur)

# internal nodes
if degree > 1:
Q = create_quadrature(ref_el, 2 * degree - 1)
qpts = Q.get_points()
Nedel = nedelec.Nedelec(ref_el, degree - 1, variant)
Nedfs = Nedel.get_nodal_basis()
Ned_at_qpts = Nedfs.tabulate(qpts)[(0,) * sd]
nodes.extend(functional.FrobeniusIntegralMoment(ref_el, Q, phi)
for phi in Ned_at_qpts)
# internal nodes
if degree > 1:
if interpolant_deg is None:
interpolant_deg = degree
Q = create_quadrature(ref_el, interpolant_deg + degree - 1)
Nedel = nedelec.Nedelec(ref_el, degree - 1, variant)
Nedfs = Nedel.get_nodal_basis()
Ned_at_qpts = Nedfs.tabulate(Q.get_points())[(0,) * sd]
nodes.extend(functional.FrobeniusIntegralMoment(ref_el, Q, phi)
for phi in Ned_at_qpts)

# sets vertices (and in 3d, edges) to have no nodes
for i in range(sd - 1):
entity_ids[i] = {}
for j in range(len(t[i])):
for j in top[i]:
entity_ids[i][j] = []

cur = 0
Expand All @@ -74,7 +70,7 @@ def __init__(self, ref_el, degree, variant, interpolant_deg):
pts_per_facet = len(pts_facet_0)

entity_ids[sd - 1] = {}
for i in range(len(t[sd - 1])):
for i in top[sd - 1]:
entity_ids[sd - 1][i] = list(range(cur, cur + pts_per_facet))
cur += pts_per_facet

Expand Down
1 change: 0 additions & 1 deletion FIAT/raviart_thomas.py
Original file line number Diff line number Diff line change
Expand Up @@ -77,7 +77,6 @@ def __init__(self, ref_el, degree, variant, interpolant_deg):
Q_ref = create_quadrature(facet, interpolant_deg + degree - 1)
Pq = polynomial_set.ONPolynomialSet(facet, degree - 1)
Pq_at_qpts = Pq.tabulate(Q_ref.get_points())[(0,)*(sd - 1)]

for f in top[sd - 1]:
cur = len(nodes)
Q = FacetQuadratureRule(ref_el, sd-1, f, Q_ref)
Expand Down

0 comments on commit 861be97

Please sign in to comment.