Skip to content

Python API

The pure-compute layer of napari-colocalization is independent of napari and Qt — you can import the metric, masking, and analysis functions from a script or notebook and use them on plain ndarrays.

Documentation index: Home · Usage · Metrics · Python API

The function reference below is generated directly from the source docstrings. The narrative and examples on this page are written by hand; the ::: blocks are filled in automatically by mkdocstrings at build time.

Module layout

napari_colocalization
├── _metrics.py     pearson, spearman, li_icq, manders, costes_threshold
├── _masking.py     shapes_to_label_mask, labels_to_label_mask, iter_regions
├── _analysis.py    analyse_pairwise, analyse_all_to_all, COLUMNS
├── _sample_data.py make_sample_data, make_sample_data_3d,
│                   make_sample_data_cbs006rbm
└── _widget.py      ColocalizationWidget (Qt-only)

The leading underscore is a convention from the napari plugin template; the symbols below are stable and intended to be imported.

Metrics — napari_colocalization._metrics

Pure-numpy colocalization metrics.

No napari / qtpy imports here — every function in this module operates on ndarrays so it can be tested headlessly and reused from scripts. PCC and Manders delegate to scikit-image; SRCC uses scipy.stats; Costes auto-threshold is implemented locally because scikit-image does not provide it.

pearson

pearson(a, b, mask=None)

Pearson correlation coefficient and two-tailed p-value.

Wraps :func:skimage.measure.pearson_corr_coeff. Returns (nan, nan) for inputs with fewer than two samples or zero variance in either channel — both edge cases for which PCC is mathematically undefined.

PARAMETER DESCRIPTION
a

Same-shape intensity arrays of any dimensionality.

TYPE: array_like

b

Same-shape intensity arrays of any dimensionality.

TYPE: array_like

mask

Boolean array with the same shape as a/b. Only pixels where mask is True are included in the calculation. None (the default) uses every pixel.

TYPE: array_like of bool DEFAULT: None

RETURNS DESCRIPTION
pcc

The Pearson correlation coefficient in [-1, 1], or nan if the inputs are degenerate.

TYPE: float

p_value

Two-tailed p-value for the null hypothesis that the coefficient is zero. nan if pcc is nan.

TYPE: float

Examples:

>>> import numpy as np
>>> from napari_colocalization._metrics import pearson
>>> rng = np.random.default_rng(0)
>>> a = rng.random((64, 64))
>>> pearson(a, a)[0]
1.0

spearman

spearman(a, b, mask=None)

Spearman rank correlation coefficient and p-value.

Uses :func:scipy.stats.spearmanr on the masked, flattened intensities. Robust to monotonic non-linearity and outliers.

PARAMETER DESCRIPTION
a

Same-shape intensity arrays of any dimensionality.

TYPE: array_like

b

Same-shape intensity arrays of any dimensionality.

TYPE: array_like

mask

Boolean array with the same shape as a/b. Only pixels where mask is True are included in the calculation.

TYPE: array_like of bool DEFAULT: None

RETURNS DESCRIPTION
rho

Spearman's rank correlation coefficient in [-1, 1], or nan if the inputs are degenerate (fewer than two samples, or zero variance in either channel).

TYPE: float

p_value

Two-tailed p-value for the rank correlation, or nan.

TYPE: float

Examples:

>>> import numpy as np
>>> from napari_colocalization._metrics import spearman
>>> a = np.linspace(0.1, 10.0, 1000)
>>> spearman(a, a ** 3)[0]    # monotonic non-linear
1.0

li_icq

li_icq(a, b, mask=None)

Li's Intensity Correlation Quotient (ICQ).

For each pixel i, form the covariance contribution P_i = (a_i - mean(a)) * (b_i - mean(b)). The ICQ is the fraction of pixels for which this product is positive, re-centred to lie on [-0.5, 0.5]:

.. math:: \mathrm{ICQ} = \frac{|{i : P_i > 0}|}{N} - 0.5

Following Li et al. (2004), values close to +0.5 indicate dependent (co-varying) staining, 0 indicates random staining, and -0.5 indicates segregated (anti-varying) staining.

PARAMETER DESCRIPTION
a

Same-shape intensity arrays of any dimensionality.

TYPE: array_like

b

Same-shape intensity arrays of any dimensionality.

TYPE: array_like

mask

Boolean array with the same shape as a / b. Only pixels where mask is True are counted. None (the default) uses every pixel.

TYPE: array_like of bool DEFAULT: None

RETURNS DESCRIPTION
icq

Value in [-0.5, 0.5], or nan for inputs with fewer than two samples or zero variance in either channel (where the means coincide with every value and the sign of P_i is degenerate).

TYPE: float

Notes

Pixels with P_i == 0 (typically because one channel equals its mean exactly at that pixel) are excluded from the fraction's numerator and denominator, mirroring the convention adopted by ImageJ's Coloc 2 plugin.

References

Li, Q. et al. (2004). A Syntaxin 1, Galpha(o), and N-type Calcium Channel Complex at a Presynaptic Nerve Terminal: Analysis by Quantitative Immunocolocalization. J. Neurosci. 24(16), 4070-4081.

Examples:

>>> import numpy as np
>>> from napari_colocalization._metrics import li_icq
>>> rng = np.random.default_rng(0)
>>> a = rng.random((128, 128))
>>> round(li_icq(a, a), 2)         # perfectly co-varying
0.5
>>> round(li_icq(a, -a), 2)        # perfectly anti-varying
-0.5

manders

manders(a, b, threshold_a, threshold_b, mask=None)

Manders' colocalization coefficients M1 and M2.

M1 is the fraction of the channel-A intensity that lies in pixels where channel B is above threshold_b. M2 is the symmetric counterpart for channel B above threshold_a. Asymmetry between the two is meaningful — it reflects the difference in how much of each channel co-occurs with the other.

Wraps :func:skimage.measure.manders_coloc_coeff, which expects a binary mask for the second image; we threshold internally and then call it twice.

PARAMETER DESCRIPTION
a

Same-shape, non-negative intensity arrays.

TYPE: array_like

b

Same-shape, non-negative intensity arrays.

TYPE: array_like

threshold_a

Per-channel thresholds. Use :func:costes_threshold to derive them automatically, or pass values you have set manually.

TYPE: float

threshold_b

Per-channel thresholds. Use :func:costes_threshold to derive them automatically, or pass values you have set manually.

TYPE: float

mask

Boolean array selecting the analysed region. None analyses every pixel.

TYPE: array_like of bool DEFAULT: None

RETURNS DESCRIPTION
m1, m2 : float

Each in [0, 1], or nan if the analysed region contains no positive intensity in either channel.

Examples:

>>> import numpy as np
>>> from napari_colocalization._metrics import manders
>>> a = np.zeros((10, 10)); a[:, :] = 1.0
>>> b = np.zeros((10, 10)); b[:5, :] = 1.0   # half of A overlaps B
>>> m1, m2 = manders(a, b, threshold_a=0.5, threshold_b=0.5)
>>> round(m1, 2), round(m2, 2)
(0.5, 1.0)

costes_threshold

costes_threshold(a, b, mask=None, n_steps=256)

Costes' iterative auto-threshold for Manders' M1 / M2.

Implements the algorithm of Costes et al. (2004):

  1. Fit a least-squares regression line b = m * a + c over the analysed region.
  2. Walk a candidate threshold T_a downward from max(a). At each step, set T_b = m * T_a + c so the threshold pair lies on the regression line.
  3. The "below-threshold" subset is every pixel with a <= T_a or b <= T_b. Compute its Pearson correlation. Stop at the first T_a where that PCC drops to zero or below — at that point the below-threshold pixels are no longer correlated, i.e. background.
PARAMETER DESCRIPTION
a

Same-shape intensity arrays.

TYPE: array_like

b

Same-shape intensity arrays.

TYPE: array_like

mask

Restrict the regression and search to this region.

TYPE: array_like of bool DEFAULT: None

n_steps

Number of candidate thresholds along [min(a), max(a)].

TYPE: int DEFAULT: 256

RETURNS DESCRIPTION
threshold_a, threshold_b : float

Per-channel thresholds suitable to feed into :func:manders. Falls back to (max(a), max(b)) when the regression slope is non-positive (no co-occurrence to threshold for) and to (min(a), min(b)) when the iteration never reaches PCC <= 0.

Notes

The Costes randomisation significance test (which scrambles pixel blocks to produce a p-value for PCC) is not implemented.

References

Costes, S.V. et al. (2004). Automatic and quantitative measurement of protein-protein colocalization in live cells. Biophys. J. 86(6), 3993-4003.

Masking — napari_colocalization._masking

Region-of-interest helpers.

These functions convert a napari Shapes or Labels layer into an integer label mask (0 = background, 1..N = regions) and iterate over the non-zero regions. The helpers duck-type the layer interface — they do not import napari — so they can be tested with simple stand-in objects.

shapes_to_label_mask

shapes_to_label_mask(shapes_layer, image_shape)

Rasterise a napari Shapes layer to an integer label mask.

Uses Shapes.to_labels(labels_shape=image_shape). The resulting mask has shape image_shape, with 0 outside any shape and 1..N inside each shape (where shape i of the layer maps to label i + 1).

PARAMETER DESCRIPTION
shapes_layer

Any object that exposes .to_labels(labels_shape=...); duck-typed so tests can pass a stand-in.

TYPE: Shapes

image_shape

The desired output shape, typically the spatial shape of the image being analysed.

TYPE: tuple of int

RETURNS DESCRIPTION
mask

Integer label image of shape image_shape.

TYPE: numpy.ndarray of int

labels_to_label_mask

labels_to_label_mask(labels_layer, image_shape)

Validate and return integer data from a napari Labels layer.

PARAMETER DESCRIPTION
labels_layer

Any object that exposes a .data attribute holding an integer ndarray.

TYPE: Labels

image_shape

The expected shape; labels_layer.data.shape must match exactly.

TYPE: tuple of int

RETURNS DESCRIPTION
mask

The label data, dtype-cast to int.

TYPE: numpy.ndarray of int

RAISES DESCRIPTION
ValueError

If the layer's data shape does not match image_shape.

iter_regions

iter_regions(label_mask)

Yield (label_id, bool_mask) for each non-zero region.

PARAMETER DESCRIPTION
label_mask

Integer label image. None is treated specially as "use the whole image".

TYPE: numpy.ndarray of int, or None

YIELDS DESCRIPTION
label_id

The integer label value (1, 2, ...). When label_mask is None, yields 0 instead.

TYPE:: int

bool_mask

Boolean mask of the same shape as label_mask, True where the label equals label_id. None when label_mask is None.

TYPE:: numpy.ndarray of bool, or None

Examples:

>>> import numpy as np
>>> from napari_colocalization._masking import iter_regions
>>> mask = np.zeros((4, 4), dtype=int)
>>> mask[:2, :2] = 1
>>> mask[2:, 2:] = 2
>>> [(label, int(m.sum())) for label, m in iter_regions(mask)]
[(1, 4), (2, 4)]

Analysis — napari_colocalization._analysis

COLUMNS is the canonical column order shared by the table, the CSV export, and the row dicts:

('region', 'channel_a', 'channel_b', 'n_pixels',
 'pcc', 'pcc_pvalue', 'srcc', 'srcc_pvalue',
 'icq',
 'm1', 'm2', 'threshold_a', 'threshold_b')

High-level orchestrators that produce a result row per region.

These accept ndarrays plus a label mask, walk every non-zero region, compute the requested metrics, and return a list of dicts suitable for direct insertion into a table or a CSV.

analyse_pairwise

analyse_pairwise(a, b, *, label_mask=None, metrics=ALL_METRICS, threshold_method='costes', threshold_a=None, threshold_b=None, channel_a='A', channel_b='B')

Compute selected metrics for each region of two grayscale arrays.

Walks every non-zero region of label_mask (or analyses the whole image when label_mask is None) and computes the requested metrics within each region. Missing metrics are written as NaN so the output schema is constant.

PARAMETER DESCRIPTION
a

Same-shape intensity arrays (2D or 3D).

TYPE: ndarray

b

Same-shape intensity arrays (2D or 3D).

TYPE: ndarray

label_mask

Integer label image (0 = background, 1..N = regions). None means analyse the whole image as one region.

TYPE: numpy.ndarray of int DEFAULT: None

metrics

Which metrics to compute. Defaults to all four.

TYPE: sequence of {'pcc', 'srcc', 'icq', 'mcc'} DEFAULT: ALL_METRICS

threshold_method

Only used when 'mcc' is in metrics. 'costes' runs :func:._metrics.costes_threshold per region; 'manual' uses the values supplied below.

TYPE: (costes, manual) DEFAULT: 'costes'

threshold_a

Required when threshold_method='manual'.

TYPE: float DEFAULT: None

threshold_b

Required when threshold_method='manual'.

TYPE: float DEFAULT: None

channel_a

Display names for the two channels — written into the result rows so they round-trip into the table and CSV.

TYPE: str DEFAULT: 'A' and 'B'

channel_b

Display names for the two channels — written into the result rows so they round-trip into the table and CSV.

TYPE: str DEFAULT: 'A' and 'B'

RETURNS DESCRIPTION
rows

One row per region. Each row has the keys listed in :data:COLUMNS.

TYPE: list of dict

RAISES DESCRIPTION
ValueError

If a and b have different shapes, or if threshold_method='manual' is given without both thresholds, or if threshold_method is unknown.

Examples:

>>> import numpy as np
>>> from napari_colocalization._analysis import analyse_pairwise
>>> rng = np.random.default_rng(0)
>>> a = rng.random((32, 32)); b = a.copy()
>>> rows = analyse_pairwise(a, b, metrics=('pcc',))
>>> rows[0]['pcc']
1.0

analyse_all_to_all

analyse_all_to_all(image, channel_axis, *, label_mask=None, metrics=ALL_METRICS, threshold_method='costes', threshold_a=None, threshold_b=None, channel_names=None)

Compute metrics for every channel pair (i, j), i < j.

Iterates over each unordered pair of channels along channel_axis and dispatches to :func:analyse_pairwise. For N channels this produces N*(N-1)/2 × R rows, where R is the number of non-zero regions in label_mask (or 1 for the whole image).

PARAMETER DESCRIPTION
image

N-dimensional array with one channel axis. The remaining axes are spatial (the plugin supports up to 3 spatial axes, i.e. 4D total).

TYPE: ndarray

channel_axis

Axis index along which channels are enumerated.

TYPE: int

label_mask

Integer label image whose shape matches image with channel_axis removed.

TYPE: numpy.ndarray of int DEFAULT: None

metrics

Forwarded to :func:analyse_pairwise.

TYPE: sequence of {'pcc', 'srcc', 'mcc'} DEFAULT: ALL_METRICS

threshold_method

Forwarded to :func:analyse_pairwise. Manual thresholds apply identically to every channel pair.

TYPE: (costes, manual) DEFAULT: 'costes'

threshold_a

Forwarded to :func:analyse_pairwise.

TYPE: float DEFAULT: None

threshold_b

Forwarded to :func:analyse_pairwise.

TYPE: float DEFAULT: None

channel_names

One name per channel along channel_axis. Defaults to ['c0', 'c1', ...].

TYPE: sequence of str DEFAULT: None

RETURNS DESCRIPTION
rows

Concatenation of the per-pair results from :func:analyse_pairwise.

TYPE: list of dict

RAISES DESCRIPTION
ValueError

If len(channel_names) does not match the channel count.

Sample data — napari_colocalization._sample_data

make_sample_data_cbs006rbm downloads the CBS006RBM benchmark image from the Colocalization Benchmark Source on first use and caches it under ~/.cache/napari-colocalization/.

Synthetic two-channel sample data for colocalization.

Synthetic samples (2D, 3D): gaussian blobs with channel B copying 60 % of channel A's blobs plus independent blobs and noise, giving PCC ~ 0.7 and Manders M1/M2 ~ 0.6.

CBS006RBM: red and blue channels with ~50 % colocalization, from the Colocalization Benchmark Source. The TIFF is downloaded once and cached in ~/.cache/napari-colocalization/.

make_sample_data

make_sample_data()

2D sample: 256x256 channels of gaussian blobs.

make_sample_data_3d

make_sample_data_3d()

3D sample: 32x128x128 (Z, Y, X) volumes of gaussian blobs.

make_sample_data_cbs006rbm

make_sample_data_cbs006rbm()

2D sample from the Colocalization Benchmark Source (CBS006RBM).

Red and blue channels with ~50 % colocalization. Downloaded once and cached under ~/.cache/napari-colocalization/.

Source: https://colocalization-benchmark.com/

Putting it together: scripted analysis

import numpy as np
import pandas as pd
from napari_colocalization._analysis import analyse_pairwise, COLUMNS

a = np.load('channel_a.npy')
b = np.load('channel_b.npy')
mask = np.load('cell_labels.npy')   # 0 = bg, 1..N = cells

rows = analyse_pairwise(
    a, b,
    label_mask=mask,
    metrics=('pcc', 'srcc', 'mcc'),
    threshold_method='costes',
    channel_a='dna', channel_b='tubulin',
)

df = pd.DataFrame(rows, columns=COLUMNS)
df.to_csv('colocalization_per_cell.csv', index=False)
print(df.describe())

Stability

  • The function signatures, return shapes and COLUMNS tuple above are intended to be stable across point releases.
  • Internal helpers prefixed with _ (private functions inside each module) may change without notice.
  • ColocalizationWidget (_widget.py) is the GUI surface and not part of the scripting API; for scripts, drive analyse_* directly.