DirichletBC initialization with ufl Expr: possible confusing/inconsistent behavior

Issue #1010 new
Francesco Ballarin created an issue

Hi, this is a very minor thing but I am opening the ticket anyway to keep track of it.

If a DOLFIN Expression is provided to a DirichletBC, it is allowed to update its user parameters and any calls to DirichletBC will reflect that. If an UFL Expr (which contains a DOLFIN Expression) is provided to a DirichletBC it will first get projected (or, as the FIXME in the relevant python file mentions, interpolated), meaning that an eventual update to the Expression user parameters will not be reflected in the DirichletBC object.

from numpy import isclose
from dolfin import *

mesh = UnitSquareMesh(3, 3)
boundary = MeshFunction("size_t", mesh, 1, 0)
boundary[0] = 1
V = FunctionSpace(mesh, "CG", 1)

g = Expression("t", t=1., element=V.ufl_element())

u1 = Function(V)
u2 = Function(V)

bc1 = DirichletBC(V, g, boundary, 1)
bc2 = DirichletBC(V, 2*g, boundary, 1)

bc1.apply(u1.vector())
bc2.apply(u2.vector())
print(u1.vector().get_local(), u2.vector().get_local())
assert isclose(2*u1.vector().get_local(), u2.vector().get_local()).all()

g.t = 5.
bc1.apply(u1.vector())
bc2.apply(u2.vector())
print(u1.vector().get_local(), u2.vector().get_local()) # note that the second array does not contain the updated value (equal to 10)
assert isclose(2*u1.vector().get_local(), u2.vector().get_local()).all() # FAILURE

Of course it very easy to work around this issue, but the inconsistent behavior may be confusing and easily overlooked.

Comments (1)

  1. Log in to comment