Skip to content

Commit

Permalink
Fix tests
Browse files Browse the repository at this point in the history
  • Loading branch information
pbrubeck committed Nov 14, 2024
2 parents b76ee5e + bb5a9bd commit 0e4c04f
Show file tree
Hide file tree
Showing 4 changed files with 30 additions and 19 deletions.
5 changes: 3 additions & 2 deletions firedrake/interpolation.py
Original file line number Diff line number Diff line change
Expand Up @@ -546,7 +546,8 @@ def __init__(
# For a VectorElement or TensorElement the correct
# VectorFunctionSpace equivalent is built from the scalar
# sub-element.
if V_dest.sub(0).value_shape != ():
ufl_scalar_element = ufl_scalar_element.sub_elements[0]
if ufl_scalar_element.reference_value_shape != ():
raise NotImplementedError(
"Can't yet cross-mesh interpolate onto function spaces made from VectorElements or TensorElements made from sub elements with value shape other than ()."
)
Expand Down Expand Up @@ -613,7 +614,7 @@ def __init__(
# I first point evaluate my expression at these locations, giving a
# P0DG function on the VOM. As described in the manual, this is an
# interpolation operation.
shape = V_dest.value_shape
shape = V_dest.ufl_function_space().value_shape
if len(shape) == 0:
fs_type = firedrake.FunctionSpace
elif len(shape) == 1:
Expand Down
40 changes: 25 additions & 15 deletions firedrake/pointquery_utils.py
Original file line number Diff line number Diff line change
Expand Up @@ -34,7 +34,7 @@ def make_wrapper(function, **kwargs):

def src_locate_cell(mesh, tolerance=None):
src = ['#include <evaluate.h>']
src.append(compile_coordinate_element(mesh.ufl_coordinate_element(), tolerance))
src.append(compile_coordinate_element(mesh, tolerance))
src.append(make_wrapper(mesh.coordinates,
forward_args=["void*", "double*", RealType_c+"*"],
kernel_name="to_reference_coords_kernel",
Expand Down Expand Up @@ -129,9 +129,10 @@ def to_reference_coords_newton_step(ufl_coordinate_element, parameters, x0_dtype
# Set up UFL form
cell = ufl_coordinate_element.cell
domain = ufl.Mesh(ufl_coordinate_element)
gdim = domain.geometric_dimension()
K = ufl.JacobianInverse(domain)
x = ufl.SpatialCoordinate(domain)
x0_element = finat.ufl.VectorElement("Real", cell, 0)
x0_element = finat.ufl.VectorElement("Real", cell, 0, dim=gdim)
x0 = ufl.Coefficient(ufl.FunctionSpace(domain, x0_element))
expr = ufl.dot(K, x - x0)

Expand Down Expand Up @@ -192,37 +193,46 @@ def predicate(index):


@PETSc.Log.EventDecorator()
def compile_coordinate_element(ufl_coordinate_element, contains_eps, parameters=None):
def compile_coordinate_element(mesh, contains_eps, parameters=None):
"""Generates C code for changing to reference coordinates.
:arg ufl_coordinate_element: UFL element of the coordinates
:returns: C code as string
Parameters
----------
mesh : MeshTopology
The mesh.
contains_eps : float
The tolerance used to verify that a point is contained by a cell.
Returns
-------
str
A string of C code.
"""
if parameters is None:
parameters = tsfc.default_parameters()
else:
_ = tsfc.default_parameters()
_.update(parameters)
parameters = _

ufl_coordinate_element = mesh.ufl_coordinate_element()
# Create FInAT element
element = tsfc.finatinterface.create_element(ufl_coordinate_element)
gdim, = ufl_coordinate_element.reference_value_shape
cell = ufl_coordinate_element.cell
extruded = isinstance(cell, ufl.TensorProductCell)

code = {
"geometric_dimension": gdim,
"topological_dimension": cell.topological_dimension(),
"geometric_dimension": mesh.geometric_dimension(),
"topological_dimension": mesh.topological_dimension(),
"celldist_l1_c_expr": celldist_l1_c_expr(element.cell, "X"),
"to_reference_coords_newton_step": to_reference_coords_newton_step(ufl_coordinate_element, parameters),
"init_X": init_X(element.cell, parameters),
"max_iteration_count": 1 if is_affine(ufl_coordinate_element) else 16,
"convergence_epsilon": 1e-12,
"dX_norm_square": dX_norm_square(cell.topological_dimension()),
"X_isub_dX": X_isub_dX(cell.topological_dimension()),
"extruded_arg": ", int const *__restrict__ layers" if extruded else "",
"extr_comment_out": "//" if extruded else "",
"non_extr_comment_out": "//" if not extruded else "",
"dX_norm_square": dX_norm_square(mesh.topological_dimension()),
"X_isub_dX": X_isub_dX(mesh.topological_dimension()),
"extruded_arg": ", int const *__restrict__ layers" if mesh.extruded else "",
"extr_comment_out": "//" if mesh.extruded else "",
"non_extr_comment_out": "//" if not mesh.extruded else "",
"IntType": as_cstr(IntType),
"ScalarType": ScalarType_c,
"RealType": RealType_c,
Expand Down
2 changes: 1 addition & 1 deletion firedrake/slate/slate.py
Original file line number Diff line number Diff line change
Expand Up @@ -198,7 +198,7 @@ def shapes(self):
"""
shapes = OrderedDict()
for i, fs in enumerate(self.arg_function_spaces):
shapes[i] = tuple(int(V.finat_element.space_dimension() * V.value_size)
shapes[i] = tuple(int(V.finat_element.space_dimension() * V.block_size)
for V in fs)
return shapes

Expand Down
2 changes: 1 addition & 1 deletion tests/regression/test_fdm.py
Original file line number Diff line number Diff line change
Expand Up @@ -228,7 +228,7 @@ def test_ipdg_direct_solver(fs):
mesh = fs.mesh()
x = SpatialCoordinate(mesh)
gdim = mesh.geometric_dimension()
ncomp = fs.ufl_element().value_size
ncomp = fs.value_size

homogenize = gdim > 2
if homogenize:
Expand Down

0 comments on commit 0e4c04f

Please sign in to comment.