Meaning of .collapse()

Imagine you have a Mixed element, lets say a Raviart Thomas space and a DG 0 space (space of cell-wise constants).

import ufl
from mpi4py import MPI

import dolfinx

mesh = dolfinx.mesh.create_unit_square(MPI.COMM_WORLD, 10, 10)
RT = ufl.FiniteElement("RT", mesh.ufl_cell(), 1)
DG0 = ufl.FiniteElement("DG", mesh.ufl_cell(), 0)
el = ufl.MixedElement(RT, DG0)
W = dolfinx.fem.FunctionSpace(mesh, el)

If we create a function in the mixed space, say up, we can check how many dofs we have in the problem:

up = dolfinx.fem.Function(W)
print(len(up.x.array))

which returns 520.
If we only want to work on the space of RT dofs, the 0th subspace of W, we call

V, V_to_W = W.sub(0).collapse()
u = dolfinx.fem.Function(V)
print(len(u.x.array), len(V_to_W))

Returning

320 320

Similarly, for the DG-0 space, we have

Q, Q_to_W = W.sub(1).collapse()
p = dolfinx.fem.Function(Q)
print(len(p.x.array), len(Q_to_W))

returning

200 200

As you can see here, collapse just gives you a function space over the ith element of the mixed element, and the second output is the map from each dof in the sub space, to the parent space.

If we keep on considering the first sub-space, the RT-space, whose dofs are defined by the facets of a cell, (with functionals being integrals over said facet), one cannot locate those geometrically, as the dof does not have a physical coordinate

dolfinx.fem.locate_dofs_geometrical(V, lambda x: x[0] < 0.5)

returns

Traceback (most recent call last):
  File "/fenics/shared/mwe123.py", line 24, in <module>
    dolfinx.fem.locate_dofs_geometrical(V, lambda x: x[0] < 0.5)
  File "/usr/local/dolfinx-real/lib/python3.10/dist-packages/dolfinx/fem/bcs.py", line 54, in locate_dofs_geometrical
    return _cpp.fem.locate_dofs_geometrical(V._cpp_object, marker)
RuntimeError: Cannot evaluate dof coordinates - this element does not have pointwise evaluation.

Thus, for the mixed-element W, one would need to work on the sub-space Q, if one want to use locate_dofs_topological, as the mixed-element cannot compute this.

As we want to have the dofs of the combined space, we send in the tuplet (W.sub(1), Q)), to being able to pass the mapping information on:

print(dolfinx.fem.locate_dofs_geometrical((W.sub(1), Q), lambda x: x[0] < 0.5))

returns

[array([ 87, 101, 119, 121, 135, 138, 156, 158, 160, 174, 177, 180, 198,
       200, 202, 204, 218, 221, 224, 227, 245, 247, 249, 251, 253, 267,
       270, 273, 276, 279, 296, 298, 300, 302, 304, 315, 318, 321, 324,
       327, 341, 343, 345, 347, 349, 358, 361, 364, 367, 370, 381, 383,
       385, 387, 389, 396, 399, 402, 405, 408, 416, 418, 420, 422, 424,
       429, 432, 435, 438, 441, 446, 448, 450, 452, 454, 457, 460, 463,
       466, 469, 471, 473, 475, 477, 480, 483, 486, 489, 491, 493, 495,
       498, 501, 504, 506, 508, 511, 514, 516, 519], dtype=int32), array([ 30,  36,  42,  43,  49,  50,  56,  57,  58,  64,  65,  66,  72,
        73,  74,  75,  81,  82,  83,  84,  90,  91,  92,  93,  94, 100,
       101, 102, 103, 104, 110, 111, 112, 113, 114, 119, 120, 121, 122,
       123, 128, 129, 130, 131, 132, 136, 137, 138, 139, 140, 144, 145,
       146, 147, 148, 151, 152, 153, 154, 155, 158, 159, 160, 161, 162,
       164, 165, 166, 167, 168, 170, 171, 172, 173, 174, 175, 176, 177,
       178, 179, 180, 181, 182, 183, 184, 185, 186, 187, 188, 189, 190,
       191, 192, 193, 194, 195, 196, 197, 198, 199], dtype=int32)]

where the first list is the dof entries in W, the second entry is the dof entries in Q.

You can verify that the second row is equivalent to the dofs in Q by calling:

print(dolfinx.fem.locate_dofs_geometrical(Q, lambda x: x[0] < 0.5))
7 Likes