Skip to content

Public API

BiodiversityObservationNetworks.AdaptiveHotspot Type
julia
AdaptiveHotspot <: BONSampler

Sampling 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)

source
BiodiversityObservationNetworks.BONSampler Type
julia
BONSampler

Abstract 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.

source
BiodiversityObservationNetworks.BalancedAcceptance Type
julia
BalancedAcceptance <: BONSampler

Balanced 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)
source
BiodiversityObservationNetworks.BiodiversityObservationNetwork Type
julia
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 keys

  • coordinates::Matrixd x n_selected positions

  • features::Union{Matrix, Missing}p x n_selected or Missing

  • inclusion::Vector — inclusion weights of the selected sites

  • sampler::BONSampler — the sampler that produced this result

source
BiodiversityObservationNetworks.CandidatePool Type
julia
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 locations

  • keys::Vector{K} — original identifiers (e.g. CartesianIndex for raster inputs)

  • coordinates::Matrix2 x n spatial coordinate matrix

  • features::Union{Matrix, Missing}p x n feature matrix (where p is the number of auxiliary variables) or Missing

  • inclusion::Vectorn-vector of relative inclusion probability summing to 1

source
BiodiversityObservationNetworks.CandidatePool Method
julia
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.

source
BiodiversityObservationNetworks.CandidatePool Method
julia
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.

source
BiodiversityObservationNetworks.CandidatePool Method
julia
CandidatePool(result::BiodiversityObservationNetwork)

Convert a BiodiversityObservationNetwork into a new CandidatePool for multi-stage sampling. Selected sites become candidates with uniform inclusion.

source
BiodiversityObservationNetworks.CubeSampling Type
julia
CubeSampling <: BONSampler

Balanced 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:

  1. Flight phase: a random walk on inclusion probabilities that preserves balance constraints while driving inclusion probabilities toward 0 or 1.

  2. 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)
source
BiodiversityObservationNetworks.DistanceToAnalogNode Type
julia
DistanceToAnalogNode <: RarityMetric

For each pixel, compute the distance in z-scored feature space to the nearest BON node.

source
BiodiversityObservationNetworks.DistanceToMedian Type
julia
DistanceToMedian

Rarity 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.

source
BiodiversityObservationNetworks.GRTS Type
julia
GRTS <: BONSampler

Generalized 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)
source
BiodiversityObservationNetworks.JensenShannon Type
julia
JensenShannon

The 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
source
BiodiversityObservationNetworks.Loarie2009 Type
julia
Loarie2009

Implements 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

source
BiodiversityObservationNetworks.MoransI Type
julia
MoransI <: SamplingMetric

Computes Moran's I on the inclusion indicator variable.

Description

Calculates spatial autocorrelation of the sample indicator δ (1 if sampled, 0 otherwise). Negative values indicate spatial inhibition (spread), which is desired for balanced sampling.

This version was proposed by (Tillé et al., 2025).

source
BiodiversityObservationNetworks.MultivariateEnvironmentalSimilarity Type
julia
MultivariateEnvironmentalSimilarity

Multivariate 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.

source
BiodiversityObservationNetworks.Pivotal Type
julia
Pivotal <: BONSampler

Local 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)
source
BiodiversityObservationNetworks.RarityMetric Type
julia
RarityMetric

Abstract type encompassing all methods for computing environmental rarity.

source
BiodiversityObservationNetworks.SamplingMetric Type
julia
SamplingMetric

Abstract supertype for metrics that evaluate a BiodiversityObservationNetwork.

Implement evaluate(metric::MyMetric, result::BiodiversityObservationNetwork) to add a new metric.

source
BiodiversityObservationNetworks.SimpleRandom Type
julia
SimpleRandom <: BONSampler

Simple random sampling (or weighted random sampling when inclusion weights are non-uniform). Selects n candidates without replacement.

source
BiodiversityObservationNetworks.SpatiallyCorrelatedPoisson Type
julia
SpatiallyCorrelatedPoisson <: BONSampler

Implements 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.

source
BiodiversityObservationNetworks.Stratified Type
julia
Stratified <: BONSampler

Stratified 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 select

  • weights::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

julia
weights = [0.1, 0.2, 0.5, 0.1, 0.1]
domain = rand(1:5, (30, 30))
sample(Stratified(weights = weights), domain)
source
BiodiversityObservationNetworks.VoronoiVariance Type
julia
VoronoiVariance <: SamplingMetric

The 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 vi as the total inclusion probability across all elements of the population in Voronoi polygon i, i.e.

vi=jiπj

then we can assess the spatial balance of a sample by measuring the distance of vi from 1 for each polygon. (Grafström, 2025) proposes the metric B, defined as

B=1ni=1n(vi1)2

to measure spatial balance, where smaller values indicate better spatial balance.

source
BiodiversityObservationNetworks.WithinRange Type
julia
WithinRange <: RarityMetric

Boolean rarity surface indicating whether each cell lies within the hyper- rectangle spanned by the per-feature minima and maxima of the BON nodes.

source
BiodiversityObservationNetworks.evaluate Method
julia
rarity(::DistanceToAnalogNode, bon, layers; mask)

For each cell, compute the distance in z-scored feature space to the nearest selected BON node in layers.

source
BiodiversityObservationNetworks.evaluate Method
julia
evaluate(::JensenShannon, bon::BiodiversityObservationNetwork, layers)
source
BiodiversityObservationNetworks.evaluate Method
julia
rarity(::WithinRange, bon, layers)
source
BiodiversityObservationNetworks.sample Function
julia
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.

source
BiodiversityObservationNetworks.spatialbalance Method
julia
spatialbalance(::MoransI, domain::Matrix, bon::BiodiversityObservationNetwork)
source

Private API

BiodiversityObservationNetworks.VelocityMetric Type
julia
VelocityMetric

Abstract type encompassing all methods for computing environmental velocity.

source
BiodiversityObservationNetworks._2d_bas Method
julia
_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.

source
BiodiversityObservationNetworks._3d_bas Method
julia
_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].

source
BiodiversityObservationNetworks._above_one_update! Method
julia
_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).

source
BiodiversityObservationNetworks._apply_update_rule! Method
julia
_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.

source
BiodiversityObservationNetworks._below_one_update! Method
julia
_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.

source
BiodiversityObservationNetworks._build_kdtree Method
julia
_build_kdtree(cpool::CandidatePool)

Build a KDTree from candidate coordinates for O(n*log(n)) nearest-neighbor queries.

source
BiodiversityObservationNetworks._cube_flight_phase Method
julia
_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.

source
BiodiversityObservationNetworks._cube_landing_phase Method
julia
_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).

source
BiodiversityObservationNetworks._get_address_length Method
julia
_get_address_length

Returns the number of digits in a GeneralizedRandomTessellatedStratified address for a specific geometry, computed based on the raster dimensions.

source
BiodiversityObservationNetworks._get_addresses Method
julia
_get_addresses(sz)

Gets an address matrix of size sz.

source
BiodiversityObservationNetworks._grts_grid_size Method
julia
_grts_grid_size

Determine grid dimensions for the GRTS tessellation. Uses metadata if available, otherwise infers from coordinate extent.

source
BiodiversityObservationNetworks._h Method
julia
_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.

source
BiodiversityObservationNetworks._histbin Method
julia
_histbin(x, edges)

Compute empirical probability mass by binning empirical values for a single feature x.

source
BiodiversityObservationNetworks._jensen_shannon Method
julia
_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.

source
BiodiversityObservationNetworks._matérn Method
julia
_matérn(d, ρ, ν)

Matérn covariance kernel evaluated at distance d, range ρ, smoothness ν. Normalized so that _matérn(0, ρ, ν) == 1.

source
BiodiversityObservationNetworks._neighbor_map Method
julia
_neighbor_map(tree, coords)

Return a dictionary with keys for all nodes i pointing to the list of indices sorted by distance to i.

source
BiodiversityObservationNetworks._neighbor_order Method
julia
_neighbor_order(tree, coords, i)

Return all neighbors of candidate i ordered by distance (excluding self).

source
BiodiversityObservationNetworks._pick_nodes Method
julia
_pick_nodes(sampler, raster, addresses)
source
BiodiversityObservationNetworks._quadrant_fill! Method
julia
_quadrant_fill!(mat)

Takes a matrix mat and splits it into quadrants randomly labeled one through four.

source
BiodiversityObservationNetworks._quadrant_split! Method
julia
_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.

source
BiodiversityObservationNetworks._sort_features_by_mahalanobis Method
julia
_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.

source
BiodiversityObservationNetworks._standardize Method
julia
_standardize

Standardizes 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: an n x d matrix, where n is the size of the population, and d is the number of predictors

  • Xsampled: an m x d matrix, where m < n is the size of the sample

source
BiodiversityObservationNetworks._update_psi Method
julia
_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].

source
BiodiversityObservationNetworks._voronoi_tesselation Method
julia
_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.

source
BiodiversityObservationNetworks.guarantees_exact_n Method

Whether the sampler guarantees exactly n selected sites.

source
BiodiversityObservationNetworks.requires_features Method

Whether the sampler requires features (auxiliary variables associated with each site) to be present.

source
BiodiversityObservationNetworks.supports_features Method

Whether the sampler can support features (auxiliary variables associated with each site).

source
BiodiversityObservationNetworks.supports_inclusion Method

Whether the sampler uses custom inclusion probabilities.

source
BiodiversityObservationNetworks.unique_permutations Method
julia
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).

source