Simplicial complexes

This module defines a minimal set of functions to generate a dimension-independent cellular complex of simplices. The name of the library was firstly used by CAD Lab at University $La Sapienza$ of Rome in years 1987/88 when started working with dimension-independent simplicial complexes (see "Dimension-independent modeling with simplicial complexes" ). The library provides combinatorial algorithms for some basic functions of geometric modelling with simplicial complexes. In particular, provides the efficient creation of simplicial complexes generated by simplicial complexes of lower dimension, the production of simplicial grids of any dimension, and the extraction of facets (i.e. of $(d-1)$-faces) of simplicial of $d$-complexes. The main aim of the simplicial functions given in this library is to provide optimal combinatorial algorithms, whose time complexity is linear in the size of the output. Such a goal is achieved by calculating each cell in the output via closed combinatorial formulas, that do not require any searching nor data structure traversal to produce their results.

Simplicial extrusion

Here we discuss the implementation of the linear extrusion of simplicial complexes according to the method discussed on papers "Dimension-independent modeling with simplicial complexes" and "Extrusion and boundary evaluation for multidimensional polyhedra". In synthesis, for each $d$-simplex in the input complex, we generate combinatorially a $(d+1)$-simplicial tube, i.e. a chain of $d+1$ simplexes of dimension $d+1$. It can be shown that if the input simplices are a simplicial complex, then the output simplices are a simplicial complex too (i.e. are correctly glued together).

In other words, if the input is a $d$complex, where all $d$-cells either intersect along a common face or are pairwise disjoints, then the output is a simplicial complex of dimension $d+1$. This method is computationally optimal, since it does not require any search or traversal of data structures. The algorithm just writes the output making a constant number $O(1)$ of operation for each one of its $n$ output $d$-cells, so that the time complexity is $\Omega(n)$, where $n = d\,m$, being $m$ the number and $d$ the dimension (and the storage size) of the input cells, represented as arrays of indices of vertices, using the Cells type.

The aim of the extrudeSimplicial function is to generate the output model vertices and cells in a multiple extrusion of a LAR model. First note that the model variable contains a pair (V, FV), where V is the array of input vertices, and FV is the array of $d$-cells (given as list of lists of vertex indices), providing the input representation of a LAR cellular complex. The pattern variable is a list of integers, whose absolute values provide the sizes of the ordered set of 1D subintervals (in local coords) specified by the pattern itself. Such subintervals are assembled in global coordinates, and each one of them is considered either solid or empty depending on the sign of the corresponding measure, which may be either positive (solid subinterval) or negative (void subinterval).

Therefore, a value pattern = [1,1,-1,1] is interpreted as the 1D simplicial complex

[0,1] $\cup$ [1,2] $\cup$ [3,4]

with five vertices W = [[0.0], [1.0], [2.0], [3.0], [4.0]] and three 1-cells [[0,1], [1,2], [3,4]].

V is the list of input $d$-vertices (each given as a list of $d$ coordinates); coords is a list of absolute translation parameters to be applied to V in order to generate the output vertices generated by the combinatorial extrusion algorithm. The cellGroups internal variable is used to select the groups of $(d+1)$-simplices corresponding to solid intervals in the input pattern.

Simplicial grids

The generation of simplicial grids of any dimension and shape using the simplexGrid is amazingly simple. The input parameter shape is either a tuple or a list of integers used to specify the shape of the created array, i.e. both the number of its dimensions (given by len(shape)) and the size of each dimension $k$ (given by the shape[k] element). The implementation starts from the LAR model of an empty simplicial model (denoted as VOID, a predefined constant) and updates the model variable extruding it iteratively according to the specs given by shape. Just notice that the returned grid modelusually has vertices with integer coordinates, that can be subsequently scaled and/or translated and/or mapped in any other way, according to the user needs.

Facet extraction from simplices

A $k$-face of a $d$-simplex is defined as the convex hull of any subset of $k$ vertices. A $(d-1)$-face of a $d$-simplex

$\sigma^d = \langle v_0, v_1, \ldots, v_d \rangle$

is also called a facet. Each of the $d+1$ facets of $\sigma^d$, obtained by removing a vertex from $\sigma^d$, is a $(d-1)$-simplex. A simplex may be oriented in two different ways according to the permutation class of its vertices. The simplex orientation is so changed by either multiplying the simplex by -1, or by executing an odd number of exchanges of its vertices.

The chain of oriented boundary facets of $\sigma^d$, usually denoted as $\partial \sigma^d$, is generated combinatorially as follows:

$\partial\, \sigma^d = \sum_{k=0}^d (-1)^d \langle v_0, \ldots, v_{k-1}, v_{k+1}, \ldots, v_d \rangle$

The larSimplexFacets function, for estraction of non-oriented $(d-1)$-facets of $d$-dimensional simplices, returns a list of $d$-tuples of integers, i.e. the input LAR representation of the topology of a cellular complex. The final steps are used to remove the duplicated facets, by transforming the sorted facets into a set of strings, so removing the duplicated elements.

Examples

Multidimensional simplicial extrusion

The algorithm for multimensional extrusion of a simplicial complex is implemented in the extrudeSimplicial function. This one can be applied to 0-, 1-, 2-, ... simplicial model, to get a 1-, 2-, 3-, .... model. A 1D pattern of linear Array type is used to specify how to decompose the added dimension.

The input and output model are a LAR model, i.e. a pair (vertices, cells), whereas pattern is an array of Int64, to be used as lateral measures of the extruded model. Note that pattern elements are assumed as either solid or empty measures, according to their (+/-) sign.

julia> Lar = LinearAlgebraicRepresentation

julia> V = [[0,0] [1,0] [2,0] [0,1] [1,1] [2,1] [0,2] [1,2] [2,2]]
2×9 Array{Int64,2}:
 0  1  2  0  1  2  0  1  2
 0  0  0  1  1  1  2  2  2

julia> FV = [[1,2,4],[2,3,5],[3,5,6],[4,5,7],[5,7,8],[6,8,9]]
6-element Array{Array{Int64,1},1}:
 [1, 2, 4]
 [2, 3, 5]
 [3, 5, 6]
 [4, 5, 7]
 [5, 7, 8]
 [6, 8, 9]

julia> pattern = repeat([1,2,-3],outer=4)
12-element Array{Int64,1}:
[1,2,-3,1,2,-3,1,2,-3,1,2,-3]

julia> model = (V,FV)
([0 1 2 0 1 2 0 1 2; 0 0 0 1 1 1 2 2 2], 
Array{Int64,1}[[1, 2, 4], [2, 3, 5], [3, 5, 6], [4, 5, 7], [5, 7, 8], [6, 8, 9]])

julia> W,FW = Lar.extrudeSimplicial(model, pattern)

julia> W
3×117 Array{Int64,2}:
 0  1  2  0  1  2  0  1  2   …   0   1   2   0   1   2   0   1   2   0   1   2
 0  0  0  1  1  1  2  2  2       2   2   2   0   0   0   1   1   1   2   2   2
 0  0  0  0  0  0  0  0  0      21  21  21  24  24  24  24  24  24  24  24  24

julia> FW
144-element Array{Array{Int64,1},1}:
 [1, 2, 4, 10]      
 [2, 4, 10, 11]     
 ⋮                  
 [96, 98, 99, 105]  
 [98, 99, 105, 107] 
 [99, 105, 107, 108]

julia> Plasm.view(W,FW)

Multidimensional grids of simplices

Generate a simplicial complex decomposition of a cubical grid of $d$-cuboids, where $d$ is the length of shape=[n_1, n_2, ..., n_d] array, so that shape defines the grid dimension $d$ and size $n_1 \times n_2 \times ... \times n_d$ as a $d$-dimensional array of cubes. Vertices (0-cells) of the grid have Int64 coordinates.

julia> Lar.simplexGrid([0]) # 0-dimensional simplicial complex
# output
([0], Array{Int64,1}[])

julia> V,EV = Lar.simplexGrid([1]) # 1-dimensional simplicial complex
# output
([0 1], Array{Int64,1}[[1, 2]])

julia> V,FV = Lar.simplexGrid([1,1]) # 2-dimensional simplicial complex
# output
([0 1 0 1; 0 0 1 1], Array{Int64,1}[[1, 2, 3], [2, 3, 4]])

julia> V,CV = Lar.simplexGrid([10,10,1]) # 3-dimensional simplicial complex
# output
([0 1 … 9 10; 0 0 … 10 10; 0 0 … 1 1], Array{Int64,1}[[1, 2, 12, 122], [2, 12, 122, 123], [12, 122, 123, 133], [2, 12, 13, 123], [12, 13, 123, 133], [13, 123, 133, 134], [2, 3, 13, 123], [3, 13, 123, 124], [13, 123, 124, 134], [3, 13, 14, 124]  …  [119, 229, 230, 240], [109, 119, 120, 230], [119, 120, 230, 240], [120, 230, 240, 241], [109, 110, 120, 230], [110, 120, 230, 231], [120, 230, 231, 241], [110, 120, 121, 231], [120, 121, 231, 241], [121, 231, 241, 242]])

julia> V
# output
3×242 Array{Int64,2}:
 0  1  2  3  4  5  6  7  8  9  10  0  1  2  3  …   1   2   3   4   5   6   7   8   9  10
 0  0  0  0  0  0  0  0  0  0   0  1  1  1  1     10  10  10  10  10  10  10  10  10  10
 0  0  0  0  0  0  0  0  0  0   0  0  0  0  0      1   1   1   1   1   1   1   1   1   1


julia> using Plasm

julia> hpc = Plasm.hpc_exploded(V,CV) # exploded visualization of the simplicial grid

julia> Plasm.view(hpc)

julia> V,HV = simplexGrid([1,1,1,1]) # 4-dim simplicial complex
# output
([0 1 … 0 1; 0 0 … 1 1; 0 0 … 1 1; 0 0 … 1 1], Array{Int64,1}[[1, 2, 3, 5, 9], [2, 3, 5, 9, 10], [3, 5, 9, 10, 11], [5, 9, 10, 11, 13], [2, 3, 5, 6, 10], [3, 5, 6, 10, 11], [5, 6, 10, 11, 13], [6, 10, 11, 13, 14], [3, 5, 6, 7, 11], [5, 6, 7, 11, 13]  …  [4, 6, 10, 11, 12], [6, 10, 11, 12, 14], [3, 4, 6, 7, 11], [4, 6, 7, 11, 12], [6, 7, 11, 12, 14], [7, 11, 12, 14, 15], [4, 6, 7, 8, 12], [6, 7, 8, 12, 14], [7, 8, 12, 14, 15], [8, 12, 14, 15, 16]])

Facets of multidimensional simplicial complexes

Compute the (d-1)-skeleton (set of facets) of a simplicial d-complex. Each of the $d+1$ facets of of a $d$-simplex $\sigma^d$, obtained by removing a vertex from $\sigma^d$, is a $(d-1)$-simplex.

julia> V,FV = Lar.simplexGrid([1,1]) # 2-dimensional complex
# output
([0 1 0 1; 0 0 1 1], Array{Int64,1}[[1, 2, 3], [2, 3, 4]])

julia> Plasm.view(V,FV)

julia> W,CW = Lar.extrudeSimplicial((V,FV), [1])
([0.0 1.0 … 0.0 1.0; 0.0 0.0 … 1.0 1.0; 0.0 0.0 … 1.0 1.0], 
Array{Int64,1}[[1,2,3,5],[2,3,5,6],[3,5,6,7],[2,3,4,6],[3,4,6,7],[4,6,7,8]])

julia> FW = Lar.simplexFacets(CW)
18-element Array{Any,1}:
[[1,3,5],[5,6,7],[3,5,7],[3,6,7],[4,6,7],[4,7,8],[4,6,8],
[6,7,8],[3,5,6],[2,3,5],[2,3,4],[3,4,7],[1,2,3],[2,4,6],[2,5,6],
[1,2,5],[2,3,6],[3,4,6]]

julia> Plasm.view(W,FW)

Main Interface

LinearAlgebraicRepresentation.simplexGridFunction
simplexGrid(shape::Array)::LAR

Generate a simplicial complex decomposition of a cubical grid of $d$-cuboids, where $d$ is the length of shape array. Vertices (0-cells) of the grid have Int64 coordinates.

Examples

julia> simplexGrid([0]) # 0-dimensional complex
# output
([0], Array{Int64,1}[])

julia> V,EV = simplexGrid([1]) # 1-dimensional complex
# output
([0 1], Array{Int64,1}[[1, 2]])

julia> V,FV = simplexGrid([1,1]) # 2-dimensional complex
# output
([0 1 0 1; 0 0 1 1], Array{Int64,1}[[1, 2, 3], [2, 3, 4]])

julia> V,CV = simplexGrid([10,10,1]) # 3-dimensional complex
# output
([0 1 … 9 10; 0 0 … 10 10; 0 0 … 1 1], Array{Int64,1}[[1, 2, 12, 122], [2, 12, 122, 123], [12, 122, 123, 133], [2, 12, 13, 123], [12, 13, 123, 133], [13, 123, 133, 134], [2, 3, 13, 123], [3, 13, 123, 124], [13, 123, 124, 134], [3, 13, 14, 124]  …  [119, 229, 230, 240], [109, 119, 120, 230], [119, 120, 230, 240], [120, 230, 240, 241], [109, 110, 120, 230], [110, 120, 230, 231], [120, 230, 231, 241], [110, 120, 121, 231], [120, 121, 231, 241], [121, 231, 241, 242]])

julia> V
# output
3×242 Array{Int64,2}:
 0  1  2  3  4  5  6  7  8  9  10  0  1  2  3  …   1   2   3   4   5   6   7   8   9  10
 0  0  0  0  0  0  0  0  0  0   0  1  1  1  1     10  10  10  10  10  10  10  10  10  10
 0  0  0  0  0  0  0  0  0  0   0  0  0  0  0      1   1   1   1   1   1   1   1   1   1

julia> using Plasm

julia> hpc = Plasm.lar2exploded_hpc(V,CV) # exploded visualization of the grid

julia> Plasm.view(hpc)

julia> V,HV = simplexGrid([1,1,1,1]) # 4-dim cellular complex from the 4D simplex
# output
([0 1 … 0 1; 0 0 … 1 1; 0 0 … 1 1; 0 0 … 1 1], Array{Int64,1}[[1, 2, 3, 5, 9], [2, 3, 5, 9, 10], [3, 5, 9, 10, 11], [5, 9, 10, 11, 13], [2, 3, 5, 6, 10], [3, 5, 6, 10, 11], [5, 6, 10, 11, 13], [6, 10, 11, 13, 14], [3, 5, 6, 7, 11], [5, 6, 7, 11, 13]  …  [4, 6, 10, 11, 12], [6, 10, 11, 12, 14], [3, 4, 6, 7, 11], [4, 6, 7, 11, 12], [6, 7, 11, 12, 14], [7, 11, 12, 14, 15], [4, 6, 7, 8, 12], [6, 7, 8, 12, 14], [7, 8, 12, 14, 15], [8, 12, 14, 15, 16]])
source
LinearAlgebraicRepresentation.simplexFacetsFunction
simplexFacets(simplices::Cells)::Cells

Compute the (d-1)-skeleton (unoriented set of facets) of a simplicial d-complex.

Example

julia> V,FV = Lar.simplexGrid([1,1]) # 2-dimensional complex
# output
([0 1 0 1; 0 0 1 1], Array{Int64,1}[[1, 2, 3], [2, 3, 4]])

julia> Plasm.view(V,FV)

julia> W,CW = Lar.extrudeSimplicial((V,FV), [1])
([0.0 1.0 … 0.0 1.0; 0.0 0.0 … 1.0 1.0; 0.0 0.0 … 1.0 1.0],
Array{Int64,1}[[1,2,3,5],[2,3,5,6],[3,5,6,7],[2,3,4,6],[3,4,6,7],[4,6,7,8]])

julia> FW = Lar.simplexFacets(CW)
18-element Array{Any,1}:
[[1,3,5],[5,6,7],[3,5,7],[3,6,7],[4,6,7],[4,7,8],[4,6,8],
[6,7,8],[3,5,6],[2,3,5],[2,3,4],[3,4,7],[1,2,3],[2,4,6],[2,5,6],
[1,2,5],[2,3,6],[3,4,6]]

julia> Plasm.view(W,FW)

Example

julia> V,(VV,EV,FV,CV) = Lar.cuboidGrid([3,3,3],true)

julia> TV = Lar.simplexFacets(CV)

julia> Plasm.view(V,TV)
source
LinearAlgebraicRepresentation.extrudeSimplicialFunction
extrudeSimplicial(model::LAR, pattern::Array)::LAR

Algorithm for multimensional extrusion of a simplicial complex. Can be applied to 0-, 1-, 2-, ... simplicial models, to get a 1-, 2-, 3-, .... model. The pattern Array is used to specify how to decompose the added dimension.

A model is a LAR model, i.e. a pair (vertices,cells) to be extruded, whereas pattern is an array of Int64, to be used as lateral measures of the extruded model. pattern elements are assumed as either solid or empty measures, according to their (+/-) sign.

Example

julia> V = [[0,0] [1,0] [2,0] [0,1] [1,1] [2,1] [0,2] [1,2] [2,2]];

julia> FV = [[1,2,4],[2,3,5],[3,5,6],[4,5,7],[5,7,8],[6,8,9]];

julia> pattern = repeat([1,2,-3],outer=4);

julia> model = (V,FV);

julia> W,FW = extrudeSimplicial(model, pattern);

julia> Plasm.view(W,FW)
source