Developers API Reference

This section lists the public and private functions exported by GraphLab.jl.

GraphLab.KdNodeType
KdNode

A node in a 2D kd-tree data structure, representing a spatial region and its recursive subdivision.

Each node contains bounding box coordinates, optional _splitting information, and links to child nodes. The node can also store associated data points and optional entry/exit traversal symbols for space-filling curve traversal.

Fields

  • bbox::Tuple{Tuple{Float64, Float64}, Tuple{Float64, Float64}}: Bounding box as (min_coords, max_coords).
  • _splitdim::Union{Int, Nothing}: _splitting dimension (1 for x-axis, 2 for y-axis), or nothing for leaf.
  • _splitval::Union{Float64, Nothing}: Coordinate value of the _splitting plane, or nothing for leaf.
  • left::Union{KdNode, Nothing}: Left child node (subtree with coordinates ≤ _splitval), or nothing.
  • right::Union{KdNode, Nothing}: Right child node (subtree with coordinates > _splitval), or nothing.
  • is_leaf::Bool: Whether the node is a leaf.
  • is_root::Bool: Whether the node is the root of the tree.
  • data::Matrix{Float64}: Data points contained in this node.
  • entry::Union{Symbol, Nothing}: Optional entry symbol for space-filling curve traversal.
  • exit::Union{Symbol, Nothing}: Optional exit symbol for space-filling curve traversal.
source
GraphLab._append_gilbert!Method
_append_gilbert!(path, grid::AbstractMatrix)

Recursively traverse the given matrix grid using a generalized Gilbert curve pattern. Appends the CartesianIndex entries to path in order.

source
GraphLab._bbox_centerMethod
_bbox_center(bbox)

Compute the center point of a bounding box.

Arguments

  • bbox: A tuple ((xmin, ymin), (xmax, ymax)) representing the bounding box corners.

Returns

  • A tuple (xcenter, ycenter) representing the center coordinates of the bounding box.

Example

julia> _bbox_center(((0.0, 0.0), (2.0, 4.0)))
(1.0, 2.0)
source
GraphLab._bounding_boxMethod
_bounding_box(coords::Matrix{Float64})

Compute the axis-aligned bounding box for a set of 2D points.

Arguments

  • coords::Matrix{Float64}: A matrix of size n × 2, where each row represents a 2D point (x, y).

Returns

  • A tuple ((xmin, ymin), (xmax, ymax)) representing the lower-left and upper-right corners of the bounding box.

Example

julia> _bounding_box([1.0 2.0; 3.0 4.0; -1.0 5.0])
((-1.0, 2.0), (3.0, 5.0))
source
GraphLab._build_treeFunction
_build_tree(coords::Matrix{Float64}, is_root::Bool=true)

Recursively build a 2D kd-tree from a set of points. Leaf nodes contain the data points inside their region.

Arguments

  • coords::Matrix{Float64}: A matrix of size n × 2, where each row is a 2D point (x, y).
  • is_root::Bool=true: Whether this node is the root of the tree (automatically set during recursion).

Returns

  • A KdNode representing the root of the kd-tree.

Example

julia> tree = _build_tree([1.0 2.0; 3.0 4.0; 0.0 5.0])
KdNode(...)
source
GraphLab._centerpointMethod

Approximate _centerpoint by coordinate-wise median of a random sample.

Arguments:

  • X::Matrix{Float64}: n × d data matrix
  • sample_size::Int: number of random rows to sample

Returns:

  • center::Vector{Float64}: estimated _centerpoint
source
GraphLab._conmapMethod

Conformal map that moves center to the origin on the unit sphere.

Arguments:

  • center::Vector{Float64}: _centerpoint on the sphere (length d+1)
  • X::Matrix{Float64}: n × (d+1) matrix of lifted points

Returns:

  • Y::Matrix{Float64}: transformed points (same size as X)
source
GraphLab._get_alignmentMethod
_get_alignment(entry::Symbol, _splitdim::Int)

Determine the alignment of an entry direction relative to the _splitting dimension.

Arguments

  • entry::Symbol: Entry direction (:L, :R, :T, or :B).
  • _splitdim::Int: _splitting dimension (1 for x-axis, 2 for y-axis).

Returns

  • :parallel if entry edge direction is parallel to the _splitting dimension.
  • :perpendicular otherwise.

Example

julia> _get_alignment(:L, 1)
:parallel
julia> _get_alignment(:T, 1)
:perpendicular
source
GraphLab._get_childMethod
_get_child(node::KdNode, side::Symbol)

Return the child node of a KdNode corresponding to the specified side.

This function retrieves either the left or right child node depending on the side and the _splitting dimension.

Arguments

  • node::KdNode: The kd-tree node.
  • side::Symbol: Side indicator (:left, :right, :bottom, or :top).

Returns

  • The child KdNode corresponding to the requested side.

Example

julia> _get_child(node, :left)
KdNode(...)
source
GraphLab._get_child_entry_exitMethod
_get_child_entry_exit(parent_entry::Symbol, parent_exit::Symbol, _splitdim::Int, side::Symbol)

Compute the entry and exit directions for a child node in a kd-tree traversal.

Given the parent node’s entry and exit directions, the _splitting dimension, and the side being traversed, this function determines the entry and exit directions for the child node according to space-filling curve traversal rules (using :cis and :tran order types and :parallel / :perpendicular alignment).

Arguments

  • parent_entry::Symbol: Entry direction at the parent node (:L, :R, :T, :B).
  • parent_exit::Symbol: Exit direction at the parent node (:L, :R, :T, :B).
  • _splitdim::Int: _splitting dimension (1 for x-axis, 2 for y-axis).
  • side::Symbol: Side of the child being traversed (:left, :right, :top, :bottom).

Returns

  • A tuple (child_entry, child_exit) indicating the entry and exit directions for the child node.

Example

julia> _get_child_entry_exit(:L, :R, 1, :left)
(:L, :R)
source
GraphLab._get_maximum_matchingMethod
_get_maximum_matching(g::Graph)

Computes a maximum weight matching of the graph g using a silent CBC solver.

Arguments

  • g::Graph: The input graph.

Returns

  • match: The matching object containing mate assignments and total weight.

Notes

  • Assumes unit weights unless otherwise specified.
source
GraphLab._get_partitionMethod
_get_partition(v::AbstractVector, k::Int)

Create a partition indicator vector by dividing indices into k contiguous groups.

This function assigns partition labels 1 to k to the indices in v, distributing them evenly into k partitions based on their position in v. Each value in v is assumed to be a node ID (1-based indexing), and the result maps each node ID to its partition.

Arguments

  • v::AbstractVector: A vector of node IDs to partition.
  • k::Int: The number of partitions.

Returns

  • A vector part of length maximum(v) where part[id] is the partition label (1 to k) assigned to node id.

Example

julia> _get_partition([1, 2, 3, 4, 5, 6], 2)
[1, 1, 1, 2, 2, 2]
source
GraphLab._gilbert_indicesMethod
_gilbert_indices(grid_dims::Tuple{Int, Int}; maj_axis=grid_dims[1] ≥ grid_dims[2] ? 1 : 2)

Generate a list of CartesianIndex{2} representing a generalized space-filling traversal of a grid with dimensions grid_dims, using a Gilbert-like curve. Traversal will favor rows (X axis) if maj_axis == 1, or columns (Y axis) if maj_axis == 2.

Arguments

  • grid_dims::Tuple{Int, Int}: a tuple (rows, cols) specifying grid size.
  • maj_axis: optional; either 1 (row-major preference) or 2 (column-major preference).
    • The curve will traverse more linearly along the maj_axis.

Returns

  • A Vector{CartesianIndex{2}} representing the traversal path.

Example

julia> order = _gilbert_indices((4, 5))
julia> println(order)
[CartesianIndex(1, 1), CartesianIndex(2, 1), CartesianIndex(2, 2), ...]
source
GraphLab._gilbert_orderMethod
_gilbert_order(grid_indices::AbstractMatrix; maj_axis=rows ≥ cols ? 1 : 2)

Return a list of CartesianIndex values from a grid, traversed recursively using Gilbert curve logic. When maj_axis == 2, the grid is transposed to preserve logic.

source
GraphLab._morton_indexMethod
_morton_index(x::Int, y::Int)

Interleave the bits of x and y to produce a Morton code (Z-order curve).

source
GraphLab._morton_indicesMethod
_morton_indices(grid_dims::Tuple{Int, Int})

Return a list of CartesianIndex{2} grid coordinates ordered by 2D Morton (Z-order) curve.

Arguments

  • grid_dims::Tuple{Int, Int}: the number of rows and columns of the grid, e.g. (64, 64).

Returns

  • A vector of CartesianIndex{2} objects representing the traversal order.

Example

julia> order = _morton_indices((4, 4))
julia> println(order)
[CartesianIndex(1, 1), CartesianIndex(1, 2), CartesianIndex(2, 1), CartesianIndex(2, 2), ...]
source
GraphLab._nested_dissectionMethod
_nested_dissection(A::SparseMatrixCSC, method::Function; coords=nothing, minsep, verbose)

Recursively computes a nested dissection ordering of a sparse matrix, handling disconnected components separately.

Arguments

  • A::SparseMatrixCSC: The adjacency matrix of the graph to partition.
  • method::Function: Function used to compute separators.
  • coords::Union{Matrix, Nothing}: Optional node coordinates if the separator method needs them.
  • minsep::Int: Minimum number of nodes to stop recursion.
  • verbose::Bool: Whether to visualize partitions and wait for user input during recursion.

Returns

  • perm::Vector{Int}: The final permutation vector.

Notes

  • Uses approximate minimum degree ordering for small components.
source
GraphLab._opposite_dirMethod
_opposite_dir(dir::Symbol)::Symbol

Return the opposite direction symbol.

Arguments

  • dir::Symbol: Direction symbol (:L, :R, :T, or :B).

Returns

  • The opposite direction symbol.

Example

julia> _opposite_dir(:L)
:R
julia> _opposite_dir(:T)
:B
source
GraphLab._order_typeMethod
_order_type(entry::Symbol, exit::Symbol)

Determine the traversal type based on entry and exit directions.

Arguments

  • entry::Symbol: Entry direction (:L, :R, :T, or :B).
  • exit::Symbol: Exit direction (:L, :R, :T, or :B).

Returns

  • :tran if entry and exit are opposite directions.
  • :cis otherwise.

Example

julia> _order_type(:L, :R)
:tran
julia> _order_type(:T, :R)
:cis
source
GraphLab._part1by1Method
_part1by1(n::Int)

Expand 16-bit integer n into 32 bits with zeros between the original bits.

source
GraphLab._partitionMethod
_partition(coords::Matrix, v::Vector)

Compute a partition based on coords using a direction vector v.

Arguments

  • coords::Matrix: Node coordinates in a 2D space.
  • v::Vector: Direction vector defining the partitioning line.

Returns

  • A tuple of two vectors: indices of nodes in each partition.

Example

julia> _partition(coords, [0, 1])
([1, 2, 3, 9, 10], [4, 5, 6, 7, 8])
source
GraphLab._partitionMethod
_partition(X, direction)

Splits a set of points into two parts by projecting them onto a given direction vector and _partitioning at the median projection value.

Arguments

  • X::Matrix{Float64}: An n × d matrix of points (each row is a point in ℝ^d or ℝ^{d+1}).
  • direction::Vector{Float64}: A direction vector along which to project the points.

Returns

  • part1::Vector{Int}: Indices of points with projection ≤ median (one side of the cut).
  • part2::Vector{Int}: Indices of points with projection > median (the other side).

Notes

  • This ensures a balanced _partition by cutting at the median of the projected values.
  • Commonly used to implement great-circle or hyperplane separators.
source
GraphLab._recursive_bisectionFunction
_recursive_bisection(method::Function, levels::Int, A::AbstractSparseMatrix, 
                     coords::Union{Matrix, Nothing}=nothing, minpoints::Int=8, 
                     vn::Vector{Int}=Int[])

Recursively partition the graph A using the given partitioning method, applying hierarchical bisection.

Arguments

  • method::Function: Partitioning method to apply (e.g., part_spectral, part_inertial).
  • levels::Int: Number of recursive partitioning levels.
  • A::AbstractSparseMatrix: Adjacency matrix of the graph.
  • coords::Union{Matrix, Nothing}=nothing: Node coordinates for spatial partitioning (if applicable).
  • minpoints::Int=8: Minimum number of nodes required to continue partitioning.
  • vn::Vector{Int}=Int[]: Vector of node indices, used for tracking original node ordering.

Returns

  • A vector of partition labels for each node, recursively refined through hierarchical bisection.

Example

julia-repl julia> _recursive_bisection(part_spectral, 3, A, coords) 1 ⋮ 4

source
GraphLab._reflectorMethod
_reflector(c::AbstractVector)

Constructs a Householder reflection matrix Q that maps the input vector c (reversed) to a multiple of the first basis vector. This is useful for aligning a given direction with a coordinate axis in geometric transformations.

Arguments

  • c::AbstractVector: A real vector of length d.

Returns

  • Q::Matrix{Float64}: A d × d orthogonal matrix representing the Householder reflection.
  • r::Float64: The leading entry of the reflected vector (i.e., norm or signed component).

Notes

  • Internally, this performs a QR decomposition of the reversed vector c[end:-1:1] and adjusts the result to restore the original ordering.
source
GraphLab._sepcircleMethod

_sepcircle(A::SparseMatrixCSC, X::Matrix{Float64}, ntrials::Int)

Try ntrials random great circle cuts on unit-sphere data.

Arguments:

  • A::SparseMatrixCSC: adjacency matrix
  • X::Matrix{Float64}: n × (d+1) unit-sphere coordinates
  • ntrials::Int: number of random directions to try

Returns:

  • bestdir::Vector{Float64}: direction vector of best cut
  • mincut::Float64: minimum edge cut found
source
GraphLab._seplineMethod

_sepline(A, xy, ntrials)

Try ntrials random straight hyperplanes in Euclidean space.

Arguments:

  • A::SparseMatrixCSC: adjacency matrix
  • coords::Matrix{Float64}: n × d original coordinates
  • ntrials::Int: number of directions to try

Returns:

  • bestdir::Vector{Float64}: best direction
  • mincut::Float64: corresponding edge cut
source
GraphLab._sepqualityMethod
_sepquality(v, A, xyz)

Evaluates the quality of a geometric separator defined by a direction vector v, by computing the number of graph edges cut by the resulting _partition.

Arguments

  • v::Vector{Float64}: A direction vector defining the separating hyperplane (e.g., great circle).
  • A::SparseMatrixCSC: The adjacency matrix of the (undirected) graph.
  • xyz::Matrix{Float64}: An n × (d+1) matrix of vertex coordinates on the unit sphere (in ℝ^{d+1}).

Returns

  • cutsize::Int: The number of edges crossing between the two sides of the _partition.

Notes

  • Vertices are _partitioned based on the sign of their projection onto v.
  • The number of crossing edges is computed using nonzero entries in A[a, b] and A[b, a]'.
source
GraphLab._splitMethod
_split(coords::Matrix{Float64}, dim::Int)

_split a set of 2D points along a specified dimension at the median.

Arguments

  • coords::Matrix{Float64}: A matrix of size n × 2, where each row is a 2D point (x, y).
  • dim::Int: _splitting dimension (1 for x-axis, 2 for y-axis).

Returns

  • left::Matrix{Float64}: Points in the left subset (≤ _split value).
  • right::Matrix{Float64}: Points in the right subset (> _split value).
  • _splitval::Float64: The coordinate value at which the _split occurs.

Example

julia> left, right, _splitval = _split([1.0 2.0; 3.0 4.0; 0.0 5.0], 1)
([0.0 5.0; 1.0 2.0], [3.0 4.0], 1.0)
source
GraphLab._stepMethod
_step(entry::Symbol, exit::Symbol, _splitdim::Int, entry_pt::Tuple{Float64, Float64}, _splitval::Float64)

Determine the traversal order of child nodes in a kd-tree based on entry and exit directions.

Arguments

  • entry::Symbol: Entry direction (:L, :R, :T, or :B).
  • exit::Symbol: Exit direction (:L, :R, :T, or :B).
  • _splitdim::Int: _splitting dimension (1 for x-axis, 2 for y-axis).
  • entry_pt::Tuple{Float64, Float64}: Coordinates of the entry point.
  • _splitval::Float64: Value of the splitting plane along `splitdim`.

Returns

  • A vector of symbols indicating the traversal order of child nodes (e.g., [:left, :right]).

Example

julia> _step(:L, :R, 1, (0.0, 0.5), 0.3)
[:left, :right]
source
GraphLab._stereodownMethod
_stereodown(xyz::AbstractMatrix)

Performs the inverse of stereographic projection, mapping points from the unit sphere in ℝ^{d} (embedded in ℝ^{d+1}) back to Euclidean space ℝ^{d}.

Arguments

  • xyz::AbstractMatrix: An n × (d + 1) matrix where each row is a point on the unit sphere in ℝ^d.

Returns

  • xy::Matrix{Float64}: An n × d matrix of projected points in Euclidean space.

Notes

  • Assumes the last coordinate in each row of xyz corresponds to the vertical axis (pole).
  • This operation is the inverse of the stereographic projection performed in _stereoup.
source
GraphLab._stereoupMethod

Stereographically lift points to the unit sphere in (d+1)-dimensional space.

coords: n × d matrix of points Returns: n × (d+1) matrix on the unit sphere

source
GraphLab._traverse_kdFunction
_traverse_kd(node::KdNode, entry::Symbol=:L, exit::Symbol=:R, last_exits::Vector{Any}=[], acc::Vector{Any}=[])

Traverse a kd-tree following an adaptive space-filling curve and collect leaf data.

This function recursively traverses a kd-tree according to entry and exit directions, determining the order of child traversal at each node using adaptive rules. At each leaf node, it collects the stored data into an accumulator.

Arguments

  • node::KdNode: The kd-tree node to traverse.
  • entry::Symbol=:L: Entry direction (:L, :R, :T, :B) at the current node.
  • exit::Symbol=:R: Exit direction at the current node.
  • last_exits::Vector{Any}=[]: List of exit points accumulated during traversal.
  • acc::Vector{Any}=[]: Accumulator for collected leaf data.

Returns

  • If node.is_root, returns acc, the accumulated leaf data in traversal order.
  • Otherwise returns the last exit point from the traversal.

Example

julia> data = _traverse_kd(tree)
[[x1, y1, id1], [x2, y2, id2], ...]
source
GraphLab._vis_graphMethod
_vis_graph(A::SparseMatrixCSC, coords::Matrix, p::Vector{Int})

Visualize the partitioning p of graph A using node coordinates coords.

Arguments

  • A: Adjacency matrix of the graph.
  • coords: Node coordinates for plotting.
  • p: Partition labels for each node.

Output

  • Displays the partitioned graph visualization.
source
GraphLab._vtxsepMethod
_vtxsep(gborder, p1border, p2border, match)

Computes a vertex separator based on alternating paths in a bipartite border graph.

Arguments

  • gborder::Graph: Bipartite graph between border nodes.
  • p1border::Vector{Int}: Border nodes from partition 1.
  • p2border::Vector{Int}: Border nodes from partition 2.
  • match: Maximum matching result on gborder.

Returns

  • (sep, new_p1border, new_p2border):
    • sep: Separator set of vertices.
    • new_p1border: Updated partition 1 border after separator removal.
    • new_p2border: Updated partition 2 border after separator removal.

Notes

  • Traverses alternating paths between matching and non-matching edges to construct the separator.
source
GraphLab.grid_graphMethod
grid_graph(n::Int, m::Int, α::Float64)

Returns the adjacency matrix A::SparseMatrixCSC and the coordinates coords::Matrix{Float64} of an n × m grid graph rotated by angle α (in radians).

Vertices are ordered row-wise: vertex i,j has index (i-1)*m + j.

source
GraphLab.wait_for_keyMethod
wait_for_key(prompt)

Prints a prompt to the standard output and waits for a single keypress from the user.

Arguments

  • prompt::String: The message to display before waiting for input.

Returns

  • nothing.
source