Skip to content
This repository has been archived by the owner on Dec 6, 2024. It is now read-only.

Commit

Permalink
tests: Update underintegration test in light of UFL changes
Browse files Browse the repository at this point in the history
The estimated quadrature degree for non-affine cells now gets a
contribution from detJ, so we need to specify the rule for the DG mass
matrix too.
  • Loading branch information
wence- committed Apr 3, 2019
1 parent ddd4b1a commit 6e5a51b
Showing 1 changed file with 15 additions and 7 deletions.
22 changes: 15 additions & 7 deletions tests/test_underintegration.py
Original file line number Diff line number Diff line change
Expand Up @@ -10,9 +10,9 @@
action, interval, quadrilateral, dot, grad)

from FIAT import ufc_cell
from FIAT.quadrature import GaussLobattoLegendreQuadratureLineRule
from FIAT.quadrature import GaussLobattoLegendreQuadratureLineRule, GaussLegendreQuadratureLineRule

from finat.point_set import GaussLobattoLegendrePointSet
from finat.point_set import GaussLobattoLegendrePointSet, GaussLegendrePointSet
from finat.quadrature import QuadratureRule, TensorProductQuadratureRule

from tsfc import compile_form
Expand All @@ -28,30 +28,38 @@ def gll_quadrature_rule(cell, elem_deg):
return finat_rule


def gl_quadrature_rule(cell, elem_deg):
fiat_cell = ufc_cell("interval")
fiat_rule = GaussLegendreQuadratureLineRule(fiat_cell, elem_deg + 1)
line_rules = [QuadratureRule(GaussLegendrePointSet(fiat_rule.get_points()),
fiat_rule.get_weights())
for _ in range(cell.topological_dimension())]
finat_rule = reduce(lambda a, b: TensorProductQuadratureRule([a, b]), line_rules)
return finat_rule


def mass_cg(cell, degree):
m = Mesh(VectorElement('Q', cell, 1))
V = FunctionSpace(m, FiniteElement('Q', cell, degree, variant='spectral'))
u = TrialFunction(V)
v = TestFunction(V)
return u*v*dx(rule=gll_quadrature_rule(cell, degree))
return u*v*dx(scheme=gll_quadrature_rule(cell, degree))


def mass_dg(cell, degree):
m = Mesh(VectorElement('Q', cell, 1))
V = FunctionSpace(m, FiniteElement('DQ', cell, degree, variant='spectral'))
u = TrialFunction(V)
v = TestFunction(V)
# In this case, the estimated quadrature degree will give the
# correct number of quadrature points by luck.
return u*v*dx
return u*v*dx(scheme=gl_quadrature_rule(cell, degree))


def laplace(cell, degree):
m = Mesh(VectorElement('Q', cell, 1))
V = FunctionSpace(m, FiniteElement('Q', cell, degree, variant='spectral'))
u = TrialFunction(V)
v = TestFunction(V)
return dot(grad(u), grad(v))*dx(rule=gll_quadrature_rule(cell, degree))
return dot(grad(u), grad(v))*dx(scheme=gll_quadrature_rule(cell, degree))


def count_flops(form):
Expand Down

0 comments on commit 6e5a51b

Please sign in to comment.