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

Possible fix for taking scalar derivatives of vector coefficients using coefficient derivative maps #124

Merged
merged 6 commits into from
Sep 12, 2023
Merged
24 changes: 24 additions & 0 deletions test/test_derivative.py
Original file line number Diff line number Diff line change
Expand Up @@ -483,6 +483,30 @@ def test_coefficient_derivatives(self):
self.assertEqual(replace(actual, fd.function_replace_map), expected)


def test_vector_coefficient_scalar_derivatives(self):
V = FiniteElement("Lagrange", triangle, 1)
VV = VectorElement("Lagrange", triangle, 1)

dv = TestFunction(V)

df = Coefficient(VV, count=0)
g = Coefficient(VV, count=1)
f = Coefficient(VV, count=2)
u = Coefficient(V, count=3)
cd = {f: df}

integrand = inner(f, g)

i0, i1, i2, i3, i4 = [Index(count=c) for c in range(5)]
expected = as_tensor(df[i1]*dv, (i1,))[i0]*g[i0]

F = integrand*dx
J = derivative(F, u, dv, cd)
fd = compute_form_data(J)
actual = fd.preprocessed_form.integrals()[0].integrand()
assert (actual*dx).signature() == (expected*dx).signature()


def test_vector_coefficient_derivatives(self):
V = VectorElement("Lagrange", triangle, 1)
VV = TensorElement("Lagrange", triangle, 1)
Expand Down
6 changes: 3 additions & 3 deletions ufl/algorithms/apply_derivatives.py
Original file line number Diff line number Diff line change
Expand Up @@ -835,9 +835,9 @@ def coefficient(self, o):
dosum = Zero(o.ufl_shape)
for do, v in zip(dos, self._v):
so, oi = as_scalar(do)
rv = len(v.ufl_shape)
oi1 = oi[:-rv]
oi2 = oi[-rv:]
rv = len(oi) - len(v.ufl_shape)
oi1 = oi[:rv]
oi2 = oi[rv:]
prod = so * v[oi2]
if oi1:
dosum += as_tensor(prod, oi1)
Expand Down