Open JHopeCollins opened 2 years ago
I have recreated the issue. Based on what @tommbendall said about last week's UFL update being the cause, I reset UFL back to 1dddf46eb5d011e7144d8780ac008dc2a87a2ce6 and ran the script without trouble. So it seems that either something that got added to FEniCS UFL since then, or indeed @rckirby 's accidental PR merge to Firedrake UFL is to blame. I think we ought to be able to find the culprit by looking at the diff for https://github.com/firedrakeproject/ufl/pull/30?
This is the culprit:
Previously the default was that a geometric quantity was considered to be is_cellwise_constant() => True
(aside: this is an unsafe default to apply to an abstract base class, but meh).
That merge brings in a fix for CellNormal
(which is only cellwise constant if the element is a linear simplex). But that means that now we have new things appearing in preprocessed forms. So one needs to figure out what the correct derivative lowering is in this case. It would be good to simplify the form. I think something like this should provoke the issue?
from tsfc.ufl_utils import compute_form_data
...
n = CellNormalI(mesh)
u = Function(V)
phi = grad(inner(u, n))[0]*dx
fd = compute_form_data(derivative(phi, u))
?
See https://github.com/FEniCS/ufl/pull/119. There was a bug in UFL which incorrectly simplified all CellNormal derivatives to 0.
~I guess TSFC knows how to compute derivatives of Jacobian,~ (not needed, below is anyway lowered to reference_grad()) The following premature CellNormal lowering could help you?
import ufl
mesh = ufl.Mesh(ufl.VectorElement("CG", ufl.triangle, 2, dim=3))
cn = ufl.algorithms.apply_geometry_lowering.apply_geometry_lowering(ufl.CellNormal(mesh))
gcn = ufl.grad(cn)
@wence-, thanks for finding that. Your simpler form does result in the same error too. Full script:
import firedrake as fd
from tsfc.ufl_utils import compute_form_data
coords_degree = 2 # works when coords_degree=1, breaks when coords_degree=2
mesh = fd.IcosahedralSphereMesh(radius=1,
refinement_level=2,
degree=coords_degree)
x = fd.SpatialCoordinate(mesh)
mesh.init_cell_orientations(x)
n = fd.CellNormal(mesh)
V = fd.FunctionSpace(mesh, "BDM", 1)
u = fd.Function(V)
phi = fd.grad(fd.inner(u, n))[0]*fd.dx
form_data = compute_form_data(fd.derivative(phi, u))
@michalhabera, thanks for the workaround, I'll see if we can use it for our form for the time being.
(That workaround appears not to help :()
Here's a minimal example that only uses UFL to demonstrate the issue:
from ufl import (Cell, CellNormal, Coefficient, FunctionSpace, Mesh,
VectorElement, derivative, dx, grad, inner)
from ufl.algorithms import compute_form_data
from ufl.algorithms.apply_geometry_lowering import apply_geometry_lowering
mesh = Mesh(VectorElement("P", Cell("triangle", 3), 2))
try_workaround = False
if try_workaround:
n = apply_geometry_lowering(CellNormal(mesh))
else:
n = CellNormal(mesh)
V = FunctionSpace(mesh, VectorElement("P", mesh.ufl_cell(), 2))
u = Coefficient(V)
phi = grad(inner(u, n))[0]*dx
fd = compute_form_data(
derivative(phi, u),
do_apply_function_pullbacks=True,
do_apply_geometry_lowering=True,
)
(That workaround appears not to help :()
Here's a minimal example that only uses UFL to demonstrate the issue:
from ufl import (Cell, CellNormal, Coefficient, FunctionSpace, Mesh, VectorElement, derivative, dx, grad, inner) from ufl.algorithms import compute_form_data from ufl.algorithms.apply_geometry_lowering import apply_geometry_lowering mesh = Mesh(VectorElement("P", Cell("triangle", 3), 2)) try_workaround = False if try_workaround: n = apply_geometry_lowering(CellNormal(mesh)) else: n = CellNormal(mesh) V = FunctionSpace(mesh, VectorElement("P", mesh.ufl_cell(), 2)) u = Coefficient(V) phi = grad(inner(u, n))[0]*dx fd = compute_form_data( derivative(phi, u), do_apply_function_pullbacks=True, do_apply_geometry_lowering=True, )
Indeed. I now remember in all manifold codes we first interpolate the normal field into a conforming VectorFunction space, this avoids the CoefficientDerivative
of ReferenceGrad
.
Alternatively, and dangerously, you could make coefficient derivatives of ReferenceGrad = 0 by monkey-patching UFL,
from ufl.algorithms.apply_derivatives import GateauxDerivativeRuleset, GenericDerivativeRuleset
def myrule(self, o):
print(f"D({o})/D({self._w}) simplified to 0")
return GenericDerivativeRuleset.independent_terminal(self, o)
GateauxDerivativeRuleset.reference_grad = myrule
The comments in apply_derivatives.py
appear to indicate that CoefficientDerivative
of ReferenceGrad
could be implemented for the simple case that causes this error, and could be allowed to still fail in more complex cases. Doing the same thing for ReferenceValue
is probably advisable.
Could I ask for an update on this issue please - I'm afraid I got rather confused trying to follow the discussions referenced above. I'm really hoping to be able to run the vector invariant form of the shallow water equations on the sphere again!
Correct fix: handle derivatives of referencegrad correctly in UFL.
Workaround:
# a quadratic approximation of the cell normal (or whatever you need)
# if you want the previous affine approximation of cell normals you could use VectorFunctionSpace(mesh, "DG", 0) instead.
cn = interpolate(CellNormal(mesh), VectorFunctionSpace(mesh, "DG", 2))
# use cn as normal (hah!)
Note that the same issue is seen for a facet normal of a higher-order cell, thus not depending on manifolds, see https://github.com/firedrakeproject/firedrake/issues/3313.
When trying to create a
LinearVariationalProblem
from the derivative of a form, UFL throws an error if:MFE:
Error message:
This seems to have only started happening since the UFL update, although I'm not sure of exactly when.