Skip to content

Unable to define a single R space when using an extruded mesh #4296

Open
@JamesJackaman

Description

@JamesJackaman

If discretising using a extruded mesh it's not currently possible to define a single R space.

Defining
FunctionSpace(mesh,'R',degree=0,vfamily='R',vdegree=0) builds a function space with DOFs equal to the number of the extruded dimension, instead of just one DOF.

If instead I define it with FunctionSpace(mesh,'R',degree=0) I run into errors in the solve command which are likely caused by the element type not being supported.

MFE:

from firedrake import *

base = IntervalMesh(100,1)
mesh = ExtrudedMesh(base,layers=100,layer_height=1/100,periodic=False)

U = FunctionSpace(mesh,TensorProductElement(FiniteElement("CG",interval,degree=1),FiniteElement("CG",interval,degree=1)))

R = FunctionSpace(mesh,'R',degree=0)

V = MixedFunctionSpace((U,R))

v = Function(V)
u, r = v.subfunctions
phi, psi = TestFunctions(V)

r.assign(100)

u, r = split(v)

F = u*phi*dx + r*psi*dx

solver_parameters = {"mat_type": "matfree",
                     "pc_type": "fieldsplit",
                     "pc_fieldsplit_type": "additive",
                     "pc_fieldsplit_0_fields": "0",
                     "pc_fieldsplit_1_fields": "1",                         
                     "fieldsplit_0": {"pc_type": "python",
                                      "pc_python_type": "firedrake.AssembledPC",
                                      "assembled_pc_type": "lu",
                                      "ksp_type": "preonly"},
                     "fieldsplit_1": {"ksp_type": "preonly",
                                      "pc_type": "python",
                                      "pc_python_type":"firedrake.AssembledPC",
                                      "assembled_pc_type":"jacobi"},
                     "ksp_monitor_true_residual": None,
                     "pc_factor_mat_solver_type": "mumps",
                     "pc_factor_shift_type":"nonzero",
                     "snes_monitor": None,
                     }

nvp = NonlinearVariationalProblem(F, v)
nvs = NonlinearVariationalSolver(nvp,solver_parameters=solver_parameters)
nvs.solve()

Error:

Traceback (most recent call last):
  File "/Users/jamesij/software/firedrake/lib/python3.13/site-packages/FIAT/expansions.py", line 253, in __new__
    expansion_set = {
                    ~
    ...<3 lines>...
        reference_element.TETRAHEDRON: TetrahedronExpansionSet,
        ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
    }[ref_el.get_shape()]
    ~^^^^^^^^^^^^^^^^^^^^
KeyError: 99

During handling of the above exception, another exception occurred:

Traceback (most recent call last):
  File "petsc4py/PETSc/petscsnes.pxi", line 335, in petsc4py.PETSc.SNES_Function
  File "/Users/jamesij/software/firedrake/lib/python3.13/site-packages/firedrake/solving_utils.py", line 432, in form_function
    ctx._assemble_residual(tensor=ctx._F, current_state=ctx._x)
    ~~~~~~~~~~~~~~~~~~~~~~^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
  File "/Users/jamesij/software/firedrake/lib/python3.13/site-packages/firedrake/assemble.py", line 1005, in assemble
    self.execute_parloops(tensor)
    ~~~~~~~~~~~~~~~~~~~~~^^^^^^^^
  File "/Users/jamesij/software/firedrake/lib/python3.13/site-packages/firedrake/assemble.py", line 1248, in execute_parloops
    for parloop in self.parloops(tensor):
                   ~~~~~~~~~~~~~^^^^^^^^
  File "/Users/jamesij/software/firedrake/lib/python3.13/site-packages/firedrake/assemble.py", line 1038, in parloops
    for local_kernel, subdomain_id in self.local_kernels:
                                      ^^^^^^^^^^^^^^^^^^
  File "/opt/homebrew/Cellar/[email protected]/3.13.3/Frameworks/Python.framework/Versions/3.13/lib/python3.13/functools.py", line 1026, in __get__
    val = self.func(instance)
  File "/Users/jamesij/software/firedrake/lib/python3.13/site-packages/firedrake/assemble.py", line 1076, in local_kernels
    kernels = tsfc_interface.compile_form(
        self._form, "form", diagonal=self.diagonal,
        parameters=self._form_compiler_params
    )
  File "petsc4py/PETSc/Log.pyx", line 188, in petsc4py.PETSc.Log.EventDecorator.decorator.wrapped_func
  File "petsc4py/PETSc/Log.pyx", line 189, in petsc4py.PETSc.Log.EventDecorator.decorator.wrapped_func
  File "petsc4py/PETSc/Log.pyx", line 188, in petsc4py.PETSc.Log.EventDecorator.decorator.wrapped_func
  File "petsc4py/PETSc/Log.pyx", line 189, in petsc4py.PETSc.Log.EventDecorator.decorator.wrapped_func
  File "/Users/jamesij/software/firedrake/lib/python3.13/site-packages/pyop2/caching.py", line 558, in wrapper
    value = func(*args, **kwargs)
  File "petsc4py/PETSc/Log.pyx", line 188, in petsc4py.PETSc.Log.EventDecorator.decorator.wrapped_func
  File "petsc4py/PETSc/Log.pyx", line 189, in petsc4py.PETSc.Log.EventDecorator.decorator.wrapped_func
  File "/Users/jamesij/software/firedrake/lib/python3.13/site-packages/pyop2/caching.py", line 558, in wrapper
    value = func(*args, **kwargs)
  File "petsc4py/PETSc/Log.pyx", line 188, in petsc4py.PETSc.Log.EventDecorator.decorator.wrapped_func
  File "petsc4py/PETSc/Log.pyx", line 189, in petsc4py.PETSc.Log.EventDecorator.decorator.wrapped_func
  File "/Users/jamesij/software/firedrake/lib/python3.13/site-packages/firedrake/tsfc_interface.py", line 220, in compile_form
    tsfc_kernel = TSFCKernel(
        f,
    ...<4 lines>...
        interface, diagonal
    )
  File "/Users/jamesij/software/firedrake/lib/python3.13/site-packages/firedrake/tsfc_interface.py", line 100, in __init__
    tree = tsfc_compile_form(form, prefix=name, parameters=parameters,
                             interface=interface,
                             diagonal=diagonal)
  File "petsc4py/PETSc/Log.pyx", line 188, in petsc4py.PETSc.Log.EventDecorator.decorator.wrapped_func
  File "petsc4py/PETSc/Log.pyx", line 189, in petsc4py.PETSc.Log.EventDecorator.decorator.wrapped_func
  File "/Users/jamesij/software/firedrake/lib/python3.13/site-packages/pyop2/caching.py", line 558, in wrapper
    value = func(*args, **kwargs)
  File "petsc4py/PETSc/Log.pyx", line 188, in petsc4py.PETSc.Log.EventDecorator.decorator.wrapped_func
  File "petsc4py/PETSc/Log.pyx", line 189, in petsc4py.PETSc.Log.EventDecorator.decorator.wrapped_func
  File "/Users/jamesij/software/firedrake/lib/python3.13/site-packages/pyop2/caching.py", line 558, in wrapper
    value = func(*args, **kwargs)
  File "/Users/jamesij/software/firedrake/lib/python3.13/site-packages/tsfc/driver.py", line 73, in compile_form
    kernel = compile_integral(integral_data, fd, prefix, parameters, interface=interface, diagonal=diagonal)
  File "/Users/jamesij/software/firedrake/lib/python3.13/site-packages/tsfc/driver.py", line 135, in compile_integral
    integrand_exprs = builder.compile_integrand(integrand, params, ctx)
  File "/Users/jamesij/software/firedrake/lib/python3.13/site-packages/tsfc/kernel_interface/common.py", line 140, in compile_integrand
    set_quad_rule(params, info.domain.ufl_cell(), info.integral_type, functions)
    ~~~~~~~~~~~~~^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
  File "/Users/jamesij/software/firedrake/lib/python3.13/site-packages/tsfc/kernel_interface/common.py", line 323, in set_quad_rule
    finat_elements = set(create_element(f.ufl_element()) for f in functions
                         if f.ufl_element().family() != "Real")
  File "/Users/jamesij/software/firedrake/lib/python3.13/site-packages/tsfc/kernel_interface/common.py", line 323, in <genexpr>
    finat_elements = set(create_element(f.ufl_element()) for f in functions
                         ~~~~~~~~~~~~~~^^^^^^^^^^^^^^^^^
  File "/Users/jamesij/software/firedrake/lib/python3.13/site-packages/finat/element_factory.py", line 324, in create_element
    finat_element, deps = _create_element(ufl_element,
                          ~~~~~~~~~~~~~~~^^^^^^^^^^^^^
                                          shape_innermost=shape_innermost,
                                          ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
                                          shift_axes=shift_axes,
                                          ^^^^^^^^^^^^^^^^^^^^^^
                                          restriction=restriction)
                                          ^^^^^^^^^^^^^^^^^^^^^^^^
  File "/Users/jamesij/software/firedrake/lib/python3.13/site-packages/finat/element_factory.py", line 354, in _create_element
    finat_element, deps = convert(ufl_element, **kwargs)
                          ~~~~~~~^^^^^^^^^^^^^^^^^^^^^^^
  File "/opt/homebrew/Cellar/[email protected]/3.13.3/Frameworks/Python.framework/Versions/3.13/lib/python3.13/functools.py", line 934, in wrapper
    return dispatch(args[0].__class__)(*args, **kw)
           ~~~~~~~~~~~~~~~~~~~~~~~~~~~^^^^^^^^^^^^^
  File "/Users/jamesij/software/firedrake/lib/python3.13/site-packages/finat/element_factory.py", line 254, in convert_mixedelement
    elements, deps = zip(*[_create_element(elem, **kwargs)
                           ~~~~~~~~~~~~~~~^^^^^^^^^^^^^^^^
  File "/Users/jamesij/software/firedrake/lib/python3.13/site-packages/finat/element_factory.py", line 354, in _create_element
    finat_element, deps = convert(ufl_element, **kwargs)
                          ~~~~~~~^^^^^^^^^^^^^^^^^^^^^^^
  File "/opt/homebrew/Cellar/[email protected]/3.13.3/Frameworks/Python.framework/Versions/3.13/lib/python3.13/functools.py", line 934, in wrapper
    return dispatch(args[0].__class__)(*args, **kw)
           ~~~~~~~~~~~~~~~~~~~~~~~~~~~^^^^^^^^^^^^^
  File "/Users/jamesij/software/firedrake/lib/python3.13/site-packages/finat/element_factory.py", line 228, in convert_finiteelement
    return lmbda(cell, element.degree(), **finat_kwargs), set()
           ~~~~~^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
  File "/Users/jamesij/software/firedrake/lib/python3.13/site-packages/finat/fiat_elements.py", line 380, in __init__
    super().__init__(FIAT.DiscontinuousLagrange(cell, degree, variant=variant))
                     ~~~~~~~~~~~~~~~~~~~~~~~~~~^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
  File "/Users/jamesij/software/firedrake/lib/python3.13/site-packages/FIAT/discontinuous_lagrange.py", line 222, in __new__
    return P0.P0(ref_el)
           ~~~~~^^^^^^^^
  File "/Users/jamesij/software/firedrake/lib/python3.13/site-packages/FIAT/P0.py", line 48, in __init__
    poly_set = polynomial_set.ONPolynomialSet(ref_el, 0)
  File "/Users/jamesij/software/firedrake/lib/python3.13/site-packages/FIAT/polynomial_set.py", line 117, in __init__
    expansion_set = expansions.ExpansionSet(ref_el, **kwargs)
  File "/Users/jamesij/software/firedrake/lib/python3.13/site-packages/FIAT/expansions.py", line 261, in __new__
    raise ValueError("Invalid reference element type.")
ValueError: Invalid reference element type.

The above exception was the direct cause of the following exception:

Traceback (most recent call last):
  File "/Users/jamesij/repos/kp/test.py", line 45, in <module>
    nvs.solve()
    ~~~~~~~~~^^
  File "petsc4py/PETSc/Log.pyx", line 188, in petsc4py.PETSc.Log.EventDecorator.decorator.wrapped_func
  File "petsc4py/PETSc/Log.pyx", line 189, in petsc4py.PETSc.Log.EventDecorator.decorator.wrapped_func
  File "/Users/jamesij/software/firedrake/lib/python3.13/site-packages/firedrake/adjoint_utils/variational_solver.py", line 108, in wrapper
    out = solve(self, **kwargs)
  File "/Users/jamesij/software/firedrake/lib/python3.13/site-packages/firedrake/variational_solver.py", line 356, in solve
    self.snes.solve(None, work)
    ~~~~~~~~~~~~~~~^^^^^^^^^^^^
  File "petsc4py/PETSc/SNES.pyx", line 1724, in petsc4py.PETSc.SNES.solve
petsc4py.PETSc.Error: error code 101
[0] SNESSolve() at /Users/jamesij/software/petsc/src/snes/interface/snes.c:4847
[0] SNESSolve_NEWTONLS() at /Users/jamesij/software/petsc/src/snes/impls/ls/ls.c:169
[0] SNESComputeFunction() at /Users/jamesij/software/petsc/src/snes/interface/snes.c:2477

Metadata

Metadata

Assignees

No one assigned

    Type

    No type

    Projects

    No projects

    Milestone

    No milestone

    Relationships

    None yet

    Development

    No branches or pull requests

    Issue actions