Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Fix for Cofunction self-assignment via interpolation #3939

Open
wants to merge 7 commits into
base: master
Choose a base branch
from
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
7 changes: 6 additions & 1 deletion firedrake/interpolation.py
Original file line number Diff line number Diff line change
Expand Up @@ -859,7 +859,12 @@ def _interpolate(self, *function, output=None, transpose=False, **kwargs):
V = self.V
result = output or firedrake.Function(V)
with function.dat.vec_ro as x, result.dat.vec_wo as out:
mul(x, out)
if x is not out:
mul(x, out)
else:
out_ = out.duplicate()
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Shouldn't a no-op be the right thing?

Copy link
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I don't think so, e.g.

u = fd.Function(space, name="u")
v = fd.Function(space, name="v")
w = fd.Cofunction(space.dual(), name="w")
...
fd.assemble(fd.action(
    fd.adjoint(fd.derivative(interpolate(v * u, space), u)),
    w), tensor=w)

Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I see, it only is a no-op in the MFE from the issue: assemble(L(w), tensor=w) with L being the identity.

Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Maybe it'd be good to add a more complicated test then

Copy link
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Done. I think there are also be issues with complex here (transpose should be Hermitian for the adjoint actions?).

mul(x, out_)
out_.copy(result=out)
return result

else:
Expand Down
20 changes: 20 additions & 0 deletions tests/firedrake/regression/test_interp_dual.py
Original file line number Diff line number Diff line change
Expand Up @@ -34,6 +34,26 @@ def f1(mesh, V1):
return Function(V1).interpolate(expr)


def test_interp_self(V1):
a = assemble(conj(TestFunction(V1)) * dx)
b = assemble(conj(TestFunction(V1)) * dx)
a.interpolate(a)
assert np.allclose(a.dat.data_ro, b.dat.data_ro)


def test_assemble_interp_adjoint_tensor(mesh, V1, f1):
a = assemble(conj(TestFunction(V1)) * dx)
# We want tensor to be a dependency of the input expression for this test
assemble(action(adjoint(Interpolate(f1 * TestFunction(V1), V1)), a),
tensor=a)

x, y = SpatialCoordinate(mesh)
f2 = Function(V1, name="f2").interpolate(
exp(x) * y)

assert np.allclose(assemble(a(f2)), assemble(Function(V1).interpolate(conj(f1 * f2)) * dx))


def test_assemble_interp_operator(V2, f1):
# Check type
If1 = Interpolate(f1, V2)
Expand Down
Loading