Simplicial Complex in Three Dimensions

We dsecribe the data structure of the simplicial complex associated to a three dimensional trianglulation give by node,elem . The node records the coordinates of vertices and elem is the pointer from local to global incices of vertices. See Basic mesh data structure.

A brief summary.

  • edge: asecond ordering, i.e. edge(:,1)<edge(:,2)

  • face: asecond ordering, i.e. face(:,1)<face(:,2)<face(:,3)

  • elem: either the positive ordering or the ascend ordering. The default one is the positive ordering and the asecond ordering is mainly used for edge and face elements.

  • Use [elem,bdFlag] = sortelem3(elem,bdFlag) to change the ordering to the ascend ordering. Note that bdFlag should be switched together.

The multigrid solvers use the original ordering of elem obtained from either uniform refinement or bisection methods. So let elemold=elem before sort.

  • Examples on the usage: Poisson3RT0; Maxwell; Maxwell2;

Outline

In [3]:
elem = [1 4 5 8; 1 4 5 7];
node = [1,0,0; 1,1,1; 1,-1,-1; 0,1,0; -2,-1,0; 1,1,-1; 0,1,1; 0,-1,-1];
showmesh3(node,elem,[],'FaceAlpha',0.25);
findelem3(node,elem);
findnode3(node,elem(:));

The basic data structure of a mesh consists of node and elem. The corresponding simplicial complex consists of vertices, edges, faces, and tetrahedron. We shall discuss three issues

  • Indexing of simplexes
  • Ordering of vertices
  • Orientation of simplexes

The indexing and ordering are related and the ordering and orientation are mixed together. However the indexing has nothing to do with the orientation. The indexing and ordering are the combinarotry structure, i.e. only elem is needed, while the orientation also depends on node, the geometry emembdding of vertices.

For indexing, ordering and orientation, there are always local and global versions. The relation between the local and global version is the most complicated issue.

Indexing of Simplexes

The indexing refers to the numbering of simplexes, e.g., which face is numbered as the first one. There are two types of the indexing: local and global. Each simplex in the simplicial complex has a unique index which is called the global index. In one tetrahedra, the four vertices and four faces have their local index from 1:4.

In the assembling procedure of finite element methods, an element-wise matrix using the local indexing is first computed and then assembled to get a big matrix using the global indexing. Thus the pointer from the local indexing to the global indexing is indispensible. For bases independent of the ordering and orientation, e.g., P1 and P2 elements, this pointer is sufficient, otherwise, the inconsistency of the local ordering/orientation and the global ordering/orientation should be taken into account.

Local indexing

The tetrahedron consists of four vertices indexed as [1 2 3 4]. Each tetrahedron contains four faces and six edges. They can be indexed as

locFace = [2 3 4; 1 3 4; 1 2 4; 1 2 3];
locEdge = [1 2; 1 3; 1 4; 2 3; 2 4; 3 4];

In locFace, the i-th face is opposite to the i-th vertices and thus this is called opposite indexing. In locEdge, it is the lexicographic indexing which is induced from the lexicographic ordering of the six edges. The ordering of vertices of each face or edge will not change the indexing. For example, the following locFacec and locEdged has the same indexing as locFace and locEdge but a different ordering of vertices.

locFacec = [2 3 4; 1 4 3; 1 2 4; 1 3 2];
locEdge = [2 1; 3 1; 4 1; 3 2; 4 2; 4 3];

Indeed any permuation of each simplex will represent the same simplex and will not change the indexing. The ordering of vertices will affect the orientation and will be discussed later.

For a face consists of three vertices [1 2 3], there are two indexing schemes of its three edges.

  • Oppoiste indexing [2 3; 3 1; 1 2]
  • Lexicographic indexing [1 2; 1 3; 2 3]

Each indexing scheme has its advantange and disadavantange and which one to chose depends on the consideration of ordering and orientation.

Global indexing and vertex pointers

Each simplex in the simplicial complex has a unqiuely index. It is represented by vertices pointer from the local index to the globa index of vertices.

The matrix elem is the pointer from local to global indices of vertices of tetrahedron, e.g. elem(t,1)=25 means the first vertex of the tetrahedron t is the 25-th vertex.

Similarly the NE x 2 matrix edge records all edges and the NF x 3 by 3 matrix face records all faces of the triangulation. These are vertices pointers. We shall discuss the elementwise pointer from the local indices to the global indices for edges and faces.

Generate index pointers for edges and faces

One can easily collect edges and faces elementwise. The issue is the duplication. For example, each interior face will be counted twice. The unique function is applied such that each edge or face has a unique global index.

Edge and Face

In [4]:
totalEdge = uint32([elem(:,[1 2]); elem(:,[1 3]); elem(:,[1 4]); ...
                    elem(:,[2 3]); elem(:,[2 4]); elem(:,[3 4])]);
sortedTotalEdge = sort(totalEdge,2);
[edge, ~, je] = unique(sortedTotalEdge,'rows');
display(edge);

totalFace = uint32([elem(:,[2 3 4]); elem(:,[1 4 3]); ...
                    elem(:,[1 2 4]); elem(:,[1 3 2])]);
sortedTotalFace = sort(totalFace,2);                
[face, i2, jf] = unique(sortedTotalFace,'rows');
display(face);
edge =

           1           4
           1           5
           1           7
           1           8
           4           5
           4           7
           4           8
           5           7
           5           8


face =

           1           4           5
           1           4           7
           1           4           8
           1           5           7
           1           5           8
           4           5           7
           4           5           8

In iFEM, N,NE,NF,NT represents the number of vertices, edges, faces and tetrahedrons, resprectively.

N = size(node,1); NT = size(elem,1); NF = size(face,1); NE = size(edge,1);

In the assembling procedure, the matrix is always computed elementwise and then assemble to a big one. A pointer from the local index of a simplex to its global index is thus indispensible.

Elementwise pointers

  • elem2node = elem
  • elem2face(1:NT, 1:4)
  • elem2edge(1:NT, 1:6)

Such information is exactly stored in the output of unique function. For example, elem2face(t,1) = 17 means the first face of t (spanned by [2 3 4]) is the 17-th element in the face matrix.

In [6]:
N = size(node,1); NT = size(elem,1); NF = size(face,1); NE = size(edge,1);
elem2edge = uint32(reshape(je,NT,6));
elem2face = uint32(reshape(jf,NT,4));
display(elem2edge);
display(elem2face);
elem2edge =

           1           2           4           5           7           9
           1           2           3           5           6           8


elem2face =

           7           5           3           1
           6           4           2           1

Face to edge Pointer

|face2edge(1:NF,1:3)| records the global indices of three edges of a face. This pointer depends on the ordering of vertices of faces and the indexing of local edges in a face. We list the following two important cases. Other combinations is possible but not attractive.

  • Ascend ordering.

All local faces and local edges are ascend ordered.

locFace = [2 3 4; 1 3 4; 1 2 4; 1 2 3];
locEdge = [1 2; 1 3; 1 4; 2 3; 2 4; 3 4];
edgeofFace = [1 2; 1 3; 2 3];
locface2edge = [4 5 6; 2 3 6; 1 3 5; 1 2 4];

  • Consistent ordering

The local face is ordered such that the corresponding orientation is consistent with the induced orientation. locFace = [2 3 4; 1 4 3; 1 2 4; 1 3 2]; locEdge = [1 2; 1 3; 1 4; 2 3; 2 4; 3 4]; edgeofFace = [2 3; 3 1; 1 2];
locface2edge = [6 5 4; 6 2 3; 5 3 1; 4 1 2];

The global one can be obtained from the composition of elem2face and locface2edge. For example, for the asecnd ordering scheme,

In [7]:
face2edge(elem2face(:,1),:) = elem2edge(:,[4 5 6]);
face2edge(elem2face(:,2),:) = elem2edge(:,[2 3 6]);
face2edge(elem2face(:,3),:) = elem2edge(:,[1 3 5]);
face2edge(elem2face(:,4),:) = elem2edge(:,[1 2 4]);

Ordering of Vertices

We discuss the ordering of vertices of simplexes. Again there are local ordering and global ordering. They may not be consistent and a sign array is used to record the inconsistency if any.

The local ordering refers to the ordering of vertices in locFace or locEdge, i.e. the ordering of the local index of vertices. For elements associated to faces or edges, the local ordering could be used in the formulation of the local basis and thus the ordering does matter.

The global ordering refers to the ordering of vertices in face or edge, i.e., the ordering of the global index of vertices. Note that that in either local or global ordering, permutation of vertices will represent the same simplex. To fix an ordering we need extra information.

elem

The local ordering is always [1 2 3 4]. Any permutation of four vertices of a tetrahedon still represents the same tetrahedron. Such freedom provide a room to record more information like:

  • global ordering of vertices
  • an orientation of element
  • refinement rules (uniform refinement or bisection)

For 2-D triangulations, three vertices of a triangle in 2-D is sorted counter-cloclwise and the first vertex is chosen as the newest vertex. Such ordering enables the efficient implementation of local refinement and coarsening in 2-D; see Bisection in Two Dimensions and Coarsening in Two Dimensions.

In 3-D, for the longest edge bisection, the newest vertex (with the highest generation) is stored as the last (4-th) vertex of a tetrahedron. For 3-D Red Refinement, the ordering determines the shape regularity of refined triangulation. Permuation of vertices in elem could deterioriate the angle condition after the refinement.

We shall reserve the ordering of elem from the mesh refinement and coarsening since they are more subtle. We switch the ordering when generating data structure for finite element basis and assemble the matrix equation. Such sorting is hidden in the subroutines when a finite element basis requiring ordering is generated.

Two types of ordering of elem is of particular importantance

  • Positive ordering
  • Ascend ordering

In the positive ordering, the four vertices are ordered such that the signed volume, the mix product of vectors (v12,v13,v14), is positive. This is the default ordering used so far. fixorder3 will switch the vertices for elements with negative volume.

v = simplexvolume(node,elem) % returns the singed volume
elem = fixorder(node,elem)   % switchs the vertices for elements with negative volume.

In the ascend ordering, the vertices of elem is sorted such that

elem(t,1) < elem(t,2) < elem(t,3) < elem(t,4). 

Such ordering will benefit the construction of local bases for high order basis or basis with orientation. This can be easily achieved by elem = sort(elem,2). Howevery, one has to rotate the boundary flag accordingly using

[elem,bdFlag] = sortelem3(elem,bdFlag);

Orientation