API Reference

Main Interface

SMLMDriftCorrection.driftcorrectFunction
driftcorrect(smld; kwargs...) -> (corrected_smld, info)
driftcorrect(smld, config::DriftConfig) -> (corrected_smld, info)

Main interface for drift correction. Uses Legendre polynomial model with entropy-based cost function and adaptive KDTree neighbor building.

Arguments

  • smld: SMLD structure containing (X, Y) or (X, Y, Z) localization coordinates (μm)
  • config: Optional DriftConfig struct (alternative to keyword arguments)

Keyword Arguments

  • quality=:singlepass: Quality tier (:fft, :singlepass, :iterative)
  • degree=2: Polynomial degree for intra-dataset drift model
  • dataset_mode=:registered: Semantic label for multi-dataset handling:
    • :registered: datasets are independent acquisitions
    • :continuous: one long acquisition split into files
  • chunk_frames=0: For continuous mode, split each dataset into chunks of this many frames
  • n_chunks=0: Alternative to chunk_frames - specify number of chunks per dataset
  • maxn=200: Maximum number of neighbors for entropy calculation
  • max_iterations=10: Maximum iterations for :iterative mode
  • convergence_tol=0.001: Convergence tolerance (μm) for :iterative mode
  • warm_start=nothing: Previous model for warm starting optimization
  • verbose=0: Verbosity level (0=quiet, 1=info, 2=debug)
  • auto_roi=false: Set to true for faster (but slightly less accurate) estimation using a dense ROI subset
  • σ_loc=0.010: Typical localization precision (μm) for ROI sizing
  • σ_target=0.001: Target drift precision (μm) for ROI sizing
  • roi_safety_factor=4.0: Safety multiplier for required localizations

Quality Tiers

  • :fft: Fast cross-correlation only (~10x faster, less accurate)
  • :singlepass: Current algorithm - parallel intra, then sequential inter (default)
  • :iterative: Full convergence - iterates intra↔inter until shift changes < tol

Returns

Tuple (corrected_smld, info) where info::DriftInfo contains:

  • model: Fitted drift model (LegendrePolynomial)
  • elapsed_s: Wall time in seconds
  • backend: Computation backend (:cpu)
  • iterations: Number of iterations completed
  • converged: Whether convergence was achieved
  • entropy: Final entropy value
  • history: Entropy per iteration (for diagnostics)

Example

# Basic usage
(smld_corrected, info) = driftcorrect(smld)

# Using DriftConfig
config = DriftConfig(; quality=:iterative, degree=3, verbose=1)
(smld_corrected, info) = driftcorrect(smld, config)

# Fast FFT-only mode
(smld_corrected, info) = driftcorrect(smld; quality=:fft)

# Warm start from previous result
(smld1, info1) = driftcorrect(smld1; degree=2)
(smld2, info2) = driftcorrect(smld2; warm_start=info1.model)

# Extract drift trajectory for plotting
traj = drift_trajectory(info.model)
source
driftcorrect(smld, info::DriftInfo; kwargs...) -> (corrected_smld, info)

Continue drift correction from a previous result using the model from info.

Keyword Arguments

  • dataset_mode=:registered: Dataset mode (:registered or :continuous)
  • max_iterations=10: Additional iterations to run
  • convergence_tol=0.001: Convergence tolerance (μm)
  • maxn=200: Maximum neighbors for entropy calculation
  • verbose=0: Verbosity level
source

See Configuration for full documentation of DriftConfig and DriftInfo.

Utility Functions

SMLMDriftCorrection.filter_emittersFunction
filter_emitters(smld, keep) -> smld

Select the emitters indexed by keep from the SMLD structure, returning a new SMLD with only those emitters. keep may be a single positive integer, a range, a vector of integers, or a BitVector mask.

Example

smld_sub = filter_emitters(smld, [1])
length(smld_sub.emitters)  # 1
source
SMLMDriftCorrection.drift_trajectoryFunction
drift_trajectory(model; dataset=nothing, frames=nothing, cumulative=false)

Extract drift trajectory from a drift model for plotting.

Arguments

  • model: Polynomial or LegendrePolynomial drift model
  • dataset: specific dataset to extract (default: all datasets)
  • frames: frame range to evaluate (default: 1:n_frames for each dataset)
  • cumulative: if true, chain datasets end-to-end showing total accumulated drift. Useful for registered acquisitions where you want to visualize what drift would look like without stage registration. Default: false.

Returns

NamedTuple with fields ready for plotting:

  • frames: frame numbers (global if multiple datasets)
  • x: x-drift values (μm)
  • y: y-drift values (μm)
  • z: z-drift values (μm) - only present for 3D models
  • dataset: dataset index for each point (useful for coloring)

Example

result = driftcorrect(smld; dataset_mode=:registered)
traj = drift_trajectory(result.model)  # Each dataset relative to ds1

# Cumulative view - chains datasets end-to-end
traj_cumul = drift_trajectory(result.model; cumulative=true)
plot(traj_cumul.frames, traj_cumul.x, label="X drift (cumulative)")
source

Drift Models

SMLMDriftCorrection.LegendrePolynomialType
LegendrePolynomial

Combined intra + inter Legendre drift model. Contains one IntraLegendre per dataset plus inter-dataset shifts. The n_frames field stores the valid frame range for drift evaluation.

source
SMLMDriftCorrection.LegendrePoly1DType
LegendrePoly1D

Univariate Legendre polynomial drift model for a single spatial dimension. Coefficients are for P1(t), P2(t), ..., P_degree(t) where t is normalized to [-1, 1].

Note: P_0 (constant) is NOT included because the inter-dataset shift already handles global offsets. This matches the convention for standard Polynomial1D which uses t^1, t^2, ... (no constant term).

The key advantage over standard polynomials: orthogonal basis functions mean each coefficient captures independent variation, leading to better-conditioned optimization.

source

Drift Application

SMLMDriftCorrection.applydriftFunction

Apply drift to coordinate using Legendre polynomial model. Uses P1 through Pdegree (no constant P_0 term).

source

Apply drift to simulated data.

source

applydrift(smld, driftmodel) -> smld

Applies a drift model to SMLM data (for simulation/testing).

source

Entropy Functions

SMLMDriftCorrection.entropy_HDFunction

H_i(D) is the entropy of the distribution p(r), where D = {d(1),...,d(L)} is the drift at all frames (1 to L). The quantity below is Gaussian Mixture Model single components summed over all localizations provided.

σx and σy are localization uncertainties.

source
SMLMDriftCorrection.ub_entropyFunction

Entropy upper bound based on maxn nearest neighbors of each localization (2D).

Arguments

  • x, y: Vectors of localization positions
  • σx, σy: Vectors of localization uncertainties
  • maxn: Maximum number of neighbors considered (default 200)
  • divmethod: Divergence method ("KL", "Symmetric", "Bhattacharyya", "Mahalanobis")
source

Entropy upper bound based on maxn nearest neighbors of each localization (3D).

source

Matrix interface for ub_entropy - extracts vectors and calls optimized version. Expects r and σ to be N×K matrices (points as rows, dimensions as columns).

source

Cross-Correlation

SMLMDriftCorrection.findshiftFunction

Perform a cross-correlation between images representing localizations in two SMLD structures and compute the shift between the two original images. histbinsize is the size of the bins in the histogram image in the same units as the localization coordinates.

source
SMLMDriftCorrection.histimage2DFunction

Produce a histogram image from the localization coordinates x and y. x and y are in arbitrary units. ROI is [xmin, xmax, ymin, ymax] of the Region Of Interest in the same units as x and y. If not provided, these values are estimated from the coordinate data. histbinsize is the size of the bins of each coordinate in the same units.

source
SMLMDriftCorrection.crosscorr2DFunction

Compute the cross-correlation between two 2D images with zero-padding.

Zero-padding to 2x size eliminates cyclic wrap-around artifacts that cause false peaks at large shifts.

source