Skip to content

Commit 65502ac

Browse files
committed
Submesh: support tuple subdomain_id
1 parent 951bcd2 commit 65502ac

2 files changed

Lines changed: 31 additions & 0 deletions

File tree

firedrake/mesh.py

Lines changed: 16 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -4876,8 +4876,24 @@ 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, list)):
4881+
# A list of subdomain ids requires us to build an internal DM label with the union
4882+
iset = PETSc.IS().createGeneral([], comm=comm or mesh.comm)
4883+
for sub in subdomain_id:
4884+
iset = iset.union(plex.getStratumIS(label_name, sub))
4885+
4886+
label_name = "temp_union"
4887+
subdomain_id = 1
4888+
plex.createLabel(label_name)
4889+
label = plex.getLabel(label_name)
4890+
label.setStratumIS(subdomain_id, iset)
4891+
48794892
subplex = dmcommon.submesh_create(plex, subdim, label_name, subdomain_id, ignore_halo, comm=comm)
48804893

4894+
if label_name == "temp_union":
4895+
plex.removeLabel(label_name)
4896+
48814897
comm = comm or mesh.comm
48824898
name = name or _generate_default_submesh_name(mesh.name)
48834899
subplex.setName(_generate_default_mesh_topology_name(name))

tests/firedrake/submesh/test_submesh_facet.py

Lines changed: 15 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -134,3 +134,18 @@ 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+
assert abs(assemble(1*dx(domain=submesh1)) - len(subdomain_id)) < 1E-12
144+
145+
V = FunctionSpace(mesh, "HDiv Trace", 0)
146+
facet_function = Function(V)
147+
DirichletBC(V, 1, subdomain_id).apply(facet_function)
148+
facet_value = 999
149+
rmesh = RelabeledMesh(mesh, [facet_function], [facet_value])
150+
submesh2 = Submesh(rmesh, mesh.topological_dimension - 1, facet_value)
151+
assert submesh2.cell_set.size == submesh1.cell_set.size

0 commit comments

Comments
 (0)