Skip to content

Commit

Permalink
Added 4d element types (#188)
Browse files Browse the repository at this point in the history
* [feature][ElementTypes] Added 4d types pentatope (4d simplex) and tesseract (4d cube). (#95)

Co-authored-by: Jack S. Hale <[email protected]>
Co-authored-by: David A. Ham <[email protected]>
Co-authored-by: Matthew Scroggs <[email protected]>

* fix merge, add test

---------

Co-authored-by: dr-robertk <[email protected]>
Co-authored-by: Jack S. Hale <[email protected]>
Co-authored-by: David A. Ham <[email protected]>
  • Loading branch information
4 people authored Aug 10, 2023
1 parent a49736c commit 9e120c6
Show file tree
Hide file tree
Showing 5 changed files with 40 additions and 6 deletions.
23 changes: 22 additions & 1 deletion test/test_cell.py
Original file line number Diff line number Diff line change
Expand Up @@ -37,6 +37,20 @@ def test_hexahedron():
assert cell.num_faces() == 6


def test_pentatope():
cell = ufl.pentatope
assert cell.num_vertices() == 5
assert cell.num_edges() == 10
assert cell.num_faces() == 10


def test_tesseract():
cell = ufl.tesseract
assert cell.num_vertices() == 16
assert cell.num_edges() == 32
assert cell.num_faces() == 24



@pytest.mark.parametrize("cell", [ufl.interval])
def test_cells_1d(cell):
Expand All @@ -53,12 +67,19 @@ def test_cells_2d(cell):


@pytest.mark.parametrize("cell", [ufl.tetrahedron, ufl.hexahedron, ufl.prism, ufl.pyramid])
def test_cells_2d(cell):
def test_cells_3d(cell):
assert cell.num_facets() == cell.num_faces()
assert cell.num_ridges() == cell.num_edges()
assert cell.num_peaks() == cell.num_vertices()


@pytest.mark.parametrize("cell", [ufl.tesseract, ufl.pentatope])
def test_cells_4d(cell):
assert cell.num_facets() == cell.num_sub_entities(3)
assert cell.num_ridges() == cell.num_faces()
assert cell.num_peaks() == cell.num_edges()


def test_tensorproductcell():
orig = ufl.TensorProductCell(ufl.interval, ufl.interval)
cell = orig.reconstruct()
Expand Down
6 changes: 4 additions & 2 deletions ufl/__init__.py
Original file line number Diff line number Diff line change
Expand Up @@ -53,6 +53,8 @@
- hexahedron
- prism
- pyramid
- pentatope
- tesseract
* Domains::
Expand Down Expand Up @@ -339,7 +341,7 @@

# Predefined convenience objects
from ufl.objects import (
vertex, interval, triangle, tetrahedron,
vertex, interval, triangle, tetrahedron, pentatope, tesseract,
quadrilateral, hexahedron, prism, pyramid, facet,
i, j, k, l, p, q, r, s,
dx, ds, dS, dP,
Expand Down Expand Up @@ -400,7 +402,7 @@
'dc', 'dC', 'dO', 'dI', 'dX',
'ds_b', 'ds_t', 'ds_tb', 'ds_v', 'dS_h', 'dS_v',
'vertex', 'interval', 'triangle', 'tetrahedron',
'prism', 'pyramid',
'prism', 'pyramid', 'pentatope', 'tesseract',
'quadrilateral', 'hexahedron', 'facet',
'i', 'j', 'k', 'l', 'p', 'q', 'r', 's',
'e', 'pi',
Expand Down
8 changes: 8 additions & 0 deletions ufl/cell.py
Original file line number Diff line number Diff line change
Expand Up @@ -196,6 +196,10 @@ def peak_types(self) -> typing.Tuple[AbstractCell, ...]:
("triangle", "quadrilateral", "quadrilateral", "quadrilateral", "triangle"), ("prism", )],
"pyramid": [tuple("vertex" for i in range(5)), tuple("interval" for i in range(8)),
("quadrilateral", "triangle", "triangle", "triangle", "triangle"), ("pyramid", )],
"pentatope": [tuple("vertex" for i in range(5)), tuple("interval" for i in range(10)),
tuple("triangle" for i in range(10)), tuple("tetrahedron" for i in range(5)), ("pentatope", )],
"tesseract": [tuple("vertex" for i in range(16)), tuple("interval" for i in range(32)),
tuple("quadrilateral" for i in range(24)), tuple("hexahedron" for i in range(8)), ("tesseract", )],
}


Expand Down Expand Up @@ -422,6 +426,8 @@ def simplex(topological_dimension: int, geometric_dimension: typing.Optional[int
return Cell("triangle", geometric_dimension)
if topological_dimension == 3:
return Cell("tetrahedron", geometric_dimension)
if topological_dimension == 4:
return Cell("pentatope", geometric_dimension)
raise ValueError(f"Unsupported topological dimension for simplex: {topological_dimension}")


Expand All @@ -435,6 +441,8 @@ def hypercube(topological_dimension, geometric_dimension=None):
return Cell("quadrilateral", geometric_dimension)
if topological_dimension == 3:
return Cell("hexahedron", geometric_dimension)
if topological_dimension == 4:
return Cell("tesseract", geometric_dimension)
raise ValueError(f"Unsupported topological dimension for hypercube: {topological_dimension}")


Expand Down
7 changes: 4 additions & 3 deletions ufl/finiteelement/elementlist.py
Original file line number Diff line number Diff line change
Expand Up @@ -12,6 +12,7 @@
# Modified by Marie E. Rognes <[email protected]>, 2010
# Modified by Lizao Li <[email protected]>, 2015, 2016
# Modified by Massimiliano Leoni, 2016
# Modified by Robert Kloefkorn, 2022

import warnings
from numpy import asarray
Expand Down Expand Up @@ -82,12 +83,12 @@ def show_elements():
# the future, add mapping name as another element property.

# Cell groups
simplices = ("interval", "triangle", "tetrahedron")
cubes = ("interval", "quadrilateral", "hexahedron")
simplices = ("interval", "triangle", "tetrahedron", "pentatope")
cubes = ("interval", "quadrilateral", "hexahedron", "tesseract")
any_cell = (None,
"vertex", "interval",
"triangle", "tetrahedron", "prism",
"pyramid", "quadrilateral", "hexahedron")
"pyramid", "quadrilateral", "hexahedron", "pentatope", "tesseract")

# Elements in the periodic table # TODO: Register these as aliases of
# periodic table element description instead of the other way around
Expand Down
2 changes: 2 additions & 0 deletions ufl/objects.py
Original file line number Diff line number Diff line change
Expand Up @@ -37,6 +37,8 @@
pyramid = Cell("pyramid", 3)
quadrilateral = Cell("quadrilateral", 2)
hexahedron = Cell("hexahedron", 3)
tesseract = Cell("tesseract", 4)
pentatope = Cell("pentatope", 4)

# Facet is just a dummy declaration for RestrictedElement
facet = "facet"

0 comments on commit 9e120c6

Please sign in to comment.