Skip to content

Commit ff69532

Browse files
committed
assemble: relax multi-domain restriction
1 parent 4dc066b commit ff69532

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
@@ -1090,17 +1090,6 @@ def local_kernels(self):
10901090
each possible combination.
10911091
10921092
"""
1093-
try:
1094-
topology, = set(d.topology.submesh_ancesters[-1] for d in self._form.ufl_domains())
1095-
except ValueError:
1096-
raise NotImplementedError("All integration domains must share a mesh topology")
1097-
1098-
for o in itertools.chain(self._form.arguments(), self._form.coefficients()):
1099-
domains = extract_domains(o)
1100-
for domain in domains:
1101-
if domain is not None and domain.topology.submesh_ancesters[-1] != topology:
1102-
raise NotImplementedError("Assembly with multiple meshes is not supported")
1103-
11041093
if isinstance(self._form, ufl.Form):
11051094
kernels = tsfc_interface.compile_form(
11061095
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
@@ -10,6 +10,37 @@
1010
cwd = abspath(dirname(__file__))
1111

1212

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

0 commit comments

Comments
 (0)