Library
SMLMData.SMLMDataSMLMData.AbstractCameraSMLMData.AbstractEmitterSMLMData.AbstractSMLDSMLMData.BasicSMLDSMLMData.BasicSMLDSMLMData.Emitter2DSMLMData.Emitter2DFitSMLMData.Emitter2DFitSMLMData.Emitter3DSMLMData.Emitter3DFitSMLMData.Emitter3DFitSMLMData.IdealCameraSMLMData.IdealCameraSMLMData.IdealCameraSMLMData.IdealCameraSMLMData.IdealCameraSMLMData.ROIBatchSMLMData.ROIBatchSMLMData.ROIBatchSMLMData.SCMOSCameraSMLMData.SCMOSCameraSMLMData.SCMOSCameraSMLMData.SMLDSMLMData.SingleROISMLMData.SmiteSMDSMLMData.SmiteSMLDAdapt.adapt_structureBase.getindexBase.iterateBase.iterateBase.lengthBase.lengthBase.showBase.showBase.showBase.showBase.sizeSMLMData.api_overviewSMLMData.cat_smldSMLMData.check_complex_fieldsSMLMData.compute_bin_edgesSMLMData.compute_bin_edgesSMLMData.compute_edges_1dSMLMData.filter_framesSMLMData.filter_roiSMLMData.format_with_commasSMLMData.get_pixel_centersSMLMData.get_valid_indicesSMLMData.has_nonzero_imagSMLMData.load_smite_2dSMLMData.load_smite_3dSMLMData.merge_smldSMLMData.physical_to_pixelSMLMData.physical_to_pixel_indexSMLMData.pixel_to_physicalSMLMData.save_smiteSMLMData.@filter
SMLMData.SMLMData — Module
SMLMDataA Julia package for working with Single Molecule Localization Microscopy (SMLM) data.
Features
- Type system for emitters, cameras, and localization data
- Physical coordinate handling (microns) with camera pixel mappings
- Filtering and ROI selection tools
- SMITE format compatibility
- Memory-efficient data structures
Basic Usage
using SMLMData
# Create a camera
cam = IdealCamera(1:512, 1:512, 0.1) # 512x512 camera with 0.1 micron pixels
# Create some emitters
emitters = [
Emitter2DFit{Float64}(1.0, 2.0, 1000.0, 10.0, 0.01, 0.01, 50.0, 2.0),
Emitter2DFit{Float64}(3.0, 4.0, 1200.0, 12.0, 0.01, 0.01, 60.0, 2.0)
]
# Create SMLD object
smld = BasicSMLD(emitters, cam, 1, 1, Dict{String,Any}())
# Filter operations
roi = filter_roi(smld, 0.0:2.0, 1.0:3.0)
bright = @filter(smld, photons > 1000)API Overview
For a comprehensive overview of the API, use the help mode on api:
?apiOr access the complete API documentation programmatically:
docs = SMLMData.api()SMLMData.AbstractCamera — Type
AbstractCameraAbstract base type for all camera implementations in single molecule localization microscopy (SMLM).
Interface Requirements
Any concrete subtype of AbstractCamera must provide:
Field Requirements:
pixel_edges_x::Vector{<:Real}: Vector of pixel edge positions in x directionpixel_edges_y::Vector{<:Real}: Vector of pixel edge positions in y direction
Units:
- All edge positions must be in physical units (microns)
- Origin (0,0) corresponds to the top-left corner of the camera
- For a camera with N×M pixels, there will be N+1 x-edges and M+1 y-edges
Coordinate Convention:
- Pixel (1,1) is centered at (pixelsizex/2, pixelsizey/2) microns
- Edge positions define the boundaries of pixels in physical space
- First edge position corresponds to the left/top edge of the first pixel
- Last edge position corresponds to the right/bottom edge of the last pixel
Notes
- Edge positions must be monotonically increasing
- The number of edges must be one more than the number of pixels in each dimension
- While pixels are typically uniform in size, this is not a requirement of the interface
SMLMData.AbstractEmitter — Type
AbstractEmitterAbstract supertype for all emitter types in single molecule localization microscopy (SMLM). All spatial coordinates are specified in physical units (microns).
SMLMData.AbstractSMLD — Type
AbstractSMLDAbstract type representing Single Molecule Localization Data (SMLD).
Interface Requirements
Any concrete subtype of AbstractSMLD must provide:
emitters::Vector{<:AbstractEmitter}: Vector of localized emitters
Additional fields may include:
- Camera information
- Acquisition parameters
- Analysis metadata
Note: All emitter coordinates must be in physical units (microns).
SMLMData.BasicSMLD — Type
BasicSMLD{T,E<:AbstractEmitter} <: AbstractSMLDBasic container for single molecule localization data.
Fields
emitters::Vector{E}: Vector of localized emitterscamera::AbstractCamera: Camera used for acquisitionn_frames::Int: Total number of frames in acquisitionn_datasets::Int: Number of datasets in the acquisitionmetadata::Dict{String,Any}: Additional dataset information
Type Parameters
T: Numeric type for coordinates (typically Float64)E: Concrete emitter type
Example
# Create camera
cam = IdealCamera(1:512, 1:512, 0.1)
# Create some emitters
emitters = [
Emitter2DFit{Float64}(1.0, 1.0, 1000.0, 10.0, 0.01, 0.01, 50.0, 2.0; frame=1),
Emitter2DFit{Float64}(5.0, 5.0, 1200.0, 12.0, 0.01, 0.01, 60.0, 2.0; frame=2)
]
# Create metadata
metadata = Dict{String,Any}(
"exposure_time" => 0.1,
"timestamp" => now(),
"sample" => "Test Sample"
)
# Create SMLD object
data = BasicSMLD(emitters, cam, 2, 1, metadata)SMLMData.BasicSMLD — Method
BasicSMLD(emitters::Vector{E}, camera::AbstractCamera,
n_frames::Int, n_datasets::Int,
metadata::Dict{String,Any}=Dict{String,Any}()) where E<:AbstractEmitterConstruct a BasicSMLD from a vector of emitters and required metadata.
Arguments
emitters::Vector{E}: Vector of localized emitterscamera::AbstractCamera: Camera used for acquisitionn_frames::Int: Total number of frames in acquisitionn_datasets::Int: Number of datasets in acquisitionmetadata::Dict{String,Any}=Dict{String,Any}(): Optional additional information
The numeric type T is inferred from the camera's pixeledgesx type.
Example
# Create with minimal metadata
data = BasicSMLD(emitters, camera, 10, 1)
# Create with additional metadata
data = BasicSMLD(emitters, camera, 10, 1, Dict(
"exposure_time" => 0.1,
"timestamp" => now()
))SMLMData.Emitter2D — Type
Emitter2D{T} <: AbstractEmitterRepresents a 2D emitter for SMLM simulations with position and brightness.
Fields
x::T: x-coordinate in micronsy::T: y-coordinate in micronsphotons::T: number of photons emitted by the fluorophore
SMLMData.Emitter2DFit — Type
Emitter2DFit{T} <: AbstractEmitterRepresents fitted 2D localization results with uncertainties and temporal/tracking information.
Fields
x::T: fitted x-coordinate in micronsy::T: fitted y-coordinate in micronsphotons::T: fitted number of photonsbg::T: fitted background in photons/pixelσ_x::T: uncertainty in x position in micronsσ_y::T: uncertainty in y position in micronsσ_photons::T: uncertainty in photon countσ_bg::T: uncertainty in background in photons/pixelframe::Int: frame number in acquisition sequencedataset::Int: identifier for specific acquisition/datasettrack_id::Int: identifier for linking localizations across frames (0 = unlinked)id::Int: unique identifier within dataset
SMLMData.Emitter2DFit — Method
Emitter2DFit{T}(x, y, photons, bg, σ_x, σ_y, σ_photons, σ_bg;
frame=0, dataset=1, track_id=0, id=0) where TConvenience constructor for 2D localization fit results with optional identification parameters.
Arguments
Required
x::T: fitted x-coordinate in micronsy::T: fitted y-coordinate in micronsphotons::T: fitted number of photonsbg::T: fitted background in photons/pixelσ_x::T: uncertainty in x position in micronsσ_y::T: uncertainty in y position in micronsσ_photons::T: uncertainty in photon countσ_bg::T: uncertainty in background level
Optional Keywords
frame::Int=1: frame number in acquisition sequencedataset::Int=1: identifier for specific acquisition/datasettrack_id::Int=0: identifier for linking localizations across framesid::Int=0: unique identifier within dataset
Example
# Create emitter with just required parameters
emitter = Emitter2DFit{Float64}(
1.0, 2.0, # x, y
1000.0, 10.0, # photons, background
0.01, 0.01, # σ_x, σ_y
50.0, 2.0 # σ_photons, σ_bg
)
# Create emitter with specific frame and dataset
emitter = Emitter2DFit{Float64}(
1.0, 2.0, 1000.0, 10.0, 0.01, 0.01, 50.0, 2.0;
frame=5, dataset=2
)SMLMData.Emitter3D — Type
Emitter3D{T} <: AbstractEmitterRepresents a 3D emitter for SMLM simulations with position and brightness.
Fields
x::T: x-coordinate in micronsy::T: y-coordinate in micronsz::T: z-coordinate in microns (axial position)photons::T: number of photons emitted by the fluorophore
SMLMData.Emitter3DFit — Type
Emitter3DFit{T} <: AbstractEmitterRepresents fitted 3D localization results with uncertainties and temporal/tracking information.
Fields
x::T: fitted x-coordinate in micronsy::T: fitted y-coordinate in micronsz::T: fitted z-coordinate in micronsphotons::T: fitted number of photonsbg::T: fitted background in photons/pixelσ_x::T: uncertainty in x position in micronsσ_y::T: uncertainty in y position in micronsσ_z::T: uncertainty in z position in micronsσ_photons::T: uncertainty in photon countσ_bg::T: uncertainty in background in photons/pixelframe::Int: frame number in acquisition sequencedataset::Int: identifier for specific acquisition/datasettrack_id::Int: identifier for linking localizations across frames (0 = unlinked)id::Int: unique identifier within dataset
SMLMData.Emitter3DFit — Method
Emitter3DFit{T}(x, y, z, photons, bg, σ_x, σ_y, σ_z, σ_photons, σ_bg;
frame=0, dataset=1, track_id=0, id=0) where TConvenience constructor for 3D localization fit results with optional identification parameters.
Arguments
Required
x::T: fitted x-coordinate in micronsy::T: fitted y-coordinate in micronsz::T: fitted z-coordinate in micronsphotons::T: fitted number of photonsbg::T: fitted background in photons/pixelσ_x::T: uncertainty in x position in micronsσ_y::T: uncertainty in y position in micronsσ_z::T: uncertainty in z position in micronsσ_photons::T: uncertainty in photon countσ_bg::T: uncertainty in background level
Optional Keywords
frame::Int=1: frame number in acquisition sequencedataset::Int=1: identifier for specific acquisition/datasettrack_id::Int=0: identifier for linking localizations across framesid::Int=0: unique identifier within dataset
Example
# Create emitter with just required parameters
emitter = Emitter3DFit{Float64}(
1.0, 2.0, -0.5, # x, y, z
1000.0, 10.0, # photons, background
0.01, 0.01, 0.02,# σ_x, σ_y, σ_z
50.0, 2.0 # σ_photons, σ_bg
)
# Create emitter with specific frame and tracking
emitter = Emitter3DFit{Float64}(
1.0, 2.0, -0.5, 1000.0, 10.0, 0.01, 0.01, 0.02, 50.0, 2.0;
frame=5, track_id=1
)SMLMData.IdealCamera — Type
IdealCamera{T} <: AbstractCameraRepresents an ideal camera with regularly spaced pixels defined by their edges in physical units (microns).
Fields
pixel_edges_x::Vector{T}: Physical positions of pixel edges in x direction (microns)pixel_edges_y::Vector{T}: Physical positions of pixel edges in y direction (microns)
The edges are computed from pixel centers, where pixel (1,1) is centered at (pixelsizex/2, pixelsizey/2) in physical coordinates.
SMLMData.IdealCamera — Method
IdealCamera(pixel_centers_x::AbstractUnitRange, pixel_centers_y::AbstractUnitRange,
pixel_size::Tuple{T, T}) where T<:RealConstruct an IdealCamera with rectangular pixels given pixel center positions and x,y pixel sizes.
Arguments
pixel_centers_x::AbstractUnitRange: Range of pixel center indices in x (typically 1:N)pixel_centers_y::AbstractUnitRange: Range of pixel center indices in y (typically 1:M)pixel_size::Tuple{T, T}: Tuple of (xsize, ysize) in microns
Returns
IdealCamera{T} where T matches the type of the pixel sizes
Type Parameters
T: Numeric type for all spatial measurements (e.g., Float64, Float32)
Example
# Create a 512x256 camera with rectangular pixels (0.1 x 0.15 microns)
cam = IdealCamera(1:512, 1:256, (0.1, 0.15))
# Create with Float32 precision
cam32 = IdealCamera(1:512, 1:256, (0.1f0, 0.15f0))Note: Pixel (1,1) is centered at (pixelsize[1]/2, pixelsize[2]/2) in physical coordinates.
SMLMData.IdealCamera — Method
IdealCamera(pixel_centers_x::AbstractUnitRange, pixel_centers_y::AbstractUnitRange, pixel_size::T) where T<:RealConstruct an IdealCamera with square pixels given pixel center positions and a scalar pixel size.
Arguments
pixel_centers_x::AbstractUnitRange: Range of pixel center indices in x (typically 1:N)pixel_centers_y::AbstractUnitRange: Range of pixel center indices in y (typically 1:M)pixel_size::Real: Size of pixels in microns
Returns
IdealCamera{T} where T matches the type of pixel_size
Type Parameters
T: Numeric type for all spatial measurements (e.g., Float64, Float32)
Example
# Create a 512x512 camera with 0.1 micron square pixels
cam = IdealCamera(1:512, 1:512, 0.1)
# Create with Float32 precision
cam32 = IdealCamera(1:512, 1:512, 0.1f0)Note: Pixel (1,1) is centered at (pixelsize/2, pixelsize/2) in physical coordinates.
SMLMData.IdealCamera — Method
IdealCamera(n_pixels_x::Integer, n_pixels_y::Integer, pixel_size::Tuple{T, T}) where T<:RealConstruct an IdealCamera with rectangular pixels directly from the number of pixels and x,y pixel sizes.
Arguments
n_pixels_x::Integer: Number of pixels in x dimensionn_pixels_y::Integer: Number of pixels in y dimensionpixel_size::Tuple{T, T}: Tuple of (xsize, ysize) in microns
Returns
IdealCamera{T} where T matches the type of the pixel sizes
Example
# Create a 512x256 camera with rectangular pixels (0.1 x 0.15 microns)
cam = IdealCamera(512, 256, (0.1, 0.15))
# Create with Float32 precision
cam32 = IdealCamera(512, 256, (0.1f0, 0.15f0))SMLMData.IdealCamera — Method
IdealCamera(n_pixels_x::Integer, n_pixels_y::Integer, pixel_size::T) where T<:RealConstruct an IdealCamera with square pixels directly from the number of pixels and pixel size.
Arguments
n_pixels_x::Integer: Number of pixels in x dimensionn_pixels_y::Integer: Number of pixels in y dimensionpixel_size::Real: Size of pixels in microns
Returns
IdealCamera{T} where T matches the type of pixel_size
Example
# Create a 512x512 camera with 0.1 micron square pixels
cam = IdealCamera(512, 512, 0.1)
# Create with Float32 precision
cam32 = IdealCamera(512, 512, 0.1f0)SMLMData.ROIBatch — Type
ROIBatch{T,N,A,C}Batch of regions of interest for efficient parallel processing with camera context.
This type serves as the standard interface for ROI-based processing across the JuliaSMLM ecosystem. ROIs are extracted by SMLMBoxer and consumed by fitting packages like GaussMLE.
Type Parameters
T: Element type of ROI data (typically Float32 or Float64)N: Dimension of data array (always 3 for ROI batches)A: Array type (e.g., Array, CuArray for GPU)C: Camera type (AbstractCamera subtype)
Fields
data::A- ROI image stack (roisize × roisize × n_rois)x_corners::Vector{Int32}- X (column) coordinates of ROI corners in camera coordinatesy_corners::Vector{Int32}- Y (row) coordinates of ROI corners in camera coordinatesframe_indices::Vector{Int32}- Frame number for each ROI (1-indexed)camera::C- Camera object (IdealCamera or SCMOSCamera) representing full imageroi_size::Int- Size of each ROI in pixels (assumed square)
Coordinate System
- Camera coordinates: 1-indexed, (1,1) = top-left of full image
- ROI corners: (x, y) = (col, row) position in camera coordinates
- ROI data: Local coordinates, (1,1) = top-left within ROI
- Frame indices: 1-indexed, matching camera frame numbering
Constructors
From arrays (main constructor)
ROIBatch(data::AbstractArray{T,3}, x_corners::Vector{Int32}, y_corners::Vector{Int32},
frame_indices::Vector{Int32}, camera::AbstractCamera)From separate x/y corner vectors
ROIBatch(data::AbstractArray{T,3}, x_corners::Vector, y_corners::Vector,
frame_indices::Vector, camera::AbstractCamera)From vector of SingleROI
ROIBatch(rois::Vector{SingleROI{T}}, camera::AbstractCamera)Validation
- ROIs must be square (data dimensions 1 == dimension 2)
- xcorners must have length nrois
- ycorners must have length nrois
- Frame indices must have length n_rois
- All arrays must have consistent n_rois
Indexing and Iteration
batch = ROIBatch(data, x_corners, y_corners, frames, camera)
# Get single ROI
roi = batch[5] # Returns SingleROI{T}
# Iterate over all ROIs
for roi in batch
process(roi.data, roi.corner, roi.frame_idx)
end
# Length
n = length(batch) # Number of ROIsGPU Support
Supports GPU transfer via Adapt.jl (KernelAbstractions.jl):
using CUDA
batch_gpu = adapt(CuArray, batch) # Transfer data to GPU
# Camera stays on host (contains metadata)Example
using SMLMData
using StaticArrays
# Create camera
camera = IdealCamera(512, 512, 0.1) # 512×512 pixels, 0.1μm/pixel
# Create ROI batch (e.g., from SMLMBoxer.getboxes)
n_rois = 100
roi_size = 11
data = rand(Float32, roi_size, roi_size, n_rois)
x_corners = rand(Int32(1):Int32(500), n_rois)
y_corners = rand(Int32(1):Int32(500), n_rois)
frames = rand(Int32(1):Int32(10), n_rois)
batch = ROIBatch(data, x_corners, y_corners, frames, camera)
# Access
println("Batch contains $(length(batch)) ROIs")
first_roi = batch[1]
println("First ROI at position: $(first_roi.corner)")See Also
SingleROI- Individual ROI typeIdealCamera,SCMOSCamera- Camera types
SMLMData.ROIBatch — Method
ROIBatch(data, x_corners, y_corners, frame_indices, camera)Construct ROIBatch from separate x and y corner vectors.
Arguments
data::AbstractArray{T,3}- ROI stack (roisize × roisize × n_rois)x_corners::Vector- X (column) coordinates of ROI cornersy_corners::Vector- Y (row) coordinates of ROI cornersframe_indices::Vector- Frame number for each ROIcamera::AbstractCamera- Camera object for full image
Example
batch = ROIBatch(data, x_corners, y_corners, frames, camera)SMLMData.ROIBatch — Method
ROIBatch(rois::Vector{SingleROI{T}}, camera)Construct ROIBatch from vector of SingleROI objects.
Arguments
rois::Vector{SingleROI{T}}- Vector of individual ROIscamera::AbstractCamera- Camera object for full image
Returns
ROIBatch{T,3,Array{T,3},typeof(camera)}
Example
rois = [SingleROI(rand(Float32, 11, 11), SVector{2,Int32}(i*10, i*10), Int32(i))
for i in 1:100]
batch = ROIBatch(rois, camera)SMLMData.SCMOSCamera — Type
SCMOSCamera{T<:Real} <: AbstractCamerasCMOS camera with pixel-dependent calibration parameters matching spec sheets.
Fields
pixel_edges_x::Vector{T}: Physical pixel edges in x (μm)pixel_edges_y::Vector{T}: Physical pixel edges in y (μm)offset::Union{T, Matrix{T}}: Dark level (ADU)gain::Union{T, Matrix{T}}: Conversion gain (e⁻/ADU)readnoise::Union{T, Matrix{T}}: Read noise (e⁻ rms)qe::Union{T, Matrix{T}}: Quantum efficiency (0-1)
Units
Calibration parameters follow camera specification sheet conventions:
offset: ADU (analog-to-digital units)
- Typical values: 100-500 ADU
- Dark level with no illumination
gain: e⁻/ADU (electrons per ADU)
- Typical values: 0.1-2.0 e⁻/ADU depending on readout mode
- Example: ORCA-Flash4.0: 0.46 e⁻/ADU (12-bit), 0.11 e⁻/ADU (16-bit)
readnoise: e⁻ rms (electrons, root-mean-square)
- Typical values: 0.3-5.0 e⁻ rms
- Example: ORCA-Flash4.0 V3: 1.6 e⁻ rms
- Example: ORCA-Quest qCMOS: 0.27 e⁻ rms
qe: dimensionless (0 to 1)
- Typical values: 0.5-0.95 at peak wavelength
- Example: ORCA-Flash4.0 V2: 0.72 at 550nm
- Example: ORCA-Fusion BT: 0.95 (back-thinned)
Physical Signal Chain
Photons → Electrons → ADU:
Incident photons (N)
↓ [× QE]
Photoelectrons (N × QE)
↓ [+ readnoise (Gaussian)]
Signal electrons (N × QE + ε), where ε ~ N(0, readnoise²)
↓ [÷ gain, + offset]
ADU readout = (N × QE + ε)/gain + offsetScalar vs Matrix Parameters
Each calibration parameter can be:
- Scalar (T): Uniform across sensor (approximation or post-calibration)
- Matrix (Matrix{T}): Per-pixel calibration map (size must match pixel grid)
Use matrices for:
- Precision SMLM (2-5% variations can affect results)
- Quantitative imaging
- Artifact correction
Use scalars for:
- Quick analysis
- Post-calibrated data
- Uniform approximation
Constructors
# Minimal - most common case (requires readnoise, others default to 0, 1, 1)
cam = SCMOSCamera(512, 512, 0.1, 1.6)
# With additional parameters
cam = SCMOSCamera(512, 512, 0.1, readnoise_map,
offset=100.0, gain=0.46, qe=0.72)
# Custom edges (advanced)
cam = SCMOSCamera(custom_edges_x, custom_edges_y,
readnoise=noise_map, gain=gain_map)Examples
# Example 1: From spec sheet (ORCA-Flash4.0 V3, 12-bit mode)
cam = SCMOSCamera(
2048, 2048, 0.065, # 2048×2048 pixels, 65nm pixel size
1.6, # From spec: 1.6 e⁻ rms readnoise
offset = 100.0, # Typical offset
gain = 0.46, # From spec: 0.46 e⁻/ADU
qe = 0.72 # 72% QE at 550nm
)
# Example 2: With calibration maps (precision SMLM)
readnoise_map = load("camera_noise.mat") # 512×512 measured values
gain_map = load("camera_gain.mat")
qe_map = load("camera_qe.mat")
cam = SCMOSCamera(
512, 512, 0.1, readnoise_map,
gain = gain_map,
qe = qe_map
)
# Example 3: Minimal (variance-only approximation)
# Common when you only have noise map, assume ideal otherwise
cam = SCMOSCamera(512, 512, 0.1, readnoise_map)
# Example 4: Ultra-low noise camera (ORCA-Quest)
cam = SCMOSCamera(
2304, 4096, 0.0044, # 4.4μm pixels
0.27, # Incredible 0.27 e⁻ rms!
offset = 100.0,
gain = 0.5,
qe = 0.85
)
# Example 5: Rectangular pixels
cam = SCMOSCamera(512, 256, (0.1, 0.15), 1.8)
# Example 6: Mixed scalar/matrix parameters
cam = SCMOSCamera(
512, 512, 0.1, readnoise_map, # Per-pixel noise
offset = 100.0, # Uniform offset
gain = 0.5, # Uniform gain
qe = qe_map # Per-pixel QE
)See Also
IdealCamera for Poisson-only noise (readnoise=0)
SMLMData.SCMOSCamera — Method
SCMOSCamera(nx, ny, pixel_size, readnoise; offset=0, gain=1, qe=1)Construct sCMOS camera from pixel dimensions and calibration parameters.
Arguments
nx::Integer: Number of pixels in xny::Integer: Number of pixels in ypixel_size::Union{T, Tuple{T,T}}: Pixel size in μm (scalar or (xsize, ysize))readnoise::Union{T, Matrix{T}}: Read noise in e⁻ rms (required)
Keywords
offset::Union{T, Matrix{T}} = 0: Dark level in ADUgain::Union{T, Matrix{T}} = 1: Conversion gain in e⁻/ADUqe::Union{T, Matrix{T}} = 1: Quantum efficiency (0-1)
Each parameter can be scalar (uniform) or Matrix{T} with size (nx, ny).
Examples
# Minimal: just readnoise (assumes calibrated data: offset=0, gain=1, qe=1)
cam = SCMOSCamera(512, 512, 0.1, 1.6)
# From spec sheet (ORCA-Flash4.0 V3)
cam = SCMOSCamera(2048, 2048, 0.065, 1.6, offset=100.0, gain=0.46, qe=0.72)
# With calibration maps
cam = SCMOSCamera(512, 512, 0.1, readnoise_map,
offset=offset_map, gain=gain_map, qe=qe_map)
# Rectangular pixels
cam = SCMOSCamera(512, 256, (0.1, 0.15), 1.8)SMLMData.SCMOSCamera — Method
SCMOSCamera(pixel_edges_x, pixel_edges_y; offset=0, gain=1, readnoise, qe=1)Construct sCMOS camera with custom pixel edge positions.
Arguments
pixel_edges_x::Vector{T}: Pixel edge positions in x (μm), length nx+1pixel_edges_y::Vector{T}: Pixel edge positions in y (μm), length ny+1
Keywords
readnoise::Union{T, Matrix{T}}: Read noise in e⁻ rms (required)offset::Union{T, Matrix{T}} = 0: Dark level in ADUgain::Union{T, Matrix{T}} = 1: Conversion gain in e⁻/ADUqe::Union{T, Matrix{T}} = 1: Quantum efficiency (0-1)
Matrix parameters must have size (nx, ny) where nx = length(pixeledgesx) - 1.
Example
# Custom non-uniform pixel grid
edges_x = [0.0, 0.1, 0.21, 0.33, 0.46] # Non-uniform spacing
edges_y = [0.0, 0.1, 0.2, 0.3]
cam = SCMOSCamera(edges_x, edges_y, readnoise=1.5, gain=0.5)SMLMData.SMLD — Type
SMLDDeprecated alias for AbstractSMLD. Use AbstractSMLD instead. This alias is provided for backward compatibility and will be removed in v1.0.
SMLMData.SingleROI — Type
SingleROI{T}Single region of interest (ROI) with location and frame context.
Fields
data::Matrix{T}- ROI image data (roisize × roisize pixels)corner::SVector{2,Int32}- ROI corner position (x, y) = (col, row) in camera pixels (1-indexed)frame_idx::Int32- Frame number in image stack (1-indexed)
Coordinate System
- Camera coordinates: 1-indexed, (1,1) is top-left pixel of full image
- Corner: (x, y) = (col, row) position where top-left of ROI starts
- ROI data: (1,1) is top-left pixel within the ROI itself
Example
# Create single ROI
roi_data = rand(Float32, 11, 11) # 11×11 pixel ROI
roi = SingleROI(roi_data, SVector{2,Int32}(100, 200), Int32(5))
# Access fields
println("ROI at camera position: ", roi.corner) # (100, 200) = (col, row)
println("From frame: ", roi.frame_idx) # Frame 5SMLMData.SmiteSMD — Type
SmiteSMDHelper structure for loading Smite SMD .mat files.
Fields
filepath::String: Path to the directory containing the .mat filefilename::String: Name of the .mat filevarname::String: Variable name in the .mat file (default: "SMD")
Example
# Load from default "SMD" variable
smd = SmiteSMD("path/to/data", "localizations.mat")
# Load from custom variable name
smd = SmiteSMD("path/to/data", "localizations.mat", "CustomSMD")SMLMData.SmiteSMLD — Type
SmiteSMLD{T,E<:AbstractEmitter} <: AbstractSMLDSMLD type compatible with the Smite SMD (Single Molecule Data) format.
Fields
emitters::Vector{E}: Vector of localized emitterscamera::AbstractCamera: Camera used for acquisitionn_frames::Int: Total number of frames in acquisitionn_datasets::Int: Number of datasets in the acquisitionmetadata::Dict{String,Any}: Additional dataset information
Type Parameters
T: Numeric type for coordinates (typically Float64)E: Concrete emitter type (typically Emitter2DFit or Emitter3DFit)
Adapt.adapt_structure — Method
Adapt.adapt_structure(to, batch::ROIBatch)Adapt ROIBatch for GPU execution via KernelAbstractions.jl/CUDA.jl.
Transfers data, xcorners, ycorners, and frame_indices to the target device. Camera remains on the host (contains metadata and variance maps).
Example
using CUDA
batch_gpu = adapt(CuArray, batch)
# Process on GPU...
batch_cpu = adapt(Array, batch_gpu) # Transfer backBase.getindex — Method
getindex(batch::ROIBatch, i::Int) -> SingleROIGet the i-th ROI from the batch.
Returns a SingleROI containing the data, corner position, and frame index.
Base.iterate — Function
iterate(batch::ROIBatch, [state]) -> Union{Nothing, Tuple{SingleROI, Int}}Iterate over ROIs in the batch.
Example
for roi in batch
println("Processing ROI at ", roi.corner)
endBase.iterate — Method
Base.iterate(smld::AbstractSMLD)
Base.iterate(smld::AbstractSMLD, state)Enable iteration over emitters in an SMLD object.
Base.length — Method
Base.length(smld::AbstractSMLD)Return the number of emitters in the SMLD object.
Base.length — Method
length(batch::ROIBatch) -> IntNumber of ROIs in the batch.
SMLMData.api_overview — Method
SMLMData.jl API Overview
This guide provides a structured overview of the SMLMData.jl package designed for Single Molecule Localization Microscopy (SMLM) data handling in Julia.
Why This Overview Exists
For Humans
- Provides a concise reference without diving into full documentation
- Offers quick-start examples for common use cases
- Shows relevant patterns more clearly than individual docstrings
- Creates an at-a-glance understanding of package capabilities
For AI Assistants
- Enables better code generation with correct API patterns
- Provides structured context about type hierarchies and relationships
- Offers consistent examples to learn from when generating code
- Helps avoid common pitfalls or misunderstandings about the API
Ecosystem Role
SMLMData is the core types package for the JuliaSMLM ecosystem. Other packages depend on SMLMData and re-export its types, so you rarely need to import SMLMData directly:
using GaussMLE # Re-exports ROIBatch, camera types, etc.
using SMLMAnalysis # Re-exports all SMLMData types for analysis workflowsDirect using SMLMData is primarily for package developers, standalone data manipulation, or learning the type system.
Key Concepts
- Emitters: Individual fluorophore localizations (2D or 3D)
- Camera: Defines pixel geometry and coordinate system
- AbstractSMLD: Container holding emitters and camera information
- Coordinates: All spatial coordinates are in microns
- Coordinate System:
- Physical space: (0,0) at top-left corner of camera
- Pixel space: (1,1) at center of top-left pixel
Type Hierarchy
AbstractEmitter # Base for all emitter types
├── Emitter2D{T} # Basic 2D emitters
├── Emitter3D{T} # Basic 3D emitters
├── Emitter2DFit{T} # 2D emitters with fit results
└── Emitter3DFit{T} # 3D emitters with fit results
AbstractCamera # Base for all camera types
├── IdealCamera{T} # Camera with regular pixel grid (Poisson noise only)
└── SCMOSCamera{T} # sCMOS camera with pixel-dependent calibration
ROI Batch Types # For batched ROI processing
├── SingleROI{T} # Single ROI with location context
└── ROIBatch{T,N,A,C} # Batch of ROIs for parallel processing
AbstractSMLD # Base for data containers
├── BasicSMLD{T,E} # General-purpose container
└── SmiteSMLD{T,E} # SMITE-compatible containerEssential Types
Emitter Types
# Basic 2D emitter
mutable struct Emitter2D{T} <: AbstractEmitter
x::T # x-coordinate in microns
y::T # y-coordinate in microns
photons::T # number of photons emitted
end
# Basic 3D emitter
mutable struct Emitter3D{T} <: AbstractEmitter
x::T # x-coordinate in microns
y::T # y-coordinate in microns
z::T # z-coordinate in microns
photons::T # number of photons emitted
end
# 2D emitter with fit results
mutable struct Emitter2DFit{T} <: AbstractEmitter
x::T # fitted x-coordinate in microns
y::T # fitted y-coordinate in microns
photons::T # fitted number of photons
bg::T # fitted background in photons/pixel
σ_x::T # uncertainty in x position in microns
σ_y::T # uncertainty in y position in microns
σ_photons::T # uncertainty in photon count
σ_bg::T # uncertainty in background level
frame::Int # frame number in acquisition sequence
dataset::Int # identifier for specific acquisition/dataset
track_id::Int # identifier for linking localizations across frames
id::Int # unique identifier within dataset
end
# 3D emitter with fit results
mutable struct Emitter3DFit{T} <: AbstractEmitter
x::T # fitted x-coordinate in microns
y::T # fitted y-coordinate in microns
z::T # fitted z-coordinate in microns
photons::T # fitted number of photons
bg::T # fitted background in photons/pixel
σ_x::T # uncertainty in x position in microns
σ_y::T # uncertainty in y position in microns
σ_z::T # uncertainty in z position in microns
σ_photons::T # uncertainty in photon count
σ_bg::T # uncertainty in background level
frame::Int # frame number in acquisition sequence
dataset::Int # identifier for specific acquisition/dataset
track_id::Int # identifier for linking localizations across frames
id::Int # unique identifier within dataset
endEmitter Constructor Examples
# Basic 2D emitter
emitter_2d = Emitter2D{Float64}(
1.5, # x-coordinate in microns
2.3, # y-coordinate in microns
1000.0 # number of photons emitted
)
# Basic 3D emitter
emitter_3d = Emitter3D{Float64}(
1.5, # x-coordinate in microns
2.3, # y-coordinate in microns
-0.5, # z-coordinate in microns (negative = below focal plane)
1000.0 # number of photons emitted
)
# 2D emitter with fit results using convenience constructor
emitter_2d_fit = Emitter2DFit{Float64}(
1.5, 2.3, # x, y coordinates in microns
1000.0, 10.0, # photons detected, background photons/pixel
0.01, 0.01, # σ_x, σ_y: position uncertainties in microns
50.0, 2.0; # σ_photons, σ_bg: photon count uncertainties
frame=5, # frame number in acquisition (1-based, default=1)
dataset=1, # dataset identifier for multi-acquisition experiments
track_id=2, # tracking ID for linked localizations (default=0 = unlinked)
id=42 # unique identifier within this dataset (default=0)
)Camera Types
# Ideal camera with uniform pixel grid (Poisson noise only)
struct IdealCamera{T} <: AbstractCamera
pixel_edges_x::Vector{T} # pixel edges in x
pixel_edges_y::Vector{T} # pixel edges in y
end
# sCMOS camera with pixel-dependent calibration parameters
struct SCMOSCamera{T} <: AbstractCamera
pixel_edges_x::Vector{T} # pixel edges in x
pixel_edges_y::Vector{T} # pixel edges in y
offset::Union{T, Matrix{T}} # dark level (ADU)
gain::Union{T, Matrix{T}} # conversion gain (e⁻/ADU)
readnoise::Union{T, Matrix{T}} # read noise (e⁻ rms)
qe::Union{T, Matrix{T}} # quantum efficiency (0-1)
endCamera Constructor Examples
# IdealCamera: Create a camera with 512x512 pixels, each 100nm (0.1μm) in size
# Convenience constructor (most common)
cam = IdealCamera(512, 512, 0.1)
# Explicit constructor using pixel center ranges
cam_explicit = IdealCamera(1:512, 1:512, 0.1)
# For non-square pixels, specify different x and y sizes
cam_rect = IdealCamera(512, 512, (0.1, 0.12))
# SCMOSCamera: Create with readnoise specification (matching spec sheets)
# Minimal (uniform readnoise, assumes offset=0, gain=1, qe=1)
cam_scmos = SCMOSCamera(512, 512, 0.1, 1.6) # 1.6 e⁻ rms readnoise
# From camera spec sheet (e.g., ORCA-Flash4.0 V3)
cam_flash = SCMOSCamera(
2048, 2048, 0.065, # 2048×2048 pixels, 65nm pixel size
1.6, # 1.6 e⁻ rms readnoise from spec
offset = 100.0, # typical dark level
gain = 0.46, # 0.46 e⁻/ADU from spec
qe = 0.72 # 72% QE at 550nm
)
# With per-pixel calibration maps (precision SMLM)
readnoise_map = load("camera_noise.mat") # 512×512 measured values
gain_map = load("camera_gain.mat")
qe_map = load("camera_qe.mat")
cam_calibrated = SCMOSCamera(512, 512, 0.1, readnoise_map,
gain=gain_map, qe=qe_map)
# Mixed scalar and matrix parameters
cam_mixed = SCMOSCamera(
512, 512, 0.1, readnoise_map, # Per-pixel noise
offset = 100.0, # Uniform offset
gain = 0.5, # Uniform gain
qe = qe_map # Per-pixel QE
)ROI Batch Types
ROI batch types provide efficient storage and processing of image regions across the JuliaSMLM ecosystem.
# Single ROI with location context
struct SingleROI{T}
data::Matrix{T} # ROI image data (roi_size × roi_size)
corner::SVector{2,Int32} # (x, y) = (col, row) corner position (1-indexed)
frame_idx::Int32 # Frame number (1-indexed)
end
# Batch of ROIs for parallel processing
struct ROIBatch{T,N,A<:AbstractArray{T,N},C<:AbstractCamera}
data::A # ROI stack (roi_size × roi_size × n_rois)
x_corners::Vector{Int32} # X (column) coordinates of ROI corners
y_corners::Vector{Int32} # Y (row) coordinates of ROI corners
frame_indices::Vector{Int32} # Frame number for each ROI
camera::C # Camera object (IdealCamera or SCMOSCamera)
roi_size::Int # Size of each ROI (square)
endROI Batch Constructor Examples
# From separate x/y corner vectors (main constructor)
camera = IdealCamera(512, 512, 0.1)
data = rand(Float32, 11, 11, 100) # 100 ROIs of 11×11 pixels
x_corners = rand(Int32(1):Int32(500), 100)
y_corners = rand(Int32(1):Int32(500), 100)
frame_indices = rand(Int32(1):Int32(50), 100)
batch = ROIBatch(data, x_corners, y_corners, frame_indices, camera)
# From vector of SingleROI
rois = [SingleROI(rand(Float32, 11, 11), SVector{2,Int32}(i*10, i*10), Int32(i))
for i in 1:100]
batch = ROIBatch(rois, camera)
# Indexing and iteration
roi = batch[5] # Get single ROI
for roi in batch
process(roi.data) # Iterate over all ROIs
end
# GPU adaptation (via Adapt.jl)
using CUDA
batch_gpu = adapt(CuArray, batch) # Transfer to GPUCoordinate System:
- Camera coordinates: 1-indexed, (1,1) = top-left of full image
- ROI corners: (x, y) = (col, row) position in camera coordinates
- ROI data: Local coordinates, (1,1) = top-left within ROI
- Frame indices: 1-indexed, matching camera frame numbering
Typical Workflow:
- SMLMBoxer extracts ROIs →
ROIBatch - GaussMLE fits ROIs →
LocalizationResult - Convert to
BasicSMLDfor analysis
SMLD Container Types
# Basic SMLD container
struct BasicSMLD{T,E<:AbstractEmitter} <: AbstractSMLD
emitters::Vector{E} # Vector of emitters
camera::AbstractCamera # Camera information
n_frames::Int # Total number of frames
n_datasets::Int # Number of datasets
metadata::Dict{String,Any} # Additional information
end
# SMITE format compatible container
struct SmiteSMLD{T,E<:AbstractEmitter} <: AbstractSMLD
emitters::Vector{E} # Vector of emitters
camera::AbstractCamera # Camera information
n_frames::Int # Total number of frames
n_datasets::Int # Number of datasets
metadata::Dict{String,Any} # Additional information
endSMLD Constructor Examples
# Create a vector of emitters
emitters = [
Emitter2DFit{Float64}(1.0, 2.0, 1000.0, 10.0, 0.01, 0.01, 50.0, 2.0),
Emitter2DFit{Float64}(3.0, 4.0, 1200.0, 12.0, 0.01, 0.01, 60.0, 2.0)
]
# Create a BasicSMLD
smld = BasicSMLD(emitters, camera, 1, 1, Dict{String,Any}())
# Add metadata
smld_with_metadata = BasicSMLD(
emitters,
camera,
10, # number of frames
1, # number of datasets
Dict{String,Any}(
"exposure_time" => 0.1,
"sample" => "Test Sample"
)
)Core Functions
Accessing the API Overview
# Get this API overview as a string programmatically
overview_text = api_overview()Coordinate Conversions
# Convert from pixel to physical coordinates (microns)
x_physical, y_physical = pixel_to_physical(px, py, pixel_size)
# Convert from physical to pixel coordinates
px, py = physical_to_pixel(x, y, pixel_size)
# Convert from physical to pixel indices (integers)
px_idx, py_idx = physical_to_pixel_index(x, y, pixel_size)
# Get physical coordinates of all pixel centers
centers_x, centers_y = get_pixel_centers(camera)Filtering Operations
# Filter by emitter properties using @filter macro
bright = @filter(smld, photons > 1000) # Select bright emitters
precise = @filter(smld, σ_x < 0.02 && σ_y < 0.02) # Select precisely localized emitters
combined = @filter(smld, photons > 1000 && σ_x < 0.02) # Combine multiple criteria
# The @filter macro supports any emitter property:
# - Basic: x, y, z (for 3D), photons
# - Fit results: bg, σ_x, σ_y, σ_z, σ_photons, σ_bg
# - Metadata: frame, dataset, track_id, id
# Select frames
frame_5 = filter_frames(smld, 5) # Single frame
early_frames = filter_frames(smld, 1:10) # Range of frames (inclusive)
specific_frames = filter_frames(smld, [1,3,5,7]) # Specific frames (uses Set for efficiency)
# Select region of interest (ROI) - coordinates in microns
# 2D ROI
roi_2d = filter_roi(smld, 1.0:5.0, 2.0:6.0) # x_range, y_range
# 3D ROI (for 3D emitters only)
roi_3d = filter_roi(smld, 1.0:5.0, 2.0:6.0, -1.0:1.0) # x, y, z rangesSMLD Operations
# Concatenate multiple SMLDs
combined = cat_smld(smld1, smld2)
combined = cat_smld([smld1, smld2, smld3])
# Merge with options to adjust frame and dataset numbering
merged = merge_smld(smld1, smld2)
merged = merge_smld([smld1, smld2, smld3])
# Merge with sequential frame numbers
sequential = merge_smld([smld1, smld2, smld3], adjust_frames=true)
# Merge with sequential dataset numbers
sequential_ds = merge_smld([smld1, smld2, smld3], adjust_datasets=true)I/O Operations
# Import from SMITE format (MATLAB)
smd = SmiteSMD("path/to/data", "localizations.mat") # Default variable name "SMD"
smd = SmiteSMD("path/to/data", "localizations.mat", "CustomSMD") # Custom variable name
# Load as 2D or 3D data
smld_2d = load_smite_2d(smd)
smld_3d = load_smite_3d(smd)
# Export to SMITE format (saved as MATLAB v7.3 format)
# Note: requires SmiteSMLD object, not BasicSMLD
smite_smld = SmiteSMLD(smld.emitters, smld.camera, smld.n_frames, smld.n_datasets, smld.metadata)
save_smite(smite_smld, "output/directory", "results.mat")Note: The SMITE loader automatically handles complex-valued fields by removing emitters with non-zero imaginary components in key fields (X, Y, Z, Photons, background, and uncertainties). Information about removed emitters is stored in the metadata as "removed_complex_emitters" => count.
Working with SMLD Objects
# Get number of emitters
n_emitters = length(smld)
# Iterate over emitters
for emitter in smld
println("Emitter at ($(emitter.x), $(emitter.y)) with $(emitter.photons) photons")
end
# Display formatted information
show(smld) # Compact view
show(stdout, MIME("text/plain"), smld) # Detailed viewCommon Workflows
Creating and Working with Emitters
# Create emitters
emitter1 = Emitter2DFit{Float64}(1.0, 2.0, 1000.0, 10.0, 0.01, 0.01, 50.0, 2.0)
emitter2 = Emitter2DFit{Float64}(3.0, 4.0, 1200.0, 12.0, 0.01, 0.01, 60.0, 2.0)
# Create camera
cam = IdealCamera(512, 512, 0.1) # 512x512 camera with 0.1 micron pixels
# Create SMLD container
emitters = [emitter1, emitter2]
smld = BasicSMLD(emitters, cam, 1, 1, Dict{String,Any}())Loading and Filtering Data
# Load from SMITE format
smd = SmiteSMD("data_directory", "localizations.mat")
smld = load_smite_2d(smd)
# Filter by quality
good_fits = @filter(smld, σ_x < 0.02 && σ_y < 0.02 && photons > 500)
# Filter by ROI
roi = filter_roi(good_fits, 10.0:20.0, 10.0:20.0)
# Filter by frames
frames_1_10 = filter_frames(roi, 1:10)Multi-Dataset Analysis
# Load multiple datasets
smd1 = SmiteSMD("experiment1", "data.mat")
smd2 = SmiteSMD("experiment2", "data.mat")
smld1 = load_smite_2d(smd1)
smld2 = load_smite_2d(smd2)
# Filter each dataset
bright1 = @filter(smld1, photons > 1000)
bright2 = @filter(smld2, photons > 1000)
# Merge datasets with sequential frame numbering
merged = merge_smld([bright1, bright2], adjust_frames=true)
# Process the merged dataset
result = @filter(merged, σ_x < 0.02 && σ_y < 0.02)
# Save the results (convert to SmiteSMLD first if needed)
result_smite = SmiteSMLD(result.emitters, result.camera, result.n_frames, result.n_datasets, result.metadata)
save_smite(result_smite, "analysis_results", "merged_filtered.mat")Complete Example
using SMLMData
# 1. Create a camera with 100nm pixels
# Camera has 512x512 pixels, each 0.1 microns (100nm) in size
cam = IdealCamera(512, 512, 0.1) # Using convenience constructor
# 2. Create emitters representing single molecule localizations
emitters = [
# Emitter at (1.0, 2.0) μm with high precision
Emitter2DFit{Float64}(1.0, 2.0, 1000.0, 10.0, 0.01, 0.01, 50.0, 2.0),
# Bright emitter at (3.0, 4.0) μm
Emitter2DFit{Float64}(3.0, 4.0, 1200.0, 12.0, 0.01, 0.01, 60.0, 2.0),
# Dimmer emitter at (5.0, 6.0) μm with lower precision
Emitter2DFit{Float64}(5.0, 6.0, 800.0, 9.0, 0.03, 0.03, 40.0, 1.5)
]
# 3. Create SMLD container to hold all data
smld = BasicSMLD(
emitters, # Vector of emitters
cam, # Camera geometry
1, # Number of frames
1, # Number of datasets
Dict{String,Any}("sample" => "Test") # Metadata
)
# 4. Filter by photons to select bright emitters
bright = @filter(smld, photons > 900) # Creates new SMLD with filtered emitters
# 5. Select region of interest (ROI)
# Select emitters in rectangular region: x ∈ [0, 4] μm, y ∈ [1, 5] μm
roi = filter_roi(bright, 0.0:4.0, 1.0:5.0)
# 6. Examine the results
println("Original dataset: $(length(smld)) emitters")
println("After filtering by brightness: $(length(bright)) emitters")
println("After ROI selection: $(length(roi)) emitters")
# 7. Access individual emitters
for (i, emitter) in enumerate(roi)
println("Emitter $i: position=($(emitter.x), $(emitter.y)) μm, photons=$(emitter.photons)")
end
# Output:
# Original dataset: 3 emitters
# After filtering by brightness: 2 emitters
# After ROI selection: 2 emitters
# Emitter 1: position=(1.0, 2.0) μm, photons=1000.0
# Emitter 2: position=(3.0, 4.0) μm, photons=1200.0Common Pitfalls and Important Notes
Coordinate System
- Physical coordinates are always in microns, not nanometers or pixels
- Pixel indices start at 1 (Julia convention), not 0
- Frame numbers start at 1 (default=1, following Julia's 1-based indexing convention)
- The origin (0,0) in physical space is at the top-left corner of the camera
Type Stability
- When creating emitters, ensure all numeric fields use the same type (e.g., all
Float64) - The
BasicSMLDconstructor automatically infers typeTfrom the camera's pixel edges - Mixing types (e.g.,
Float32andFloat64) can lead to performance issues
Filtering
- The
@filtermacro creates a new SMLD object; it doesn't modify the original - Filtering by frames with a vector uses
Setinternally for O(1) lookup performance - Applying an ROI filter to incompatible emitter types will throw an error
SMITE Format
- Complex-valued fields in SMITE files are automatically handled by removing affected emitters
- The loader adds metadata about removed emitters:
"removed_complex_emitters" => count - SMITE files are saved in MATLAB v7.3 format (HDF5-based)
Memory Considerations
- Large datasets benefit from using appropriate numeric types (e.g.,
Float32vsFloat64) - The
filter_framesfunction with specific frame lists is optimized for sparse selections - Iterating over emitters is memory-efficient (doesn't create intermediate arrays)
Common Mistakes
# WRONG: Using pixel units instead of microns
emitter = Emitter2D{Float64}(100, 200, 1000.0) # ❌ Likely pixel coordinates
# CORRECT: Using micron coordinates
emitter = Emitter2D{Float64}(10.0, 20.0, 1000.0) # ✓ Physical coordinates
# WRONG: Modifying original SMLD
bright = @filter(smld, photons > 1000)
# smld is unchanged!
# CORRECT: Working with the filtered result
bright = @filter(smld, photons > 1000)
# Use 'bright' for further analysisapi_overview() returns this documentation as a plain String.
SMLMData.cat_smld — Method
cat_smld(smlds::Vector{<:AbstractSMLD})
cat_smld(smlds::AbstractSMLD...)Concatenate multiple SMLD objects into a single SMLD.
Arguments
smlds: Vector of AbstractSMLD objects or multiple AbstractSMLD arguments
Returns
New AbstractSMLD containing all emitters from inputs
Notes
- Camera must be identical across all SMLDs
- n_frames is set to maximum frame number across all inputs
- n_datasets is set to maximum dataset number across all inputs
- Metadata from first SMLD is used, with conflicts noted in metadata
Examples
# Concatenate two SMLDs
combined = cat_smld(smld1, smld2)
# Concatenate multiple SMLDs
combined = cat_smld(smld1, smld2, smld3)
# Concatenate vector of SMLDs
combined = cat_smld([smld1, smld2, smld3])SMLMData.check_complex_fields — Method
check_complex_fields(s, fields)Check if any of the given fields in s are complex and have non-zero imaginary components. Returns a tuple with:
- Boolean indicating if any fields are complex with non-zero imaginary parts
- Dict mapping field names to arrays of indices with non-zero imaginary parts
SMLMData.compute_bin_edges — Method
compute_bin_edges(centers_x::AbstractUnitRange, centers_y::AbstractUnitRange, pixel_size::Tuple{Real, Real})Compute pixel edges in both dimensions for rectangular pixels.
Arguments
centers_x::AbstractUnitRange: Range of pixel center indices in xcenters_y::AbstractUnitRange: Range of pixel center indices in ypixel_size::Tuple{Real, Real}: Tuple of (xsize, ysize) in microns
Returns
Tuple{Vector{Float64}, Vector{Float64}}: (edgesx, edgesy) in physical units (microns)
SMLMData.compute_bin_edges — Method
compute_bin_edges(centers_x::AbstractUnitRange, centers_y::AbstractUnitRange, pixel_size::T) where TCompute pixel edges in both dimensions. Returns vectors with same type as pixel_size.
SMLMData.compute_edges_1d — Method
compute_edges_1d(centers::AbstractUnitRange, pixel_size::T) where T<:RealCompute pixel edges in one dimension. Maintains the numeric type of pixelsize. The first edge starts at 0 and each pixel has width pixelsize.
Arguments
centers::AbstractUnitRange: Range of pixel center indicespixel_size::T: Size of pixels in microns
Returns
Vector{T}: Edge positions in physical units (microns), starting at 0
SMLMData.filter_frames — Method
filter_frames(smld::AbstractSMLD, frame::Integer)
filter_frames(smld::AbstractSMLD, frames::Union{AbstractVector,AbstractRange})Efficiently select emitters from specified frames.
Arguments
smld::AbstractSMLD: Input SMLD structureframes: Single frame number, vector of frame numbers, or range of frames
Returns
New AbstractSMLD containing only emitters from specified frames
Examples
# Single frame
frame_5 = filter_frames(smld, 5)
# Range of frames
early = filter_frames(smld, 1:10)
# Multiple specific frames
selected = filter_frames(smld, [1,3,5,7])SMLMData.filter_roi — Method
filter_roi(smld::AbstractSMLD, x_range, y_range)
filter_roi(smld::AbstractSMLD, x_range, y_range, z_range)Efficiently select emitters within a region of interest.
Arguments
smld::AbstractSMLD: Input SMLD structurex_range: Range or tuple for x coordinates (microns)y_range: Range or tuple for y coordinates (microns)z_range: Optional range or tuple for z coordinates (microns)
Returns
New AbstractSMLD containing only emitters within the specified ROI
Examples
# 2D ROI
region = filter_roi(smld, 1.0:5.0, 2.0:6.0)
region = filter_roi(smld, (1.0, 5.0), (2.0, 6.0))
# 3D ROI
volume = filter_roi(smld, 1.0:5.0, 2.0:6.0, -1.0:1.0)SMLMData.format_with_commas — Method
format_with_commas(n::Integer)Format an integer with thousands separators for better readability.
SMLMData.get_pixel_centers — Method
get_pixel_centers(cam::AbstractCamera)Calculate the physical coordinates of all pixel centers for any camera type.
For each dimension, the center positions are computed as the midpoint between consecutive edge positions. This works for both regular (uniform pixel size) and irregular (varying pixel size) cameras.
Arguments
cam::AbstractCamera: Any camera type that implements the AbstractCamera interface with pixeledgesx and pixeledgesy fields in physical units (microns)
Returns
Tuple{Vector, Vector}: (centersx, centersy) where each vector contains the physical coordinates (in microns) of pixel centers along that dimension
Example
# For a 512x512 camera with 0.1 micron pixels
cam = IdealCamera(1:512, 1:512, 0.1)
centers_x, centers_y = get_pixel_centers(cam)
# First pixel center should be at (0.05, 0.05) microns
@assert centers_x[1] ≈ 0.05
@assert centers_y[1] ≈ 0.05SMLMData.get_valid_indices — Method
get_valid_indices(s, complex_indices)Get indices of elements that don't have complex values with non-zero imaginary parts in any field.
SMLMData.has_nonzero_imag — Method
has_nonzero_imag(value)Check if a value has a non-zero imaginary component. Works for both scalar values and arrays.
SMLMData.load_smite_2d — Method
load_smite_2d(smd::SmiteSMD)Load a 2D Smite SMD .mat file and convert it to SmiteSMLD format. Checks for complex fields and removes emitters with non-zero imaginary components.
Arguments
smd::SmiteSMD: SmiteSMD object specifying the file to load
Returns
SmiteSMLD containing 2D localizations
Notes
- All spatial coordinates are converted to microns
- If PixelSize is not specified in the file, defaults to 0.1 microns
- Emitters with non-zero imaginary components will be excluded with a warning
- Fields are converted from Float32 to Float64 as needed
SMLMData.load_smite_3d — Method
load_smite_3d(smd::SmiteSMD)Load a 3D Smite SMD .mat file and convert it to SmiteSMLD format. Checks for complex fields and removes emitters with non-zero imaginary components.
Arguments
smd::SmiteSMD: SmiteSMD object specifying the file to load
Returns
SmiteSMLD containing 3D localizations
Notes
- All spatial coordinates are converted to microns
- If PixelSize is not specified in the file, defaults to 0.1 microns
- Emitters with non-zero imaginary components will be excluded with a warning
- Fields are converted from Float32 to Float64 as needed
SMLMData.merge_smld — Method
merge_smld(smlds::Vector{<:AbstractSMLD}; adjust_frames=false, adjust_datasets=false)
merge_smld(smlds::AbstractSMLD...; adjust_frames=false, adjust_datasets=false)Merge multiple SMLD objects with options to adjust frame and dataset numbering.
Arguments
smlds: Vector of AbstractSMLD objects or multiple AbstractSMLD argumentsadjust_frames: If true, adjusts frame numbers to be sequentialadjust_datasets: If true, adjusts dataset numbers to be sequential
Returns
New AbstractSMLD containing all emitters with adjusted numbering if requested
Notes
- Camera must be identical across all SMLDs
- When adjust_frames=true, frame numbers are made sequential across all inputs
- When adjust_datasets=true, dataset numbers are made sequential
- Metadata includes information about the merge operation
Examples
# Simple merge
merged = merge_smld(smld1, smld2)
# Merge with frame number adjustment
merged = merge_smld(smld1, smld2, adjust_frames=true)
# Merge multiple with both adjustments
merged = merge_smld([smld1, smld2, smld3],
adjust_frames=true,
adjust_datasets=true)SMLMData.physical_to_pixel — Method
physical_to_pixel(x::Real, y::Real, pixel_size::Real)Convert physical coordinates (in microns) to pixel coordinates.
Arguments
x::Real: x coordinate in microns (0,0 is top-left of image)y::Real: y coordinate in microns (0,0 is top-left of image)pixel_size::Real: size of a pixel in microns
Returns
Tuple{Float64, Float64}: (px,py) pixel coordinates where (1,1) is center of top-left pixel
Example
# For a camera with 0.1 micron pixels
px, py = physical_to_pixel(0.05, 0.05, 0.1) # Point 0.05,0.05 microns from origin
# Returns (1.0, 1.0) - center of first pixelSMLMData.physical_to_pixel_index — Method
physical_to_pixel_index(x::Real, y::Real, pixel_size::Real)Convert physical coordinates (in microns) to integer pixel indices. Returns the pixel that contains the given physical coordinate.
Arguments
x::Real: x coordinate in microns (0,0 is top-left of image)y::Real: y coordinate in microns (0,0 is top-left of image)pixel_size::Real: size of a pixel in microns
Returns
Tuple{Int, Int}: (px,py) pixel indices where (1,1) is top-left pixel
Example
# For a camera with 0.1 micron pixels
px, py = physical_to_pixel_index(0.05, 0.05, 0.1) # Point at center of first pixel
# Returns (1, 1)SMLMData.pixel_to_physical — Method
pixel_to_physical(px::Real, py::Real, pixel_size::T) where TConvert pixel coordinates to physical coordinates (in microns). Returns coordinates with the same type as pixel_size.
SMLMData.save_smite — Method
save_smite(smld::SmiteSMLD, filepath::String, filename::String)Save SmiteSMLD data back to SMITE's SMD .mat format.
Arguments
smld::SmiteSMLD: SMLD object to savefilepath::String: Directory path where to save the filefilename::String: Name of the output .mat file
Notes
- Saves in MATLAB v7.3 format
- Preserves all metadata fields
SMLMData.@filter — Macro
@filter(smld, condition)Filter SMLD emitters using a natural condition syntax. Transforms expressions at compile time into efficient filtering operations.
Examples
# Simple conditions
bright = @filter(smld, photons > 1000)
early = @filter(smld, frame < 10)
# Compound conditions
good_fits = @filter(smld, σ_x < 0.02 && σ_y < 0.02)
roi = @filter(smld, 1.0 <= x <= 5.0 && 1.0 <= y <= 5.0)