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

FEniCS backend: assignment and interpolation updates #613

Merged
merged 5 commits into from
Jan 7, 2025
Merged
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
30 changes: 21 additions & 9 deletions tests/fenics/test_equations.py
Original file line number Diff line number Diff line change
Expand Up @@ -490,22 +490,34 @@ def forward_J(x):

dJ = compute_gradient(J, x)

min_order = taylor_test(forward_J, x, J_val=J_val, dJ=dJ, seed=1.0e-3)
assert min_order > 1.99
dMs = tuple(var_new(x) for _ in range(2))
for i, dm in enumerate(dMs):
if complex_mode:
interpolate_expression(
dm, cos((i + 1) * pi * X[0]) * exp(1.0j * X[0] + X[1]))
else:
interpolate_expression(
dm, cos((i + 1) * pi * X[0]) * exp(X[1]))

min_order = taylor_test(forward_J, x, J_val=J_val, dJ=dJ, seed=1.0e-2,
dM=dMs[0])
assert min_order > 2.00

ddJ = Hessian(forward_J)
min_order = taylor_test(forward_J, x, J_val=J_val, ddJ=ddJ, seed=1.0e-3)
assert min_order > 2.98
min_order = taylor_test(forward_J, x, J_val=J_val, ddJ=ddJ, seed=1.0e-2,
dM=dMs[0])
assert min_order > 3.00

min_order = taylor_test_tlm(forward_J, x, tlm_order=1, seed=1.0e-3)
assert min_order > 1.99
min_order = taylor_test_tlm(forward_J, x, tlm_order=1, seed=1.0e-3,
dMs=dMs[:1])
assert min_order > 2.00

min_order = taylor_test_tlm_adjoint(forward_J, x, adjoint_order=1,
seed=1.0e-3)
assert min_order > 1.99
seed=1.0e-3, dMs=dMs[:1])
assert min_order > 2.00

min_order = taylor_test_tlm_adjoint(forward_J, x, adjoint_order=2,
seed=1.0e-3)
seed=1.0e-3, dMs=dMs)
assert min_order > 1.99

h.close()
Expand Down
28 changes: 8 additions & 20 deletions tests/fenics/test_patches.py
Original file line number Diff line number Diff line change
Expand Up @@ -282,28 +282,16 @@ def forward(m):
u.assign(-2.0)
u.assign(u + 2.0 * m)

m_ = Function(space, name="m")
assign_fn(m_, m)
m = m_
del m_

u_ = Function(space, name="u")
assign_fn(u_, u)
assign_fn(u_, m)
u = u_
del u_

one = Function(space, name="one")
assign_fn(one, Constant(1.0))

v = Function(space, name="v")
assign_fn(v, u)
v.assign(u + one)
assign_fn(v, Constant(0.0))
v.assign(u + v + one)
v.assign(2.5 * u + 3.6 * v + 4.7 * m)

J = Functional(name="J")
J.assign(((v - 1.0) ** 4) * dx)
J.assign(((v - 0.5) ** 4) * dx)
return J

m = Constant(2.0, name="m")
Expand All @@ -313,29 +301,29 @@ def forward(m):
stop_manager()

J_val = J.value
assert abs(J_val - 342974.2096) < 1.0e-9
assert abs(J_val - 1.5 ** 4) < 1.0e-14

dJ = compute_gradient(J, m)

dm = Constant(1.0)

min_order = taylor_test(forward, m, J_val=J_val, dJ=dJ, dM=dm)
assert min_order > 1.99
assert min_order > 2.00

ddJ = Hessian(forward)
min_order = taylor_test(forward, m, J_val=J_val, ddJ=ddJ, dM=dm)
assert min_order > 2.99
assert min_order > 3.00

min_order = taylor_test_tlm(forward, m, tlm_order=1, dMs=(dm,))
assert min_order > 1.99
assert min_order > 2.00

min_order = taylor_test_tlm_adjoint(forward, m, adjoint_order=1,
dMs=(dm,))
assert min_order > 1.99
assert min_order > 2.00

min_order = taylor_test_tlm_adjoint(forward, m, adjoint_order=2,
dMs=(dm, dm))
assert min_order > 1.99
assert min_order > 2.00


@pytest.mark.fenics
Expand Down
34 changes: 23 additions & 11 deletions tests/firedrake/test_equations.py
Original file line number Diff line number Diff line change
Expand Up @@ -813,23 +813,35 @@ def forward_J(x):

dJ = compute_gradient(J, x)

min_order = taylor_test(forward_J, x, J_val=J_val, dJ=dJ, seed=1.0e-3)
assert min_order > 1.99
dMs = tuple(var_new(x) for _ in range(2))
for i, dm in enumerate(dMs):
if complex_mode:
interpolate_expression(
dm, cos((i + 1) * pi * X[0]) * exp(1.0j * X[0] + X[1]))
else:
interpolate_expression(
dm, cos((i + 1) * pi * X[0]) * exp(X[1]))

min_order = taylor_test(forward_J, x, J_val=J_val, dJ=dJ, seed=1.0e-2,
dM=dMs[0])
assert min_order > 2.00

ddJ = Hessian(forward_J)
min_order = taylor_test(forward_J, x, J_val=J_val, ddJ=ddJ, seed=1.0e-3)
assert min_order > 2.98
min_order = taylor_test(forward_J, x, J_val=J_val, ddJ=ddJ, seed=1.0e-2,
dM=dMs[0])
assert min_order > 2.99

min_order = taylor_test_tlm(forward_J, x, tlm_order=1, seed=1.0e-3)
assert min_order > 1.99
min_order = taylor_test_tlm(forward_J, x, tlm_order=1, seed=1.0e-3,
dMs=dMs[:1])
assert min_order > 2.00

min_order = taylor_test_tlm_adjoint(forward_J, x, adjoint_order=1,
seed=1.0e-3)
assert min_order > 1.99
seed=1.0e-3, dMs=dMs[:1])
assert min_order > 2.00

min_order = taylor_test_tlm_adjoint(forward_J, x, adjoint_order=2,
seed=1.0e-3)
assert min_order > 1.98
seed=1.0e-3, dMs=dMs)
assert min_order > 2.00

h.close()

Expand Down Expand Up @@ -1153,7 +1165,7 @@ def forward(m):
ddJ = Hessian(forward)
min_order = taylor_test(forward, m, J_val=J_val, ddJ=ddJ, seed=1.0e-3,
dM=dm)
assert min_order > 2.97
assert min_order > 2.96

min_order = taylor_test_tlm(forward, m, tlm_order=1, seed=1.0e-4,
dMs=(dm,))
Expand Down
4 changes: 2 additions & 2 deletions tests/firedrake/test_patches.py
Original file line number Diff line number Diff line change
Expand Up @@ -889,8 +889,8 @@ def forward(b):
u_v.setValue(n, 1)
u_v.assemblyBegin()
u_v.assemblyEnd()
error = abs(assemble(y(u) - b(Function(space_1).interpolate(u))))
assert error == 0
error_norm = abs(assemble(y(u) - b(Function(space_1).interpolate(u))))
assert error_norm == 0

J_val = J.value

Expand Down
29 changes: 12 additions & 17 deletions tlm_adjoint/fenics/backend_patches.py
Original file line number Diff line number Diff line change
Expand Up @@ -467,32 +467,27 @@ def Function_assign(self, orig, orig_args, rhs):
annotate = annotation_enabled()
tlm = tlm_enabled()
if isinstance(rhs, backend_Function):
# Prevent a new vector being created
if not space_eq(rhs.function_space(), self.function_space()):
raise ValueError("Invalid assignment")

if space_eq(rhs.function_space(), self.function_space()):
if rhs is not self:
var_assign(self, rhs)

if annotate or tlm:
eq = Assignment(self, rhs)
assert not eq._pre_process_required
eq._post_process()
else:
value = var_new(self)
orig(value, rhs)
var_assign(self, value)
if rhs is not self:
# Prevent a new vector being created
var_assign(self, rhs)

if annotate or tlm:
eq = ExprInterpolation(self, rhs)
eq = Assignment(self, rhs)
assert not eq._pre_process_required
eq._post_process()
else:
orig_args()

if annotate or tlm:
eq = ExprInterpolation(self, expr_new_x(rhs, self))
assert not eq._pre_process_required
eq._post_process()
if isinstance(rhs, backend_Constant):
eq = ExprInterpolation(self, rhs)
assert not eq._pre_process_required
eq._post_process()
else:
raise NotImplementedError("Case not implemented")


@manager_method(backend_Function, "copy", patch_without_manager=True)
Expand Down
5 changes: 1 addition & 4 deletions tlm_adjoint/fenics/interpolation.py
Original file line number Diff line number Diff line change
Expand Up @@ -56,10 +56,7 @@ def value_shape(self):
value, = value
var_assign(x, value)
elif isinstance(x, backend_Function):
try:
x.assign(expr)
except RuntimeError:
x.interpolate(Expr())
x.interpolate(Expr())
else:
raise TypeError(f"Unexpected type: {type(x)}")
else:
Expand Down
Loading