Skip to content

Commit

Permalink
Fix BaseForm.__add__ simplification of Zero
Browse files Browse the repository at this point in the history
  • Loading branch information
nbouziani committed Mar 26, 2024
1 parent fbd288e commit 68f1c26
Show file tree
Hide file tree
Showing 2 changed files with 9 additions and 1 deletion.
8 changes: 8 additions & 0 deletions test/test_duals.py
Original file line number Diff line number Diff line change
Expand Up @@ -101,6 +101,9 @@ def test_addition():
V = FunctionSpace(domain_2d, f_2d)
V_dual = V.dual()

fvector_2d = FiniteElement("Lagrange", triangle, 1, (2, ), identity_pullback, H1)
W = FunctionSpace(domain_2d, fvector_2d)

u = TrialFunction(V)
v = TestFunction(V)

Expand Down Expand Up @@ -132,6 +135,11 @@ def test_addition():
res -= ZeroBaseForm((v,))
assert res == L

# Simplification with respect to ufl.Zero
a_W = Matrix(W, W)
res = a_W + Zero(W.ufl_element().value_shape)
assert res == a_W


def test_scalar_mult():
domain_2d = Mesh(FiniteElement("Lagrange", triangle, 1, (2, ), identity_pullback, H1))
Expand Down
2 changes: 1 addition & 1 deletion ufl/form.py
Original file line number Diff line number Diff line change
Expand Up @@ -137,7 +137,7 @@ def __add__(self, other):
if isinstance(other, (int, float)) and other == 0:
# Allow adding 0 or 0.0 as a no-op, needed for sum([a,b])
return self
elif isinstance(other, Zero) and not (other.ufl_shape or other.ufl_free_indices):
elif isinstance(other, Zero):
# Allow adding ufl Zero as a no-op, needed for sum([a,b])
return self

Expand Down

0 comments on commit 68f1c26

Please sign in to comment.