fissa.core module#

Main FISSA user interface.

Authors:
class fissa.core.Experiment(images, rois, folder=None, nRegions=None, expansion=None, method=None, alpha=None, max_iter=None, tol=None, max_tries=None, ncores_preparation=- 1, ncores_separation=- 1, lowmemory_mode=False, datahandler=None, verbosity=1)[source]#

Bases: object

Class-based interface for running FISSA on experiments.

Uses the methodology described in FISSA: A neuropil decontamination toolbox for calcium imaging signals.

Parameters
imagesstr or list

The raw imaging data. Should be one of:

  • the path to a directory containing TIFF files (string),

  • a list of paths to TIFF files (list of strings),

  • a list of array_like data already loaded into memory, each shaped (n_frames, height, width).

Note that each TIFF or array is considered a single trial.

roisstr or list

The region of interest (ROI) definitions. Should be one of:

  • the path to a directory containing ImageJ ZIP files (string),

  • the path of a single ImageJ ZIP file (string),

  • a list of ImageJ ZIP files (list of strings),

  • a list of arrays, each encoding a ROI polygons,

  • a list of lists of binary arrays, each representing a ROI mask.

This can either be a single roiset for all trials, or a different roiset for each trial.

folderstr, optional

Path to a cache directory from which pre-extracted data will be loaded if present, and saved to otherwise. If folder is unset, the experiment data will not be saved.

nRegionsint, default=4

Number of neuropil regions and signals to use. Default is 4. Use a higher number for densely labelled tissue.

expansionfloat, default=1

Expansion factor for each neuropil region, relative to the ROI area. Default is 1. The total neuropil area will be nRegions * expansion * area(ROI).

method“nmf” or “ica”, default=”nmf”

Which blind source-separation method to use. Either "nmf" for non-negative matrix factorization, or "ica" for independent component analysis. Default is "nmf" (recommended).

alphafloat, default=0.1

Sparsity regularizaton weight for NMF algorithm. Set to zero to remove regularization. Default is 0.1.

max_iterint, default=20000

Maximum number of iterations before timing out on an attempt.

New in version 1.0.0.

tolfloat, default=1e-4

Tolerance of the stopping condition.

New in version 1.0.0.

max_triesint, default=1

Maximum number of random initial states to try. Each random state will be optimized for max_iter iterations before timing out.

New in version 1.0.0.

ncores_preparationint or None, default=-1

The number of parallel subprocesses to use during the data preparation steps of separation_prep(). These steps are ROI and neuropil subregion definitions, and extracting raw signals from TIFFs.

If set to None or -1 (default), the number of processes used will equal the number of threads on the machine. If this is set to -2, the number of processes used will be one less than the number of threads on the machine; etc.

Note that the preparation process can be quite memory-intensive and it may be necessary to reduce the number of processes from the default.

ncores_separationint or None, default=-1

The number of parallel subprocesses to use during the signal separation steps of separate().

If set to None or -1 (default), the number of processes used will equal the number of threads on the machine. If this is set to -2, the number of processes used will be one less than the number of threads on the machine; etc.

The separation routine requires less memory per process than the preparation routine, and so ncores_separation be often be set higher than ncores_preparation.

lowmemory_modebool, optional

If True, FISSA will load TIFF files into memory frame-by-frame instead of holding the entire TIFF in memory at once. This option reduces the memory load, and may be necessary for very large inputs. Default is False.

datahandlerfissa.extraction.DataHandlerAbstract, optional

A custom datahandler object for handling ROIs and calcium data can be given here. See fissa.extraction for example datahandler classes. The default datahandler is DataHandlerTifffile. If datahandler is set, the lowmemory_mode parameter is ignored.

verbosityint, default=1

How verbose the processing will be. Increase for more output messages. Processing is silent if verbosity=0.

New in version 1.0.0.

Attributes
resultnumpy.ndarray

A numpy.ndarray of shape (n_rois, n_trials), each element of which is itself a numpy.ndarray shaped (n_signals, n_timepoints).

The final output of FISSA, with separated signals ranked in order of their weighting toward the raw cell ROI signal relative to their weighting toward other mixed raw signals. The ordering is such that experiment.result[roi, trial][0, :] is the signal with highest score in its contribution to the raw neuronal signal. Subsequent signals are sorted in order of diminishing score. The units are same as raw (candelas per unit area).

This field is only populated after separate() has been run; until then, it is set to None.

roi_polysnumpy.ndarray

A numpy.ndarray of shape (n_rois, n_trials), each element of which is itself a list of length nRegions + 1, each element of which is a list of length n_contour containing a numpy.ndarray of shape (n_nodes, 2).

Polygon contours describing the outline of each region.

For contiguous ROIs, the outline of the i_roi-th ROI used in the i_trial-th trial is described by the array at experiment.roi_polys[i_roi, i_trial][0][0]. This array consists of n_nodes rows, each representing the coordinate of a node in (y, x) format. For non-contiguous ROIs, a contour is needed for each disconnected polygon making up the total aggregate ROI. These contours are found at experiment.roi_polys[i_roi, i_trial][0][i_contour].

Similarly, the nRegions neuropil regions are each described by the polygons experiment.roi_polys[i_roi, i_trial][i_neurpil + 1][i_contour], respectively.

meanslist of n_trials numpy.ndarray, each shaped (height, width)

The temporal-mean image for each trial (i.e. for each TIFF file, the average image over all of its frames).

rawnumpy.ndarray

A numpy.ndarray of shape (n_rois, n_trials), each element of which is itself a numpy.ndarray shaped (n_signals, n_timepoints).

For each ROI and trial (raw[i_roi, i_trial]) we extract a temporal trace of the average value within the spatial area of each of the nRegions + 1 regions. The 0-th region is the i_roi-th ROI (raw[i_roi, i_trial][0]). The subsequent nRegions vectors are the traces for each of the neuropil regions.

The units are the same as the supplied imagery (candelas per unit area).

sepnumpy.ndarray

A numpy.ndarray of shape (n_rois, n_trials), each element of which is itself a numpy.ndarray shaped (n_signals, n_timepoints).

The separated signals, before output signals are ranked according to their matching against the raw signal from within the ROI. Separated signal i for a specific ROI and trial can be found at experiment.sep[roi, trial][i, :].

This field is only populated after separate() has been run; until then, it is set to None.

mixmatnumpy.ndarray

A numpy.ndarray of shape (n_rois, n_trials), each element of which is itself a numpy.ndarray shaped (n_rois, n_signals).

The mixing matrix, which maps from experiment.raw to experiment.sep. Because we use the collate the traces from all trials to determine separate the signals, the mixing matrices for a given ROI are the same across all trials. This means all n_trials elements in mixmat[i_roi, :] are identical.

This field is only populated after separate() has been run; until then, it is set to None.

infonumpy.ndarray shaped (n_rois, n_trials) of dicts

Information about the separation routine.

Each dictionary in the array has the following fields:

convergedbool

Whether the separation model converged, or if it ended due to reaching the maximum number of iterations.

iterationsint

The number of iterations which were needed for the separation model to converge.

max_iterationsint

Maximum number of iterations to use when fitting the separation model.

random_stateint or None

Random seed used to initialise the separation model.

This field is only populated after separate() has been run; until then, it is set to None.

deltaf_rawnumpy.ndarray

A numpy.ndarray of shape (n_rois, n_trials), each element of which is itself a numpy.ndarray shaped (1, n_timepoint).

The amount of change in fluorence relative to the baseline fluorence (Δf/f0).

This field is only populated after calc_deltaf() has been run; until then, it is set to None.

Changed in version 1.0.0: The shape of the interior arrays changed from (n_timepoint, ) to (1, n_timepoint).

deltaf_resultnumpy.ndarray

A numpy.ndarray of shape (n_rois, n_trials), each element of which is itself a numpy.ndarray shaped (n_signals, n_timepoints).

The amount of change in fluorence relative to the baseline fluorence (Δf/f0). By default, the baseline is taken from raw because the minimum values in result are typically zero. See calc_deltaf() for details.

This field is only populated after calc_deltaf() has been run; until then, it is set to None.

Methods

calc_deltaf(freq[, use_raw_f0, across_trials])

Calculate Δf/f0 for raw and result traces.

clear([verbosity])

Clear prepared data, and all data downstream of prepared data.

clear_separated([verbosity])

Clear separated data, and all data downstream of separated data.

load([path, force, skip_clear])

Load data from cache file in npz format.

save_prep([destination])

Save prepared raw signals, extracted from images, to an npz file.

save_separated([destination])

Save separated signals to an npz file.

save_to_matlab([fname])

Save the results to a MATLAB file.

separate([redo_prep, redo_sep])

Separate all the trials with FISSA algorithm.

separation_prep([redo])

Prepare and extract the data to be separated.

to_matfile([fname, legacy])

Save the results to a MATLAB file.

calc_deltaf(freq, use_raw_f0=True, across_trials=True)[source]#

Calculate Δf/f0 for raw and result traces.

The outputs are found in the deltaf_raw and deltaf_result attributes, which can be accessed at experiment.deltaf_raw and experiment.deltaf_result.

Parameters
freqfloat

Imaging frequency, in Hz.

use_raw_f0bool, optional

If True (default), use an f0 estimate from the raw ROI trace for both raw and result traces. If False, use individual f0 estimates for each of the traces.

across_trialsbool, optional

If True, we estimate a single baseline f0 value across all trials. If False, each trial will have their own baseline f0, and Δf/f0 value will be relative to the trial-specific f0. Default is True.

clear(verbosity=None)[source]#

Clear prepared data, and all data downstream of prepared data.

New in version 1.0.0.

Parameters
verbosityint, optional

Whether to show the data fields which were cleared. By default, the object’s verbosity attribute is used.

clear_separated(verbosity=None)[source]#

Clear separated data, and all data downstream of separated data.

New in version 1.0.0.

Parameters
verbosityint, optional

Whether to show the data fields which were cleared. By default, the object’s verbosity attribute is used.

load(path=None, force=False, skip_clear=False)[source]#

Load data from cache file in npz format.

New in version 1.0.0.

Parameters
pathstr, optional

Path to cache file (.npz format) or a directory containing "prepared.npz" and/or "separated.npz" files. Default behaviour is to use the folder parameter which was provided when the object was initialised is used (experiment.folder).

forcebool, optional

Whether to load the cache even if its experiment parameters differ from the properties of this experiment. Default is False.

skip_clearbool, optional

Whether to skip clearing values before loading. Default is False.

property nCell#
property nTrials#
save_prep(destination=None)[source]#

Save prepared raw signals, extracted from images, to an npz file.

New in version 1.0.0.

Parameters
destinationstr, optional

Path to output file. The default destination is "prepared.npz" within the cache directory experiment.folder.

save_separated(destination=None)[source]#

Save separated signals to an npz file.

New in version 1.0.0.

Parameters
destinationstr, optional

Path to output file. The default destination is "separated.npz" within the cache directory experiment.folder.

save_to_matlab(fname=None)[source]#

Save the results to a MATLAB file.

Deprecated since version 1.0.0: Use experiment.to_matfile(legacy=True) instead.

This will generate a .mat file which can be loaded into MATLAB to provide structs: ROIs, result, raw.

If Δf/f0 was calculated, these will also be stored as df_result and df_raw, which will have the same format as result and raw.

These can be interfaced with as follows, for ROI 0, trial 0:

ROIs.cell0.trial0{1}

Polygon outlining the ROI.

ROIs.cell0.trial0{2}

Polygon outlining the first (of nRegions) neuropil region.

result.cell0.trial0(1, :)

Final extracted neuronal signal.

result.cell0.trial0(2, :)

Contaminating signal.

raw.cell0.trial0(1, :)

Raw measured cell signal, average over the ROI.

raw.cell0.trial0(2, :)

Raw signal from first (of nRegions) neuropil region.

Parameters
fnamestr, optional

Destination for output file. Default is a file named "matlab.mat" within the cache save directory for the experiment (the folder argument when the Experiment instance was created).

separate(redo_prep=False, redo_sep=False)[source]#

Separate all the trials with FISSA algorithm.

After running separate, data can be found as follows:

experiment.sep

Raw separation output, without being matched. Signal i for a specific ROI and trial can be found as experiment.sep[roi, trial][i, :].

experiment.result

Final output, in order of presence in the ROI. Signal i for a specific ROI and trial can be found at experiment.result[roi, trial][i, :]. Note that the ordering is such that i = 0 is the signal most strongly present in the ROI, and subsequent entries are in diminishing order.

experiment.mixmat

The mixing matrix, which maps from experiment.raw to experiment.sep.

experiment.info

Information about separation routine, iterations needed, etc.

Parameters
redo_prepbool, optional

Whether to redo the preparation. Default is False. Note that if this is true, we set redo_sep = True as well.

redo_sepbool, optional

Whether to redo the separation. Default is False. Note that this parameter is ignored if redo_prep is set to True.

separation_prep(redo=False)[source]#

Prepare and extract the data to be separated.

For each trial, performs the following steps:

  • Load in data as arrays.

  • Load in ROIs as masks.

  • Grow and seaparate ROIs to define neuropil regions.

  • Using neuropil and original ROI regions, extract traces from data.

After running this you can access the raw data (i.e. pre-separation) as experiment.raw and experiment.rois. experiment.raw is a list of arrays. experiment.raw[roi, trial] gives you the traces of a specific ROI and trial, across the ROI and neuropil regions. experiment.roi_polys is a list of lists of arrays. experiment.roi_polys[roi, trial][region][0] gives you the polygon for the region for a specific ROI, trial and region. region=0 is the ROI itself (i.e. the outline of the neuron cell), and region>0 gives the different neuropil regions. For separable masks, it is possible multiple outlines are found, which can be accessed as experiment.roi_polys[roi, trial][region][i], where i is the outline index.

Parameters
redobool, optional

If False, we load previously prepared data when possible. If True, we re-run the preparation, even if it has previously been run. Default is False.

to_matfile(fname=None, legacy=False)[source]#

Save the results to a MATLAB file.

New in version 1.0.0.

This will generate a MAT-file (.mat) which can be loaded into MATLAB. The MAT-file contains structs for all the experiment output attributes (roi_polys, result, raw, etc.) and analysis parameters (expansion, nRegions, alpha, etc.). If Δf/f0 was calculated with calc_deltaf(), deltaf_result and deltaf_raw are also included.

These can be interfaced with as illustrated below.

result{1, 1}(1, :)

The separated signal for the first ROI and first trial. This is equivalent to experiment.result[0, 0][0, :] when interacting with the Experiment object in Python.

result{roi, trial}(1, :)

The separated signal for the roi-th ROI and trial-th trial. This is equivalent to experiment.result[roi - 1, trial - 1][0, :] when interacting with the Experiment object in Python.

result{roi, trial}(2, :)

A contaminating signal.

raw{roi, trial}(1, :)

Raw measured neuronal signal, averaged over the ROI. This is equivalent to experiment.raw[roi - 1, trial - 1][0, :] when interacting with the Experiment object in Python.

raw{roi, trial}(2, :)

Raw signal from first neuropil region (of nRegions).

roi_polys{roi, trial}{1}

Polygon outlining the ROI, as an n-by-2 array of coordinates.

roi_polys{roi, trial}{2}

Polygon outlining the first neuropil region (of nRegions), as an n-by-2 array of coordinates.

Parameters
fnamestr, optional

Destination for output file. The default is a file named "separated.mat" within the cache save directory for the experiment (the folder argument when the Experiment instance was created).

legacybool, default=False

Whether to use the legacy format of save_to_matlab(). This also changes the default output name to "matlab.mat".

Examples

Here are some example MATLAB plots.

Plotting raw and decontaminated traces:

% Load the FISSA output data
S = load('separated.mat')
% Separated signal for the third ROI, second trial
roi = 3; trial = 2;
% Plot the raw and result traces for the ROI signal
figure; hold on;
plot(S.raw{roi, trial}(1, :));
plot(S.result{roi, trial}(1, :));
title(sprintf('ROI %d, Trial %d', roi, trial));
xlabel('Time (frame number)');
ylabel('Signal intensity (candela per unit area)');
legend({'Raw', 'Result'});

If all ROIs are contiguous and described by a single contour, the the mean image and ROI locations for one trial can be plotted as follows:

% Load the FISSA output data
S = load('separated.mat')
trial = 1;
figure; hold on;
% Plot the mean image
imagesc(squeeze(S.means(trial, :, :)));
colormap('gray');
% Plot ROI locations
for i_roi = 1:size(S.result, 1);
    contour = S.roi_polys{i_roi, trial}{1};
    plot(contour(:, 2), contour(:, 1));
end
set(gca, 'YDir', 'reverse');
fissa.core.extract(image, rois, nRegions=4, expansion=1, datahandler=None, verbosity=1, label=None, total=None)[source]#

Extract data for all ROIs in a single 3d array or TIFF file.

New in version 1.0.0.

Parameters
imagestr or array_like shaped (time, height, width)

The imaging data. Either a path to a multipage TIFF file, or 3d array_like data.

roisstr or list of array_like

The regions-of-interest, specified by either a string containing a path to an ImageJ roi zip file, or a list of arrays encoding polygons, or list of binary arrays representing masks.

nRegionsint, default=4

Number of neuropil regions to draw. Use a higher number for densely labelled tissue. Default is 4.

expansionfloat, default=1

Expansion factor for the neuropil region, relative to the ROI area. Default is 1. The total neuropil area will be nRegions * expansion * area(ROI).

datahandlerfissa.extraction.DataHandlerAbstract, optional

A datahandler object for handling ROIs and calcium data. The default is DataHandlerTifffile.

verbosityint, default=1

Level of verbosity. The options are:

  • 0: No outputs.

  • 1: Print extraction start.

  • 2: Print extraction end.

  • 3: Print start of each step within the extraction process.

labelstr or int, optional

The label for the current trial. Only used for reporting progress.

totalint, optional

Total number of trials. Only used for reporting progress.

Returns
tracesnumpy.ndarray shaped (n_rois, nRegions + 1, n_frames)

The raw signal, determined as the average fluorence trace extracted from each ROI and neuropil region.

Each vector traces[i_roi, 0, :] contains the traces for the i_roi-th ROI. The following nRegions arrays in traces[i_roi, 1 : nRegions + 1, :] contain the traces from the nRegions grown neuropil regions surrounding the i_roi-th ROI.

polyslist of list of list of numpy.ndarray shaped (n_nodes, 2)

Polygon contours describing the outline of each region.

For contiguous ROIs, the outline of the i_roi-th ROI is described by the array at polys[i_roi][0][0]. This array is n_nodes rows, each representing the coordinate of a node in (y, x) format. For non-contiguous ROIs, a contour is needed for each disconnected polygon making up the total aggregate ROI. These contours are found at polys[i_roi][0][i_contour].

Similarly, the nRegions neuropil regions are each described by the polygons polys[i_roi][i_neurpil + 1][i_contour] respectively.

meannumpy.ndarray shaped (height, width)

Mean image.

fissa.core.run_fissa(images, rois, folder=None, freq=None, return_deltaf=False, deltaf_across_trials=True, export_to_matfile=False, **kwargs)[source]#

Function-based interface to run FISSA on an experiment.

New in version 1.0.0.

Uses the methodology described in FISSA: A neuropil decontamination toolbox for calcium imaging signals.

Parameters
imagesstr or list

The raw recording data. Should be one of:

  • the path to a directory containing TIFF files (string),

  • a list of paths to TIFF files (list of strings),

  • a list of array_like data already loaded into memory, each shaped (n_frames, height, width).

Note that each TIFF/array is considered a single trial.

roisstr or list

The roi definitions. Should be one of:

  • the path to a directory containing ImageJ ZIP files (string),

  • the path of a single ImageJ ZIP file (string),

  • a list of ImageJ ZIP files (list of strings),

  • a list of arrays, each encoding a ROI polygons,

  • a list of lists of binary arrays, each representing a ROI mask.

This can either be a single roiset for all trials, or a different roiset for each trial.

folderstr, optional

Path to a cache directory from which pre-extracted data will be loaded if present, and saved to otherwise. If folder is unset, the experiment data will not be saved.

freqfloat, optional

Imaging frequency, in Hz. Required if return_deltaf=True.

return_deltafbool, optional

Whether to return Δf/f0. Otherwise, the decontaminated signal is returned scaled against the raw recording. Default is False.

deltaf_across_trialsbool, default=True

If True, we estimate a single baseline f0 value across all trials when computing Δf/f0. If False, each trial will have their own baseline f0, and Δf/f0 value will be relative to the trial-specific f0. Default is True.

export_to_matfilebool or str or None, default=False

Whether to export the data to a MATLAB-compatible .mat file. If export_to_matfile is a string, it is used as the path to the output file. If export_to_matfile=True, the matfile is saved to the default path of "separated.mat" within the folder directory, and folder must be set. If this is None, the matfile is exported to the default path if folder is set, and otherwise is not exported. Default is False.

**kwargs

Additional keyword arguments as per Experiment.

Returns
result2d numpy.ndarray of 2d numpy.ndarrays of np.float64

The vector result[roi, trial][0, :] is the trace from ROI roi in trial trial. If return_deltaf=True, this is Δf/f0; otherwise, it is the decontaminated signal scaled as per the raw signal. f0 is the baseline as calculated from the raw signal.

raw2d numpy.ndarray of 2d numpy.ndarrays of np.float64

The raw traces without separation. The vector raw[c, t][0, :] is the ROI trace from cell c in trial t. The vector raw[c, t][i, :] for i>=1 the trace from cell c in trial t, from neuropil region i-1. If return_deltaf=True, this is Δf/f0; otherwise it’s the raw extracted signal.

fissa.core.separate_trials(raw, alpha=0.1, max_iter=20000, tol=0.0001, max_tries=1, method='nmf', verbosity=1, label=None, total=None)[source]#

Separate signals within a set of 2d arrays.

New in version 1.0.0.

Parameters
rawlist of n_trials array_like, each shaped (nRegions + 1, observations)

Raw signals. A list of 2-d arrays, each of which contains observations of mixed signals, mixed in the same way across all trials. The nRegions signals must be the same for each trial, and the 0-th region, raw[trial][0], should be from the region of interest for which a matching source signal should be identified.

alphafloat, default=0.1

Sparsity regularizaton weight for NMF algorithm. Set to zero to remove regularization. Default is 0.1. (Only used for method="nmf".)

max_iterint, default=20000

Maximum number of iterations before timing out on an attempt.

tolfloat, default=1e-4

Tolerance of the stopping condition.

max_triesint, default=1

Maximum number of random initial states to try. Each random state will be optimized for max_iter iterations before timing out.

method{“nmf”, “ica”}, default=”nmf”

Which blind source-separation method to use. Either "nmf" for non-negative matrix factorization, or "ica" for independent component analysis. Default is "nmf".

verbosityint, default=1

Level of verbosity. The options are:

  • 0: No outputs.

  • 1: Print separation start.

  • 2: Print separation end.

  • 3: Print progress details during separation.

labelstr or int, optional

Label/name or index of the ROI currently being processed. Only used for progress messages.

totalint, optional

Total number of ROIs. Only used for reporting progress.

Returns
Xseplist of n_trials numpy.ndarray, each shaped (nRegions + 1, observations)

The separated signals, unordered.

Xmatchlist of n_trials numpy.ndarray, each shaped (nRegions + 1, observations)

The separated traces, ordered by matching score against the raw ROI signal.

Xmixmatnumpy.ndarray, shaped (nRegions + 1, nRegions + 1)

Mixing matrix.

convergencedict

Metadata for the convergence result, with the following keys and values:

convergedbool

Whether the separation model converged, or if it ended due to reaching the maximum number of iterations.

iterationsint

The number of iterations which were needed for the separation model to converge.

max_iterationsint

Maximum number of iterations to use when fitting the separation model.

random_stateint or None

Random seed used to initialise the separation model.