SMLMFrameConnection
Documentation for SMLMFrameConnection.
SMLMFrameConnection.ParamStructSMLMFrameConnection.combinelocalizationsSMLMFrameConnection.compress_connectIDSMLMFrameConnection.compress_connectID!SMLMFrameConnection.computeclusterinfoSMLMFrameConnection.connect1DS!SMLMFrameConnection.connectlocalizationsSMLMFrameConnection.connectlocalizations!SMLMFrameConnection.create_costmatrixSMLMFrameConnection.defineidealFC!SMLMFrameConnection.estimatedensitiesSMLMFrameConnection.estimateparamsSMLMFrameConnection.frameconnectSMLMFrameConnection.frameconnect!SMLMFrameConnection.linkclustersSMLMFrameConnection.linkclusters!SMLMFrameConnection.organizeclustersSMLMFrameConnection.preclusterSMLMFrameConnection.solveLAP
SMLMFrameConnection.ParamStruct — TypeParamStruct(initialdensity::Vector{Float64},
k_on::Float64, k_off::Float64, k_bleach::Float64, p_miss::Float64,
nsigmadev::Float64, maxframegap::Int, nnearestclusters::Int)Structure of parameters needed for frame-connection.
Fields
initialdensity: Density of emitters at the start of the experiment. (see estimatedensities()) (emitters/pixel^2)nnearestclusters: 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. (see estimateparams())(1/frame)k_off: Rate at which visible emitters are converted to the reversible dark state. (see estimateparams())(1/frame)k_bleach: Rate at which visible emitters are irreversibly photobleached. (see estimateparams())(1/frame)p_miss: Probability of missing a localization of a visible emitter.nsigmadev: Multiplier of localization errors that defines a pre-clustering distance threshold. (default = 5)(see precluster())(pixels)maxframegap: Maximum frame gap between temporally adjacent localizations in a precluster. (default = 5)(see precluster())(frames)nmaxnn: 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())
SMLMFrameConnection.combinelocalizations — Methodsmld_combined = combinelocalizations(smld::SMLMData.SMLD2D)Combine clustered localizations in smld into higher precision localizations.
Description
This function combines localizations in smld that share the same value of smld.connectID. Localizations are combined assuming they arose from independent measurements of the same position with Gaussian errors.
SMLMFrameConnection.compress_connectID! — Methodcompress_connectID!(connectID::Vector{Int32})Make connectID consists of the set of integers 1:length(unique(connectID))`
SMLMFrameConnection.compress_connectID — MethodconnectID = compress_connectID(connectID::Vector{Int32})Make connectID consists of the set of integers 1:length(unique(connectID))`
SMLMFrameConnection.computeclusterinfo — Methodclusterdurations, 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.
SMLMFrameConnection.connect1DS! — Methodconnect1DS!(smld::SMLMData.SMLD2D, dataset::Int; maxframegap::Int = 5)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: SMLMData.SMLD2D with smld.connectID populated to indicate emitter ID.dataset: Dataset number to be connected.maxframegap: Maximum frame gap allowed between localizations connected in the "ideal" result.
SMLMFrameConnection.connectlocalizations! — Methodconnectlocalizations!(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.
SMLMFrameConnection.connectlocalizations — MethodconnectID = 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.
SMLMFrameConnection.create_costmatrix — Methodcostmatrix = create_costmatrix(clusterdata::Vector{Matrix{Float32}},
params::ParamStruct, nframes::Int64)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.
SMLMFrameConnection.defineidealFC! — Methodsmld, smld_combined = defineidealFC!(smld::SMLMData.SMLD2D;
maxframegap::Int = 5)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 SMLD2D structure smld with connectID 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 maxframegap 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 maxframegap frames of its previous blink). Note that localizations are not allowed to be connected across datasets.
Inputs
-smld: SMLMData.SMLD2D with smld.connectID populated to indicate emitter ID. -maxframegap: Maximum frame gap allowed between localizations connected in the "ideal" result.
Outputs
-smld_out: Ideal frame-connection result on input smld with localizations combined as appropriate. -smld_connected: Copy of the input smld with smldconnected.connectID modified to reflect blinking event ID (i.e., `smldout` is generated as smldout = combinelocalizations(smldconnected))
SMLMFrameConnection.estimatedensities — Methodinitialdensity = estimatedensities(smld::SMLMData.SMLD2D,
clusterdata::Vector{Matrix{Float32}}, params::ParamStruct)Estimate local emitter densities for clusters in smld and clusterdata.
Description
The initial local densities initialdensity 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.
SMLMFrameConnection.estimateparams — Methodk_on, k_off, k_bleach, p_miss, nemitters =
estimateparams(smld::SMLMData.SMLD2D, clusterdata::Vector{Matrix{Float32}})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
SMLMFrameConnection.frameconnect! — Methodsmld, smld_preclustered, smld_combined, params = frameconnect!(smld::SMLMData.SMLD2D;
nnearestclusters::Int=2, nsigmadev::Float64=5.0, maxframegap::Int=5, nmaxnn::Int=2)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
- 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.
Inputs
smld: SMLD2D structure containing the localizations that should be connected (see SMLMData.SMLD2D for organization/fields).nnearestclusters: Number of nearest preclusters used for local density estimates. (default = 2)(see estimatedensities())nsigmadev: Multiplier of localization errors that defines a pre-clustering distance threshold. (default = 5)(see precluster())(pixels)maxframegap: Maximum frame gap between temporally adjacent localizations in a precluster. (default = 5)(see precluster())(frames)nmaxnn: 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())
Outputs
smld: Inputsmldwith field connectID updated to reflect connected localizations (however localizations remain uncombined).smld_preclustered: Copy of the inputsmldwith the field connectID populated to reflect the results of pre-clustering.smld_combined: Final frame-connection result (i.e.,smldwith localizations that seem to be from the same blinking event combined into higher precision localizations).params: Structure of parameters used in the algorithm, with some copied directly from the option kwargs to this function, and others calculated internally (see SMLMFrameConnection.ParamStruct).
SMLMFrameConnection.frameconnect — Methodsmld, smld_preclustered, smld_combined, params = frameconnect(smld::SMLMData.SMLD2D;
nnearestclusters::Int=2, nsigmadev::Float64=5.0, maxframegap::Int=5, nmaxnn::Int=2)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
- 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.
Inputs
-smld: SMLD2D structure containing the localizations that should be connected (see SMLMData.SMLD2D for organization/fields). -nnearestclusters: Number of nearest preclusters used for local density estimates. (default = 2)(see estimatedensities()) -nsigmadev: Multiplier of localization errors that defines a pre-clustering distance threshold. (default = 5)(see precluster())(pixels) -maxframegap: Maximum frame gap between temporally adjacent localizations in a precluster. (default = 5)(see precluster())(frames) -nmaxnn: 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())
Outputs
-smld: Input smld with field connectID updated to reflect connected localizations (however localizations remain uncombined). -smld_preclustered: Copy of the input smld with the field connectID populated to reflect the results of pre-clustering. -smld_combined: Final frame-connection result (i.e., smld with localizations that seem to be from the same blinking event combined into higher precision localizations). -params: Structure of parameters used in the algorithm, with some copied directly from the option kwargs to this function, and others calculated internally (see SMLMFrameConnection.ParamStruct).
SMLMFrameConnection.linkclusters! — MethodconnectID, maxconnectID = linkclusters!(connectID::Int64,
maxconnectID::Int64, updateind::Vector{Int64},
assignment::Vector{Int64})Force localizations linked by assignment to share same the same connectID
SMLMFrameConnection.linkclusters — MethodconnectID, maxconnectID = linkclusters(
connectID::Int64, maxconnectID::Int64,
updateind::Vector{Int64}, assignment::Vector{Int64})Force localizations linked by assignment to share same the same connectID
SMLMFrameConnection.organizeclusters — Methodclusterdata = organizeclusters(smld::SMLMData.SMLD2D)Organize pre-clusters into a vector indexing distinct clusters.
Description
Pre-clusters in the input smld, as related by their shared integer value k_off smld.connectID, are organized into a vector of matrices clusterdata. Each index of clusterdata corresponds to a distinct cluster.
SMLMFrameConnection.precluster — Functionsmld_preclustered = precluster(smld::SMLMData.SMLD2D,
params::ParamStruct = ParamStruct())Cluster localizations in smld based on distance and time thresholds in params.
Description
Localizations in the input structure smd are clustered together based on their spatiotemporal separations. All localizations within a spatial threshold of params.nsigmadev*mean([smld.σ_x smld.σ_y) and a temporal threshold of params.maxframegap of one another will be clustered together, meaning that these localizations now share the same unique integer value for smld.connectID.
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.
SMLMFrameConnection.solveLAP — Methodassignment, 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.