Frame-connection for 2D 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.
# Kwargs form (most common)
(combined, info) = frameconnect(smld::BasicSMLD; kwargs...)
# Config form (reusable settings)
(combined, info) = frameconnect(smld::BasicSMLD, config::FrameConnectConfig)Main entry point. Connects repeated localizations and combines them.
Returns tuple (combined, info):
combined::BasicSMLD: Main output - combined high-precision localizationsinfo::FrameConnectInfo: Track assignments and algorithm metadata
combinelocalizations(smld::BasicSMLD) -> BasicSMLDCombines emitters sharing the same track_id using MLE weighted mean. Use when track_id is already populated.
(smld_connected, smld_combined) = defineidealFC(smld::BasicSMLD; max_frame_gap::Int=5)For simulated data where track_id indicates ground-truth emitter ID. Useful for validation/benchmarking.
@kwdef struct FrameConnectConfig <: AbstractSMLMConfig
n_density_neighbors::Int = 2 # Nearest preclusters for local density estimation
max_sigma_dist::Float64 = 5.0 # Sigma multiplier for preclustering distance threshold
max_frame_gap::Int = 5 # Maximum frame gap for temporal adjacency
max_neighbors::Int = 2 # Maximum nearest-neighbors for precluster membership
endConfiguration parameters for frame connection. Use with frameconnect(smld, config) or pass fields as kwargs to frameconnect(smld; kwargs...).
Example:
# Create config with custom settings
config = FrameConnectConfig(max_frame_gap=10, max_sigma_dist=3.0)
(combined, info) = frameconnect(smld, config)
# Equivalent kwargs form
(combined, info) = frameconnect(smld; max_frame_gap=10, max_sigma_dist=3.0)struct FrameConnectInfo{T} <: AbstractSMLMInfo
connected::BasicSMLD{T} # Input with track_id assigned (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 # Rate: dark → visible (1/frame)
k_off::Float64 # Rate: visible → dark (1/frame)
k_bleach::Float64 # Rate: photobleaching (1/frame)
p_miss::Float64 # Probability of missed detection
initial_density::Vector{Float64} # Emitter density per cluster (emitters/μm²)
elapsed_s::Float64 # Wall time in seconds
algorithm::Symbol # Algorithm used (:lap)
n_preclusters::Int # Number of preclusters formed
endAccess track assignments:
(combined, info) = frameconnect(smld)
track_ids = [e.track_id for e in info.connected.emitters]mutable struct ParamStruct
initial_density::Vector{Float64} # Emitter density per cluster (emitters/μm²)
n_density_neighbors::Int # Clusters for density estimation
k_on::Float64 # Rate: dark → visible (1/frame)
k_off::Float64 # Rate: visible → dark (1/frame)
k_bleach::Float64 # Rate: photobleaching (1/frame)
p_miss::Float64 # Probability of missing localization when on
max_sigma_dist::Float64 # Preclustering threshold multiplier
max_frame_gap::Int # Max frame gap in preclusters
max_neighbors::Int # Max nearest-neighbors for preclustering
endPhotophysics terms:
k_on: Transition rate from dark (non-fluorescent) to visible statek_off: Transition rate from visible to dark statek_bleach: Irreversible photobleaching ratep_miss: Probability localization algorithm fails to detect an "on" fluorophore
Input BasicSMLD must contain emitters with position uncertainties.
Required fields:
x,y: Position (microns)σ_x,σ_y: Position uncertainties (microns) - must be > 0 for MLEframe: Frame number (1-based)
Optional fields:
photons,σ_photons: Photon count (summed in output)bg,σ_bg: Backgrounddataset: Dataset identifier (default: 1)track_id: Set to 0 for input (populated by algorithm)
Output combined contains emitters with:
- Combined position via MLE weighted mean:
x = Σ(x/σ²) / Σ(1/σ²) - Reduced uncertainties:
σ = √(1/Σ(1/σ²)) - Summed photons
track_idindicating connection group
- SMLMData.jl (v0.7): BasicSMLD, Emitter types, AbstractSMLMConfig, AbstractSMLMInfo
- Hungarian.jl: Linear assignment problem solver
- NearestNeighbors.jl: Spatial clustering
- Optim.jl: Parameter estimation
Schodt & Lidke (2021), "Spatiotemporal Clustering of Repeated Super-Resolution Localizations via Linear Assignment Problem", Front. Bioinform. https://doi.org/10.3389/fbinf.2021.724325