Skip to content

Commit

Permalink
Fix 3D manifolds from vector and tensorfunctionspace (#2655)
Browse files Browse the repository at this point in the history
* Fix 3D manifolds from vector and tensorfunctionspace

* Flake8

* Flake8 again
  • Loading branch information
jorgensd committed May 15, 2023
1 parent 4cb7eb8 commit 45a606c
Show file tree
Hide file tree
Showing 2 changed files with 27 additions and 9 deletions.
10 changes: 7 additions & 3 deletions python/dolfinx/fem/function.py
Original file line number Diff line number Diff line change
Expand Up @@ -632,11 +632,13 @@ def VectorFunctionSpace(mesh: Mesh,
ufl_e = basix.ufl.element(element.family(), element.cell_type, # type: ignore
element.degree(), element.lagrange_variant, # type: ignore
element.dpc_variant, element.discontinuous, # type: ignore
shape=None if dim is None else (dim, ), gdim=mesh.geometry.dim, rank=1)
shape=(mesh.geometry.dim, ) if dim is None else (dim, ),
gdim=mesh.geometry.dim, rank=1)
except AttributeError:
ed = ElementMetaData(*element)
ufl_e = basix.ufl.element(ed.family, mesh.basix_cell(), ed.degree,
shape=None if dim is None else (dim, ), gdim=mesh.geometry.dim, rank=1)
shape=(mesh.geometry.dim, ) if dim is None else (dim, ),
gdim=mesh.geometry.dim, rank=1)
return FunctionSpace(mesh, ufl_e)


Expand All @@ -647,7 +649,9 @@ def TensorFunctionSpace(mesh: Mesh, element: typing.Union[ElementMetaData, typin
raise ValueError("Cannot create tensor element containing a non-scalar.")

e = ElementMetaData(*element)
gdim = mesh.geometry.dim
shape_ = (gdim, gdim) if shape is None else shape
ufl_element = basix.ufl.element(e.family, mesh.basix_cell(),
e.degree, shape=shape, symmetry=symmetry,
e.degree, shape=shape_, symmetry=symmetry,
gdim=mesh.geometry.dim, rank=2)
return FunctionSpace(mesh, ufl_element)
26 changes: 20 additions & 6 deletions python/test/unit/fem/test_function_space.py
Original file line number Diff line number Diff line change
Expand Up @@ -4,16 +4,15 @@
#
# SPDX-License-Identifier: LGPL-3.0-or-later
"""Unit tests for the FunctionSpace class"""
import basix
import numpy as np
import pytest

import basix
from basix.ufl import mixed_element, element
from dolfinx.fem import Function, FunctionSpace, VectorFunctionSpace
from basix.ufl import element, mixed_element
from dolfinx.fem import (Function, FunctionSpace, TensorFunctionSpace,
VectorFunctionSpace)
from dolfinx.mesh import create_mesh, create_unit_cube
from ufl import Cell, Mesh, TestFunction, TrialFunction, grad

from mpi4py import MPI
from ufl import Cell, Mesh, TestFunction, TrialFunction, grad


@pytest.fixture
Expand Down Expand Up @@ -252,3 +251,18 @@ def test_vector_function_space_cell_type():
# is correct
V = VectorFunctionSpace(mesh, ('Lagrange', 1))
assert V.ufl_element().cell() == cell


@pytest.mark.skip_in_parallel
def test_manifold_spaces():
vertices = [(0.0, 0.0, 1.0), (1.0, 1.0, 1.0), (1.0, 0.0, 0.0), (0.0, 1.0,
0.0)]
cells = [(0, 1, 2), (0, 1, 3)]
domain = Mesh(element("Lagrange", "triangle", 1, gdim=3, rank=1))
mesh = create_mesh(MPI.COMM_WORLD, cells, vertices, domain)
QV = VectorFunctionSpace(mesh, ("Lagrange", 1))
QT = TensorFunctionSpace(mesh, ("Lagrange", 1))
u = Function(QV)
v = Function(QT)
assert u.ufl_shape == (3,)
assert v.ufl_shape == (3, 3)

0 comments on commit 45a606c

Please sign in to comment.