Conley Theory Functions
Multivector Fields
ConleyDynamics.mvf_information
— Functionmvf_information(lc::LefschetzComplex, mvf::CellSubsets)
Extract basic information about a multivector field.
The input argument lc
contains the Lefschetz complex, and mvf
describes the multivector field. 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:
"N mv"
: Number of multivectors"N critical"
: Number of critcal multivectors"N regular"
: Number of regular multivectors"Lengths critical"
: Length distribution of critical multivectors"Lengths regular"
: Length distribution of regular multivectors
In the last two cases, the dictionary entries are vectors of pairs (length,frequency)
, where each pair indicates that there are frequency
multivectors of length length
.
ConleyDynamics.create_mvf_hull
— Functioncreate_mvf_hull(lc::LefschetzComplex, mvfbase::Vector{Vector{Int}})
Create the smallest multivector field containing the given sets.
The resulting multivector field has the property that every set of the form mvfbase[k]
is contained in a minimal multivector. Notice that these sets do not have to be disjoint, and that not even their locally closed hulls have to be disjoint. In the latter case, this leads to two such sets having to be contained in the same multivector. If the sets in mvfbase
are poorly chosen, one might end up with extremely large multivectors due to the above potential merging of locally closed hulls.
create_mvf_hull(lc::LefschetzComplex, mvfbase::Vector{Vector{String}})
Create the smallest multivector field containing the given sets.
The resulting multivector field has the property that every set of the form mvfbase[k]
is contained in a minimal multivector. Notice that these sets do not have to be disjoint, and that not even their locally closed hulls have to be disjoint. In the latter case, this leads to two such sets having to be contained in the same multivector. If the sets in mvfbase
are poorly chosen, one might end up with extremely large multivectors due to the above potential merging of locally closed hulls.
ConleyDynamics.create_planar_mvf
— Functioncreate_planar_mvf(lc::LefschetzComplex, coords::Vector{Vector{Float64}}, vf)
Create a planar multivector field from a regular vector field.
The function expects a planar Lefschetz complex lc
and a coordinate vector coords
of coordinates for all the 0-dimensional cells in the complex. Moreover, the underlying vector field is specified by the function vf(z::Vector{Float64})::Vector{Float64}
, where both the input and output vectors have length two. The function create_planar_mvf
returns a multivector field mvf
on lc
, which can then be further analyzed using for example the function connection_matrix
.
The input data lc
and coords
can be generated using one of the following methods:
create_cubical_rectangle
create_simplicial_rectangle
create_simplicial_delaunay
In each case, the provided coordinate vector can be transformed to the correct bounding box values using the conversion function convert_planar_coordinates
.
Example 1
Suppose we define a sample vector field using the commands
function samplevf(x::Vector{Float64})
#
# Sample vector field with nontrivial Morse decomposition
#
x1, x2 = x
y1 = x1 * (1.0 - x1*x1 - 3.0*x2*x2)
y2 = x2 * (1.0 - 3.0*x1*x1 - x2*x2)
return [y1, y2]
end
One first creates a triangulation of the enclosing box, which in this case is given by [-2,2] x [-2,2]
using the commands
n = 21
lc, coords = create_simplicial_rectangle(n,n);
coordsN = convert_planar_coordinates(coords,[-2.0,-2.0],[2.0,2.0]);
The multivector field is then generated using
mvf = create_planar_mvf(lc,coordsN,samplevf);
and the commands
cm = connection_matrix(lc, mvf);
cm.conley
full_from_sparse(cm.matrix)
finally show that this vector field gives rise to a Morse decomposition with nine Morse sets, and twelve connecting orbits. Using the commands
fname = "morse_test.pdf"
plot_planar_simplicial_morse(lc, coordsN, fname, cm.morse, pv=true)
these Morse sets can be visualized. The image will be saved in fname
.
Example 2
An example with periodic orbits can be generated using the vector field
function samplevf2(x::Vector{Float64})
#
# Sample vector field with nontrivial Morse decomposition
#
x1, x2 = x
c0 = x1*x1 + x2*x2
c1 = (c0 - 4.0) * (c0 - 1.0)
y1 = -x2 + x1 * c1
y2 = x1 + x2 * c1
return [-y1, -y2]
end
The Morse decomposition can now be computed via
n2 = 51
lc2, coords2 = create_cubical_rectangle(n2,n2);
coords2N = convert_planar_coordinates(coords2,[-4.0,-4.0],[4.0,4.0]);
mvf2 = create_planar_mvf(lc2,coords2N,samplevf2);
cm2 = connection_matrix(lc2, mvf2);
cm2.conley
cm2.poset
full_from_sparse(cm2.matrix)
fname2 = "morse_test2.pdf"
plot_planar_cubical_morse(lc2, fname2, cm2.morse, pv=true)
In this case, one obtains three Morse sets: One is a stable equilibrium, one is an unstable periodic orbit, and the last is a stable periodic orbit.
ConleyDynamics.create_spatial_mvf
— Functioncreate_spatial_mvf(lc::LefschetzComplex, coords::Vector{Vector{Float64}}, vf)
Create a spatial multivector field from a regular vector field.
The function expects a three-dimensional Lefschetz complex lc
and a coordinate vector coords
of coordinates for all the 0-dimensional cells in the complex. Moreover, the underlying vector field is specified by the function vf(z::Vector{Float64})::Vector{Float64}
, where both the input and output vectors have length three. The function create_spatial_mvf
returns a multivector field mvf
on lc
, which can then be further analyzed using for example the function connection_matrix
.
The input data lc
and coords
has to be of one of the following two types:
lc
is a tetrahedral mesh of a region in three dimensions. In other words, the underlying Lefschetz complex is in fact a simplicial complex, and the vectorcoords
contains the vertex coordinates.lc
is a three-dimensional cubical complex, i.e., it is the closure of a collection of three-dimensional cubes in space. The vertex coordinates can be slightly perturbed from the original position in the cubical lattice, as long as the overall structure of the complex stays intact. In that case, the faces are interpreted as Bezier surfaces with straight edges.
Example 1
Suppose we define a sample vector field using the commands
function samplevf(x::Vector{Float64})
#
# Sample vector field with nontrivial Morse decomposition
#
x1, x2, x3 = x
y1 = x1 * (1.0 - x1*x1)
y2 = -x2
y3 = -x3
return [y1, y2, y3]
end
One first creates a cubical complex covering the interesting dynamics, say the trapping region [-1.5,1.5] x [-1,1] x [-1,1]
, using the commands
lc, coords = create_cubical_box(3,3,3);
coordsN = convert_spatial_coordinates(coords,[-1.5,-1.0,-1.0],[1.5,1.0,1.0]);
The multivector field is then generated using
mvf = create_spatial_mvf(lc,coordsN,samplevf);
and the commands
cm = connection_matrix(lc, mvf);
cm.conley
full_from_sparse(cm.matrix)
finally show that this vector field gives rise to a Morse decomposition with three Morse sets, and two connecting orbits.
ConleyDynamics.extract_multivectors
— Functionextract_multivectors(lc::LefschetzComplex, mvf::Vector{Vector{Int}},
scells::Vector{Int})
Extract all multivectors containing a provided selection of cells.
The function returns all multivectors which contain at least one of the cells in the input vector scells
. The return argument has type Vector{Vector{Int}}
.
extract_multivectors(lc::LefschetzComplex, mvf::Vector{Vector{String}},
scells::Vector{String})
Extract all multivectors containing a provided selection of cells.
The function returns all multivectors which contain at least one of the cells in the input vector scells
. The return argument has type Vector{Vector{String}}
.
ConleyDynamics.planar_nontransverse_edges
— Functionplanar_nontransverse_edges(lc::LefschetzComplex, coords::Vector{Vector{Float64}}, vf;
npts::Int=100)
Find all edges of a planar Lefschetz complex which are not flow transverse.
The Lefschetz complex is given in lc
, the coordinates of all vertices of the complex in coords
, and the vector field is specified in vf
. The optional parameter npts
determines how many points along an edge are evaluated for the transversality check. The function returns a list of nontransverse edges as Vector{Int}
, which contains the edge indices.
Conley Index Computations
ConleyDynamics.isoinvset_information
— Functionisoinvset_information(lc::LefschetzComplex, mvf::CellSubsets, iis::Cells)
Compute basic information about an isolated invariant set.
The input argument lc
contains the Lefschetz complex, and mvf
describes the multivector field. The isolated invariant set is specified in the argument iis
. The function returns the information in the form of a Dict{String,Any}
. The command keys
can be used to see the keyset of the return dictionary. These describe the following information:
"Conley index"
contains the Conley index of the isolated invariant set."N multivectors"
contains the number of multivectors in the isolated invariant set.
ConleyDynamics.conley_index
— Functionconley_index(lc::LefschetzComplex, subcomp::Vector{String})
Determine the Conley index of a Lefschetz complex subset.
The function raises an error if the subset subcomp
is not locally closed. The computations are performed over the field associated with the Lefschetz complex boundary matrix.
conley_index(lc::LefschetzComplex, subcomp::Vector{Int})
Determine the Conley index of a Lefschetz complex subset.
The function raises an error if the subset subcomp
is not locally closed. The computations are performed over the field associated with the Lefschetz complex boundary matrix.
ConleyDynamics.morse_sets
— Functionmorse_sets(lc::LefschetzComplex, mvf::CellSubsets; poset::Bool=false)
Find the nontrivial Morse sets of a multivector field on a Lefschetz complex.
The input argument lc
contains the Lefschetz complex, and mvf
describes the multivector field. The function returns the nontrivial Morse sets as a Vector{Vector{Int}}
. If the optional argument poset=true
is added, then the function returns both the Morse sets and the adjacency matrix of the Hasse diagram of the underlying poset.
ConleyDynamics.morse_interval
— Functionmorse_interval(lc::LefschetzComplex, mvf::CellSubsets,
ms::CellSubsets)
Find the isolated invariant set for a Morse set interval.
The input argument lc
contains the Lefschetz complex, and mvf
describes the multivector field. The collection of Morse sets are contained inms
. All of these sets should be Morse sets in the sense of being strongly connected components of the flow graph. (Nevertheless, this will be enforced in the function!) In other words, the sets in ms
should be determined using the function morse_sets
!
The function returns the smallest isolated invariant set which contains the Morse sets and their connections as a Vector{Int}
.
ConleyDynamics.restrict_dynamics
— Functionrestrict_dynamics(lc::LefschetzComplex, mvf::CellSubsets, lcsub::Cells)
Restrict a multivector field to a Lefschetz subcomplex.
For a given multivector field mvf
on a Lefschetz complex lc
, and a subcomplex which is given by the locally closed set represented by lcsub
, create the associated Lefschetz subcomplex lcreduced
and the induced multivector field mvfreduced
on the subcomplex. The multivectors of the new multivector field are the intersections of the original multivectors and the subcomplex.
ConleyDynamics.remove_exit_set
— Functionremove_exit_set(lc::LefschetzComplex, mvf::CellSubsets)
Exit set removal for a multivector field on a Lefschetz subcomplex.
It is assumed that the Lefschetz complex lc
is a topological manifold and that mvf
contains a multivector field that is created via either create_planar_mvf
or create_spatial_mvf
. The function identifies cells on the boundary at which the flows exits the region covered by the Lefschetz complex. If this exit set is closed, we have found an isolated invariant set and the function returns a Lefschetz complex lcr
restricted to it, as well as the restricted multivector field mvfr
. If the exit set is not closed, a warning is displayed and the function returns the restricted Lefschetz complex and multivector field obtained by removing the closure of the exit set. In the latter case, unexpected results might be obtained.
Connection Matrix Computation
ConleyDynamics.connection_matrix
— Functionconnection_matrix(lc::LefschetzComplex, mvf::CellSubsets;
[returnbasis::Bool])
Compute a connection matrix for the multivector field mvf
on the Lefschetz complex lc
over the field associated with the Lefschetz complex boundary matrix.
The function returns an object of type ConleyMorseCM
. If the optional argument returnbasis::Bool=true
is given, then the function also returns a dictionary which gives the basis for the connection matrix columns in terms of the original labels.
ConleyDynamics.cm_reduce!
— Functioncm_reduce!(matrix::SparseMatrix, psetvec::Vector{Int};
[returnbasis::Bool],[returntm::Bool])
Compute the connection matrix.
Assumes that matrix
is upper triangular and filtered according to psetvec
. Modifies the argument matrix
.
Return values:
cmatrix
: Connection matrixcmatrix_cols
: Columns of the connection matrix in the boundarybasisvecs
(optional): If the argumentreturnbasis=true
is given, this returns information about the computed basis. The k-th entry ofbasisvecs
is a vector containing the columns making up the k-th basis vector, which corresponds to columncmatrix_cols[k]
.tmatrix
(optional): If the argumentreturntm=true
is given in addition toreturnbasis=true
, then instead ofbasisvecs
the function returns the complete transformation matrix. In this case,basicvecs
is not returned.