Lefschetz Complex Functions
Simplicial Complexes
ConleyDynamics.create_simplicial_complex
— Functioncreate_simplicial_complex(labels::Vector{String},
simplices::Vector{Vector{Int}};
p::Int=2)
Initialize a Lefschetz complex from a simplicial complex. The complex is over the rationals if p=0
, and over GF(p)
if p>0
.
The vector labels
contains a label for every vertex, while simplices
contains all the highest-dimensional simplices necessary to define the simplicial complex. Every simplex is represented as a vector of Int
, with entries corresponding to the vertex indices.
Note that the labels all have to have the same character length!
create_simplicial_complex(labels::Vector{String},
simplices::Vector{Vector{String}};
p::Int=2)
Initialize a Lefschetz complex from a simplicial complex. The complex is over the rationals if p=0
, and over GF(p)
if p>0
.
The vector labels
contains a label for every vertex, while simplices
contains all the highest-dimensional simplices necessary to define the simplicial complex.
ConleyDynamics.create_simplicial_rectangle
— Functioncreate_simplicial_rectangle(nx::Int, ny::Int; p::Int=2)
Create a simplicial complex covering a rectangle in the plane. The complex is over the rationals if p=0
, and over GF(p)
if p>0
.
The rectangle is given by the subset [0,nx] x [0,ny]
of the plane. Each unit square is represented by four triangles, which meet in the center point of the square. Labels have the following meaning:
- The label
XXXYYYb
corresponds to the point(XXX, YYY)
. - The label
XXXYYYc
corresponds to(XXX + 1/2, YYY + 1/2)
.
The number of characters in XXX
and YYY
matches the number of digits of the larger number of nx
and ny
. The function returns the following objects:
- A simplicial complex
sc::LefschetzComplex
. - A vector
coords::Vector{Vector{Float64}}
of vertex coordinates.
ConleyDynamics.create_simplicial_delaunay
— Functioncreate_simplicial_delaunay(boxw::Real, boxh::Real, pdist::Real, attmpt::Int;
p::Int=2)
Create a planar Delaunay triangulation inside a box. The complex is over the rationals if p=0
, and over GF(p)
if p>0
.
The function selects a random sample of points inside the rectangular box [0,boxw] x [0,boxh]
, while trying to maintain a minimum distance of pdist
between the points. The argument attmpt
specifies the number of attempts when trying to add points. A standard value is 20, and larger values tend to fill holes better, but at the expense of runtime. From the random sample, the function then creates a Delaunay triangulation, and returns the following objects:
- A simplicial complex
sc::LefschetzComplex
. - A vector
coords::Vector{Vector{Float64}}
of vertex coordinates.
Note that the function does not provide a full triangulation of the given rectangle. Close to the boundary there will be gaps.
create_simplicial_delaunay(boxw::Real, boxh::Real, npoints::Int;
p::Int=2)
Create a planar Delaunay triangulation inside a box. The complex is over the rationals if p=0
, and over GF(p)
if p>0
.
The function selects a random sample of npoints
points inside the rectangular box [0,boxw] x [0,boxh]
. From the random sample, the function then creates a Delaunay triangulation, and returns the following objects:
- A simplicial complex
sc::LefschetzComplex
. - A vector
coords::Vector{Vector{Float64}}
of vertex coordinates.
Note that the function does not provide a full triangulation of the given rectangle. Close to the boundary there will be gaps.
Cubical Complexes
ConleyDynamics.create_cubical_complex
— Functioncreate_cubical_complex(cubes::Vector{String}; p::Int=2)
Initialize a Lefschetz complex from a cubical complex. The complex is over the rationals if p=0
, and over GF(p)
if p>0
.
The vector cubes
contains a list of all the highest-dimensional cubes necessary to define the cubical complex. Every cube is represented as a string as follows:
- d integers, which correspond to the coordinates of a point in d-dimensional Euclidean space
- a point
.
- d integers 0 or 1, which give the interval length in the respective dimension
The first d integers all have to occupy the same number of characters. In addition, if the occupied space is L characters for each coordinate, the coordinates only can take values from 0 to 10^L - 2. This is due to the fact that the boundary operator will add one to certain coordinates, and they still need to be representable withing the same L digits.
For example, the string 030600.101
corresponds to the point (3,6,0)
in three dimensions. The dimensions are 1, 0, and 1, and therefore this string corresponds to the cube [3,4] x [6] x [0,1]
. The same cube could have also been represented by 360.101
or by 003006000.101
.
Note that the labels all have to have the same format!
Example
julia> cubes = ["00.11", "01.01", "02.10", "11.10", "11.01", "22.00"];
julia> lc = create_cubical_complex(cubes);
julia> lc.ncells
17
julia> homology(lc)
3-element Vector{Int64}:
2
1
0
ConleyDynamics.create_cubical_rectangle
— Functioncreate_cubical_rectangle(nx::Int, ny::Int;
p::Int=2, randomize::Real=0.0)
Create a cubical complex covering a rectangle in the plane. The complex is over the rationals if p=0
, and over GF(p)
if p>0
.
The rectangle is given by the subset [0,nx] x [0,ny]
of the plane, and each unit square gives a two-dimensional cube in the resulting cubical complex. The function returns the following objects:
- A cubical complex
cc::LefschetzComplex
- A vector
coords::Vector{Vector{Float64}}
of vertex coordinates
If the optional parameter randomize
is assigned a positive real fraction r
less that 0.5, then the actual coordinates will be randomized. They are chosen uniformly from discs of radius r
centered at each vertex.
ConleyDynamics.create_cubical_box
— Functioncreate_cubical_box(nx::Int, ny::Int, nz::Int;
p::Int=2, randomize::Real=0.0)
Create a cubical complex covering a box in space. The complex is over the rationals if p=0
, and over GF(p)
if p>0
.
The box is given by the subset [0,nx] x [0,ny] x [0,nz]
of space, and each unit cube gives a three-dimensional cube in the resulting cubical complex. The function returns the following objects:
- A cubical complex
cc::LefschetzComplex
- A vector
coords::Vector{Vector{Float64}}
of vertex coordinates
If the optional parameter randomize
is assigned a positive real fraction r
less that 0.5, then the actual coordinates will be randomized. They are chosen uniformly from balls of radius r
centered at each vertex.
ConleyDynamics.cube_field_size
— Functioncube_field_size(cube::String)
Determine the field sizes of a given cube label.
The function returns the dimension of the ambient space in the first output parameter pointdim
, and the length of the individual coordinate fields in the second return variable pointlen
.
Example
julia> cube_field_size("011654003020.0110")
(4, 3)
ConleyDynamics.cube_information
— Functioncube_information(cube::String)
Compute a cube's coordinate information.
The function returns an integer vector with the cubes coordinate information. The return vector intinfo
contains in its components the following data:
1:pointdim
: Coordinates of the anchor point1+pointdim:2*pointdim
: Interval length in each dimension1+2*pointdim
: Dimension of the cube
Note that pointdim
equals the dimension of the points specifying the cube.
Example
julia> cube_information("011654003.011")
7-element Vector{Int64}:
11
654
3
0
1
1
2
ConleyDynamics.cube_label
— Functioncube_label(pointdim::Int, pointlen::Int, pointinfo::Vector{Int})
Create a label from a cube's coordinate information.
The dimension of the ambient Eucliden space is pointdim
, while the field length for each coordinate is pointlen
. The vector pointinfo
has to be of length at least two times pointdim
. The first pointdim
entries contain the coordinates of the anchor point, while the next pointdim
entries are either 0 or 1 depending on the size of the interval. For example, if poindim = 3
and pointinfo = [1,2,3,1,0,1]
, then we represent the cube in three-dimensional space given by [1,2] x [2] x [3 4]
.
Example
julia> cube_label(3,2,[10,23,5,1,1,0])
"102305.110"
ConleyDynamics.get_cubical_coords
— Functionget_cubical_coords(cc::LefschetzComplex)
Compute the vertex coordinates for a cubical complex.
The variable cc
has to contain a cubical complex, and the function returns a vector of coordinates for the vertices of the complex, that can then be used for plotting. ```
Lefschetz Complex Queries
ConleyDynamics.lefschetz_field
— Functionfieldstr = lefschetz_field(lc::LefschetzComplex)
Returns the Lefschetz complex coefficient field.
ConleyDynamics.lefschetz_is_closed
— Functionlefschetz_is_closed(lc::LefschetzComplex, subcomp::Vector{Int})
Determine whether a Lefschetz complex subset is closed.
lefschetz_is_closed(lc::LefschetzComplex, subcomp::Vector{String})
Determine whether a Lefschetz complex subset is closed.
ConleyDynamics.lefschetz_is_locally_closed
— Functionlefschetz_is_locally_closed(lc::LefschetzComplex, subcomp::Vector{Int})
Determine whether a Lefschetz complex subset is locally closed.
lefschetz_is_locally_closed(lc::LefschetzComplex, subcomp::Vector{String})
Determine whether a Lefschetz complex subset is locally closed.
Topological Features
ConleyDynamics.lefschetz_boundary
— Functionlefschetz_boundary(lc::LefschetzComplex, cellI::Int)
Compute the support of the boundary of a Lefschetz complex cell.
This method returns the boundary support as a Vector{Int}
.
lefschetz_boundary(lc::LefschetzComplex, cellS::String)
Compute the support of the boundary of a Lefschetz complex cell.
This method returns the boundary support as a Vector{String}
.
ConleyDynamics.lefschetz_coboundary
— Functionlefschetz_coboundary(lc::LefschetzComplex, cellI::Int)
Compute the support of the coboundary of a Lefschetz complex cell.
This method returns the boundary support as a Vector{Int}
.
lefschetz_coboundary(lc::LefschetzComplex, cellS::String)
Compute the support of the coboundary of a Lefschetz complex cell.
This method returns the boundary support as a Vector{String}
.
ConleyDynamics.lefschetz_closure
— Functionlefschetz_closure(lc::LefschetzComplex, subcomp::Vector{Int})
Compute the closure of a Lefschetz complex subset.
lefschetz_closure(lc::LefschetzComplex, subcomp::Vector{String})
Compute the closure of a Lefschetz complex subset.
ConleyDynamics.lefschetz_openhull
— Functionlefschetz_openhull(lc::LefschetzComplex, subcomp::Vector{Int})
Compute the open hull of a Lefschetz complex subset.
lefschetz_openhull(lc::LefschetzComplex, subcomp::Vector{String})
Compute the open hull of a Lefschetz complex subset.
ConleyDynamics.lefschetz_lchull
— Functionlefschetz_lchull(lc::LefschetzComplex, subcomp::Vector{Int})
Compute the locally closed hull of a Lefschetz complex subset.
The locally closed hull is the smallest locally closed set which contains the given cells. It is the intersection of the closure and the open hull.
lefschetz_lchull(lc::LefschetzComplex, subcomp::Vector{String})
Compute the locally closed hull of a Lefschetz complex subset.
The locally closed hull is the smallest locally closed set which contains the given cells. It is the intersection of the closure and the open hull.
ConleyDynamics.lefschetz_clomo_pair
— Functionlefschetz_clomopair(lc::LefschetzComplex, subcomp::Vector{Int})
Determine the closure-mouth-pair associated with a Lefschetz complex subset.
The function returns the pair (closure,mouth)
.
lefschetz_clomopair(lc::LefschetzComplex, subcomp::Vector{String})
Determine the closure-mouth-pair associated with a Lefschetz complex subset.
The function returns the pair (closure,mouth)
.
ConleyDynamics.lefschetz_skeleton
— Functionlefschetz_skeleton(lc::LefschetzComplex, subcomp::Vector{Int}, skdim::Int)
Compute the skdim
-dimensional skeleton of a Lefschetz complex subset.
The computed skeleton is for the closure of the subcomplex given by subcomp
.
lefschetz_skeleton(lc::LefschetzComplex, subcomp::Vector{String}, skdim::Int)
Compute the skdim
-dimensional skeleton of a Lefschetz complex subset.
The computed skeleton is for the closure of the subcomplex given by subcomp
.
lefschetz_skeleton(lc::LefschetzComplex, skdim::Int)
Compute the skdim
-dimensional skeleton of a Lefschetz complex.
The computed skeleton is for the full Lefschetz complex.
ConleyDynamics.manifold_boundary
— Functionmanifold_boundary(lc::LefschetzComplex)
Extract the manifold boundary from a Lefschetz complex.
The function expects a Lefschetz complex which represents a compact d-dimensional manifold with boundary. It returns a list of all cells which lie on the topological boundary of the manifold, in the form of a Vector{Int}
.
Lefschetz Subcomplexes
ConleyDynamics.lefschetz_subcomplex
— Functionlefschetz_subcomplex(lc::LefschetzComplex, subcomp::Vector{Int})
Extract a subcomplex from a Lefschetz complex. The subcomplex has to be locally closed, and is given by the collection of cells in subcomp
.
lefschetz_subcomplex(lc::LefschetzComplex, subcomp::Vector{String})
Extract a subcomplex from a Lefschetz complex. The subcomplex has to be locally closed, and is given by the collection of cells in subcomp
.
ConleyDynamics.lefschetz_closed_subcomplex
— Functionlefschetz_closed_subcomplex(lc::LefschetzComplex, subcomp::Vector{Int})
Extract a closed subcomplex from a Lefschetz complex. The subcomplex is the closure of the collection of cells given in subcomp
.
lefschetz_closed_subcomplex(lc::LefschetzComplex, subcomp::Vector{String})
Extract a closed subcomplex from a Lefschetz complex. The subcomplex is the closure of the collection of cells given in subcomp
.
ConleyDynamics.permute_lefschetz_complex
— Functionpermute_lefschetz_complex(lc::LefschetzComplex,
permutation::Vector{Int})
Permute the indices of a Lefschetz complex.
The vector permutation
contains a permutation of the indices for the given Lefschetz complex lc
. If no permutation is specified, or if the length of the vector is not correct, then a randomly generated one will be used. Note that the permutation has to respect the ordering of the cells by dimension, otherwise an error is raised. In other words, the permutation has to decompose into permutations within each dimension. This is automatically done if no permutation is explicitly specified.
Lefschetz Helper Functions
ConleyDynamics.lefschetz_gfp_conversion
— Functionlcgfp = lefschetz_gfp_conversion(lc::LefschetzComplex, p::Int)
Convert a Lefschetz complex to the same complex over a finite field.
It is expected that the boundary matrix of the given Lefschetz complex lc
is defined over the rationals, and that the target characteristic p
is a prime.
ConleyDynamics.lefschetz_filtration
— Functionlefschetz_filtration(lc::LefschetzComplex, fvalues::Vector{Int})
Compute a filtration on a Lefschetz subset.
The considered Lefschetz complex is given in lc
. The vector fvalues
assigns an integer between 0 and N to every cell in lc
. For every k the complex L_k
is given by the closure of all cells with values between 1 and k. The function returns the following variables:
lcsub
: The subcomplexL_N
fvalsub
: The filtration on the subcomplex with values 1,...,N
Example
julia> labels = ["A","B","C","D","E","F","G"];
julia> simplices = [["A","B","D"],["B","D","E"],["B","C","E"],["C","E","F"],["F","G"]];
julia> sc = create_simplicial_complex(labels,simplices);
julia> filtration = [0,0,0,0,0,0,0,1,1,0,1,2,0,4,2,4,0,5,3,7,6];
julia> lcsub, fvalsub = lefschetz_filtration(sc,filtration);
julia> phinf, phint = persistent_homology(lcsub,fvalsub);
julia> phinf
3-element Vector{Vector{Int64}}:
[1]
[]
[]
julia> phint
3-element Vector{Vector{Tuple{Int64, Int64}}}:
[]
[(1, 5), (2, 7), (4, 6)]
[]
lefschetz_filtration(lc::LefschetzComplex, strfilt::Vector{Vector{String}})
Compute a filtration on a Lefschetz subset.
The considered Lefschetz complex is given in lc
. The vector of string vectors strfilt
contains the necessary simplices to build the filtration. The list strfilt[k]
contains the simplices that are added at the k-th step, together with their closures. Thus, for every k the complex L_k
is given by the closure of all cells listed in strfilt[i]
for i
between 1 and k. The function returns the following variables:
lcsub
: The subcomplexL_N
, whereN = length(strfilt)
fvalsub
: The filtration on the subcomplex with values 1,...,N
Example
julia> labels = ["A","B","C","D","E","F","G"];
julia> simplices = [["A","B","D"],["B","D","E"],["B","C","E"],["C","E","F"],["F","G"]];
julia> sc = create_simplicial_complex(labels,simplices);
julia> strfiltration = [["AB","AD","BD"],["BE","DE"],["BCE"],["CF","EF"],["ABD"],["CEF"],["BDE"]];
julia> lcsub, fvalsub = lefschetz_filtration(sc, strfiltration);
julia> phinf, phint = persistent_homology(lcsub,fvalsub);
julia> phinf
3-element Vector{Vector{Int64}}:
[1]
[]
[]
julia> phint
3-element Vector{Vector{Tuple{Int64, Int64}}}:
[]
[(1, 5), (2, 7), (4, 6)]
[]
Cell Subset Helper Functions
ConleyDynamics.convert_cells
— Functionconvert_cells(lc::LefschetzComplex, cl::Vector{Int})
Convert cell list cl
in the Lefschetz complex lc
from index form to label form.
convert_cells(lc::LefschetzComplex, cl::Vector{String})
Convert cell list cl
in the Lefschetz complex lc
from label form to index form.
ConleyDynamics.convert_cellsubsets
— Functionconvert_cellsubsets(lc::LefschetzComplex, clsub::Vector{Vector{Int}})
Convert CellSubsets clsub
in the Lefschetz complex lc
from index form to label form.
convert_cellsubsets(lc::LefschetzComplex, clsub::Vector{Vector{String}})
Convert CellSubsets clsub
in the Lefschetz complex lc
from label form to index form.
Coordinate Helper Functions
ConleyDynamics.convert_planar_coordinates
— Functionconvert_planar_coordinates(coords::Vector{Vector{Float64}},
p0::Vector{Float64},
p1::Vector{Float64})
Convert a given collection of planar coordinates.
The vector coords
contains pairs of coordinates, which are then transformed to fit into the box with vertices p0 = (p0x,p0y)
and p1 = (p1x,p1y)
. It is assumed that p0
denotes the lower left box corner, while p1
is the upper right corner. The function shifts and scales the coordinates in such a way that every side of the box contains at least one point. Upon completion, it returns a new coordinate vector coordsNew
.
More precisely, if the x-coordinates are spanning the interval [xmin,xmax]
and the y-coordinates span [ymin,ymax]
, then the point (x,y)
is transformed to (xn,yn)
with:
xn = p0x + (p1x-p0x) * (x-cxmin) / (cxmax-cxmin)
yn = p0y + (p1y-p0y) * (y-cymin) / (cymax-cymin)
ConleyDynamics.convert_spatial_coordinates
— Functionconvert_spatial_coordinates(coords::Vector{Vector{Float64}},
p0::Vector{Float64},
p1::Vector{Float64})
Convert a given collection of spatial coordinates.
The vector coords
contains triples of coordinates, which are then transformed to fit into the box with vertices p0 = (p0x,p0y,p0z)
and p1 = (p1x,p1y,p1z)
. It is assumed that each coordinate of p0
is strictly smaller than the corresponding coordinate of p1
. The function shifts and scales the coordinates in such a way that every side of the box contains at least one point. Upon completion, it returns a new coordinate vector coordsNew
.
More precisely, if the x-coordinates are spanning the interval [xmin,xmax]
, the y-coordinates span [ymin,ymax]
, and the z-coordinates span [zmin,zmax]
, then the point (x,y,z)
is transformed to (xn,yn,zn)
with:
xn = p0x + (p1x-p0x) * (x-cxmin) / (cxmax-cxmin)
yn = p0y + (p1y-p0y) * (y-cymin) / (cymax-cymin)
zn = p0z + (p1z-p0z) * (z-czmin) / (czmax-czmin)