Plotting Functions
Visualizing Simplicial Complexes
ConleyDynamics.plot_planar_simplicial
— Functionplot_planar_simplicial(sc::LefschetzComplex,
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.
ConleyDynamics.plot_planar_simplicial_morse
— Functionplot_planar_simplicial_morse(sc::LefschetzComplex,
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])
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. Finally if one passes the argument pv=true
, then in addition to saving the file a preview is displayed.
Visualizing Cubical Complexes
ConleyDynamics.plot_planar_cubical
— Functionplot_planar_cubical(cc::LefschetzComplex,
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::LefschetzComplex,
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.
ConleyDynamics.plot_planar_cubical_morse
— Functionplot_planar_cubical_morse(cc::LefschetzComplex,
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::LefschetzComplex,
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.