SMLMFrameConnection

Frame-connection for 2D single molecule localization microscopy (SMLM) data: linking localizations from the same fluorophore blinking event across consecutive frames into single, higher-precision localizations. Uses spatiotemporal LAP assignment to optimally connect temporally adjacent detections based on spatial proximity and estimated blinking kinetics. See Schodt & Lidke 2021.

Installation

using Pkg
Pkg.add("SMLMFrameConnection")

Quick Start

using SMLMFrameConnection

# Run frame connection on your BasicSMLD with Emitter2DFit emitters
(combined, info) = frameconnect(smld)

# combined is the main output - higher precision localizations
# info contains track assignments and algorithm metadata

Input Requirements

Input BasicSMLD must contain Emitter2DFit emitters with:

Required:

  • x, y: Position coordinates (microns)
  • σ_x, σ_y: Position uncertainties (microns) - must be > 0 for MLE weighting
  • frame: Frame number (1-based)

Optional:

  • photons, bg: Photometry (summed in output)
  • dataset: Dataset identifier (default: 1)

Algorithm Overview

  1. Preclustering: Groups spatially and temporally adjacent localizations into candidate clusters
  2. Parameter Estimation: Estimates fluorophore blinking kinetics (k_on, k_off, k_bleach) and emitter density from the data
  3. Frame Connection: Uses Linear Assignment Problem (LAP) to optimally assign localizations to emitters based on spatial proximity and estimated photophysics
  4. Combination: Combines connected localizations using MLE weighted mean for improved precision

Outputs

OutputDescription
combinedMain output - combined high-precision localizations
info.connectedOriginal localizations with track_id assigned
info.n_tracksNumber of tracks formed
info.elapsed_sWall time in seconds

Parameters

frameconnect(smld;
    n_density_neighbors = 2,   # Clusters for density estimation
    max_sigma_dist = 5.0,        # Distance threshold multiplier
    max_frame_gap = 5,        # Max frame gap for connections
    max_neighbors = 2              # Nearest neighbors for preclustering
)
  • max_sigma_dist: Higher values allow connections over larger distances
  • max_frame_gap: Increase for dyes with long dark states (dSTORM: 10-20)

API Reference

SMLMFrameConnection.CalibrationConfigType
CalibrationConfig

Configuration for optional uncertainty calibration within frame connection.

Calibration analyzes frame-to-frame jitter within connected tracks to estimate a motion variance (σ_motion²) and CRLB scale factor (), then applies corrected uncertainties before track combination: Σ_corrected = σ_motion² I + k² Σ_CRLB.

Fields

  • clamp_k_to_one::Bool=true: Clamp k_scale to minimum of 1.0 (CRLB is a lower bound)
  • filter_high_chi2::Bool=false: Filter tracks with high chi² pairs before fitting
  • chi2_filter_threshold::Float64=6.0: Chi² threshold for track filtering (per pair)
source
SMLMFrameConnection.CalibrationResultType
CalibrationResult

Output diagnostics from uncertainty calibration. Returned in FrameConnectInfo.calibration when calibration is enabled.

The calibration model fits: observed_var = A + B * CRLB_var where A = σ_motion² (additive motion/vibration variance) and B = k² (CRLB scale factor).

Fields

  • sigma_motion_nm::Float64: Estimated motion std dev in nm (√A × 1000)
  • k_scale::Float64: CRLB scale factor (√B, or clamped √B if clampkto_one)
  • A::Float64: Fit intercept (σ_motion² in μm²)
  • B::Float64: Fit slope (k²)
  • A_sigma::Float64: Standard error of A
  • B_sigma::Float64: Standard error of B
  • r_squared::Float64: R² of the WLS fit
  • mean_chi2::Float64: Mean chi² across all frame-to-frame pairs
  • n_pairs::Int: Number of frame-to-frame pairs used in fit
  • n_tracks_used::Int: Number of tracks contributing pairs
  • n_tracks_filtered::Int: Number of tracks removed by chi² filter
  • bin_centers::Vector{Float64}: Bin centers (CRLB variance) for diagnostic plots
  • bin_observed::Vector{Float64}: Bin observed variance for diagnostic plots
  • frame_shifts::Dict{Int, Vector{NTuple{2,Float64}}}: Per-dataset (dx,dy) shifts for jitter plots
  • calibration_applied::Bool: Whether calibration was actually applied (false on fallback)
  • warning::String: Warning message if calibration fell back (empty if OK)
source
SMLMFrameConnection.FrameConnectConfigType
FrameConnectConfig

Configuration parameters for frame connection algorithm.

Fields

  • n_density_neighbors::Int=2: Number of nearest preclusters used for local density estimates (see estimatedensities)
  • max_sigma_dist::Float64=5.0: Multiplier of localization errors that defines a pre-clustering distance threshold (see precluster)
  • max_frame_gap::Int=5: Maximum frame gap between temporally adjacent localizations in a precluster (see precluster)
  • max_neighbors::Int=2: Maximum number of nearest-neighbors inspected for precluster membership (see precluster)
  • calibration::Union{CalibrationConfig, Nothing}=nothing: Optional uncertainty calibration

Example

# Without calibration (default)
(combined, info) = frameconnect(smld)

# With calibration
config = FrameConnectConfig(calibration=CalibrationConfig())
(combined, info) = frameconnect(smld, config)

# With calibration + chi2 filtering
config = FrameConnectConfig(
    calibration=CalibrationConfig(filter_high_chi2=true, chi2_filter_threshold=4.0)
)
(combined, info) = frameconnect(smld, config)
source
SMLMFrameConnection.FrameConnectInfoType
FrameConnectInfo{T}

Secondary output from frameconnect() containing track assignments and algorithm metadata.

Fields

  • connected::BasicSMLD{T}: Input SMLD with track_id assigned (localizations uncombined)
  • n_input::Int: Number of input localizations
  • n_tracks::Int: Number of tracks formed
  • n_combined::Int: Number of output localizations
  • k_on::Float64: Estimated on rate (1/frame)
  • k_off::Float64: Estimated off rate (1/frame)
  • k_bleach::Float64: Estimated bleach rate (1/frame)
  • p_miss::Float64: Probability of missed detection
  • initial_density::Vector{Float64}: Initial density estimate per cluster (emitters/μm²)
  • elapsed_s::Float64: Wall time in seconds
  • algorithm::Symbol: Algorithm used (:lap)
  • n_preclusters::Int: Number of preclusters formed
  • calibration::Union{CalibrationResult, Nothing}: Calibration diagnostics (nothing if disabled)

Rate Parameter Interpretation

The rate parameters describe the photophysics of blinking fluorophores:

  • k_on: Rate at which dark emitters convert to visible state
  • k_off: Rate at which visible emitters convert to reversible dark state
  • k_bleach: Rate at which visible emitters are irreversibly photobleached
  • Duty cycle = kon / (kon + k_off)
  • For typical dSTORM: kon << koff (low duty cycle, mostly dark with brief blinks)

Example

(combined, info) = frameconnect(smld)
println("Connected $(info.n_input) → $(info.n_combined) localizations")
println("Formed $(info.n_tracks) tracks in $(info.elapsed_s)s")
# Access track assignments via info.connected
# Check calibration results if enabled
if info.calibration !== nothing
    println("k_scale = $(info.calibration.k_scale)")
    println("σ_motion = $(info.calibration.sigma_motion_nm) nm")
end
source
SMLMFrameConnection.ParamStructType
ParamStruct(initial_density::Vector{Float64},
    k_on::Float64, k_off::Float64, k_bleach::Float64, p_miss::Float64,
    max_sigma_dist::Float64, max_frame_gap::Int, n_density_neighbors::Int)

Structure of parameters needed for frame-connection.

Fields

  • initial_density: Density of emitters at the start of the experiment. (see estimatedensities()) (emitters/μm²)
  • n_density_neighbors: Number of nearest preclusters used for local density estimates. (default = 2)(see estimatedensities())
  • k_on: Rate at which dark emitters convert to the visible state (1/frame).
  • k_off: Rate at which visible emitters convert to the reversible dark state (1/frame).
  • k_bleach: Rate at which visible emitters are irreversibly photobleached (1/frame).
  • p_miss: Probability of missing a localization of a visible emitter.
  • max_sigma_dist: Multiplier of localization errors that defines a pre-clustering distance threshold. (default = 5)(see precluster())(unitless)
  • max_frame_gap: Maximum frame gap between temporally adjacent localizations in a precluster. (default = 5)(see precluster())(frames)
  • max_neighbors: Maximum number of nearest-neighbors inspected for precluster membership. Ideally, this would be set to inf, but that's not feasible for most data. (default = 2)(see precluster())

Note: kon, koff, and kbleach are transition rates, not duty cycle fractions. The duty cycle (fraction of time in ON state) is kon/(kon + koff). For typical dSTORM, kon << koff gives low duty cycle (mostly dark, brief blinks).

source
SMLMFrameConnection.analyze_calibrationMethod
analyze_calibration(smld::BasicSMLD, config::CalibrationConfig) -> CalibrationResult

Analyze frame-to-frame jitter in connected tracks to estimate uncertainty calibration parameters.

Collects consecutive-frame pairs within each track, computes observed positional variance vs CRLB variance, and fits the model: observed_var = A + B * CRLB_var using weighted least squares on binned data.

Falls back gracefully (returns calibration_applied=false) if fewer than 100 frame-to-frame pairs are available or if the fit quality is poor (R² < 0.1).

source
SMLMFrameConnection.apply_calibrationMethod
apply_calibration(smld::BasicSMLD, result::CalibrationResult) -> BasicSMLD

Apply calibrated uncertainties to emitters. Returns a new SMLD with corrected covariance: Σ_corrected = σ_motion² I + k² Σ_CRLB.

Only modifies σx, σy, σ_xy. Positions and other fields are unchanged.

source
SMLMFrameConnection.combinelocalizationsMethod
smld_combined = combinelocalizations(smld::BasicSMLD{T,E}) where {T, E<:SMLMData.AbstractEmitter}

Combine clustered localizations in smld into higher precision localizations.

Description

This function combines localizations in smld that share the same value of track_id. Localizations are combined assuming they arose from independent measurements of the same position with Gaussian errors.

source
SMLMFrameConnection.computeclusterinfoMethod
clusterdurations, nobservations = 
    computeclusterinfo(clusterdata::Vector{Matrix{Float32}})

Compute the durations of and number of member localizations in clusters.

Description

This method computes the duration of each cluster and the number of localizations in each cluster present in the input clusterdata.

source
SMLMFrameConnection.connect1DSMethod
connect1DS(smld::BasicSMLD{T,E}, dataset::Int,
           connectID::Vector{Int}, maxID::Int;
           max_frame_gap::Int = 5) where {T, E<:SMLMData.AbstractEmitter}

Define the "ideal" frame-connection result for simulated smld with one dataset.

Description

This is a helper function which defines the ideal FC result for a single dataset. The user should probably not bother using this, and instead should call defineidealFC().

Inputs

  • smld: BasicSMLD with track_id populated to indicate emitter ID.
  • dataset: Dataset number to be connected.
  • connectID: Current connection IDs for all emitters.
  • maxID: Current maximum ID value.
  • max_frame_gap: Maximum frame gap allowed between localizations connected in the "ideal" result.

Returns

  • Updated connectID array and new maxID.
source
SMLMFrameConnection.connectlocalizations!Method
connectlocalizations!(connectID::Vector{Int64},
    clusterdata::Vector{Matrix{Float32}},
    params::ParamStruct, nframes::Int64)

Connect localizations in clusterdata by solving a linear assignment problem.

Description

Connect localizations in clusterdata by giving "connected" localizations the same integer value for their field smd.connectID. Associations are made by solving a linear assignment problem which designates the ideal connections between localizations.

source
SMLMFrameConnection.connectlocalizationsMethod
connectID = connectlocalizations(connectID::Vector{Int64}, 
    clusterdata::Vector{Matrix{Float32}}, 
    params::ParamStruct, nframes::Int64)

Connect localizations in clusterdata by solving a linear assignment problem.

Description

Connect localizations in clusterdata by giving "connected" localizations the same integer value for their field smd.connectID. Associations are made by solving a linear assignment problem which designates the ideal connections between localizations.

source
SMLMFrameConnection.create_costmatrixMethod
costmatrix = create_costmatrix(clusterdata::Vector{Matrix{Float32}},
    params::ParamStruct, clusterind::Int64, nframes::Int64,
    base_rho_on::Vector{Float64})

Create a cost matrix for connections between localizations in clusterdata.

Description

The costs associated with connecting/not connecting the localizations in the input set of data in clusterdata are placed in entries of the output matrix costmatrix. For N localizations in clusterdata, the output costmatrix will be an 2Nx2N block matrix. The NxN blocks are attributed the following meanings: the upper-left "connection" block corresponds to adding localizations to existing clusters of other localizations, the bottom-left "birth" block corresponds to introduction of a new cluster, and the upper-right "death" block corresponds to preventing admission of any more localizations to a cluster. The bottom-right "auxillary" block is the transpose of the "connection" block. Forbidden connections are indicated by missing.

base_rho_on is a precomputed unit-density on-state density curve of length nframes, shared across all clusters to avoid per-cluster allocation.

source
SMLMFrameConnection.defineidealFCMethod
smld_connected, smld_combined = defineidealFC(
    smld::BasicSMLD{T,E};
    max_frame_gap::Int = 5) where {T, E<:SMLMData.AbstractEmitter}

Define the "ideal" frame-connection result for a simulated smld.

Description

This function defines the "ideal" frame connection result from a simulation. That is to say, for a simulated BasicSMLD structure smld with track_id field populated to indicate emitter membership of localizations, this function will generate an "ideal" FC result which combines all blinking events that appeared with frame gaps less than max_frame_gap of one another. Note that for very high duty cycles, multiple blinking events might be mistakingly combined by this method (i.e., if the emitter blinks back on within max_frame_gap frames of its previous blink). Note that localizations are not allowed to be connected across datasets.

Inputs

  • smld: BasicSMLD with track_id populated to indicate emitter ID.
  • max_frame_gap: Maximum frame gap allowed between localizations connected in the "ideal" result.

Outputs

  • smld_connected: Copy of the input smld with track_id modified to reflect blinking event ID.
  • smld_combined: Ideal frame-connection result with localizations combined.
source
SMLMFrameConnection.estimatedensitiesMethod
initial_density = estimatedensities(smld::BasicSMLD{T,E},
    clusterdata::Vector{Matrix{Float32}}, params::ParamStruct) where {T, E<:SMLMData.AbstractEmitter}

Estimate local emitter densities for clusters in smld and clusterdata.

Description

The initial local densities initial_density around each pre-cluster present in smld/clusterdata are estimated based on the local density of pre-clusters throughout the entire set of data as well as some of the rate parameters provided in params.

source
SMLMFrameConnection.estimateparamsMethod
k_on, k_off, k_bleach, p_miss, nemitters =
    estimateparams(smld::BasicSMLD{T,E},
                   clusterdata::Vector{Matrix{Float32}}) where {T, E<:SMLMData.AbstractEmitter}

Estimate rate parameters from the clusters in smld and clusterdata.

Description

The rate parameters k_on, k_off, and k_bleach, the probability of missing a localization for a visible emitter p_miss, and the underlying number of emitters nemitters are estimated from the pre-clusters in smld/clusterdata. This is done by assuming that pre-clusters are, on average, representative of a single blinking event of a single emitter.

Citation

Schodt David J., Lidke Keith A., "Spatiotemporal Clustering of Repeated Super-Resolution Localizations via Linear Assignment Problem", Front. Bioinform., 20 October 2021 https://doi.org/10.3389/fbinf.2021.724325

source
SMLMFrameConnection.frameconnectMethod
(combined, info) = frameconnect(smld::BasicSMLD, config::FrameConnectConfig)
(combined, info) = frameconnect(smld::BasicSMLD; kwargs...)

Connect repeated localizations of the same emitter in smld.

Description

Repeated localizations of the same emitter present in smld are connected and combined into higher precision localizations of that emitter. This is done by

  1. forming pre-clusters of localizations, 2) estimating rate parameters from

the pre-clusters, 3) solving a linear assignment problem for connecting localizations in each pre-cluster, and 4) combining the connected localizations using their MLE position estimate assuming Gaussian noise.

Arguments

  • smld::BasicSMLD: Localizations to connect. Must contain emitters with valid position uncertainties (σx, σy).
  • config::FrameConnectConfig: Configuration parameters (optional, can use kwargs instead)

Keyword Arguments (equivalent to FrameConnectConfig fields)

  • n_density_neighbors::Int=2: Number of nearest preclusters used for local density estimates (see estimatedensities)
  • max_sigma_dist::Float64=5.0: Multiplier of localization errors that defines a pre-clustering distance threshold (see precluster)
  • max_frame_gap::Int=5: Maximum frame gap between temporally adjacent localizations in a precluster (see precluster)
  • max_neighbors::Int=2: Maximum number of nearest-neighbors inspected for precluster membership (see precluster)

Returns

A tuple (combined, info):

  • combined::BasicSMLD: Connected localizations combined into higher precision results
  • info::FrameConnectInfo: Track assignments and algorithm metadata (see FrameConnectInfo)

Example

# Using kwargs (most common)
(combined, info) = frameconnect(smld)
(combined, info) = frameconnect(smld; max_frame_gap=10)

# Using config struct
config = FrameConnectConfig(max_frame_gap=10, max_sigma_dist=3.0)
(combined, info) = frameconnect(smld, config)

println("Connected $(info.n_input) → $(info.n_combined) localizations")
println("Formed $(info.n_tracks) tracks from $(info.n_preclusters) preclusters")

# Access track assignments for downstream analysis
track_ids = [e.track_id for e in info.connected.emitters]
source
SMLMFrameConnection.linkclusters!Method
connectID, maxconnectID = linkclusters!(connectID::Int64, 
    maxconnectID::Int64, updateind::Vector{Int64}, 
    assignment::Vector{Int64})

Force localizations linked by assignment to share same the same connectID

source
SMLMFrameConnection.linkclustersMethod
connectID, maxconnectID = linkclusters(
    connectID::Int64, maxconnectID::Int64, 
    updateind::Vector{Int64}, assignment::Vector{Int64})

Force localizations linked by assignment to share same the same connectID

source
SMLMFrameConnection.organizeclustersMethod
clusterdata = organizeclusters(smld::BasicSMLD{T,E}) where {T, E<:SMLMData.AbstractEmitter}

Organize pre-clusters into a vector indexing distinct clusters.

Description

Pre-clusters in the input smld, as related by their shared integer value of track_id, are organized into a vector of matrices clusterdata. Each index of clusterdata corresponds to a distinct cluster.

source
SMLMFrameConnection.preclusterMethod
smld_preclustered = precluster(smld::BasicSMLD{T,E},
    params::ParamStruct = ParamStruct()) where {T, E<:SMLMData.AbstractEmitter}

Cluster localizations in smld based on distance and time thresholds in params.

Description

Localizations in the input structure smld are clustered together based on their spatiotemporal separations. All localizations within a spatial threshold of params.max_sigma_dist*mean([σ_x σ_y]) and a temporal threshold of params.max_frame_gap of one another will be clustered together, meaning that these localizations now share the same unique integer value for their track_id field.

Notes

Pre-clustering allows localizations observed in the same frame to be in the same cluster. This is done to prevent exclusion of the "correct" localization from its ideal cluster due to a previously included "incorrect" localization into that cluster.

source
SMLMFrameConnection.solveLAPMethod
assignment, cost = solveLAP(costmatrix::Matrix{Float32})

Solve the linear assignment problem designated by costmatrix.

Description

Solve the linear assignment problem (LAP) defined by costmatrix using the package Hungarian.jl. For now, this method is just a wrapper around Hungarian.jl, however the intention is to add more checks and options to this method in the future.

source
SMLMFrameConnection.to_emitter2dfitMethod
to_emitter2dfit(e::AbstractEmitter, track_id::Int) -> Emitter2DFit

Convert any AbstractEmitter to Emitter2DFit with specified track_id. Used internally to support different input emitter types while producing consistent Emitter2DFit output.

source

Citation

David J. Schodt and Keith A. Lidke, "Spatiotemporal Clustering of Repeated Super-Resolution Localizations via Linear Assignment Problem", Frontiers in Bioinformatics, 2021. DOI: 10.3389/fbinf.2021.724325