Skip to content
Open
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
12 changes: 8 additions & 4 deletions firedrake/formmanipulation.py
Original file line number Diff line number Diff line change
Expand Up @@ -3,9 +3,9 @@
import collections

from ufl import as_tensor, as_vector, split
from ufl.classes import Zero, FixedIndex, ListTensor, ZeroBaseForm
from ufl.classes import FixedIndex, Form, ListTensor, Zero, ZeroBaseForm
from ufl.algorithms.map_integrands import map_integrand_dags
from ufl.algorithms import expand_derivatives
from ufl.algorithms import expand_derivatives, expand_indices
from ufl.corealg.map_dag import MultiFunction, map_expr_dags

from pyop2 import MixedDat
Expand Down Expand Up @@ -78,9 +78,13 @@ def split(self, form, argument_indices):
assert (len(idx) == 1 for idx in self.blocks.values())
assert (idx[0] == 0 for idx in self.blocks.values())
return form
# TODO find a way to distinguish empty Forms avoiding expand_derivatives
f = map_integrand_dags(self, form)
if expand_derivatives(f).empty():

# TODO find a better way to distinguish empty Forms
f_expanded = expand_derivatives(f)
if isinstance(f_expanded, Form):
f_expanded = expand_indices(f_expanded)
if f_expanded.empty():
# Get ZeroBaseForm with the right shape
f = ZeroBaseForm(tuple(map(self._subspace_argument, form.arguments())))
return f
Expand Down
12 changes: 12 additions & 0 deletions tests/firedrake/regression/test_split.py
Original file line number Diff line number Diff line change
Expand Up @@ -125,3 +125,15 @@ def test_split_coefficient_not_argument():
as_vector([TestFunction(V), 0])),
w, wr)
assert J00.signature() == expect.signature()


def test_split_zero_block():
mesh = UnitSquareMesh(1, 1)
V = FunctionSpace(mesh, "DG", 0)
Z = V * V * V * V
J = inner(TrialFunction(Z), TestFunction(Z))*dx
splitter = ExtractSubBlock()

J00 = splitter.split(J, (0, 1))

assert isinstance(J00, ZeroBaseForm)
Loading