API

These are the full docstrings for UlamMethod.jl.

Index

Functions

UlamMethod.P_openMethod
P_open(ulam_result)

Return ulam_result.O2O, i.e. the transition matrix without nirvana.

source
UlamMethod.bins_disMethod
bins_dis(ulam_result)

Return the bins that contained data but were removed when the largest strongly connected component was taken.

source
UlamMethod.countsMethod
counts(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).

source
UlamMethod.membershipMethod
membership(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).

source
UlamMethod.ulam_methodMethod
ulam_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: A Trajectories object, holding the short-range trajectory data.
  • binner: A BinningAlgorithm that specifies the algorithm used to partition the boundary into bins.

Optional Arguments


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.

source
UlamMethod.vizFunction
viz(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 geometries
  • alpha - transparency in [0,1]
  • colormap - color scheme/map from ColorSchemes.jl
  • colorrange - minimum and maximum color values
  • showsegments - visualize segments
  • segmentcolor - color of segments
  • segmentsize - width of segments
  • showpoints - visualize points
  • pointmarker - marker of points
  • pointcolor - color of points
  • pointsize - 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.
source
UlamMethod.viz!Function
viz!(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.

source

Types

UlamMethod.BinningAlgorithmType
abstract 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 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.
  • 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.

source
UlamMethod.BinsType
struct Bins{Dim, M, CRS}

A container type for bins of dimension Dim embedded in manifold M with coordinate reference system CRS.

Fields

  • bins: A vector of Polytope{Dim, M, CRS} objects.

Methods

points(bins)

Return a vector of raw vertices for each bin.

source
UlamMethod.BoundaryType
struct 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: A Polytope defining 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.

source
UlamMethod.DataReinjectionType
struct DataReinjection

Reinject 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()
source
UlamMethod.HexagonBinnerType
struct 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: A Boundary object.
  • bins: A Bins object.
  • idx2pos: A vector such that idx2pos[i] gives the position (in bins) of the bin initially (before removing dataless and disconnected bins) labelled i. idx2pos[i] == nothing if this bin was removed.

Constructor

HexagonBinner(nbins, boundary; hardclip = true)
source
UlamMethod.HyperRectangleType
struct HyperRectangle{Dim, M, CRS}

An extention of Meshes to dimensions greater than 3.

The extension is minimal and implements no methods and does not interface with CRS.

Fields

  • min: The coordinates of the minimum of the HyperRectangle as an NTuple.
  • max: The coordinates of the maximum of the HyperRectangle as an NTuple.
  • vertices: [min, max], required by Meshes.

Constructor

HyperRectangle(min, max)
source
UlamMethod.HyperRectangleBinnerType
struct 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: A Boundary object.
  • bins: A Bins object.
  • idx2pos: A vector such that idx2pos[i] gives the position (in bins) of the bin initially (before removing dataless and disconnected bins) labelled i. idx2pos[i] == nothing if this bin was removed.
  • ranges: A vector such that ranges[i] is the AbstractRange giving the gridpoints of the bins along dimension i.

Constructor

HyperRectangleBinner(nbins, boundary)

nbins can be:

  • An Integer, in which case a HyperRectangleBinner will be constructed attempting to distribute boxes proportionally across dimensions such that the total number of boxes is as close as possible to nbins.
  • An NTuple{Dim, Integer}, in which case dimension i receives nbins[i] bins.
source
UlamMethod.LineBinnerType
struct LineBinner{M, CRS}

Bin a one dimensional Segment (line segement) with nbins equally-spaced bins.

Fields

  • boundary: A Boundary object.
  • bins: A Bins object.
  • idx2pos: A vector such that idx2pos[i] gives the position (in bins) of the bin initially (before removing dataless and disconnected bins) labelled i. idx2pos[i] == nothing if this bin was removed.

Constructor

LineBinner(nbins, boundary; hardclip = true)
source
UlamMethod.RectangleBinnerType
struct 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: A Boundary object.
  • bins: A Bins object.
  • idx2pos: A vector such that idx2pos[i] gives the position (in bins) of the bin initially (before removing dataless and disconnected bins) labelled i. idx2pos[i] == nothing if this bin was removed.

Constructor

RectangleBinner(nbins, boundary; hardclip = true)

nbins can be:

  • An Integer, in which case a RectangleBinner will be constructed attempting to distribute boxes proportionally across x and ysuch that the total number of boxes is as close as possible tonbins`.
  • An NTuple{2, Integer}, in which case dimension i receives nbins[i] bins.
  • An NTuple{2, AbstractVector}, in which case nbins[i] defines the edges of the bins in dimension i.
source
UlamMethod.ReinjectionAlgorithmType
abstract type ReinjectionAlgorithm

The 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 matrix Pijin 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.

source
UlamMethod.SourceReinjectionType
struct SourceReinjection{Dim}

Reinject the data at particular locations.

Fields

  • points: A vector of the form [p1, p2, ...] where pi is a point such as (x) or (x, y). Data are reinjected uniformly across all bins that contain at least one member of points.
  • fallback: The ReinjectionAlgorithm to apply in case no bins contain any points or there is no data to reinject.

Constructor

SourceReinjection(points; fallback = DataReinjection())
source
UlamMethod.StationaryReinjectionType
struct StationaryReinjection

Reinject the data weighted by the stationary distribution of P_open.

Fields

  • fallback: The ReinjectionAlgorithm to apply in case the stationary distribution can not be computed.

Constructor

StationaryReinjection(; fallback = DataReinjection())
source
UlamMethod.TrajectoriesType
struct Trajectories{Dim}

A container for trajectory data of dimension Dim.

Fields

  • x0: Initial coordinates, should be a matrix of size Dim x n_traj.
  • xT: Final coordinates, should be a matrix of size Dim x n_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

source
UlamMethod.TriangleBinnerType
struct 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: A Boundary object.
  • bins: A Bins object.
  • idx2pos: A vector such that idx2pos[i] gives the position (in bins) of the bin initially (before removing dataless and disconnected bins) labelled i. idx2pos[i] == nothing if this bin was removed.

Constructor

TriangleBinner(nbins, boundary; hardclip = true)
source
UlamMethod.UlamResultType
struct 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: The BinningAlgorithm used 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.

source
UlamMethod.VoronoiBinnerType
struct 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: A Boundary object.
  • bins: A Bins object.
  • idx2pos: A vector such that idx2pos[i] gives the position (in bins) of the bin initially (before removing dataless and disconnected bins) labelled i. idx2pos[i] == nothing if this bin was removed.

Constructor

VoronoiBinner(nbins, boundary, traj; hardclip = true)
source

Earth Polygons

UlamMethod.EarthPolygons.North_Atlantic_clippedConstant
North_Atlantic_clipped

Vertices defining a polygon contained in the rectangle with lower left corner (-100, -9) and upper right corner (15, 39). This polygon clips away the land, leaving only the ocean in the rectangle.

source