Facet permutation meaning

If I run this example with 2 processors:

from dolfinx import mesh
from mpi4py import MPI
import numpy as np

comm = MPI.COMM_WORLD
rank = comm.Get_rank()

def main():
    # Meant to be ran with 2 procs
    els_side = 4
    domain  = mesh.create_unit_square(MPI.COMM_WORLD, els_side, els_side, mesh.CellType.quadrilateral)
    domain.topology.create_connectivity(1,0)
    normal_facet = 8
    upside_down_facet = 10
    owner_cell_normal_facet = 3
    owner_cell_upside_down_facet = 2
    local_index_normal_facet = 3
    local_index_upside_down_facet = 3
    num_facets_per_cell = 4
    con_facet_nodes = domain.topology.connectivity(1, 0)
    domain.topology.create_entity_permutations()
    fperms = domain.topology.get_facet_permutations()
    if rank==0:
        l = con_facet_nodes.links(normal_facet)
        print(f"Rank {rank}, normal facet goes from {domain.geometry.x[l[0]]} to \
                {domain.geometry.x[l[1]]}", flush=True)
        l = con_facet_nodes.links(upside_down_facet)
        print(f"Rank {rank}, upside down facet goes from {domain.geometry.x[l[0]]} to \
                {domain.geometry.x[l[1]]}", flush=True)
        print(f"Facet permutation info for normal facet = {np.unpackbits(fperms[owner_cell_normal_facet*num_facets_per_cell+local_index_normal_facet])}")
        print(f"Facet permutation info for upside down facet = {np.unpackbits(fperms[owner_cell_upside_down_facet*num_facets_per_cell+local_index_upside_down_facet])}")

if __name__=="__main__":
    main()

I get:

Rank 0, normal facet goes from [0.5  0.75 0.  ] to                 [0.5 1.  0. ]
Rank 0, upside down facet goes from [0.5  0.75 0.  ] to                 [0.5 0.5 0. ]
Facet permutation info for normal facet = [0 0 0 0 0 0 0 0]
Facet permutation info for upside down facet = [0 0 0 0 0 0 0 1]

So one facet is going from bottom node to upper node and the other one is going in the opposite direction. I need to retrieve this facet orientation / permutation information but it is not so clear to me how to do so from topology.get_facet_permutations. The documentation says

The data can be unpacked with numpy.unpack_bits.

I get 8 bits for each facet. How to interpret them?

You then need to read the paper describing facet permutations and transformations:

or the explanation in DOLFINx on when they are required (section 8.4)

1 Like

Thanks for the references! I have yet to give it an in depth look but I understand that I am supposed to deal with permutations using Basix. I had a quick look and I saw no function in Basix that takes uint8s. Any tips?

The functions in question are: basix/python/basix/finite_element.py at 92c54d2dbd3f10771453e1f859c20b5f6ca0da42 · FEniCS/basix · GitHub

which are also available through DOLFINx: dolfinx.fem — DOLFINx 0.10.0.0 documentation

Please note that you only need facet permutations if you are considering hand written interior facet kernels. They are not needed for cell integrals or exterior facet integrals (ffcx/ffcx/ir/integral.py at 7164d8cfec00eec4dda0a8752819beeb7c0009dc · FEniCS/ffcx · GitHub). It might also be needed for co-dim 1 assembly by hand (i.e. integration of data from a mesh and a submesh of facets in a single form).

1 Like

Thank you! I am doing assembly by hand on internal facets with data interpolated from another mesh.

Thank you for the help so far. I’ve been looking into these functions, they allow me to work with 32 bit integers that represent cell permutations and reorder data. The facet permutations are 8 bit integers. Is there a way to work with those?

I’ve been playing around with T_apply and I’m not sure it does anything on the Python side. I modified this test as follows:

data1 = np.array(list(range(size**2)), dtype=np.float32)
data0 = data1.copy()
e.T_apply(data1, size, cell_info)
assert np.allclose(data0, data1)

and it still passes. It should fail since I’m comparing data pre and post permutation.

What element type and degree are you using? For low degree elements, these permutations are often the identity

1 Like

I was playing around with quad elements of degree 1 with made up cell permutations and nothing changed. Then I looked into the test. I ran it with with pytest parametrization and all the 305 cases pass this assertion. Here is a fork with the change, it passes the CI up to linting.

EDIT: I tried making the same change in another test and it does fail there (commit, CI).

EDIT 2: Thanks to your comment @mscroggs I am starting to understand that it is not a Basix DOF transformation that I need since the type of element I am using has dof_transformations_are_identity = True.

In my case, I have uint8 that describes the facet permutation and I need to reoder its Gauss points accordingly. So far the most helpful comment I found in the code was in permutationcomputation.h:

/// 1. The facet rotation and reflection data is encoded so that:
///
/// - n % 2 gives the number of reflections to apply
/// - n // 2 gives the number of rotations to apply
/// …

So for N Gauss points for a given facet I have to figure out the corresponding permutation.

I am reading this as recommended by @dokken

since I want to understand what are exactly these rotations and reflections. The most relevant paragraph I found is:

For sub-entities that are faces, we use two base transformations: one of these represents moving the vertices one position clockwise (i.e., a 120° rotation for a equilateral triangle or a 90° rotation for a square); the other represents a reflection of the face. […]

I’ve looked a bit but I can’t find anywhere where reflection is explained textually although there are many schematics:

EDIT: I think I found the right functions in ffcx: permute_quadrature_interval, permute_quadrature_triangle and permute_quadrature_quadrilateral