Skip to content

Commit debd46e

Browse files
authored
Fix regression in FiniteElement Python wrapper (#3552)
* Fix documentation in finite element wrapper and inverse permutation * Further typing fixes * Add dof inverse permutation check
1 parent fcc7719 commit debd46e

File tree

2 files changed

+48
-24
lines changed

2 files changed

+48
-24
lines changed

python/dolfinx/fem/element.py

Lines changed: 18 additions & 22 deletions
Original file line numberDiff line numberDiff line change
@@ -263,13 +263,15 @@ def signature(self) -> str:
263263
"""String identifying the finite element."""
264264
return self._cpp_object.signature
265265

266-
def T_apply(self, x: npt.NDArray[np.floating], cell_permutations: np.int32, dim: int) -> None:
266+
def T_apply(
267+
self, x: npt.NDArray[np.floating], cell_permutations: npt.NDArray[np.uint32], dim: int
268+
) -> None:
267269
"""Transform basis functions from the reference element ordering and orientation to the
268270
globally consistent physical element ordering and orientation.
269271
270272
Args:
271-
x: Data to transform (in place). The shape is ``(m, n)``, where `m` is the number of
272-
dgerees-of-freedom and the storage is row-major.
273+
x: Data to transform (in place). The shape is ``(num_cells, n, dim)``, where ``n`` is
274+
the number degrees-of-freedom and the data is flattened (row-major).
273275
cell_permutations: Permutation data for the cell.
274276
dim: Number of columns in ``data``.
275277
@@ -279,37 +281,31 @@ def T_apply(self, x: npt.NDArray[np.floating], cell_permutations: np.int32, dim:
279281
"""
280282
self._cpp_object.T_apply(x, cell_permutations, dim)
281283

282-
def Tt_apply(self, x: npt.NDArray[np.floating], cell_permutations: np.int32, dim: int) -> None:
284+
def Tt_apply(
285+
self, x: npt.NDArray[np.floating], cell_permutations: npt.NDArray[np.uint32], dim: int
286+
) -> None:
283287
"""Apply the transpose of the operator applied by T_apply().
284288
285289
Args:
286-
x: Data to transform (in place). The shape is ``(m, n)``, where `m` is the number of
287-
dgerees-of-freedom and the storage is row-major.
288-
cell_permutations: Permutation data for the cell.
289-
dim: Number of columns in `data`.
290-
291-
Note:
292-
Exposed for testing. Function is not vectorised across multiple cells. Please see
293-
`basix.numba_helpers` for performant versions.
290+
x: Data to transform (in place). The shape is ``(num_cells, n, dim)``, where ``n`` is
291+
the number degrees-of-freedom and the data is flattened (row-major).
292+
cell_permutations: Permutation data for the cells
293+
dim: Number of columns in ``data``.
294294
"""
295295
self._cpp_object.Tt_apply(x, cell_permutations, dim)
296296

297297
def Tt_inv_apply(
298-
self, x: npt.NDArray[np.floating], cell_permutations: np.int32, dim: int
298+
self, x: npt.NDArray[np.floating], cell_permutations: npt.NDArray[np.uint32], dim: int
299299
) -> None:
300300
"""Apply the inverse transpose of the operator applied by T_apply().
301301
302302
Args:
303-
x: Data to transform (in place). The shape is ``(m, n)``, where ``m`` is the number of
304-
dgerees-of-freedom and the storage is row-major.
305-
cell_permutations: Permutation data for the cell.
306-
dim: Number of columns in `data`.
307-
308-
Note:
309-
Exposed for testing. Function is not vectorised across multiple cells. Please see
310-
``basix.numba_helpers`` for performant versions.
303+
x: Data to transform (in place). The shape is ``(num_cells, n, dim)``, where ``n`` is
304+
the number degrees-of-freedom and the data is flattened (row-major).
305+
cell_permutations: Permutation data for the cells
306+
dim: Number of columns in ``data``.
311307
"""
312-
self._cpp_object.Tt_apply(x, cell_permutations, dim)
308+
self._cpp_object.Tt_inv_apply(x, cell_permutations, dim)
313309

314310

315311
def finiteelement(

python/test/unit/fem/test_dof_permuting.py

Lines changed: 30 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -1,4 +1,4 @@
1-
# Copyright (C) 2009-2020 Garth N. Wells, Matthew W. Scroggs and Jorgen S. Dokken
1+
# Copyright (C) 2009-2024 Garth N. Wells, Matthew W. Scroggs and Jorgen S. Dokken
22
#
33
# This file is part of DOLFINx (https://www.fenicsproject.org)
44
#
@@ -16,7 +16,7 @@
1616
from basix.ufl import element
1717
from dolfinx import default_real_type
1818
from dolfinx.fem import Function, assemble_scalar, form, functionspace
19-
from dolfinx.mesh import create_mesh
19+
from dolfinx.mesh import create_mesh, create_unit_cube
2020

2121

2222
def randomly_ordered_mesh(cell_type):
@@ -392,3 +392,31 @@ def tangent2(x):
392392

393393
value = assemble_scalar(form(_form))
394394
assert np.isclose(value, 0.0, rtol=1.0e-6, atol=1.0e-6)
395+
396+
397+
@pytest.mark.parametrize(
398+
"data_types",
399+
[
400+
(np.float32, np.float32),
401+
(np.float64, np.float64),
402+
(np.float32, np.complex64),
403+
(np.float64, np.complex128),
404+
],
405+
)
406+
@pytest.mark.parametrize("space_order", range(3, 5))
407+
def test_permutation_wrappers(space_order, data_types):
408+
s_type, d_type = data_types
409+
domain = create_unit_cube(MPI.COMM_WORLD, 5, 3, 4, dtype=s_type)
410+
V = functionspace(domain, ("N1curl", space_order))
411+
u = Function(V, name="u", dtype=d_type)
412+
u.interpolate(lambda x: np.array([x[0], x[1], x[2]]))
413+
414+
arr = u.x.array[V.dofmap.list]
415+
416+
domain.topology.create_entity_permutations()
417+
cell_perm = domain.topology.get_cell_permutation_info()
418+
org_data = arr.copy()
419+
V.element.Tt_apply(arr.reshape(-1), cell_perm, 1)
420+
V.element.Tt_inv_apply(arr.reshape(-1), cell_perm, 1)
421+
eps = 100 * np.finfo(s_type).eps
422+
np.testing.assert_allclose(org_data.reshape(-1), arr.reshape(-1), atol=eps)

0 commit comments

Comments
 (0)