Skip to content

Commit

Permalink
Merge pull request FEniCS#41 from firedrakeproject/vertexcell
Browse files Browse the repository at this point in the history
Add Point/Vertex Cell
  • Loading branch information
wence- authored Jun 26, 2020
2 parents ccc94c3 + 09c34de commit e74e070
Show file tree
Hide file tree
Showing 6 changed files with 65 additions and 7 deletions.
3 changes: 3 additions & 0 deletions AUTHORS
Original file line number Diff line number Diff line change
Expand Up @@ -73,3 +73,6 @@ Contributors:

Cyrus Cheng
email: [email protected]

Reuben W. Hill
email: [email protected]
5 changes: 4 additions & 1 deletion FIAT/P0.py
Original file line number Diff line number Diff line change
Expand Up @@ -19,7 +19,10 @@ def __init__(self, ref_el):
entity_ids = {}
nodes = []
vs = numpy.array(ref_el.get_vertices())
bary = tuple(numpy.average(vs, 0))
if ref_el.get_dimension() == 0:
bary = ()
else:
bary = tuple(numpy.average(vs, 0))

nodes = [functional.PointEvaluation(ref_el, bary)]
entity_ids = {}
Expand Down
40 changes: 37 additions & 3 deletions FIAT/expansions.py
Original file line number Diff line number Diff line change
Expand Up @@ -109,6 +109,34 @@ def xi_tetrahedron(eta):
return xi1, xi2, xi3


class PointExpansionSet(object):
"""Evaluates the point basis on a point reference element."""

def __init__(self, ref_el):
if ref_el.get_spatial_dimension() != 0:
raise ValueError("Must have a point")
self.ref_el = ref_el
self.base_ref_el = reference_element.Point()

def get_num_members(self, n):
return 1

def tabulate(self, n, pts):
"""Returns a numpy array A[i,j] = phi_i(pts[j]) = 1.0."""
assert n == 0
return numpy.ones((1, len(pts)))

def tabulate_derivatives(self, n, pts):
"""Returns a numpy array of size A where A[i,j] = phi_i(pts[j])
but where each element is an empty tuple (). This maintains
compatibility with the interfaces of the interval, triangle and
tetrahedron expansions."""
deriv_vals = numpy.empty_like(self.tabulate(n, pts), dtype=tuple)
deriv_vals.fill(())

return deriv_vals


class LineExpansionSet(object):
"""Evaluates the Legendre basis on a line reference element."""

Expand Down Expand Up @@ -377,7 +405,9 @@ def tabulate_jet(self, n, pts, order=1):
def get_expansion_set(ref_el):
"""Returns an ExpansionSet instance appopriate for the given
reference element."""
if ref_el.get_shape() == reference_element.LINE:
if ref_el.get_shape() == reference_element.POINT:
return PointExpansionSet(ref_el)
elif ref_el.get_shape() == reference_element.LINE:
return LineExpansionSet(ref_el)
elif ref_el.get_shape() == reference_element.TRIANGLE:
return TriangleExpansionSet(ref_el)
Expand All @@ -390,11 +420,15 @@ def get_expansion_set(ref_el):
def polynomial_dimension(ref_el, degree):
"""Returns the dimension of the space of polynomials of degree no
greater than degree on the reference element."""
if ref_el.get_shape() == reference_element.LINE:
if ref_el.get_shape() == reference_element.POINT:
if degree > 0:
raise ValueError("Only degree zero polynomials supported on point elements.")
return 1
elif ref_el.get_shape() == reference_element.LINE:
return max(0, degree + 1)
elif ref_el.get_shape() == reference_element.TRIANGLE:
return max((degree + 1) * (degree + 2) // 2, 0)
elif ref_el.get_shape() == reference_element.TETRAHEDRON:
return max(0, (degree + 1) * (degree + 2) * (degree + 3) // 6)
else:
raise Exception("Unknown reference element type.")
raise ValueError("Unknown reference element type.")
7 changes: 6 additions & 1 deletion FIAT/polynomial_set.py
Original file line number Diff line number Diff line change
Expand Up @@ -75,7 +75,12 @@ def tabulate(self, pts, jet_order=0):
for i in range(jet_order + 1):
alphas = mis(self.ref_el.get_spatial_dimension(), i)
for alpha in alphas:
D = form_matrix_product(self.dmats, alpha)
if len(self.dmats) > 0:
D = form_matrix_product(self.dmats, alpha)
else:
# special for vertex without defined point location
assert pts == [()]
D = numpy.eye(1)
result[alpha] = numpy.dot(self.coeffs,
numpy.dot(numpy.transpose(D),
base_vals))
Expand Down
12 changes: 12 additions & 0 deletions FIAT/reference_element.py
Original file line number Diff line number Diff line change
Expand Up @@ -6,6 +6,7 @@
#
# Modified by David A. Ham ([email protected]), 2014
# Modified by Lizao Li ([email protected]), 2016

"""
Abstract class and particular implementations of finite element
reference simplex geometry/topology.
Expand Down Expand Up @@ -482,6 +483,15 @@ def __init__(self):
topology = {0: {0: (0,)}}
super(Point, self).__init__(POINT, verts, topology)

def construct_subelement(self, dimension):
"""Constructs the reference element of a cell subentity
specified by subelement dimension.
:arg dimension: subentity dimension (integer). Must be zero.
"""
assert dimension == 0
return self


class DefaultLine(Simplex):
"""This is the reference line with vertices (-1.0,) and (1.0,)."""
Expand Down Expand Up @@ -968,6 +978,8 @@ def ufc_cell(cell):
return UFCQuadrilateral()
elif celltype == "hexahedron":
return UFCHexahedron()
elif celltype == "vertex":
return ufc_simplex(0)
elif celltype == "interval":
return ufc_simplex(1)
elif celltype == "triangle":
Expand Down
5 changes: 3 additions & 2 deletions test/unit/test_fiat.py
Original file line number Diff line number Diff line change
Expand Up @@ -20,7 +20,7 @@
import pytest

from FIAT.reference_element import LINE, ReferenceElement
from FIAT.reference_element import UFCInterval, UFCTriangle, UFCTetrahedron
from FIAT.reference_element import Point, UFCInterval, UFCTriangle, UFCTetrahedron
from FIAT.lagrange import Lagrange
from FIAT.discontinuous_lagrange import DiscontinuousLagrange # noqa: F401
from FIAT.discontinuous_taylor import DiscontinuousTaylor # noqa: F401
Expand Down Expand Up @@ -49,7 +49,7 @@
from FIAT.enriched import EnrichedElement # noqa: F401
from FIAT.nodal_enriched import NodalEnrichedElement


P = Point()
I = UFCInterval() # noqa: E741
T = UFCTriangle()
S = UFCTetrahedron()
Expand Down Expand Up @@ -110,6 +110,7 @@ def __init__(self, a, b):
"P0(I)",
"P0(T)",
"P0(S)",
"DiscontinuousLagrange(P, 0)",
"DiscontinuousLagrange(I, 0)",
"DiscontinuousLagrange(I, 1)",
"DiscontinuousLagrange(I, 2)",
Expand Down

0 comments on commit e74e070

Please sign in to comment.