Multi-color data analysis#
This module provides functions to analyze multi-color data. In particular, there is support for
registration of color channels, i.e., finding transforms for mapping coordinates between channels (
Registrator
) as a prerequisite for finding colocalizationsfinding single-molecule colocalizations (
find_colocalizations()
)detecting and plotting single-molecule codiffusion (
find_codiffusion()
,plot_codiffusion()
)merging single-molecule data from different channels (
merge_channels()
)selection of images and single-molecule data according to excitation type (
FrameSelector
)
Examples
To demonstrate how to find colocalizations, create fake data first:
>>> loc1 = pandas.DataFrame([[1, 1], [11, 11]], columns=["x", "y"])
>>> loc2 = pandas.DataFrame([[1, 2], [22, 22]], columns=["x", "y"])
>>> loc1["frame"] = loc2["frame"] = 0
Now, find the colocalizations using find_colocalizations()
:
>>> coloc = find_colocalizations(loc1, loc2, max_dist=2)
>>> coloc
channel1 channel2
x y frame x y frame
0 1 1 0 1 2 0
Detecting codiffusing particles works similar using
find_codiffusion()
:
>>> trc1 = sdt.io.load("tracking1.h5") # load data
>>> trc2 = sdt.io.load("tracking2.h5")
>>> codiff = find_codiffusion(trc1, trc2)
A single codiffusing pair can be plot with help of plot_codiffusion()
using the result of find_codiffusion()
:
>>> plot_codiffusion(codiff, particle=0)
Data from two color channels can be merged so that colocalizing features only
appear once using merge_channels()
:
>>> merged = merge_channels(loc1, loc2, max_dist=2)
>>> merged
x y frame
0 1 1 0
1 11 11 0
2 22 22 0
Programming reference#
- class sdt.multicolor.Registrator(feat1=None, feat2=None, columns={}, channel_names=['channel1', 'channel2'])[source]#
Registration of two fluorescence microscopy emission channels
This class provides an easy-to-use interface to determine maps between the channels’ coordinates using localization data from fiducial markers. It is based on the algorithm published by Preibisch et al. [Preibisch2010].
Examples
Let’s assume that multiple images/sequences of fluorescent beads have been acquired, which are visible in both emission channels. First, the beads need to be localized (e.g. using the
sdt.gui.locator
application). These need to be loaded:>>> bead_loc = [sdt.io.load(b) for b in glob.glob("beads*.h5")]
Next, we define ROIs for the two channels using
sdt.roi.ROI
and choose the bead localizations with respect to the ROIs:>>> r1 = sdt.roi.ROI((0, 0), size=(200, 100)) >>> r2 = sdt.roi.ROI((0, 200), size=(200, 100)) >>> beads_r1 = [r1(b) for b in bead_loc] >>> beads_r2 = [r2(b) for b in bead_loc]
Now, calculate the transform that overlays the two channels using
Registrator
:>>> corr = Registrator(beads_r1, beads_r2) >>> corr.determine_parameters() >>> corr.test() # Plot results
This can now be used to transform i.e. image data from channel 1 so that it can be overlaid with channel 1:
>>> with io.ImageSequence("image.tif") as s: ... img = s[0] # Load first frame (number 0) >>> img_r1 = r1(img) # Get channel 1 part of the image >>> img_r2 = r2(img) # Get channel 2 part of the image >>> img_r1_corr = corr(img_r1, channel=1) # Transform channel 1 image >>> overlay = numpy.dstack([img_r1, img_r2, np.zeros_like(img_r1)]) >>> matplotlib.pyplot.imshow(overlay) # Plot overlay
Similarly, coordinates of single molecule data can be transformed:
>>> loc = sdt.io.load("loc.h5") # Load data >>> loc_r1 = r1(loc) # Get data from channel 1 >>> loc_r2 = r2(loc) # Get data from channel 2 >>> loc_r1_corr = corr(loc_r1, channel=1) # Transform channel 1 data >>> matplotlib.pyplot.scatter(loc_r1_corr["x"], loc_r1_corr["y"], ... marker="+") >>> matplotlib.pyplot.scatter(loc_r2["x"], loc_r2["y"], marker="x")
There is also support for saving and loading a
Registrator
instance to/from YAML:>>> with open("output.yaml", "w") as f: >>> sdt.io.yaml.safe_dump(corr, f) >>> with open("output.yaml", "r") as f: >>> corr_loaded = sdt.io.yaml.safe_load(f)
References
[Preibisch2010]Preibisch, S.; Saalfeld, S.; Schindelin, J. & Tomancak, P.: “Software for bead-based registration of selective plane illumination microscopy data”, Nature Methods, Springer Science and Business Media LLC, 2010, 7, 418–419
- Parameters:
feat1 (Sequence[DataFrame]) – Set the feat1 and feat2 attribute (turning it into a list if it is a single DataFrame). Can also be None, but in this case
find_pairs()
anddetermine_parameters()
will not work. Defaults to None.feat2 (Sequence[DataFrame]) – Set the feat1 and feat2 attribute (turning it into a list if it is a single DataFrame). Can also be None, but in this case
find_pairs()
anddetermine_parameters()
will not work. Defaults to None.channel_names (Sequence[str]) – Set the channel_names attribute.
columns (Dict) – Override default column names as defined in
config.columns
. Relevant name are coords and time. That means, if the DataFrames have coordinate columns “x” and “z”, and a time column “alt_frame”, setcolumns={"coords": ["x", "z"], "time": "alt_frame"}
. This is used to set thecolumns
attribute.
- feat1: Sequence[DataFrame]#
Positions of beads (as found by a localization algorithm) in the first channel. Each DataFrame corresponds to one image (sequence), thus multiple bead images can be used to increase the accuracy.
- channel_names: Sequence[str]#
Channel names
- pairs: DataFrame#
Pairs found by
determine_parameters()
.
- parameters1: ndarray#
Array describing the affine transformation of data from channel 1 to channel 2.
- parameters2: ndarray#
Array describing the affine transformation of data from channel 2 to channel 1.
- determine_parameters(n_neighbors=3, ambiguity_factor=5.0, max_error=1.0)[source]#
Determine the parameters for the affine transformation
This takes the localizations of
feat1
and tries to match them with those offeat2
. Then a fit is used to determine the affine transformation between the channels.This is a convenience function that calls
find_pairs()
andfit_parameters()
; see the documentation of the methods for details.- Parameters:
n_neighbors (int) – Number of neighboring beads to consider to find matching features across channels.
ambiguity_factor (float) – A low value (around 1) will accept pairs even if there are similar possible matches. The higher this value, the less are ambigous results accepted.
max_error (float) – Maximum error (i.e., distance between transformed position from channel 1 and matched position in channel 2) to consider a feature pair not an outlier and thus remove it from the transformation fit.
- __call__(data: pd.DataFrame, channel: int | str, inplace: Literal[True], mode: str, cval: float | Callable[[np.ndarray], float], columns: Mapping)[source]#
- __call__(data: pd.DataFrame, channel: int | str, inplace: Literal[False], mode: str, cval: float | Callable[[np.ndarray], float], columns: Mapping) pd.DataFrame
- __call__(data: helper.Slicerator | helper.Pipeline, channel: int | str, inplace: Literal[False], mode: str, cval: float | Callable[[np.ndarray], float], columns: Mapping) helper.Pipeline
Apply transformation to data
This can be done either on coordinates (e. g. resulting from feature localization) or directly on raw images.
- Parameters:
data – Data to be processed. If a pandas.Dataframe, the feature coordinates will be corrected. Otherwise,
sdt.helper.Pipeline
is used to correct image data using an affine transformation.channel – If features are in the first channel (corresponding to the feat1 arg of the constructor), set to 1 or the first channel’s name. If features are in the second channel, set to 2 or the second channel’s name. Depending on this, a transformation will be applied to the coordinates of features to match the other channel (mathematically speaking depending on this parameter either the “original” transformation or its inverse are applied.)
inplace – Only has an effect if data is a DataFrame. If True, the feature coordinates will be corrected in place.
mode – How to fill points outside of the uncorrected image boundaries. Possibilities are “constant”, “nearest”, “reflect” or “wrap”.
cval – What value to use for mode=”constant”. If this is callable, it should take a single argument (the uncorrected image) and return a scalar, which will be used as the fill value.
columns – Override default column names in case data is a
pandas.DataFrame
. The only relevant name is coords. That means, if the DataFrame has coordinate columns “x” and “z”, setcolumns={"coords": ["x", "z"]}
.
- test(ax=None)[source]#
Test validity of the correction parameters
This plots the affine transformation functions and the coordinates of the pairs that were matched in the channels. If everything went well, the dots (i. e. pair coordinates) should lie on the line (transformation function).
- Parameters:
ax (Sequence | None) – Two matplotib axes instances for plotting. If None, allocate new axes using
matplotlib.pyplot.subplots()
.
- save(file, fmt='npz', key=('chromatic_param1', 'chromatic_param2'))[source]#
Save transformation parameters to file
- Parameters:
file (BinaryIO | str | Path) – File name or an open file-like object to save data to.
fmt (str) – Format to save data. Either numpy (“npz”) or MATLAB (“mat”). Defaults to “npz”.
key (Tuple[str, str]) – Names of two transformations in the saved file.
- classmethod load(file, fmt='npz', key=('chromatic_param1', 'chromatic_param2'))[source]#
Read paramaters from a file and construct a Registrator
- Parameters:
file (BinaryIO | str | Path) – File name or an open file-like object to load data from.
fmt (str) – Format to save data. Either numpy (“npz”) or MATLAB (“mat”) or determine_shiftstretch’s wrp (“wrp”). Defaults to “npz”.
key (Tuple[str, str]) – Name of the variables in the saved file (does not apply to “wrp”). Defaults to (“chromatic_param1”, “chromatic_param2”).
- Return type:
A class instance with the parameters read from the file.
- find_pairs(n_neighbors=3, ambiguity_factor=5.0)[source]#
Match features of feat1 with features of feat2
Find the geomtric signature for each feature in each channel. Those with best matching signatures are taken to be the same feature.
- Parameters:
n_neighbors (int) – Number of neighboring beads to consider for signature calculation.
ambiguity_factor (float) – Accept only matches where the distance (with respect to the geometric signature) between the beads of the second best match is at least
ambiguity_factor
times as large as the distance for the best match.
- fit_parameters(max_error=1.0)[source]#
Determine parameters for the affine transformation
An affine transformation is used to map x and y coordinates of feat1 to to those of feat2. This requires
pairs
to be set, e.g. by runningfind_pairs()
. The result is saved asparameters1
(transform of channel 1 coordinates to channel 2) andparameters2
(transform of channel 2 coordinates to channel 1) attributes.parameters1
is calculated by determining the affine transformation between thepairs
entries using a RANSAC algorithm.parameters2
is its inverse. Therefore, results may be slightly different depending on what is channel1 and what is channel2.- Parameters:
max_error (float) – Maximum error (i.e., distance between transformed position from channel 1 and matched position in channel 2) to consider a feature pair not an outlier.
- sdt.multicolor.find_colocalizations(features1, features2, max_dist=2.0, keep_unmatched=False, channel_names=['channel1', 'channel2'], columns={}, **kwargs)[source]#
Match localizations in one channel to localizations in another
For every localization in features1 find localizations in features2 (in the same frame) that are in a circle with radius max_dist around it, then pick the closest one.
- Parameters:
features1 (DataFrame) – Localization data
features2 (DataFrame) – Localization data
max_dist (float) – Maximum distance between features to still be considered colocalizing.
keep_unmatched (bool) – If True, also keep non-colocalized features in the result DataFrame. Non-colocalized features have NaNs in the columns of the channel they don’t appear in.
channel_names (Sequence[str]) – Names of the two channels.
columns (Mapping) – Override default column names as defined in
config.columns
. Relevant names are coords and time. This means, if your DataFrame has coordinate columns “x” and “z” and the time column “alt_frame”, setcolumns={"coords": ["x", "z"], "time": "alt_frame"}
.keep_non_coloc – Deprecated alias for keep_unmatched
- Returns:
The DataFrame has a multi-index for columns with the top level
given by the channel_names parameter. Each line of DataFrame
corresponds to one pair of colocalizing particles.
- Return type:
DataFrame
- sdt.multicolor.calc_pair_distance(data, channel_names=None, columns={})[source]#
Calculate distances between colocalized features
- Parameters:
data (pandas.DataFrame) – Colocalized feature data, e.g. output of
find_colocalizations()
.channel_names (list of str or None, optional) – Names of the channels. If None, use the first two entries of teh top level of data’s MultiIndex. Defaults to None.
columns (dict, optional) – Override default column names as defined in
config.columns
. The only relevant name is coords. This means, if your DataFrame has coordinate columns “x” and “z”, set ``columns={“coords”: [“x”, “z”]}`.
- Returns:
Distances between colocalized features. The series’ index is the same as in the input DataFrame.
- Return type:
pandas.Series
- sdt.multicolor.find_codiffusion(tracks1, tracks2, abs_threshold=3, rel_threshold=0.75, return_data='data', feature_pairs=None, max_dist=2.0, keep_unmatched='gaps', channel_names=['channel1', 'channel2'], columns={})[source]#
Find codiffusion in tracking data
First, find pairs of localizations, the look up to which tracks they belong to and thus match tracks.
- Parameters:
tracks1 (DataFrame) – Tracking data. In addition to what is required by
find_colocalizations()
, a “particle” column has to be present.tracks2 – Tracking data. In addition to what is required by
find_colocalizations()
, a “particle” column has to be present.abs_threshold – Minimum number of matched localizations per track pair.
rel_threshold (float, optional) – Minimum fraction of matched localizations that have to belong to the same pair of tracks. E. g., if localizations of a particle in the first channel match five localizations of a particle in the second channel, and there are eight frames between the first and the last match, that fraction would be 5/8.
return_data – Whether to return the full data of the codiffusing particles (
"data"
), only their particle numbers ("numbers"
), or both ("both"
).feature_pairs – If
find_colocalizations()
has already been called on the data, the result can be passed to save re-computation. If None,find_colocalizations()
is called in this function.max_dist – max_dist parameter for
find_colocalizations()
call.keep_unmatched – If
"all"
, also keep all non-colocalized features in the result DataFrame. Non-colocalized features have NaNs in the columns of the channel they don’t appear in.If"gaps"
, only keep non-colocalized features within codiffusing parts tracks. If"none"
, remove all non-colocalized entries.channel_names – Names of the two channels.
columns – Override default column names as defined in
config.columns
. Relevant names are coords, particle, and time. This means, if your DataFrame has coordinate columns “x” and “z” and the time column “alt_frame”, setcolumns={"coords": ["x", "z"], "time": "alt_frame"}
.
- Returns:
If return_data is
"data"
or"both"
, a DataFrame with fulldata (from tracks1 and tracks2) of the codiffusing particles is
returned. The DataFrame has a multi-index for columns with the top level
being the two channels. Each line of DataFrame corresponds to one pair of
colocalizing particles.
If return_data is
"numbers"
or"both"
, an array with fourcolumns is returned.
Each row’s first entry is a particle number in the first channel. The
second entry is the matching particle number in the second channel.
Third and fourth columns are start and end frame, respectively.
- sdt.multicolor.plot_codiffusion(data, particle, ax=None, cmap=None, show_legend=True, legend_loc=0, linestyles=['-', '--', ':', '-.'], channel_names=None, columns={})[source]#
Plot trajectories of codiffusing particles
Each step is colored differently so that by comparing colors one can figure out which steps in one channel correspond to which steps in the other channel.
- Parameters:
data (pandas.DataFrame or tuple of pandas.DataFrames) – Tracking data of codiffusing particles. This can be a
pandas.DataFrame
with a MultiIndex for columns as e. g. returned byfind_codiffusion()
(i. e. matching indices in the DataFrames correspond to matching localizations) or a tuple of DataFrames, one for each channel.particle (int or tuple of int) – Specify particle ID. In case data is a list of DataFrames and the particles have different IDs, one can pass the tuple of IDs.
ax (matplotlib.axes.Axes) – To be used for plotting. If None,
matplotlib.pyplot.gca()
will be used. Defaults to None.cmap (matplotlib.colors.Colormap, optional) – To be used for coloring steps. Defaults to the “Paired” map of matplotlib.
show_legend (bool, optional) – Whether to print a legend or not. Defaults to True.
legend_loc (int or str) – Is passed as the loc parameter to matplotlib’s axes’ legend method. Defaults to 0.
channel_names (list of str or None, optional) – Names of the channels. If None, use the first two entries of teh top level of data’s MultiIndex if it has one (i. e. if it is a DataFrame as created by
find_codiffusion()
), otherwise use [“channel1”, “channel2”]. Defaults to None.columns (dict, optional) – Override default column names as defined in
config.columns
. Relevant names are coords, particle, and time. This means, if your DataFrame has coordinate columns “x” and “z” and the time column “alt_frame”, setcolumns={"coords": ["x", "z"], "time": "alt_frame"}
.
- sdt.multicolor.merge_channels(features1, features2, max_dist=2.0, mean_pos=False, return_data='data', columns={})[source]#
Merge features of features1 and features2
Concatenate all of features1 and those entries of features2 that do not colocalize with any of features1 (in the same frame).
- Parameters:
features1 (pandas.DataFrame) – Localization data
features2 (pandas.DataFrame) – Localization data
max_dist (float, optional) – Maximum distance between features to still be considered colocalizing. Defaults to 2.
mean_pos (bool, optional) – When entries are merged (i. e., if an entry in features1 is close enough to one in features2), calculate the center of mass of the two channel’s entries. If False, use the position in features1. All other (non-coordinate) columns are taken from features1 in any case. Defaults to False.
return_data ({"data", "index", "both"}, optional) – Whether to return the full data of merged features, only indices (that is, the DatatFrame’s indices) of features in features2 that have no counterpart in features1, or both. Defaults to “data”.
columns (dict, optional) – Override default column names as defined in
config.columns
. Relevant names are coords and time. This means, if your DataFrame has coordinate columns “x” and “z” and the time column “alt_frame”, setcolumns={"coords": ["x", "z"], "time": "alt_frame"}
.
- Returns:
data (pandas.DataFrame) – DataFrame of merged features. Returned if return_data is “data” or “both”.
feat2_index (pandas.Index) – Indices of features out of features2 that have no counterpart in features1. One can construct the DataFrame (as returned if return_data is “data” or “both”) by
pandas.concat([features1, features2.loc[feat2_index]], keys=channel_names, ignore_index=True)
. Returned if return_data is “index” or “both”.
- class sdt.multicolor.FrameSelector(excitation_seq)[source]#
Select images and datapoints of a certain excitation type
For instance, if employing alternating excitation for a FRET experiment, this can be used to select only the donor or only the acceptor frames. For details on how the excitation sequence is specified, see
excitation_seq
.Examples
>>> # Sequence of 8 "images" >>> img_seq = numpy.array([numpy.full((3, 3), i) for i in range(8)]) >>> sel = FrameSelector("odda") >>> sel.select(img_seq, "o") array([[[0, 0, 0], [0, 0, 0], [0, 0, 0]], [[4, 4, 4], [4, 4, 4], [4, 4, 4]]]) >>> sel.select(img_seq, "d") array([[[1, 1, 1], [1, 1, 1], [1, 1, 1]], [[2, 2, 2], [2, 2, 2], [2, 2, 2]], [[5, 5, 5], [5, 5, 5], [5, 5, 5]], [[6, 6, 6], [6, 6, 6], [6, 6, 6]]])
A flexible sequence specification uses
?
as a multiplicator. This causes repetition of a subsequence such that the excitation sequence length matches the image sequence length. In the following example,"da"
is repeated three times.>>> sel = FrameSelector("o + da * ? + c") >>> sel.select(img_seq, "d") array([[[1, 1, 1], [1, 1, 1], [1, 1, 1]], [[3, 3, 3], [3, 3, 3], [3, 3, 3]], [[5, 5, 5], [5, 5, 5], [5, 5, 5]]])
- Parameters:
excitation_seq (str) – Set
excitation_seq
- excitation_seq: str#
Excitation sequence. Use different letters to describe different excitation types, e.g., “d” for donor excitation, “a” for acceptor, etc.
One needs only specify the shortest sequence that is repeated, i. e.
"ddddaddddadddda"
is equivalent to"dddda"
. It is possible to specifiy repetition via the multiplication operator."c + da * 3 + c"
is equivalent to"cdadadac"
. Also, flexible sequences can be specified using?
as a multiplicator. In"c + da * ? + c"
,"da"
is repeated an appropriate number of times to generate a sequence of given length. This length is either derived from the image sequence the class instance is applied to or given explicitly (seeselect()
).If the sequence is empty, the FrameSelector does nothing. E.g.,
select()
andrenumber_frames()
just return their argument unaltered.
- eval_seq(n_frames=None)[source]#
Evaluate excitation sequence
- Parameters:
n_frames (int | None) – Number of frames. This needs to be given if sequence is flexible, i.e., contains multiplication by
?
. In this case,?
evaluates to a number such that the total sequence length is n_frames. If n_frames is equal to -1,?
is replaced by 1. This is useful for testing and finding e.g. the first occurence of a frame without depending on actual data.- Return type:
Array of characters representing evaluated sequence
- renumber_frames(frame_nos, which, restore=False, n_frames=None)[source]#
Renumber a sequence of frame numbers
The numbers can be with respect to the original image sequence. In this case, set
restore=False
and this function will return frame numbers with respect to an image sequence which was returned byselect()
. Any frames that don’t belong to the excitation type given by which are mapped to -1. Ifrestore=True
is set, this works the opposite other way.- Parameters:
frame_nos (ndarray) – Array of frame numbers to be mapped
which (str | Iterable[str]) – Excitation type(s), one or multiple characters out of
excitation_seq
.restore (bool) – If False, renumber the frames (as in
select()
). If True, restore original frame numbers (as inrestore_frame_numbers()
).n_frames (int | None) –
- Return type:
New frame numbers
- find_other_frames(data, which, other, how='nearest-up')[source]#
For given excitation type, find frames of other type
Typically, this is used to find the closest frame of a certain type, e.g., find the next acceptor excitation frame for donor excitation frames in FRET ALEX experiments.
Examples
>>> fs = FrameSelector("dddda") >>> d_img = fs.select(image_sequence, "d") >>> a_img = fs.find_other_frames(image_sequence, "d", "a) >>> for di, ai in zip(d_img, a_img): ... # `ai` is the frame of type "a" closest to `di` ... pass
- Parameters:
data (int | Sequence) – If this is a sequence, return a subsequence of the same length as
select(data, which)
. Each entry of the returned sequence is the corresponding entry to whatselect(data, which)
returns. Passing an int is equivalent to usingnumpy.arange(data)
.which (str | Iterable[str]) – Excitation types for which to find other excitation type’s frames
other (str | Iterable[str]) – Other frames’ excitation type(s)
how (str) – How to do the interpolation. “nearest-up” will return the nearest frame. In case of a tie, the “other” frame after the current one is used. “nearest” will return the “other” frame before the current one in case of a tie. “previous” and “next” return the “other” frame before and after the current one, respectively.
- Returns:
Subsequence of data with entries being the “other” entries to
select(data, which)
.
- Return type:
Sequence
- select(data, which, renumber=False, n_frames=None, columns={})[source]#
Get only data corresponding to a certain excitation type
- Parameters:
data (DataFrame | Iterable[ndarray]) – Localization data or sequence of images
which (str | Iterable[str]) – Excitation type(s), one or multiple characters out of
excitation_seq
.renumber (bool) –
Renumber frames so that only frames for excitation types corresponding to which are counted.
After selecting only data corresponding to an excitation type (e.g. “a”) in a
pandas.DataFrame
, all other types are missing and frame number are not consecutive any more. This can be a problem for tracking or diffusion analysis. Settingrenumber=True
will work around this.Only applicable if data is a
pandas.DataFrame
.n_frames (int | None) – Number of frames. This needs to be given if sequence is flexible, i.e., contains multiplication by
?
. In this case,?
evaluates to a number such that the total sequence length is n_frames. If data is an image sequence, n_frames can be inferred fromlen(data)
.columns (Mapping) – Override default column names in case data is a
pandas.DataFrame
. The only relevant name is time. That means, if the DataFrame has frame number columns “frame2”, setcolumns={"time": "frame2"}
.
- Returns:
Selected frames. If data is a DataFrame, return only lines
corresponding to excitation types given by which. If
renumber=True
, a modified copy of data is returned.If data is an image sequence, first an attempt is made at
indexing the images directly. If data supports this (e.g. if it
is a
numpy.ndarray
), the indexed version is returned.Otherwise, img_seq is converted to
helper.Slicerator
, indexed, and returned.
- Return type:
DataFrame | ndarray | Slicerator
- restore_frame_numbers(data, which, n_frames=None, columns={})[source]#
Undo frame renumbering from
select()
data is modified in place for this purpose.
- Parameters:
data (DataFrame) – Localization data
which (str | Iterable[str]) – Excitation type(s), one or multiple characters out of
excitation_seq
.n_frames (int | None) – Number of frames. This needs to be given if sequence is flexible, i.e., contains multiplication by
?
. In this case,?
evaluates to a number such that the total sequence length is n_frames.columns (Mapping) – Override default column names in case data is a
pandas.DataFrame
. The only relevant name is time. That means, if the DataFrame has frame number columns “frame2”, setcolumns={"time": "frame2"}
.
Low level helper functions#
- sdt.multicolor.find_closest_pairs(coords1, coords2, max_dist)[source]#
Find best matches between coordinates
Given two coordinate arrays (coords1 and coords2), find pairs of points with the minimum distance between them. Each point will be in at most one pair (in contrast to usual KD-tree queries).
- Parameters:
coords1 (array-like, shape(n, m)) – Arrays of m-dimensional coordinate tuples. The dimension (m) has to be the same in both arrays while the number of points (n) can differ.
coords2 (array-like, shape(n, m)) – Arrays of m-dimensional coordinate tuples. The dimension (m) has to be the same in both arrays while the number of points (n) can differ.
max_dist (float) – Maximum distance for two points to be considered a pair
- Returns:
Each row describes one match. The first entry is the index of a point in coords1. The second entry is the index of its match in coords2.
- Return type:
numpy.ndarray, shape(n, 2), dtype(int)