Skip to content

Commit c33d1de

Browse files
committed
assemble: relax multi-domain restriction
1 parent d5391b3 commit c33d1de

File tree

2 files changed

+31
-11
lines changed

2 files changed

+31
-11
lines changed

firedrake/assemble.py

Lines changed: 0 additions & 11 deletions
Original file line numberDiff line numberDiff line change
@@ -1057,17 +1057,6 @@ def local_kernels(self):
10571057
each possible combination.
10581058
10591059
"""
1060-
try:
1061-
topology, = set(d.topology.submesh_ancesters[-1] for d in self._form.ufl_domains())
1062-
except ValueError:
1063-
raise NotImplementedError("All integration domains must share a mesh topology")
1064-
1065-
for o in itertools.chain(self._form.arguments(), self._form.coefficients()):
1066-
domains = extract_domains(o)
1067-
for domain in domains:
1068-
if domain is not None and domain.topology.submesh_ancesters[-1] != topology:
1069-
raise NotImplementedError("Assembly with multiple meshes is not supported")
1070-
10711060
if isinstance(self._form, ufl.Form):
10721061
kernels = tsfc_interface.compile_form(
10731062
self._form, "form", diagonal=self.diagonal,

tests/firedrake/submesh/test_submesh_solve.py

Lines changed: 31 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -7,6 +7,37 @@
77
cwd = abspath(dirname(__file__))
88

99

10+
# Solve an uncoupled system of two equations in one shot.
11+
def test_submesh_solve_uncoupled():
12+
mesh0 = UnitTriangleMesh()
13+
mesh1 = UnitSquareMesh(1, 1, quadrilateral=True)
14+
V0 = FunctionSpace(mesh0, "P", 1)
15+
V1 = FunctionSpace(mesh1, "Q", 1)
16+
V = V0 * V1
17+
u = TrialFunction(V)
18+
v = TestFunction(V)
19+
u0, u1 = split(u)
20+
v0, v1 = split(v)
21+
dx0 = Measure(
22+
"dx", domain=mesh0,
23+
)
24+
dx1 = Measure(
25+
"dx", domain=mesh1,
26+
)
27+
a = (
28+
inner(u0, v0) * dx0
29+
+ inner(u1, v1) * dx1
30+
)
31+
L = (
32+
inner(Constant(3.), v0) * dx0
33+
+ inner(Constant(4.), v1) * dx1
34+
)
35+
solution = Function(V)
36+
solve(a == L, solution)
37+
assert np.allclose(solution.subfunctions[0].dat.data_ro, 3.)
38+
assert np.allclose(solution.subfunctions[1].dat.data_ro, 4.)
39+
40+
1041
def _solve_helmholtz(mesh):
1142
V = FunctionSpace(mesh, "CG", 1)
1243
u = TrialFunction(V)

0 commit comments

Comments
 (0)