Skip to content

Commit

Permalink
add form_compiler_parameters to Interpolator
Browse files Browse the repository at this point in the history
  • Loading branch information
Pablo Brubeck Martinez committed Jan 18, 2023
1 parent 1d333bd commit 7692be6
Show file tree
Hide file tree
Showing 2 changed files with 13 additions and 7 deletions.
15 changes: 9 additions & 6 deletions firedrake/interpolation.py
Original file line number Diff line number Diff line change
Expand Up @@ -82,9 +82,11 @@ class Interpolator(object):
:class:`Interpolator` is also collected).
"""
def __init__(self, expr, V, subset=None, freeze_expr=False, access=op2.WRITE, bcs=None):
def __init__(self, expr, V, subset=None, freeze_expr=False, access=op2.WRITE, bcs=None, form_compiler_parameters=None):
try:
self.callable, arguments = make_interpolator(expr, V, subset, access, bcs=bcs)
self.callable, arguments = make_interpolator(expr, V, subset, access,
bcs=bcs,
parameters=form_compiler_parameters)
except FIAT.hdiv_trace.TraceError:
raise NotImplementedError("Can't interpolate onto traces sorry")
self.arguments = arguments
Expand Down Expand Up @@ -155,7 +157,7 @@ def interpolate(self, *function, output=None, transpose=False):


@PETSc.Log.EventDecorator()
def make_interpolator(expr, V, subset, access, bcs=None):
def make_interpolator(expr, V, subset, access, bcs=None, parameters=None):
assert isinstance(expr, ufl.classes.Expr)

arguments = extract_arguments(expr)
Expand Down Expand Up @@ -216,7 +218,7 @@ def make_interpolator(expr, V, subset, access, bcs=None):
if len(V) > 1:
raise NotImplementedError(
"UFL expressions for mixed functions are not yet supported.")
loops.extend(_interpolator(V, tensor, expr, subset, arguments, access, bcs=bcs))
loops.extend(_interpolator(V, tensor, expr, subset, arguments, access, bcs=bcs, parameters=parameters))

if bcs and len(arguments) == 0:
loops.extend([partial(bc.apply, f) for bc in bcs])
Expand All @@ -230,7 +232,7 @@ def callable(loops, f):


@utils.known_pyop2_safe
def _interpolator(V, tensor, expr, subset, arguments, access, bcs=None):
def _interpolator(V, tensor, expr, subset, arguments, access, bcs=None, parameters=None):
try:
expr = ufl.as_ufl(expr)
except ufl.UFLException:
Expand Down Expand Up @@ -283,7 +285,8 @@ def _interpolator(V, tensor, expr, subset, arguments, access, bcs=None):
assert subset.superset == cell_set
cell_set = subset

parameters = {}
if parameters is None:
parameters = {}
parameters['scalar_type'] = utils.ScalarType

# We need to pass both the ufl element and the finat element
Expand Down
5 changes: 4 additions & 1 deletion firedrake/preconditioners/hiptmair.py
Original file line number Diff line number Diff line change
Expand Up @@ -180,7 +180,10 @@ def coarsen(self, pc):

if G_callback is None:
from firedrake.preconditioners.hypre_ams import chop
interp_petscmat = chop(Interpolator(dminus(test), V, bcs=bcs + coarse_space_bcs).callable().handle)
fcp = {"mode": "vanilla",}
interp_petscmat = chop(Interpolator(dminus(test), V,
bcs=bcs + coarse_space_bcs,
form_compiler_parameters=fcp).callable().handle)
else:
interp_petscmat = G_callback(V, coarse_space, bcs, coarse_space_bcs)

Expand Down

0 comments on commit 7692be6

Please sign in to comment.