Skip to content

Commit 8f6011e

Browse files
committed
Submesh support tuple subdomain_id
1 parent 951bcd2 commit 8f6011e

2 files changed

Lines changed: 26 additions & 0 deletions

File tree

firedrake/mesh.py

Lines changed: 12 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -4876,8 +4876,20 @@ def Submesh(mesh, subdim=None, subdomain_id=None, label_name=None, name=None, ig
48764876
label_name = dmcommon.CELL_SETS_LABEL
48774877
elif subdim == dim - 1:
48784878
label_name = dmcommon.FACE_SETS_LABEL
4879+
4880+
if isinstance(subdomain_id, tuple):
4881+
points = np.concatenate([plex.getStratumIS(label_name, sub).indices for sub in subdomain_id])
4882+
label_name = "temp_union"
4883+
subdomain_id = 1
4884+
plex.createLabel(label_name)
4885+
for p in points:
4886+
plex.setLabelValue(label_name, p, subdomain_id)
4887+
48794888
subplex = dmcommon.submesh_create(plex, subdim, label_name, subdomain_id, ignore_halo, comm=comm)
48804889

4890+
if label_name == "temp_union":
4891+
plex.removeLabel(label_name)
4892+
48814893
comm = comm or mesh.comm
48824894
name = name or _generate_default_submesh_name(mesh.name)
48834895
subplex.setName(_generate_default_mesh_topology_name(name))

tests/firedrake/submesh/test_submesh_facet.py

Lines changed: 14 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -134,3 +134,17 @@ def test_submesh_facet_all_facets():
134134
rmesh = RelabeledMesh(mesh, [facet_function], [facet_value])
135135
submesh2 = Submesh(rmesh, mesh.topological_dimension - 1, facet_value)
136136
assert submesh2.cell_set.size == submesh1.cell_set.size
137+
138+
139+
def test_submesh_facet_subdomain_id_tuple():
140+
mesh = UnitCubeMesh(2, 2, 2)
141+
subdomain_id = (1, 3, 6)
142+
submesh1 = Submesh(mesh, mesh.topological_dimension - 1, subdomain_id=subdomain_id)
143+
144+
V = FunctionSpace(mesh, "HDiv Trace", 0)
145+
facet_function = Function(V)
146+
DirichletBC(V, 1, subdomain_id).apply(facet_function)
147+
facet_value = 999
148+
rmesh = RelabeledMesh(mesh, [facet_function], [facet_value])
149+
submesh2 = Submesh(rmesh, mesh.topological_dimension - 1, facet_value)
150+
assert submesh2.cell_set.size == submesh1.cell_set.size

0 commit comments

Comments
 (0)