You signed in with another tab or window. Reload to refresh your session.You signed out in another tab or window. Reload to refresh your session.You switched accounts on another tab or window. Reload to refresh your session.Dismiss alert
Describe the bug
The derivative of an Action may not be aBaseForm subclass.
Steps to Reproduce
from firedrake import *
import ufl
mesh = UnitIntervalMesh(2)
space = FunctionSpace(mesh, "Lagrange", 1)
test = TestFunction(space)
u = Function(space, name="u")
v = Cofunction(space.dual(), name="v")
w = Cofunction(space.dual(), name="w")
F = v(u)
dF = ufl.algorithms.expand_derivatives(
ufl.derivative(F, v))
print(f"{type(dF)=}") # Function, should be Action
F = w(u)
dF = ufl.algorithms.expand_derivatives(
ufl.derivative(F, v))
print(f"{type(dF)=}") # int, should be ZeroBaseForm
Expected behavior
The result of ufl.algorithms.expand_derivative(ufl.derivative(form, ...)) where form is a BaseForm should always be a BaseForm.
Error message
This leads to difficulties when constructing tangent-linear right-hand-sides, e.g.
from firedrake import *
import ufl
mesh = UnitIntervalMesh(2)
space = FunctionSpace(mesh, "Lagrange", 1)
test = TestFunction(space)
u = Function(space, name="u")
v = Cofunction(space.dual(), name="v")
form = v(u)
tlm_v = Cofunction(space.dual(), name="tlm_v")
ufl.action(ufl.algorithms.expand_derivatives(ufl.derivative(form, v)), tlm_v)
leads to the error
Traceback (most recent call last):
File "/home/maddison/test.py", line 15, in <module>
ufl.action(ufl.algorithms.expand_derivatives(ufl.derivative(form, v)), tlm_v)
File "/home/maddison/build/firedrake/firedrake/src/ufl/ufl/formoperators.py", line 109, in action
form = as_form(form)
File "/home/maddison/build/firedrake/firedrake/src/ufl/ufl/form.py", line 688, in as_form
raise ValueError(f"Unable to convert object to a UFL form: {ufl_err_str(form)}")
ValueError: Unable to convert object to a UFL form: <Coefficient id=127176171859712>
Additional Info
The natural solution would be for Functions to be BaseForms with one Coargument argument, but I can imagine that might cause a lot of problems.
The text was updated successfully, but these errors were encountered:
Describe the bug
The derivative of an
Action
may not be aBaseForm
subclass.Steps to Reproduce
Expected behavior
The result of
ufl.algorithms.expand_derivative(ufl.derivative(form, ...))
whereform
is aBaseForm
should always be aBaseForm
.Error message
This leads to difficulties when constructing tangent-linear right-hand-sides, e.g.
leads to the error
Additional Info
The natural solution would be for
Function
s to beBaseForm
s with oneCoargument
argument, but I can imagine that might cause a lot of problems.The text was updated successfully, but these errors were encountered: