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

Expunge Expr.ufl_domain Part II #3259

Merged
merged 12 commits into from
Nov 29, 2023
6 changes: 3 additions & 3 deletions firedrake/function.py
Original file line number Diff line number Diff line change
Expand Up @@ -616,16 +616,16 @@ def at(self, arg, *args, **kwargs):
dont_raise = kwargs.get('dont_raise', False)

tolerance = kwargs.get('tolerance', None)
mesh = self.function_space().mesh()
if tolerance is None:
tolerance = self.ufl_domain().tolerance
tolerance = mesh.tolerance
else:
self.ufl_domain().tolerance = tolerance
mesh.tolerance = tolerance

# Handle f.at(0.3)
if not arg.shape:
arg = arg.reshape(-1)

mesh = self.function_space().mesh()
if mesh.variable_layers:
raise NotImplementedError("Point evaluation not implemented for variable layers")

Expand Down
8 changes: 4 additions & 4 deletions firedrake/interpolation.py
Original file line number Diff line number Diff line change
Expand Up @@ -9,7 +9,7 @@
import finat.ufl
from ufl.algorithms import extract_arguments, extract_coefficients
from ufl.algorithms.signature import compute_expression_signature
from ufl.domain import extract_unique_domain
from ufl.domain import as_domain, extract_unique_domain

from pyop2 import op2
from pyop2.caching import disk_cached
Expand Down Expand Up @@ -191,7 +191,7 @@ class Interpolator(abc.ABC):
"""

def __new__(cls, expr, V, **kwargs):
target_mesh = V.ufl_domain()
target_mesh = as_domain(V)
source_mesh = extract_unique_domain(expr) or target_mesh
if target_mesh is not source_mesh:
if isinstance(target_mesh.topology, firedrake.mesh.VertexOnlyMeshTopology):
Expand Down Expand Up @@ -323,7 +323,7 @@ def __init__(
# setup
V_dest = V
src_mesh = extract_unique_domain(expr)
dest_mesh = V_dest.ufl_domain()
dest_mesh = as_domain(V_dest)
src_mesh_gdim = src_mesh.geometric_dimension()
dest_mesh_gdim = dest_mesh.geometric_dimension()
if src_mesh_gdim != dest_mesh_gdim:
Expand Down Expand Up @@ -702,7 +702,7 @@ def make_interpolator(expr, V, subset, access, bcs=None):
assert isinstance(expr, ufl.classes.Expr)

arguments = extract_arguments(expr)
target_mesh = V.ufl_domain()
target_mesh = as_domain(V)
if len(arguments) == 0:
source_mesh = extract_unique_domain(expr) or target_mesh
vom_onto_other_vom = (
Expand Down
3 changes: 2 additions & 1 deletion firedrake/ufl_expr.py
Original file line number Diff line number Diff line change
Expand Up @@ -3,6 +3,7 @@
from ufl.duals import is_dual
from ufl.split_functions import split
from ufl.algorithms import extract_arguments, extract_coefficients
from ufl.domain import as_domain

import firedrake
from firedrake import utils, function, cofunction
Expand Down Expand Up @@ -230,7 +231,7 @@ def derivative(form, u, du=None, coefficient_derivatives=None):
raise ValueError("Taking derivative of form wrt u, but form contains coefficients from u.subfunctions."
"\nYou probably meant to write split(u) when defining your form.")

mesh = form.ufl_domain()
mesh = as_domain(form)
if not mesh:
raise ValueError("Expression to be differentiated has no ufl domain."
"\nDo you need to add a domain to your Constant?")
Expand Down
2 changes: 1 addition & 1 deletion tests/regression/test_subdomain_integrals.py
Original file line number Diff line number Diff line change
Expand Up @@ -44,7 +44,7 @@ def test_solve_cell_subdomains(form, u):
solve(form == 0, u)
expect = Function(u.function_space())

mesh = u.ufl_domain()
mesh = u.function_space().mesh()
expect.interpolate(Constant(1.0), subset=mesh.cell_subset(1))
expect.interpolate(Constant(0.5), subset=mesh.cell_subset(2))

Expand Down
2 changes: 1 addition & 1 deletion tests/regression/test_vfs_component_bcs.py
Original file line number Diff line number Diff line change
Expand Up @@ -144,7 +144,7 @@ def test_cant_integrate_subscripted_VFS(V):
f = Function(V)
f.assign(Constant([2, 1]))
assert np.allclose(assemble(f.sub(0)*dx),
assemble(Constant(2)*dx(domain=f.ufl_domain())))
assemble(Constant(2)*dx(domain=V.mesh())))


@pytest.mark.parametrize("cmpt",
Expand Down