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
from firedrake import *
mesh = UnitSquareMesh(10, 10)
V = FunctionSpace(mesh, "CG", 1)
x, y = SpatialCoordinate(mesh)
u = Function(V)
u.interpolate(x)
print(u.at([0.5, 0.5]))
mesh.coordinates.dat.data[:] **= 2
print(u.at((0.5, 0.5)))
The second point interpolation at (0.5, 0.5) fails with a PointNotInDomainError. This is because it uses a spatialindex that got created and cached on the first point interpolation call and thus returns candidate elements in which the point would have originally been located, but these elements have since moved elsewhere. Firedrake then tries these elements, and finds none of them contain the queried point. This can be fixed by call mesh.clear_spatialindex() after each move.
Note that this specifically only addresses interpolating in individual points via the at() (or __call__()) method. I believe function-to-function interpolate now uses VertexOnlyMesh. I don't currently have any evidence that is not doing the right thing, but we should be really careful as AFAIK this stuff isn't really supported (if the warning here https://www.firedrakeproject.org/mesh-coordinates.html is still valid) but we should add some tests to verify it is working as we expect, and to ensure it isn't broken in the future. The same goes for when we use cross-mesh project() as this probably involves a spatial index that is stored in yet another place (super mesh).
The text was updated successfully, but these errors were encountered:
Call clear_spatial_index() on the moved mesh when updating its
coordinates in MA mover. Otherwise if a point was queried (e.g. with
Function.at()) before the update, the same spatial index will still be
used when the query is repeated after the update but at that time the
point may be located in a completely different cell resulting in
PointNotInDomainErrors.
Consider the following code:
The second point interpolation at (0.5, 0.5) fails with a PointNotInDomainError. This is because it uses a spatialindex that got created and cached on the first point interpolation call and thus returns candidate elements in which the point would have originally been located, but these elements have since moved elsewhere. Firedrake then tries these elements, and finds none of them contain the queried point. This can be fixed by call
mesh.clear_spatialindex()
after each move.Note that this specifically only addresses interpolating in individual points via the
at()
(or__call__()
) method. I believe function-to-function interpolate now uses VertexOnlyMesh. I don't currently have any evidence that is not doing the right thing, but we should be really careful as AFAIK this stuff isn't really supported (if the warning here https://www.firedrakeproject.org/mesh-coordinates.html is still valid) but we should add some tests to verify it is working as we expect, and to ensure it isn't broken in the future. The same goes for when we use cross-mesh project() as this probably involves a spatial index that is stored in yet another place (super mesh).The text was updated successfully, but these errors were encountered: