Advanced Usage
- Advanced Usage
- High level overview
- Loading trajectories
- Defining Boundaries
- Binning Algorithms
- Reinjection algorithms
- Computing the main results
- Working with results
- Visualization
- References
High level overview
Use of this package proceeds as follows.
- Create
Trajectoriesusing your trajectory data - Select the computational domain (
Boundary). Data inside the boundary are considered "interior" and data outside the boundary are considered "exterior" or "in nirvana"[1] [2]. Use theAutoBoundaryfunction to attempt to select a reasonable boundary automatically. In 2D, theAutoBoundary2Dfunction is also available, which may find a tighter boundary thanAutoBoundary. - Select a
BinningAlgorithmand set any parameters it has. - Optionally, select a
ReinjectionAlgorithmto handle the behavior of trajectories points from nirvana to the interior. - Run
ulam_method(traj, binner; reinj_algo). - Inspect and process the result with
P_closed,bins,membershipand other functions. - Visualize the result with
vizandviz!.
Loading trajectories
UlamMethod.Trajectories — Typestruct Trajectories{Dim}A container for trajectory data of dimension Dim.
Fields
x0: Initial coordinates, should be a matrix of sizeDimxn_traj.xT: Final coordinates, should be a matrix of sizeDimxn_traj.
Constructors
Trajectories(x0, xT)Construct Trajectories directly from data.
Trajectories(dim, n_traj; corners = nothing, mu_sigma = nothing)Construct n_traj random Trajectories of dimension dim for testing.
The x0 will lie inside the box defined by corners = ((xmin, ymin, zmin, ...), (xmax, ymax, zmax, ...)) which defaults to the unit box when not provided.
The xT are displaced from the x0 by normal distributions along each dimension with means and variances defined by mu_sigma = [(mu1, sigma1), (mu2, sigma2), ...]. If not provided, default to mu = 1.0, sigma = 1.0
Defining Boundaries
UlamMethod.Boundary — Typestruct Boundary{Dim, M, CRS}A container type for the computational domain within which Ulam's method is applied.
The boundary has dimension Dim embedded in manifold M with coordinate reference system CRS.
The boundary is partitioned according to an BinningAlgorithm.
Fields
boundary: APolytopedefining the boundary.
Methods
points(boundary)Calculate a Dim x N matrix of boundary vertices.
1D Constructors
Boundary(x_start, x_end)The boundary is a continuous line segment Segment between x_start and x_end.
2D Constructors
Boundary(verts)The boundary is a closed polygon Ngon with vertices verts, which should be be a vector of (x, y) coordinates or a 2 x N matrix.
Boundary(xmin, xmax, ymin, ymax)A convenience constructor is also provided for the case of a rectangular boundary.
≥3D Constructors
Boundary(corner_min, corner_max)The boundary is a HyperRectangle with minimum and maximum vertices corner_min, corner_max where the corners are NTuples, (x_min, y_min, z_min, ...), (x_max, y_max, z_max, ...).
Automatric Boundary Construction
AutoBoundary(traj; nirvana = 0.10)Construct a rectangular boundary in arbitrary dimensions based on the Trajectories in traj.
The boundary is placed such that a fraction nirvana of the data is in nirvana. For example, with the default value of 0.10, roughly 10% of all of the datapoints (split between x0 and xT) will be outside the boundary with a roughly equal amount on each side.
The shape of the boundary can be further controlled by providing nirvana as a vector of length Dim of tuples of the form (min, max) such that a fraction min (respectively max) will be "below" (respectively "above") the boundary along each dimension.
AutoBoundary2D(traj; nirvana = 0.10, tol = 0.001)Construct a tight boundary in two dimensions based on the Trajectories in traj.
The boundary in question is formed by scaling the convex hull of the trajectory data until a fraction nirvana of the data are bouside the boundary to within tol percent.
Binning Algorithms
UlamMethod.BinningAlgorithm — Typeabstract type BinningAlgorithm{Dim}An abstract type for binning algorithms of dimension Dim.
Implementation
Each subtype should have a fields:
boundary::Boundary{Dim, CRS}`bins::Bins{Dim, CRS}idx2pos::Vector{Union{Int64, Nothing}}
Each subtype should implement the following:
- A function
membership(data::Matrix, binner)which takes aDim x N_pointsmatrix of points and returns a vectormembwherememb[i] = jifdata[:,i]is insidebinner.bins[j]andmemb[i] = nothingifdata[:,i]is not inside any bin. - A function
membership(traj::Trajectories, binner)which returns(memb_x0, memb_xT). In many cases this will be given by(membership(traj.x0, binner), membership(traj.xT, binner)).
Methods
points(binner)Return a vector of raw vertices for each bin.
1D Algorithms
UlamMethod.LineBinner — Typestruct LineBinner{M, CRS}Bin a one dimensional Segment (line segement) with nbins equally-spaced bins.
Fields
boundary: ABoundaryobject.bins: ABinsobject.idx2pos: A vector such thatidx2pos[i]gives the position (in bins) of the bin initially (before removing dataless and disconnected bins) labelledi.idx2pos[i] == nothingif this bin was removed.
Constructor
LineBinner(nbins, boundary; hardclip = true)2D Algorithms
UlamMethod.RectangleBinner — Typestruct RectangleBinner{M, CRS}Bin a two dimensional Polygon with a tight covering of rectangles. The rectangles are chosen to be as close as possible to squares and so the final number of bins may be slightly different than the number requested.
Fields
boundary: ABoundaryobject.bins: ABinsobject.idx2pos: A vector such thatidx2pos[i]gives the position (in bins) of the bin initially (before removing dataless and disconnected bins) labelledi.idx2pos[i] == nothingif this bin was removed.
Constructor
RectangleBinner(nbins, boundary; hardclip = true)nbins can be:
- An
Integer, in which case aRectangleBinnerwill be constructed attempting to distribute boxes proportionally acrossxand ysuch that the total number of boxes is as close as possible tonbins`. - An
NTuple{2, Integer}, in which case dimensionireceivesnbins[i]bins. - An
NTuple{2, AbstractVector}, in which casenbins[i]defines the edges of the bins in dimensioni.
UlamMethod.TriangleBinner — Typestruct TriangleBinner{M, CRS}Bin a two dimensional Polygon with a covering of equilateral triangles.
The final number of bins may be slightly different than the number requested.
Fields
boundary: ABoundaryobject.bins: ABinsobject.idx2pos: A vector such thatidx2pos[i]gives the position (in bins) of the bin initially (before removing dataless and disconnected bins) labelledi.idx2pos[i] == nothingif this bin was removed.
Constructor
TriangleBinner(nbins, boundary; hardclip = true)UlamMethod.HexagonBinner — Typestruct HexagonBinner{M, CRS}Bin a two dimensional Polygon with a covering of regular hexagons.
The final number of bins may be slightly different than the number requested.
Fields
boundary: ABoundaryobject.bins: ABinsobject.idx2pos: A vector such thatidx2pos[i]gives the position (in bins) of the bin initially (before removing dataless and disconnected bins) labelledi.idx2pos[i] == nothingif this bin was removed.
Constructor
HexagonBinner(nbins, boundary; hardclip = true)UlamMethod.VoronoiBinner — Typestruct VoronoiBinner{Dim, M, CRS}Bin Dim-dimensional (where Dim ∈ (1, 2)) dataset based on a Voronoi tesellation of points generated by a k-means clustering of initial trajectory data (i.e. Trajectories.x0).
Fields
boundary: ABoundaryobject.bins: ABinsobject.idx2pos: A vector such thatidx2pos[i]gives the position (in bins) of the bin initially (before removing dataless and disconnected bins) labelledi.idx2pos[i] == nothingif this bin was removed.
Constructor
VoronoiBinner(nbins, boundary, traj; hardclip = true)≥3D Algorithms
UlamMethod.HyperRectangleBinner — Typestruct HyperRectangleBinner{Dim, M, CRS}Bin Dim-dimensional (where Dim ≥ 3) dataset based on a tight covering of identical HyperRectangles.
The hyperrectangles are chosen to be as close as possible to hypercubes.
The HyperRectangleBinner does not have hard clipping functionality. However, bins with no data are still removed in the main Ulam's method calculation.
The final number of bins may be different (larger) than the number requested.
Fields
boundary: ABoundaryobject.bins: ABinsobject.idx2pos: A vector such thatidx2pos[i]gives the position (in bins) of the bin initially (before removing dataless and disconnected bins) labelledi.idx2pos[i] == nothingif this bin was removed.ranges: A vector such thatranges[i]is theAbstractRangegiving the gridpoints of the bins along dimensioni.
Constructor
HyperRectangleBinner(nbins, boundary)nbins can be:
- An
Integer, in which case aHyperRectangleBinnerwill be constructed attempting to distribute boxes proportionally across dimensions such that the total number of boxes is as close as possible tonbins. - An
NTuple{Dim, Integer}, in which case dimensionireceivesnbins[i]bins.
Reinjection algorithms
UlamMethod.ReinjectionAlgorithm — Typeabstract type ReinjectionAlgorithmThe abstract type for the algorithm controlling how trajectories pointing from nirvana to the interior are redistributed.
Each subtype reinj_algo should implement the following:
- a function
reinject!(binner, Pij, reinj_algo)that modifies the transition matrixPijin place given
the binning algorithm binner. This function is applied after the largest strongly connected component of the original Pij is taken, but before the matrix is row-stochasticized.
UlamMethod.DataReinjection — Typestruct DataReinjectionReinject the data according to which bins are actually hit by transitions from the exterior to the interior.
Since this is applied by default in the main Ulam's method algorithm, the associated reinject! function for DataReinjection does nothing.
Constructor
DataReinjection()UlamMethod.UniformReinjection — Typestruct UniformReinjectionReinject the data uniformly across all available bins.
Constructor
UniformReinjection()UlamMethod.SourceReinjection — Typestruct SourceReinjection{Dim}Reinject the data at particular locations.
Fields
points: A vector of the form[p1, p2, ...]wherepiis a point such as(x)or(x, y). Data are reinjected uniformly across all bins that contain at least one member ofpoints.fallback: TheReinjectionAlgorithmto apply in case no bins contain any points or there is no data to reinject.
Constructor
SourceReinjection(points; fallback = DataReinjection())UlamMethod.StationaryReinjection — Typestruct StationaryReinjectionReinject the data weighted by the stationary distribution of P_open.
Fields
fallback: TheReinjectionAlgorithmto apply in case the stationary distribution can not be computed.
Constructor
StationaryReinjection(; fallback = DataReinjection())Computing the main results
UlamMethod.ulam_method — Functionulam_method(traj, binner; reinj_algo)Run the main Ulam's method calculation and return an UlamResult.
This is the lower level method with all options available.
Arguments
traj: ATrajectoriesobject, holding the short-range trajectory data.binner: ABinningAlgorithmthat specifies the algorithm used to partition the boundary into bins.
Optional Arguments
reinj_algo: AReinjectionAlgorithmthat specifies how trajectories pointing from nirvana to the interior should be reinjected. DefaultDataReinjection.
ulam_method(traj, nbins; nirvana = 0.10)Run the main Ulam's method calculation and return an UlamResult.
This is the higher level "convenience" method with most options selected automatically.
Nirvana Optional Argument
The boundary is placed such that a fraction nirvana of the data is in nirvana. For example, with the default value of 0.10, roughly 10% of all of the datapoints (split between x0 and xT) will be outside the boundary with a roughly equal amount on each side.
The shape of the boundary can be further controlled by providing nirvana as a vector of length Dim of tuples of the form (min, max) such that a fraction min (respectively max) will be "below" (respectively "above") the boundary along each dimension.
Working with results
UlamMethod.UlamResult — Typestruct UlamResult{Dim, M, CRS}A container for the result of the main Ulam's method calculation.
This consists of three main components, the transition probability matrix, the associated bins and the bins that were disconnected.
The transition probability matrix is of the form
| P_O2O | P_O2ω |
| P_ω2O | 0 |
Where O represents the "open" (interior) region and ω reprsents the "nirvana" (exterior) region. The union of O and ω form a closed region.
Fields
P_O2O: As in diagram.P_O2ω: As in diagram.P_ω2O: As in diagram.binner: TheBinningAlgorithmused to bin the data.bins_dis: The bins that contained data but were removed when the largest strongly connected component was taken.
Methods
P_closed(UlamResult)Access the full matrix.
P_open(UlamResult)Access P_O2O.
bins(UlamResult)Access bins.
bins_dis(UlamResult)Access bins_dis.
membership(data_or_traj, ulam_result)Compute the membership of data_or_traj for the binned domain. See the membership function for more information.
counts(data_or_traj, ulam_result)Compute the bin counts of data_or_traj for the binned domain. See the counts function for more information.
UlamMethod.P_closed — FunctionP_closed(ulam_result)Return the full transition matrix.
UlamMethod.bins — Functionbins(ulam_result)Return the bins associated to the transition matrix.
UlamMethod.bins_dis — Functionbins_dis(ulam_result)Return the bins that contained data but were removed when the largest strongly connected component was taken.
UlamMethod.membership — Functionmembership(data, binner)Takes a Dim x N_points matrix of points and returns a vector memb where memb[i] = j if data[:,i] is inside binner.bins[j] and memb[i] = nothing if data[:,i] is not inside any bin.
membership(data, ulam_result)Compute membership(data, ulam_result.binner).
membership(traj, binner)Compute (membership(traj.x0, binner), membership(traj.xT, binner)).
membership(traj, ulam_result)Compute membership(traj, ulam_result.binner).
UlamMethod.points — Functionpoints(bins)Return a vector of raw vertices for each bin. Can also be applied to a BinningAlgorithm.
UlamMethod.counts — Functioncounts(data, ulam_result)Compute a vector counts such that counts[i] is equal to the number of points inside the ith bin and counts[end] is the number of points outside the boundary (in nirvana).
counts(traj, ulam_result)Compute (counts(traj.x0, ulam_result), counts(traj.xT, ulam_result).
Visualization
UlamMethod.viz — Functionviz(object; kwargs...)This function wraps Meshes.viz for UlamMethod objects for plotting.
Methods
viz(boundary; kwargs...)Plot a boundary.
viz(bins; kwargs...)Plot the bins.
viz(binner; kwargs...)Plot binner.bins.
viz(data; kwargs...)Plot a 2 x N matrix of points.
viz(traj; kwargs...)Plot traj.x0 and traj.xT on the same graph.
viz(ulam; bins_labels = true, bins_labels_size = 16)Plot the bins, boundary and a heatmap of the stationary distribution of P_closed(ulam).
Optionally, superimpose bin labels with bin_labels and adjust their size with bins_labels_size.
viz(traj, ulam; bins_labels = true, bins_labels_size = 16)Plot a heatmap of the counts in each bin of both traj.x0 and traj.xT.
Optionally, superimpose bin labels with bin_labels and adjust their size with bins_labels_size.
The docstring for Meshes.viz is provided below.
viz(object; [options])Visualize Meshes.jl object with various options.
Available options
color- color of geometriesalpha- transparency in [0,1]colormap- color scheme/map from ColorSchemes.jlcolorrange- minimum and maximum color valuesshowsegments- visualize segmentssegmentcolor- color of segmentssegmentsize- width of segmentsshowpoints- visualize pointspointmarker- marker of pointspointcolor- color of pointspointsize- size of points
The option color can be a single scalar or a vector of scalars. For Mesh subtypes, the length of the vector of colors determines if the colors should be assigned to vertices or to elements.
Examples
Different coloring methods (vertex vs. element):
# vertex coloring (i.e. linear interpolation)
viz(mesh, color = 1:nvertices(mesh))
# element coloring (i.e. discrete colors)
viz(mesh, color = 1:nelements(mesh))Different strategies to show the boundary of geometries (showsegments vs. boundary):
# visualize boundary with showsegments
viz(polygon, showsegments = true)
# visualize boundary with separate call
viz(polygon)
viz!(boundary(polygon))Notes
This function will only work in the presence of a Makie.jl backend via package extensions in Julia v1.9 or later versions of the language.
UlamMethod.viz! — Functionviz!(object; kwargs...)This function wraps Meshes.viz! for UlamMethod objects for plotting.
Methods
viz!(boundary; kwargs...)Plot a boundary.
viz!(bins; kwargs...)Plot the bins.
viz!(binner; kwargs...)Plot binner.bins.
viz!(data; kwargs...)Plot a 2 x N matrix of points.
The docstring for Meshes.viz! is provided below.
viz!(object; [options])Visualize Meshes.jl object in an existing scene with options forwarded to viz.
References
- 1Miron, Philippe, et al. "Transition paths of marine debris and the stability of the garbage patches." Chaos: An Interdisciplinary Journal of Nonlinear Science 31.3 (2021): 033101.
- 2In brief, nirvana is an extra state appended to an open system to close it; trajectories which point from inside the domain to the outisde of the domain transition to this nirvana state. Trajectories which point from outside the domain to the inside are transitions "from" nirvana - how exactly these data are reinjected is controlled by the
ReinjectionAlgorithm.