Public API
BiodiversityObservationNetworks.AdaptiveHotspot Type
AdaptiveHotspot <: BONSamplerSampling for hotspots. Takes an auxiliary variable (e.g. uncertainty) to target samples toward. Proposed by (Andrade-Pacheco et al., 2022).
Starts at the global maximum of the target/uncertainty surface. Subsequent points are chosen to maximize a trade-off between the target value and spatial diversity (measured via the determinant of a kernel matrix).
Requires a CandidatePool with features. The first feature row is used as the target/uncertainty surface; subsequent rows are ignored.
Fields
n::Int: number of sites to select (default 50)scale: Matérn kernel range parameter ρ (default 1.0)smoothness: Matérn kernel smoothness ν (default 0.5, equivalent to an exponential kernel)
BiodiversityObservationNetworks.BONSampler Type
BONSamplerAbstract supertype for all spatial sampling algorithms. Each concrete sampler is a struct with at minimum an n::Int field specifying the desired sample size.
Implement _sample(rng, sampler, candidatepool) to add a new algorithm.
BiodiversityObservationNetworks.BalancedAcceptance Type
BalancedAcceptance <: BONSamplerBalanced Acceptance Sampling (BAS) using Halton sequences, proposed by (Robertson et al., 2025).
Generates spatially balanced samples by mapping Halton sequences to the candidate coordinate space. When inclusion weights are non-uniform, a third Halton dimension acts as a threshold for acceptance, preferentially selecting higher-weighted candidates.
Fields
n::Int: number of sites to select (default 50)
BiodiversityObservationNetworks.BiodiversityObservationNetwork Type
BiodiversityObservationNetwork{K}The output of sample. Holds selected sites together with the sampler and domain that produced them, enabling downstream evaluation without passing the domain separately.
Fields
sites::Vector{K}: selected keyscoordinates::Matrix—d x n_selectedpositionsfeatures::Union{Matrix, Missing}—p x n_selectedorMissinginclusion::Vector— inclusion weights of the selected sitessampler::BONSampler— the sampler that produced this result
BiodiversityObservationNetworks.CandidatePool Type
CandidatePool{K}The internal representation of a sampling domain: a set of n candidate locations, each with coordinates. Each candidate coordinate can potentially have features (auxiliary variables) and inclusion weights, which may be required/optional for some sampling methods.
Fields
n::Int— the number of candidate locationskeys::Vector{K}— original identifiers (e.g.CartesianIndexfor raster inputs)coordinates::Matrix—2 x nspatial coordinate matrixfeatures::Union{Matrix, Missing}—p x nfeature matrix (where p is the number of auxiliary variables) orMissinginclusion::Vector—n-vector of relative inclusion probability summing to 1
BiodiversityObservationNetworks.CandidatePool Method
CandidatePool(mat::AbstractMatrix; mask = missing, inclusion = missing)Convert a matrix into a CandidatePool. Each cell becomes a candidate with grid-index coordinates (row, col). An optional boolean mask restricts which cells are included. An optional inclusion matrix provides per-cell weights.
BiodiversityObservationNetworks.CandidatePool Method
CandidatePool(matrices::AbstractVector{<:AbstractMatrix}; mask=missing, inclusion=missing)Convert a vector of same-sized matrices into a CandidatePool with features. Each matrix contributes one feature row.
BiodiversityObservationNetworks.CandidatePool Method
CandidatePool(result::BiodiversityObservationNetwork)Convert a BiodiversityObservationNetwork into a new CandidatePool for multi-stage sampling. Selected sites become candidates with uniform inclusion.
BiodiversityObservationNetworks.CubeSampling Type
CubeSampling <: BONSamplerBalanced sampling with respect to auxiliary variables, using the Cube method from (Deville and Tillé, 2022).
Selects a sample whose Horvitz–Thompson estimates of auxiliary variables match the population as closely as possible. The algorithm runs in two phases:
Flight phase: a random walk on inclusion probabilities that preserves balance constraints while driving inclusion probabilities toward 0 or 1.
Landing phase: uses integer programming (via the JuMP + HiGHS packages) to resolve any remaining fractional probabilities once the flight phase terminates.
Requires a CandidatePool with features; pass a vector of matrices or construct a pool with explicit features.
Fields
n::Int: number of sites to select (default 50)
BiodiversityObservationNetworks.DistanceToAnalogNode Type
DistanceToAnalogNode <: RarityMetricFor each pixel, compute the distance in z-scored feature space to the nearest BON node.
sourceBiodiversityObservationNetworks.DistanceToMedian Type
DistanceToMedianRarity score defined as Euclidean distance in feature space to the per-feature median across the raster stack. Optionally, features can be PCA-transformed prior to z-scoring.
sourceBiodiversityObservationNetworks.GRTS Type
GRTS <: BONSamplerGeneralized Random Tessellation Stratified (GRTS) sampling.
GRTS produces spatially balanced samples by recursively partitioning the domain into quadrants, assigning each cell a hierarchical address. The recursive tessellation ensures spatial spread without requiring distance computations.
Originally proposed by (Stevens and Olsen, 2004).
Fields
n::Int: number of sites to select (default 50)
BiodiversityObservationNetworks.JensenShannon Type
JensenShannonThe JensenShannon method for evaluating BiodiversityObservationNetwork design computes how representative the selected BON sites are of the environmental variables across the domain using Jensen-Shannon divergence, a method for measuring the distance between two distributions.
Fields
nbins::Int: number of bins to use when computing empirical probability mass
BiodiversityObservationNetworks.Loarie2009 Type
Loarie2009Implements the climate velocity metric from: Loarie, S., Duffy, P., Hamilton, H. et al. The velocity of climate change. Nature 462, 1052–1055 (2009). https://doi.org/10.1038/nature08649
sourceBiodiversityObservationNetworks.MoransI Type
MoransI <: SamplingMetricComputes Moran's I on the inclusion indicator variable.
Description
Calculates spatial autocorrelation of the sample indicator
This version was proposed by (Tillé et al., 2025).
sourceBiodiversityObservationNetworks.MultivariateEnvironmentalSimilarity Type
MultivariateEnvironmentalSimilarityMultivariate Environmental Similarity Surface (MESS). For each cell, compute the minimum over features of a per-feature similarity score derived from the ECDF of the training distribution, following the standard MESS definition.
sourceBiodiversityObservationNetworks.Pivotal Type
Pivotal <: BONSamplerLocal Pivotal Method (LPM) for spatially balanced sampling.
Iteratively pairs nearby candidates and updates their inclusion probabilities so that one tends toward selection and the other toward exclusion. Repeating over the domain produces a spatially balanced sample that respects per-unit inclusion probabilities.
Proposed in (Grafström et al., 2022).
Fields
n::Int: number of sites to select (default 50)
BiodiversityObservationNetworks.RarityMetric Type
RarityMetricAbstract type encompassing all methods for computing environmental rarity.
sourceBiodiversityObservationNetworks.SamplingMetric Type
SamplingMetricAbstract supertype for metrics that evaluate a BiodiversityObservationNetwork.
Implement evaluate(metric::MyMetric, result::BiodiversityObservationNetwork) to add a new metric.
BiodiversityObservationNetworks.SimpleRandom Type
SimpleRandom <: BONSamplerSimple random sampling (or weighted random sampling when inclusion weights are non-uniform). Selects n candidates without replacement.
BiodiversityObservationNetworks.SpatiallyCorrelatedPoisson Type
SpatiallyCorrelatedPoisson <: BONSamplerImplements Spatially Correlated Poisson Sampling (SCPS) from (Grafström, 2025).
Fields
n::Int: expected number of sites to select (default 50)
Description
Iterates through units, selecting them based on inclusion probabilities, and dynamically adjusting the probabilities of neighboring units to maintain spatial balance.
sourceBiodiversityObservationNetworks.Stratified Type
Stratified <: BONSamplerStratified sampling. The first feature row of the CandidatePool is treated as a discrete stratum label. Draws are allocated across strata via a Multinomial draw on the stratum weights; within each stratum, sites are drawn without replacement.
Recommended usage is to pass a single matrix/SDMLayer of discrete integer stratum labels.
Fields
n::Int: number of sites to selectweights::Union{Vector, Missing}: per-stratum weights. Index corresponds to label.missing(default) uses weights proportional to area.
Notes
Does not guarantee exactly n selected sites when a Multinomial draw allocates more draws to a stratum than it has candidates. In that case the stratum is exhausted and the total falls short.
Usage
weights = [0.1, 0.2, 0.5, 0.1, 0.1]
domain = rand(1:5, (30, 30))
sample(Stratified(weights = weights), domain)BiodiversityObservationNetworks.VoronoiVariance Type
VoronoiVariance <: SamplingMetricThe VoronoiVariance method for characterizing the spatial balance of a sample is based on the initial method proposed by (Stevens and Olsen, 2004), and then extended by (Grafström, 2025).
For a given BiodiversityObservationNetwork bon, the Voronoi tesselation splits the plane into a series of polygons, where the i-th polygon consists of all points in the plane whose nearest node in bon is the i-th node.
These polygons can then be used to assess the spatial balance of a sample.
In an ideally balanced sample, the sum of the inclusion probabilities across each polygon i would equal 1, because in expectation exactly one unit would be sampled in that region.
If we define i, i.e.
then we can assess the spatial balance of a sample by measuring the distance of B, defined as
to measure spatial balance, where smaller values indicate better spatial balance.
sourceBiodiversityObservationNetworks.WithinRange Type
WithinRange <: RarityMetricBoolean rarity surface indicating whether each cell lies within the hyper- rectangle spanned by the per-feature minima and maxima of the BON nodes.
sourceBiodiversityObservationNetworks.evaluate Method
rarity(::DistanceToAnalogNode, bon, layers; mask)For each cell, compute the distance in z-scored feature space to the nearest selected BON node in layers.
BiodiversityObservationNetworks.evaluate Method
evaluate(::JensenShannon, bon::BiodiversityObservationNetwork, layers)BiodiversityObservationNetworks.sample Function
sample([rng::AbstractRNG,] sampler::BONSampler, domain; mask=missing, inclusion=missing)Select sites from domain using sampler. Pass an explicit rng for reproducibility.
domain can be any type accepted by CandidatePool — a Matrix, a Vector of matrices, a CandidatePool, or a BiodiversityObservationNetwork (for multi-stage sampling).
CandidatePool can also accept types from SpeciesDistributionToolkit, including [SDMLayer]s from the [SimpleSDMLayers] subpackage, which represent rasters with geospatial metadata, and vector data types from the [SimpleSDMPolygons] subpackage. The support for these types are included as extensions, meaning their functionality is only loaded once [SpeciesDistributionToolkit] (or one of the corresponding subpackages) is loaded.
BiodiversityObservationNetworks.spatialbalance Method
spatialbalance(::MoransI, domain::Matrix, bon::BiodiversityObservationNetwork)Private API
BiodiversityObservationNetworks.VelocityMetric Type
VelocityMetricAbstract type encompassing all methods for computing environmental velocity.
sourceBiodiversityObservationNetworks._2d_bas Method
_2d_bas(sampler, domain)2D BAS using Halton bases [2,3] to generate spatially spread candidate cells, accepting those that fall on unmasked locations until num_nodes are selected.
BiodiversityObservationNetworks._3d_bas Method
_3d_bas(rng, sampler, cpool)3D BAS using Halton bases [2,3,5]. A candidate (i,j,z) is accepted if the cell is unmasked and z < inclusion[i,j].
BiodiversityObservationNetworks._above_one_update! Method
_above_one_update!(inclusion, pool, i, j, complete_flags, selected_flag)Apply the LPM update when the paired units have πᵢ + πⱼ ≥ 1.
One of the two units is set to 1 (selected) and the other is reduced to πᵢ + πⱼ - 1. The unit whose probability reaches 1 is marked both as included (selected_flag) and complete (complete_flags).
BiodiversityObservationNetworks._apply_update_rule! Method
_apply_update_rule!(inclusion, pool, i, j, selected_flag, complete_flags)Dispatch to the appropriate LPM update based on whether πᵢ + πⱼ is below or at/above one.
BiodiversityObservationNetworks._below_one_update! Method
_below_one_update!(inclusion, pool, i, j, complete_flags)Apply the LPM update when the paired units have πᵢ + πⱼ < 1.
One of the two units is set to 0 (not selected) and the other is increased to πᵢ + πⱼ. The unit whose probability reaches 0 is marked complete.
BiodiversityObservationNetworks._build_kdtree Method
_build_kdtree(cpool::CandidatePool)Build a KDTree from candidate coordinates for O(n*log(n)) nearest-neighbor queries.
BiodiversityObservationNetworks._cube_flight_phase Method
_cube_flight_phase(πₖ, x)Run the flight phase of the cube method.
πₖ are current inclusion probabilities; x (auxiliary matrix) is augmented with a first row of πₖ' so sample size fixed. The method repeatedly finds a direction in the null space of the constraint matrix and pushes a small subset of probabilities to 0 or 1 while preserving balances in expectation.
BiodiversityObservationNetworks._cube_landing_phase Method
_cube_landing_phase(pikstar, πₖ, x)Run the landing phase if some probabilities are still fractional after flight. Formulate a small linear program over the fractional units to select a 0/1 sample whose auxiliary totals match pikstar in expectation with minimal cost C(s) (per Deville & Tillé, 2004).
BiodiversityObservationNetworks._get_address_length Method
_get_address_lengthReturns the number of digits in a GeneralizedRandomTessellatedStratified address for a specific geometry, computed based on the raster dimensions.
BiodiversityObservationNetworks._get_addresses Method
_get_addresses(sz)Gets an address matrix of size sz.
BiodiversityObservationNetworks._grts_grid_size Method
_grts_grid_sizeDetermine grid dimensions for the GRTS tessellation. Uses metadata if available, otherwise infers from coordinate extent.
sourceBiodiversityObservationNetworks._h Method
_h(K)Entropy-like diversity score based on the log-determinant of the kernel matrix K. Larger values encourage selecting points that are diverse with respect to previously chosen sites.
BiodiversityObservationNetworks._histbin Method
_histbin(x, edges)Compute empirical probability mass by binning empirical values for a single feature x.
BiodiversityObservationNetworks._jensen_shannon Method
_jensen_shannon(Xfull, Xsampled; nbins=20)Compute Jensen-Shannon divergence between two empirical distributions of different size: Xfull, the features across the entire domain, and Xsampled, the features at sampled sites.
BiodiversityObservationNetworks._matérn Method
_matérn(d, ρ, ν)Matérn covariance kernel evaluated at distance d, range ρ, smoothness ν. Normalized so that _matérn(0, ρ, ν) == 1.
BiodiversityObservationNetworks._neighbor_map Method
_neighbor_map(tree, coords)Return a dictionary with keys for all nodes i pointing to the list of indices sorted by distance to i.
BiodiversityObservationNetworks._neighbor_order Method
_neighbor_order(tree, coords, i)Return all neighbors of candidate i ordered by distance (excluding self).
BiodiversityObservationNetworks._pick_nodes Method
_pick_nodes(sampler, raster, addresses)BiodiversityObservationNetworks._quadrant_fill! Method
_quadrant_fill!(mat)Takes a matrix mat and splits it into quadrants randomly labeled one through four.
BiodiversityObservationNetworks._quadrant_split! Method
_quadrant_split!(mat, grid_size)Splits a matrix mat into nested quadrants, where the side-length of a submatrix to be split into quadrants is given by grid_size.
BiodiversityObservationNetworks._sort_features_by_mahalanobis Method
_sort_features_by_mahalanobis(features, inclusion)Order units to stabilize the flight phase by spreading early decisions across feature space. Units are sorted by Mahalanobis distance from the inclusion- weighted mean feature vector.
sourceBiodiversityObservationNetworks._standardize Method
_standardizeStandardizes the values of a matrix of predictors across the entire population Xfull, and a set of predictors associated with the sampled sites, Xsampled, by scaling each predictor to [0,1].
Xsampled is standardized based on the minimum and maximum values of each predictor across the population, so both matrices return on the same scale.
Arguments:
Xfull: annxdmatrix, wherenis the size of the population, anddis the number of predictorsXsampled: anmxdmatrix, wherem<nis the size of the sample
BiodiversityObservationNetworks._update_psi Method
_update_psi(Ψ, B)Given current working vector Ψ and constraint block B, compute a feasible update along a null-space direction u by maximizing step sizes λ₁, λ₂ that keep probabilities within [0,1].
BiodiversityObservationNetworks._voronoi_tesselation Method
_voronoi_tesselation(bon, domain)Construct Voronoi polygons for the nodes in a BiodiversityObservationNetwork within the given domain. Output polygons are clipped to the domain extent.
BiodiversityObservationNetworks.guarantees_exact_n Method
Whether the sampler guarantees exactly n selected sites.
BiodiversityObservationNetworks.requires_features Method
Whether the sampler requires features (auxiliary variables associated with each site) to be present.
sourceBiodiversityObservationNetworks.supports_features Method
Whether the sampler can support features (auxiliary variables associated with each site).
sourceBiodiversityObservationNetworks.supports_inclusion Method
Whether the sampler uses custom inclusion probabilities.
sourceBiodiversityObservationNetworks.unique_permutations Method
unique_permutations(x::T, prefix = T()) where {T}Generate all unique permutations for a multiset x without repetition of duplicates. Based on StackOverflow (https://stackoverflow.com/questions/65051953/julia-generate-all-non-repeating-permutations-in-set-with-duplicates).