Plotting Functions
Visualizing Simplicial Complexes (Luxor)
ConleyDynamics.plot_planar_simplicial — Function
plot_planar_simplicial(sc::AbstractComplex,
coords::Vector{<:Vector{<:Real}},
fname::String;
[mvf::CellSubsets=Vector{Vector{Int}}([]),]
[labeldir::Vector{<:Real}=Vector{Int}([]),]
[labeldis::Real=8,]
[hfac::Real=1.2,]
[vfac::Real=1.2,]
[sfac::Real=0,]
[pdim::Vector{Bool}=[true,true,true],]
[pv::Bool=false])Create an image of a planar simplicial complex, and if specified, a Forman vector field on it.
The vector coords contains coordinates for every one of the vertices of the simplicial complex sc. The image will be saved in the file with name fname, and the ending determines the image type. Accepted are .pdf, .svg, .png, and .eps.
If the optional mvf is specified and is a Forman vector field, then this Forman vector field is drawn as well. The optional vector labeldir contains directions for the vertex labels, and labeldis the distance from the vertex. The directions have to be reals between 0 and 4, with 0,1,2,3 corresponding to E,N,W,S. The optional constants hfac and vfac contain the horizontal and vertical scale vectors, while sfac describes a uniform scale. If sfac=0 the latter is automatically determined. The vector pdim specifies in which dimensions cells are drawn; the default shows vertices, edges, and triangles. Finally if one passes the argument pv=true, then in addition to saving the file a preview is displayed.
Examples
Suppose we have created a simplicial complex using the commands
sc, coords = create_simplicial_delaunay(300, 300, 30, 20)
fname = "sc_plot_test.pdf"Then the following code creates an image of the simplicial complex without labels, but with a preview:
plot_planar_simplicial(sc, coords, fname, pv=true)If we want to see the labels, we can use
ldir = fill(0.5, sc.ncells);
plot_planar_simplicial(sc, coords, fname, labeldir=ldir, labeldis=10, pv=true)This command puts all labels in the North-East direction at a distance of 10.
plot_planar_simplicial(ec::EuclideanComplex,
fname::String;
kwargs...)Create an image of a planar simplicial complex from an EuclideanComplex.
The vertex positions for every cell are taken directly from ec.coords[k], which stores the complete vertex coordinates for each cell k. This means the function works correctly even when some vertices are no longer present as cells in the complex (e.g. after subcomplex creation), as long as the higher-dimensional cells still carry their coordinates.
ConleyDynamics.plot_planar_simplicial_morse — Function
plot_planar_simplicial_morse(sc::AbstractComplex,
coords::Vector{<:Vector{<:Real}},
fname::String,
morsesets::CellSubsets;
[hfac::Real=1.2,]
[vfac::Real=1.2,]
[sfac::Real=0,]
[pdim::Vector{Bool}=[false,true,true],]
[pv::Bool=false,]
[ci::Bool=false])Create an image of a planar simplicial complex, together with Morse sets, or also selected multivectors.
The vector coords contains coordinates for every one of the vertices of the simplicial complex sc. The image will be saved in the file with name fname, and the ending determines the image type. Accepted are .pdf, .svg, .png, and .eps.
The vector morsesets contains a list of Morse sets, or more general, subsets of the simplicial complex. For every k, the set described by morsesets[k] will be shown in a distinct color.
The optional constants hfac and vfac contain the horizontal and vertical scale vectors for the margins, while sfac describes a uniform scale. If sfac=0 the latter is automatically determined. The vector pdim specifies in which dimensions cells are drawn; the default only shows edges and triangles. If one passes the argument pv=true, then in addition to saving the file a preview is displayed. Lastly, passing the argument ci=true will color code the Morse sets according to their respective Conley indices.
plot_planar_simplicial_morse(ec::EuclideanComplex,
fname::String,
morsesets::CellSubsets;
kwargs...)Create an image of a planar simplicial complex with Morse sets from an EuclideanComplex.
The vertex positions for every cell are taken directly from ec.coords[k], which stores the complete vertex coordinates for each cell k. This means the function works correctly even when some vertices are no longer present as cells in the complex (e.g. after subcomplex creation), as long as the higher-dimensional cells still carry their coordinates.
Visualizing Cubical Complexes (Luxor)
ConleyDynamics.plot_planar_cubical — Function
plot_planar_cubical(cc::AbstractComplex,
coords::Vector{<:Vector{<:Real}},
fname::String;
[hfac::Real=1.2,]
[vfac::Real=1.2,]
[cubefac::Real=0,]
[pdim::Vector{Bool}=[true,true,true],]
[pv::Bool=false])Create an image of a planar cubical complex.
The vector coords contains coordinates for every one of the vertices of the cubical complex cc. The image will be saved in the file with name fname, and the ending determines the image type. Accepted are .pdf, .svg, .png, and .eps. The optional constants hfac and vfac contain the horizontal and vertical scale vectors. The optional argument cubefac specifies the side length of an elementary cube for plotting, and it will be automatically determined otherwise. The vector pdim specifies which cell dimensions should be plotted, with pdim[k] representing dimension k-1. Finally if one passes the argument pv=true, then in addition to saving the file a preview is displayed.
Examples
Suppose we have created a cubical complex using the commands
cubes = ["00.11", "01.01", "02.10", "11.10", "11.01", "22.00"]
coords = [[0,0],[0,1],[0,2],[1,0],[1,1],[1,2],[2,1],[2,2]]
cc = create_cubical_complex(cubes)
fname = "cc_plot_test.pdf"Then the following code creates an image of the simplicial complex without labels, but with a preview:
plot_planar_cubical(cc, coords, fname, pv=true)If one only wants to plot the edges in the complex, but not the vertices or rectangles, then one can use:
plot_planar_cubical(cc, coords, fname, pv=true, pdim=[false,true,false])plot_planar_cubical(cc::AbstractComplex,
fname::String;
[hfac::Real=1.2,]
[vfac::Real=1.2,]
[cubefac::Real=0,]
[pdim::Vector{Bool}=[true,true,true],]
[pv::Bool=false])Create an image of a planar cubical complex.
This is an alternative method which does not require the specification of the vertex coordinates. They will be taken from the cube vertex labels.
plot_planar_cubical(ec::EuclideanComplex,
fname::String;
kwargs...)Create an image of a planar cubical complex from an EuclideanComplex.
The vertex positions for every cell are taken directly from ec.coords[k], which stores the complete vertex coordinates for each cell k. This means the function works correctly even when some vertices are no longer present as cells in the complex (e.g. after subcomplex creation), as long as the higher-dimensional cells still carry their coordinates.
ConleyDynamics.plot_planar_cubical_morse — Function
plot_planar_cubical_morse(cc::AbstractComplex,
coords::Vector{<:Vector{<:Real}},
fname::String,
morsesets::CellSubsets;
[hfac::Real=1.2,]
[vfac::Real=1.2,]
[cubefac::Real=0,]
[pdim::Vector{Bool}=[false,true,true],]
[pv::Bool=false])Create an image of a planar cubical complex, together with Morse sets, or also selected multivectors.
The vector coords contains coordinates for every one of the vertices of the cubical complex cc. The image will be saved in the file with name fname, and the ending determines the image type. Accepted are .pdf, .svg, .png, and .eps.
The vector morsesets contains a list of Morse sets, or more general, subsets of the cubical complex. For every k, the set described by morsesets[k] will be shown in a distinct color.
The optional constants hfac and vfac contain the horizontal and vertical scale vectors for the margins, while cubefac describes a uniform scale. If cubefac=0 the latter is automatically determined. The vector pdim specifies in which dimensions cells are drawn; the default only shows edges and squares. Finally if one passes the argument pv=true, then in addition to saving the file a preview is displayed.
plot_planar_cubical_morse(cc::AbstractComplex,
fname::String,
morsesets::CellSubsets;
[hfac::Real=1.2,]
[vfac::Real=1.2,]
[cubefac::Real=0,]
[pdim::Vector{Bool}=[false,true,true],]
[pv::Bool=false])Create an image of a planar cubical complex, together with Morse sets, or also selected multivectors.
This is an alternative method which does not require the specification of the vertex coordinates. They will be taken from the cube vertex labels.
plot_planar_cubical_morse(ec::EuclideanComplex,
fname::String,
morsesets::CellSubsets;
kwargs...)Create an image of a planar cubical complex with Morse sets from an EuclideanComplex.
The vertex positions for every cell are taken directly from ec.coords[k], which stores the complete vertex coordinates for each cell k. This means the function works correctly even when some vertices are no longer present as cells in the complex (e.g. after subcomplex creation), as long as the higher-dimensional cells still carry their coordinates.
Plot Data Types (Plots.jl)
ConleyDynamics.SimplicialComplexPlot — Type
SimplicialComplexPlotWrapper type for plotting a planar simplicial complex via Plots.jl. Constructed automatically by plot_simplicial; can also be passed directly to plot(::SimplicialComplexPlot).
ConleyDynamics.CubicalComplexPlot — Type
CubicalComplexPlotWrapper type for plotting a planar cubical complex via Plots.jl.
ConleyDynamics.SimplicialMorsePlot — Type
SimplicialMorsePlotWrapper type for plotting a planar simplicial complex with Morse sets via Plots.jl.
ConleyDynamics.CubicalMorsePlot — Type
CubicalMorsePlotWrapper type for plotting a planar cubical complex with Morse sets via Plots.jl.
ConleyDynamics.MVFPlot — Type
MVFPlotWrapper type for plotting a multivector field on a planar complex via Plots.jl. Each multivector is rendered as a stadium (pill) shaped region.
ConleyDynamics.MVRegionPlot — Type
MVRegionPlotWrapper type for plotting a single multivector on a planar complex via Plots.jl. The multivector is rendered as an inflated convex hull per maximal cell.
Visualizing Simplicial Complexes (Plots.jl)
ConleyDynamics.plot_simplicial — Function
plot_simplicial(ec::EuclideanComplex; kwargs...) -> Plots.PlotPlot a planar simplicial complex. Requires using Plots before using ConleyDynamics.
Optional keyword arguments:
mvf::CellSubsets=[]: Forman vector field to overlay (arrows + critical cells)labeldir::Vector{<:Real}=[]: label directions per cell (0=E,1=N,2=W,3=S)labeldis::Real=0.05: label offset in data-space unitspdim::Vector{Bool}=[true,true,true]: which dimensions (0,1,2) to draw
ConleyDynamics.plot_simplicial_morse — Function
plot_simplicial_morse(ec::EuclideanComplex, morsesets::CellSubsets;
kwargs...) -> Plots.PlotPlot a planar simplicial complex with Morse sets highlighted. Requires using Plots before using ConleyDynamics.
Optional keyword arguments:
pdim::Vector{Bool}=[false,true,true]: which dimensions (0,1,2) to drawci::Bool=false: color Morse sets by their Conley indexaddcritical::Bool=false: when true, cells absent from morsesets are shown as implicit singletons
ConleyDynamics.plot_simplicial_mvf — Function
plot_simplicial_mvf(ec::EuclideanComplex, mvf::CellSubsets;
kwargs...) -> Plots.PlotPlot a planar simplicial complex with multivector field regions shaded. Each multivector is rendered as a single colored region, inset from cell boundaries so adjacent multivectors are visually distinct. Requires using Plots before using ConleyDynamics.
Optional keyword arguments:
pdim::Vector{Bool}=[true,true,true]: which dimensions to show in backgroundtubefac::Real=0.05: tube half-width as fraction of average edge lengthmvfcolor::Union{String,Vector{String}}="darkorange": color(s) for multivector regions; a single string applies to all, a vector cycles through multivectors in ordermvfalpha::Real=0.3: fill opacity (0.0 = transparent, 1.0 = opaque)addcritical::Bool=true: when true, cells absent from the MVF are shown as implicit singletonsci::Bool=false: color each multivector by its Conley index instead of usingmvfcolor
ConleyDynamics.plot_simplicial_mv — Function
plot_simplicial_mv(ec::EuclideanComplex, mv; kwargs...) -> Plots.PlotPlot a single multivector on a planar simplicial complex background. The multivector is rendered as an inflated convex hull of barycenters, one hull per maximal cell of the multivector. Requires using Plots before using ConleyDynamics.
mv may be a Vector{Int} (cell indices) or Vector{String} (cell labels).
Optional keyword arguments:
pdim::Vector{Bool}=[true,true,true]: which background dimensions to drawtubefac::Real=0.05: inflation radius as fraction of average edge lengthmvfcolor::String="darkorange": X11 color name for the regionmvfalpha::Real=0.3: fill opacity (0.0 = transparent, 1.0 = opaque)
Visualizing Cubical Complexes (Plots.jl)
ConleyDynamics.plot_cubical — Function
plot_cubical(ec::EuclideanComplex; kwargs...) -> Plots.PlotPlot a planar cubical complex. Requires using Plots before using ConleyDynamics.
Optional keyword arguments:
mvf::CellSubsets=[]: multivector field to overlay (arrows + critical-cell dots)pdim::Vector{Bool}=[true,true,true]: which dimensions (0,1,2) to draw
ConleyDynamics.plot_cubical_morse — Function
plot_cubical_morse(ec::EuclideanComplex, morsesets::CellSubsets;
kwargs...) -> Plots.PlotPlot a planar cubical complex with Morse sets highlighted. Requires using Plots before using ConleyDynamics.
Optional keyword arguments:
pdim::Vector{Bool}=[false,true,true]: which dimensions (0,1,2) to drawci::Bool=false: color Morse sets by their Conley indexaddcritical::Bool=false: when true, cells absent from morsesets are shown as implicit singletons
ConleyDynamics.plot_cubical_mvf — Function
plot_cubical_mvf(ec::EuclideanComplex, mvf::CellSubsets;
kwargs...) -> Plots.PlotPlot a planar cubical complex with multivector field regions shaded. Requires using Plots before using ConleyDynamics.
Optional keyword arguments:
pdim::Vector{Bool}=[true,true,true]: which dimensions to show in backgroundtubefac::Real=0.05: tube half-width as fraction of average edge lengthmvfcolor::Union{String,Vector{String}}="darkorange": color(s) for multivector regionsmvfalpha::Real=0.3: fill opacity (0.0 = transparent, 1.0 = opaque)addcritical::Bool=true: when true, cells absent from the MVF are shown as implicit singletonsci::Bool=false: color each multivector by its Conley index instead of usingmvfcolor
ConleyDynamics.plot_cubical_mv — Function
plot_cubical_mv(ec::EuclideanComplex, mv; kwargs...) -> Plots.PlotPlot a single multivector on a planar cubical complex background. See plot_simplicial_mv for details.