Dofmap of fem.VectorFunctionSpace vs. fem.FunctionSpace(mesh, MixedElement)

Hi all,

please consider the following code:

import dolfinx
import ufl

mesh = dolfinx.mesh.create_unit_square(MPI.COMM_WORLD, 1, 1, dolfinx.mesh.CellType.quadrilateral)
V = dolfinx.fem.VectorFunctionSpace(mesh, ("CG", 1))
print(V.dofmap.list)
print(V.sub(0).dofmap.list)
print(V.sub(1).dofmap.list)

fe = ufl.FiniteElement("CG", ufl.quadrilateral, 1)
element = ufl.MixedElement([fe, fe])
W = dolfinx.fem.FunctionSpace(mesh, element)

print(W.dofmap.list)
print(W.sub(0).dofmap.list)
print(W.sub(1).dofmap.list)

The output is (copied from interactive session)

In [49]: V.dofmap.list
Out[49]: 
<AdjacencyList> with 1 nodes
  0: [0 1 2 3 ]

In [50]: V.sub(0).dofmap.list
Out[50]: 
<AdjacencyList> with 1 nodes
  0: [0 2 4 6 ]

In [51]: V.sub(1).dofmap.list
Out[51]: 
<AdjacencyList> with 1 nodes
  0: [1 3 5 7 ]

In [52]: W.dofmap.list
Out[52]: 
<AdjacencyList> with 1 nodes
  0: [0 1 2 3 4 5 6 7 ]

In [53]: W.sub(0).dofmap.list
Out[53]: 
<AdjacencyList> with 1 nodes
  0: [0 1 2 3 ]

In [54]: W.sub(1).dofmap.list
Out[54]: 
<AdjacencyList> with 1 nodes
  0: [4 5 6 7 ]

Remarks/Questions:

  1. I would have expected an array of length 8 for both V.dofmap.list and W.dofmap.list
  2. I would have expected the same (local) layout of the dofs for V.sub(0).dofmap.list and W.sub(0).dofmap.list
  3. Under the hood basix defines the finite element and layout of the dofs independent of how the space is constructed? (or does the local ordering of dofs come from ufl in the case of the space W?)
  4. Could someone please explain the above findings?

Thanks,
Philipp

PS: I have fenicsx installed via conda-forge:

fenics-basix              0.5.0           py310hbf28c38_0    conda-forge
fenics-dolfinx            0.5.1           py310h117620d_100    conda-forge
fenics-ffcx               0.5.0              pyhb871ab6_1    conda-forge
fenics-libbasix           0.5.0                heb3b609_0    conda-forge
fenics-libdolfinx         0.5.1              hcde13a9_100    conda-forge
fenics-ufcx               0.5.0                hb871ab6_1    conda-forge
fenics-ufl                2022.2.0           pyhd8ed1ab_0    conda-forge

The difference is that V.dofmap.bs=2 while W.dofmap.bs=1.
This is used internally to avoid doing extra computations with basis function when using vector spaces.

For a mixed element, the dofs are most usually ordered (dofs_element_1, dofs_element_2),
while for a vector element, you have that they are ordered (dof_0_x, dof_0_y, dof_0_z, dof_1_x, … , dof_N_x, dof_N_y, dof_N_z).

If you are using the same finite element in the MixedElement and the VectorElement, basix does the same thing (basix doesn’t really have a notion of vector spaces:

it only uses it to interface with ufl.

Thank you for the explanation @dokken !

I observed the same thing when I increased the size of the mesh, e.g.

from dolfinx import mesh, fem
msh  = mesh.create_unit_square(MPI.COMM_WORLD, 2, 2)
FE   = FiniteElement("Lagrange", msh.ufl_cell(), 1)
ME   = MixedElement([FE, FE])
V_ME = fem.FunctionSpace(msh, ME)
_, dofmaps = V_ME.sub(0).collapse()
print(dofmaps)

It returns [0, 1, 2, 6, 8, 10, 12, 14, 16]. Except for the first triangle element, other DOFs are ordered.

I am wondering if the feature of not-entirely ordered dofmap in a mixed element could cause interpolation issue especially across two non-matching meshes.

It should not be a problem, as the ordering assumption is not used anywhere in the dolfinx source code