Lefschetz Complex Functions

Simplicial Complexes

ConleyDynamics.create_simplicial_complexFunction
create_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.

Warning

Note that the labels all have to have the same character length!

source
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.

source
ConleyDynamics.create_simplicial_rectangleFunction
create_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.
source
ConleyDynamics.create_simplicial_delaunayFunction
create_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.

source
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.

source
ConleyDynamics.simplicial_torusFunction
sc = simplicial_torus(p::Int)

Create a triangulation of the two-dimensional torus.

The function returns a simplicial complex which represents a two-dimensional torus. The argument p specifies the characteristic of the underlying field. This triangulation is taken from Figure 6.4 in Munkres' book on Algebraic Topology. The boundary vertices are labeled as letters as in the book, the five center vertices are labeled by 1 through 5.

Examples

julia> println(homology(simplicial_torus(0)))
[1, 2, 1]

julia> println(homology(simplicial_torus(2)))
[1, 2, 1]

julia> println(homology(simplicial_torus(3)))
[1, 2, 1]
source
ConleyDynamics.simplicial_klein_bottleFunction
sc = simplicial_klein_bottle(p::Int)

Create a triangulation of the two-dimensional Klein bottle.

The function returns a simplicial complex which represents the two-dimensional Klein bottle. The argument p specifies the characteristic of the underlying field. This triangulation is taken from Figure 6.6 in Munkres' book on Algebraic Topology. The boundary vertices are labeled as letters as in the book, the five center vertices are labeled by 1 through 5.

Examples

julia> println(homology(simplicial_klein_bottle(0)))
[1, 1, 0]

julia> println(homology(simplicial_klein_bottle(2)))
[1, 2, 1]

julia> println(homology(simplicial_klein_bottle(3)))
[1, 1, 0]
source
ConleyDynamics.simplicial_projective_planeFunction
sc = simplicial_projective_plane(p::Int)

Create a triangulation of the projective plane.

The function returns a simplicial complex which represents the projective plane. The argument p specifies the characteristic of the underlying field. This triangulation is taken from Figure 6.6 in Munkres' book on Algebraic Topology. The boundary vertices are labeled as letters as in the book, the five center vertices are labeled by 1 through 5.

Examples

julia> println(homology(simplicial_projective_plane(0)))
[1, 0, 0]

julia> println(homology(simplicial_projective_plane(2)))
[1, 1, 1]

julia> println(homology(simplicial_projective_plane(3)))
[1, 0, 0]
source
ConleyDynamics.simplicial_torsion_spaceFunction
sc = simplicial_torsion_space(n::Int, p::Int)

Create a triangulation of a space with 1-dimensional torsion.

The function returns a simplicial complex which has the following integer homology groups:

  • In dimension 0 it is the group of integers.
  • In dimension 1 it is the integers modulo n.
  • In dimension 2 it is the trivial group.

In other words, the simplicial complex has nontrivial torsion in dimension 1. It is a triangulation of an n-gon, in which all boundary edges are oriented counterclockwise, and all of these edges are identified. The parameter p specifies the characteristic of the underlying field.

Examples

julia> println(homology(simplicial_torsion_space(6,2)))
[1, 1, 1]

julia> println(homology(simplicial_torsion_space(6,3)))
[1, 1, 1]

julia> println(homology(simplicial_torsion_space(6,5)))
[1, 0, 0]
source

Cubical Complexes

ConleyDynamics.create_cubical_complexFunction
create_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.

Warning

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
source
ConleyDynamics.create_cubical_rectangleFunction
create_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.

source
ConleyDynamics.create_cubical_boxFunction
create_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.

source
ConleyDynamics.cube_field_sizeFunction
cube_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)
source
ConleyDynamics.cube_informationFunction
cube_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 point
  • 1+pointdim:2*pointdim: Interval length in each dimension
  • 1+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
source
ConleyDynamics.cube_labelFunction
cube_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"
source
ConleyDynamics.get_cubical_coordsFunction
get_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. ```

source

Lefschetz Complex Creation

ConleyDynamics.create_lefschetz_gf2Function
create_lefschetz_gf2(defcellbnd)

Create a Lefschetz complex over GF(2) by specifying its essential cells and boundaries.

The input argument defcellbnd has to be a vector of vectors. Each entry defcellbnd[k] has to be of one of the following two forms:

  • [String, Int, String, String, ...]: The first String contains the label for the cell k, followed by its dimension in the second entry. The remaining entries are for the labels of the cells which make up the boundary.
  • [String, Int]: This shorther form is for cells with empty boundary. The first entry denotes the cell label, and the second its dimension.

The cells of the resulting Lefschetz complex correspond to the union of all occurring labels. Cell labels that only occur in the boundary specification are assumed to have empty boundary, and they do not have to be specified separately in the second form above. However, if their boundary is not empty, they have to be listed via the above first form as well.

Examples

julia> defcellbnd = [["A",0], ["a",1,"B","C"], ["b",1,"B","C"]];

julia> push!(defcellbnd, ["c",1,"B","C"]);

julia> push!(defcellbnd, ["alpha",2,"b","c"]);

julia> lc = create_lefschetz_gf2(defcellbnd);

julia> lc.labels
7-element Vector{String}:
 "A"
 "B"
 "C"
 "a"
 "b"
 "c"
 "alpha"

julia> homology(lc)
3-element Vector{Int64}:
 2
 1
 0
source
ConleyDynamics.lefschetz_subcomplexFunction
lefschetz_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.

source
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.

source
ConleyDynamics.lefschetz_closed_subcomplexFunction
lefschetz_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.

source
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.

source
ConleyDynamics.lefschetz_reductionFunction
lefschetz_reduction(lc::LefschetzComplex, redpairs::Vector{Vector{Int}})

Apply a sequence of elementary reductions to a Lefschetz complex.

The reduction pairs have to be specified in the argument redpairs. Each entry has to be a vector of length two which contains an elementary reduction pair in index form. In particular, the dimensions of the two cells in the pair have to differ by one, and once the pair is reached in the reduction sequence, one cell has to be a face of the other. The function returns a new Lefschetz complex, where all cells in redpairs have been removed.

source
lefschetz_reduction(lc::LefschetzComplex, redpairs::Vector{Vector{String}})

Apply a sequence of elementary reductions to a Lefschetz complex.

The reduction pairs have to be specified in the argument redpairs. Each entry has to be a vector of length two which contains an elementary reduction pair in label form. In particular, the dimensions of the two cells in the pair have to differ by one, and once the pair is reached in the reduction sequence, one cell has to be a face of the other. The function returns a new Lefschetz complex, where all cells in redpairs have been removed.

source
lefschetz_reduction(lc::LefschetzComplex, r1::Int, r2::Int)

Apply a single elementary reduction to a Lefschetz complex.

This method expects that the two cells r1 and r2 which form the reduction pair are given in index form. The function returns the reduced Lefschetz complex.

source
lefschetz_reduction(lc::LefschetzComplex, r1::String, r2::String)

Apply a single elementary reduction to a Lefschetz complex.

This method expects that the two cells r1 and r2 which form the reduction pair are given in label form. The function returns the reduced Lefschetz complex.

source
ConleyDynamics.permute_lefschetz_complexFunction
permute_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.

source

Lefschetz Complex Queries

ConleyDynamics.lefschetz_informationFunction
lefschetz_information(lc::LefschetzComplex)

Extract basic information about a Lefschetz complex.

The input argument lc contains the Lefschetz complex. The function returns the information in the form of a Dict{String,Any}. You can use the command keys to see the keyset of the return dictionary:

  • "Dimension": Dimension of the Lefschetz complex
  • "Coefficient field": Underlying coefficient field
  • "Euler characteristic": Euler characteristic of the complex
  • "Homology": Betti numbers of the Lefschetz complex
  • "Boundary sparsity": Sparsity percentage of the boundary matrix
  • "Number of cells": Total number of cells in the complex
  • "Cell counts by dim": Cell counts by dimension

In the last case, the dictionary entry is a vector of pairs (dimension, cell count).

source
ConleyDynamics.lefschetz_is_closedFunction
lefschetz_is_closed(lc::LefschetzComplex, subcomp::Vector{Int})

Determine whether a Lefschetz complex subset is closed.

source
lefschetz_is_closed(lc::LefschetzComplex, subcomp::Vector{String})

Determine whether a Lefschetz complex subset is closed.

source
ConleyDynamics.lefschetz_is_locally_closedFunction
lefschetz_is_locally_closed(lc::LefschetzComplex, subcomp::Vector{Int})

Determine whether a Lefschetz complex subset is locally closed.

source
lefschetz_is_locally_closed(lc::LefschetzComplex, subcomp::Vector{String})

Determine whether a Lefschetz complex subset is locally closed.

source

Topological Features

ConleyDynamics.lefschetz_boundaryFunction
lefschetz_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}.

source
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}.

source
ConleyDynamics.lefschetz_coboundaryFunction
lefschetz_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}.

source
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}.

source
ConleyDynamics.lefschetz_closureFunction
lefschetz_closure(lc::LefschetzComplex, subcomp::Vector{Int})

Compute the closure of a Lefschetz complex subset.

source
lefschetz_closure(lc::LefschetzComplex, subcomp::Vector{String})

Compute the closure of a Lefschetz complex subset.

source
ConleyDynamics.lefschetz_interiorFunction
lefschetz_interior(lc::LefschetzComplex, subcomp::Vector{Int})

Compute the interior of a Lefschetz complex subset.

source
lefschetz_interior(lc::LefschetzComplex, subcomp::Vector{String})

Compute the interior of a Lefschetz complex subset.

source
ConleyDynamics.lefschetz_topboundaryFunction
lefschetz_topboundary(lc::LefschetzComplex, subcomp::Vector{Int})

Compute the topological boundary of a Lefschetz complex subset.

In contrast to the algebraic boundary defined via the boundary operator, this function computes the boundary of the Lefschetz complex subset specified in subcomp if the Lefschetz complex is interpreted as a finite topological space. In other words, the topological boundary is the set difference of the closure and the interior of the subset.

source
lefschetz_topboundary(lc::LefschetzComplex, subcomp::Vector{String})

Compute the topological boundary of a Lefschetz complex subset.

In contrast to the algebraic boundary defined via the boundary operator, this function computes the boundary of the Lefschetz complex subset specified in subcomp if the Lefschetz complex is interpreted as a finite topological space. In other words, the topological boundary is the set difference of the closure and the interior of the subset.

source
ConleyDynamics.lefschetz_openhullFunction
lefschetz_openhull(lc::LefschetzComplex, subcomp::Vector{Int})

Compute the open hull of a Lefschetz complex subset.

source
lefschetz_openhull(lc::LefschetzComplex, subcomp::Vector{String})

Compute the open hull of a Lefschetz complex subset.

source
ConleyDynamics.lefschetz_lchullFunction
lefschetz_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.

source
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.

source
ConleyDynamics.lefschetz_clomo_pairFunction
lefschetz_clomopair(lc::LefschetzComplex, subcomp::Vector{Int})

Determine the closure-mouth-pair associated with a Lefschetz complex subset.

The function returns the pair (closure,mouth).

source
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).

source
ConleyDynamics.lefschetz_skeletonFunction
lefschetz_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.

source
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.

source
lefschetz_skeleton(lc::LefschetzComplex, skdim::Int)

Compute the skdim-dimensional skeleton of a Lefschetz complex.

The computed skeleton is for the full Lefschetz complex.

source
ConleyDynamics.manifold_boundaryFunction
manifold_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}.

source

Filters on Lefschetz Complexes

ConleyDynamics.create_random_filterFunction
create_random_filter(lc::LefschetzComplex)

Create a random injective filter on a Lefschetz complex.

The function creates a random injective filter on a Lefschetz complex. The filter is created by assigning integers to cell groups, increasing with dimension. Within each dimension the assignment is random, but all filter values of cells of dimension k are less than all filter values of cells with dimension k+1. The function returns the filter as Vector{Int}, with indices corresponding to the cell indices in the Lefschetz complex.

source
ConleyDynamics.filter_shallow_pairsFunction
filter_shallow_pairs(lc::LefschetzComplex, phi)

Find all shallow pairs for a filter.

This function finds all shallow pairs for the filter phi. These are face-coface pairs (x,y) whose dimensions differ by one, and such that y has the smallest filter value on the coboundary of x, and x has the largest filter value on the boundary of y.

If the filter is injective, these pairs give rise to a Forman vector field on the underlying Lefschetz complex. For noninjective filters this is not true in general.

source
ConleyDynamics.filter_induced_mvfFunction
filter_induced_mvf(lc::LefschetzComplex, phi)

Compute the multivector field induced by a filter.

This function returns the smallest multivector field which has the property that every shallow pair is contained in a multivector. For injective filters this is a Forman vector field, but in the noninjective case it can be a general multivector field.

source

Lefschetz Helper Functions

ConleyDynamics.lefschetz_gfp_conversionFunction
lcgfp = 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.

source
ConleyDynamics.lefschetz_filtrationFunction
lefschetz_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 subcomplex L_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)]
 []
source
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 subcomplex L_N, where N = 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)]
 []
source

Cell Subset Helper Functions

ConleyDynamics.convert_cellsFunction
convert_cells(lc::LefschetzComplex, cl::Vector{Int})

Convert cell list cl in the Lefschetz complex lc from index form to label form.

source
convert_cells(lc::LefschetzComplex, cl::Vector{String})

Convert cell list cl in the Lefschetz complex lc from label form to index form.

source
ConleyDynamics.convert_cellsubsetsFunction
convert_cellsubsets(lc::LefschetzComplex, clsub::Vector{Vector{Int}})

Convert CellSubsets clsub in the Lefschetz complex lc from index form to label form.

source
convert_cellsubsets(lc::LefschetzComplex, clsub::Vector{Vector{String}})

Convert CellSubsets clsub in the Lefschetz complex lc from label form to index form.

source

Coordinate Helper Functions

ConleyDynamics.convert_planar_coordinatesFunction
convert_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)
source
ConleyDynamics.convert_spatial_coordinatesFunction
convert_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)
source