Point evaluation

Firedrake can evaluate Functions at arbitrary physical points. This feature can be useful for the evaluation of the result of a simulation. Two APIs are offered to this feature: a Firedrake-specific one, and one from UFL.

Firedrake API

Firedrake offers a convenient API for evaluating functions at arbitrary points via at():

# evaluate f at a 1-dimensional point
f.at(0.3)

# evaluate f at two 1-dimensional points, or at one 2-dimensional point
# (depending on f's geometric dimension)
f.at(0.2, 0.4)

# evaluate f at one 2-dimensional point
f.at([0.2, 0.4])

# evaluate f at two 2-dimensional point
f.at([0.2, 0.4], [1.2, 0.5])

# evaluate f at two 2-dimensional point (same as above)
f.at([[0.2, 0.4], [1.2, 0.5]])

While in these examples we have only shown lists, other iterables such as tuples and numpy arrays are also accepted. The following are equivalent:

f.at(0.2, 0.4)
f.at((0.2, 0.4))
f.at([0.2, 0.4])
f.at(numpy.array([0.2, 0.4]))

For a single point, the result is a numpy array, or a tuple of numpy arrays in case of mixed functions. When evaluating multiple points, the result is a list of values for each point. To summarise:

  • Single point, non-mixed: numpy array

  • Single point, mixed: tuple of numpy arrays

  • Multiple points, non-mixed: list of numpy arrays

  • Multiple points, mixed: list of tuples of numpy arrays

Points outside the domain

When any point is outside the domain of the function, PointNotInDomainError exception is raised. If dont_raise=True is passed to at(), the result is None for those points which fall outside the domain.

mesh = UnitIntervalMesh(8)
f = mesh.coordinates

f.at(1.2)                   # raises exception
f.at(1.2, dont_raise=True)  # returns None

f.at(0.5, 1.2)                   # raises exception
f.at(0.5, 1.2, dont_raise=True)  # returns [0.5, None]

Warning

Point evaluation on immersed manifolds is not supported yet, due to the difficulty of specifying a physical point on the manifold.

Evaluation on a moving mesh

If you move the mesh, by changing the mesh coordinates, then the bounding box tree that Firedrake maintains to ensure fast point evaluation must be rebuilt. To do this, after moving the mesh, call clear_spatial_index() on the mesh you have just moved.

Evaluation with a distributed mesh

There is limited support for point evaluation when running Firedrake in parallel. There is no special API, but there are some restrictions:

  • Point evaluation is a collective operation.

  • Each process must ask for the same list of points.

  • Each process will get the same values.

UFL API

UFL reserves the function call operator for evaluation:

f([0.2, 0.4])

will evaluate \(f\) at \((0.2, 0.4)\). UFL does not accept multiple points at once, and cannot configure what to do with a point which is not in the domain. The advantage of this syntax is that it works on any Expr, for example:

(f*sin(f)([0.2, 0.4])

will evaluate \(f \cdot \sin(f)\) at \((0.2, 0.4)\).

Note

The expression itself is not translated into C code. While the evaluation of a function uses the same infrastructure as the Firedrake API, which uses generated C code, the expression tree is evaluated by UFL in Python.