Skip to content

Commit 99e29fb

Browse files
authored
Vertex submesh support (#3455)
* Add tests for submesh of vertices * Redirect basix branch * Ruff formatting * Add fixes and simplifications for submesh of vertices * Tidy up test * Ignore when edim becomes negative by setting to zero in parametrized test * Revert branches * Remove formatting error
1 parent 8261439 commit 99e29fb

File tree

3 files changed

+16
-10
lines changed

3 files changed

+16
-10
lines changed

cpp/dolfinx/mesh/utils.h

Lines changed: 6 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -1250,7 +1250,12 @@ create_subgeometry(const Mesh<T>& mesh, int dim,
12501250
// Create sub-geometry coordinate element
12511251
CellType sub_coord_cell
12521252
= cell_entity_type(geometry.cmap().cell_shape(), dim, 0);
1253-
fem::CoordinateElement<T> sub_cmap(sub_coord_cell, geometry.cmap().degree(),
1253+
// Special handling if point meshes, as they only support constant basis
1254+
// functions
1255+
int degree = geometry.cmap().degree();
1256+
if (sub_coord_cell == CellType::point)
1257+
degree = 0;
1258+
fem::CoordinateElement<T> sub_cmap(sub_coord_cell, degree,
12541259
geometry.cmap().variant());
12551260

12561261
// Sub-geometry input_global_indices

python/dolfinx/mesh.py

Lines changed: 1 addition & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -675,11 +675,10 @@ def create_submesh(
675675
submsh, entity_map, vertex_map, geom_map = _cpp.mesh.create_submesh(
676676
msh._cpp_object, dim, entities
677677
)
678-
submsh_ufl_cell = ufl.Cell(to_string(submsh.topology.cell_type))
679678
submsh_domain = ufl.Mesh(
680679
basix.ufl.element(
681680
"Lagrange",
682-
submsh_ufl_cell.cellname(),
681+
to_string(submsh.topology.cell_type),
683682
submsh.geometry.cmap.degree,
684683
basix.LagrangeVariant(submsh.geometry.cmap.variant),
685684
shape=(submsh.geometry.dim,),

python/test/unit/mesh/test_mesh.py

Lines changed: 9 additions & 7 deletions
Original file line numberDiff line numberDiff line change
@@ -503,35 +503,37 @@ def boundary_2(x):
503503

504504

505505
# TODO Test that submesh of full mesh is a copy of the mesh
506-
@pytest.mark.parametrize("d", [2, 3])
506+
@pytest.mark.parametrize("d", [1, 2, 3])
507507
@pytest.mark.parametrize("n", [3, 6])
508508
@pytest.mark.parametrize("codim", [0, 1, 2])
509509
@pytest.mark.parametrize("marker", [lambda x: x[0] >= 0.5, lambda x: x[0] >= -1])
510510
@pytest.mark.parametrize("ghost_mode", [GhostMode.none, GhostMode.shared_facet])
511511
@pytest.mark.parametrize("simplex", [True, False])
512512
def test_submesh_full(d, n, codim, marker, ghost_mode, simplex):
513-
if d == codim:
514-
pytest.xfail("Cannot create vertex submesh")
515-
if d == 2:
513+
if d == 1:
514+
mesh = create_unit_interval(MPI.COMM_WORLD, n, ghost_mode=ghost_mode)
515+
elif d == 2:
516516
ct = CellType.triangle if simplex else CellType.quadrilateral
517517
mesh = create_unit_square(MPI.COMM_WORLD, n, n, ghost_mode=ghost_mode, cell_type=ct)
518518
else:
519519
ct = CellType.tetrahedron if simplex else CellType.hexahedron
520520
mesh = create_unit_cube(MPI.COMM_WORLD, n, n, n, ghost_mode=ghost_mode, cell_type=ct)
521521

522-
edim = mesh.topology.dim - codim
522+
edim = max(mesh.topology.dim - codim, 0)
523523
entities = locate_entities(mesh, edim, marker)
524524
submesh, entity_map, vertex_map, geom_map = create_submesh(mesh, edim, entities)
525525
submesh_topology_test(mesh, submesh, entity_map, vertex_map, edim)
526526
submesh_geometry_test(mesh, submesh, entity_map, geom_map, edim)
527527

528528

529-
@pytest.mark.parametrize("d", [2, 3])
529+
@pytest.mark.parametrize("d", [1, 2, 3])
530530
@pytest.mark.parametrize("n", [3, 6])
531531
@pytest.mark.parametrize("boundary", [boundary_0, boundary_1, boundary_2])
532532
@pytest.mark.parametrize("ghost_mode", [GhostMode.none, GhostMode.shared_facet])
533533
def test_submesh_boundary(d, n, boundary, ghost_mode):
534-
if d == 2:
534+
if d == 1:
535+
mesh = create_unit_interval(MPI.COMM_WORLD, n, ghost_mode=ghost_mode)
536+
elif d == 2:
535537
mesh = create_unit_square(MPI.COMM_WORLD, n, n, ghost_mode=ghost_mode)
536538
else:
537539
mesh = create_unit_cube(MPI.COMM_WORLD, n, n, n, ghost_mode=ghost_mode)

0 commit comments

Comments
 (0)