scikit-image

Форк
0
/
_colocalization.py 
304 строки · 12.0 Кб
1
import numpy as np
2
from scipy.stats import pearsonr
3

4
from .._shared.utils import check_shape_equality, as_binary_ndarray
5

6
__all__ = [
7
    'pearson_corr_coeff',
8
    'manders_coloc_coeff',
9
    'manders_overlap_coeff',
10
    'intersection_coeff',
11
]
12

13

14
def pearson_corr_coeff(image0, image1, mask=None):
15
    r"""Calculate Pearson's Correlation Coefficient between pixel intensities
16
    in channels.
17

18
    Parameters
19
    ----------
20
    image0 : (M, N) ndarray
21
        Image of channel A.
22
    image1 : (M, N) ndarray
23
        Image of channel 2 to be correlated with channel B.
24
        Must have same dimensions as `image0`.
25
    mask : (M, N) ndarray of dtype bool, optional
26
        Only `image0` and `image1` pixels within this region of interest mask
27
        are included in the calculation. Must have same dimensions as `image0`.
28

29
    Returns
30
    -------
31
    pcc : float
32
        Pearson's correlation coefficient of the pixel intensities between
33
        the two images, within the mask if provided.
34
    p-value : float
35
        Two-tailed p-value.
36

37
    Notes
38
    -----
39
    Pearson's Correlation Coefficient (PCC) measures the linear correlation
40
    between the pixel intensities of the two images. Its value ranges from -1
41
    for perfect linear anti-correlation to +1 for perfect linear correlation.
42
    The calculation of the p-value assumes that the intensities of pixels in
43
    each input image are normally distributed.
44

45
    Scipy's implementation of Pearson's correlation coefficient is used. Please
46
    refer to it for further information and caveats [1]_.
47

48
    .. math::
49
        r = \frac{\sum (A_i - m_A_i) (B_i - m_B_i)}
50
        {\sqrt{\sum (A_i - m_A_i)^2 \sum (B_i - m_B_i)^2}}
51

52
    where
53
        :math:`A_i` is the value of the :math:`i^{th}` pixel in `image0`
54
        :math:`B_i` is the value of the :math:`i^{th}` pixel in `image1`,
55
        :math:`m_A_i` is the mean of the pixel values in `image0`
56
        :math:`m_B_i` is the mean of the pixel values in `image1`
57

58
    A low PCC value does not necessarily mean that there is no correlation
59
    between the two channel intensities, just that there is no linear
60
    correlation. You may wish to plot the pixel intensities of each of the two
61
    channels in a 2D scatterplot and use Spearman's rank correlation if a
62
    non-linear correlation is visually identified [2]_. Also consider if you
63
    are interested in correlation or co-occurence, in which case a method
64
    involving segmentation masks (e.g. MCC or intersection coefficient) may be
65
    more suitable [3]_ [4]_.
66

67
    Providing the mask of only relevant sections of the image (e.g., cells, or
68
    particular cellular compartments) and removing noise is important as the
69
    PCC is sensitive to these measures [3]_ [4]_.
70

71
    References
72
    ----------
73
    .. [1] https://docs.scipy.org/doc/scipy/reference/generated/scipy.stats.pearsonr.html  # noqa
74
    .. [2] https://docs.scipy.org/doc/scipy/reference/generated/scipy.stats.spearmanr.html  # noqa
75
    .. [3] Dunn, K. W., Kamocka, M. M., & McDonald, J. H. (2011). A practical
76
           guide to evaluating colocalization in biological microscopy.
77
           American journal of physiology. Cell physiology, 300(4), C723–C742.
78
           https://doi.org/10.1152/ajpcell.00462.2010
79
    .. [4] Bolte, S. and Cordelières, F.P. (2006), A guided tour into
80
           subcellular colocalization analysis in light microscopy. Journal of
81
           Microscopy, 224: 213-232.
82
           https://doi.org/10.1111/j.1365-2818.2006.01706.x
83
    """
84
    image0 = np.asarray(image0)
85
    image1 = np.asarray(image1)
86
    if mask is not None:
87
        mask = as_binary_ndarray(mask, variable_name="mask")
88
        check_shape_equality(image0, image1, mask)
89
        image0 = image0[mask]
90
        image1 = image1[mask]
91
    else:
92
        check_shape_equality(image0, image1)
93
        # scipy pearsonr function only takes flattened arrays
94
        image0 = image0.reshape(-1)
95
        image1 = image1.reshape(-1)
96

97
    return tuple(float(v) for v in pearsonr(image0, image1))
98

99

100
def manders_coloc_coeff(image0, image1_mask, mask=None):
101
    r"""Manders' colocalization coefficient between two channels.
102

103
    Parameters
104
    ----------
105
    image0 : (M, N) ndarray
106
        Image of channel A. All pixel values should be non-negative.
107
    image1_mask : (M, N) ndarray of dtype bool
108
        Binary mask with segmented regions of interest in channel B.
109
        Must have same dimensions as `image0`.
110
    mask : (M, N) ndarray of dtype bool, optional
111
        Only `image0` pixel values within this region of interest mask are
112
        included in the calculation.
113
        Must have same dimensions as `image0`.
114

115
    Returns
116
    -------
117
    mcc : float
118
        Manders' colocalization coefficient.
119

120
    Notes
121
    -----
122
    Manders' Colocalization Coefficient (MCC) is the fraction of total
123
    intensity of a certain channel (channel A) that is within the segmented
124
    region of a second channel (channel B) [1]_. It ranges from 0 for no
125
    colocalisation to 1 for complete colocalization. It is also referred to
126
    as M1 and M2.
127

128
    MCC is commonly used to measure the colocalization of a particular protein
129
    in a subceullar compartment. Typically a segmentation mask for channel B
130
    is generated by setting a threshold that the pixel values must be above
131
    to be included in the MCC calculation. In this implementation,
132
    the channel B mask is provided as the argument `image1_mask`, allowing
133
    the exact segmentation method to be decided by the user beforehand.
134

135
    The implemented equation is:
136

137
    .. math::
138
        r = \frac{\sum A_{i,coloc}}{\sum A_i}
139

140
    where
141
        :math:`A_i` is the value of the :math:`i^{th}` pixel in `image0`
142
        :math:`A_{i,coloc} = A_i` if :math:`Bmask_i > 0`
143
        :math:`Bmask_i` is the value of the :math:`i^{th}` pixel in
144
        `mask`
145

146
    MCC is sensitive to noise, with diffuse signal in the first channel
147
    inflating its value. Images should be processed to remove out of focus and
148
    background light before the MCC is calculated [2]_.
149

150
    References
151
    ----------
152
    .. [1] Manders, E.M.M., Verbeek, F.J. and Aten, J.A. (1993), Measurement of
153
           co-localization of objects in dual-colour confocal images. Journal
154
           of Microscopy, 169: 375-382.
155
           https://doi.org/10.1111/j.1365-2818.1993.tb03313.x
156
           https://imagej.net/media/manders.pdf
157
    .. [2] Dunn, K. W., Kamocka, M. M., & McDonald, J. H. (2011). A practical
158
           guide to evaluating colocalization in biological microscopy.
159
           American journal of physiology. Cell physiology, 300(4), C723–C742.
160
           https://doi.org/10.1152/ajpcell.00462.2010
161

162
    """
163
    image0 = np.asarray(image0)
164
    image1_mask = as_binary_ndarray(image1_mask, variable_name="image1_mask")
165
    if mask is not None:
166
        mask = as_binary_ndarray(mask, variable_name="mask")
167
        check_shape_equality(image0, image1_mask, mask)
168
        image0 = image0[mask]
169
        image1_mask = image1_mask[mask]
170
    else:
171
        check_shape_equality(image0, image1_mask)
172
    # check non-negative image
173
    if image0.min() < 0:
174
        raise ValueError("image contains negative values")
175

176
    sum = np.sum(image0)
177
    if sum == 0:
178
        return 0
179
    return np.sum(image0 * image1_mask) / sum
180

181

182
def manders_overlap_coeff(image0, image1, mask=None):
183
    r"""Manders' overlap coefficient
184

185
    Parameters
186
    ----------
187
    image0 : (M, N) ndarray
188
        Image of channel A. All pixel values should be non-negative.
189
    image1 : (M, N) ndarray
190
        Image of channel B. All pixel values should be non-negative.
191
        Must have same dimensions as `image0`
192
    mask : (M, N) ndarray of dtype bool, optional
193
        Only `image0` and `image1` pixel values within this region of interest
194
        mask are included in the calculation.
195
        Must have ♣same dimensions as `image0`.
196

197
    Returns
198
    -------
199
    moc: float
200
        Manders' Overlap Coefficient of pixel intensities between the two
201
        images.
202

203
    Notes
204
    -----
205
    Manders' Overlap Coefficient (MOC) is given by the equation [1]_:
206

207
    .. math::
208
        r = \frac{\sum A_i B_i}{\sqrt{\sum A_i^2 \sum B_i^2}}
209

210
    where
211
        :math:`A_i` is the value of the :math:`i^{th}` pixel in `image0`
212
        :math:`B_i` is the value of the :math:`i^{th}` pixel in `image1`
213

214
    It ranges between 0 for no colocalization and 1 for complete colocalization
215
    of all pixels.
216

217
    MOC does not take into account pixel intensities, just the fraction of
218
    pixels that have positive values for both channels[2]_ [3]_. Its usefulness
219
    has been criticized as it changes in response to differences in both
220
    co-occurence and correlation and so a particular MOC value could indicate
221
    a wide range of colocalization patterns [4]_ [5]_.
222

223
    References
224
    ----------
225
    .. [1] Manders, E.M.M., Verbeek, F.J. and Aten, J.A. (1993), Measurement of
226
           co-localization of objects in dual-colour confocal images. Journal
227
           of Microscopy, 169: 375-382.
228
           https://doi.org/10.1111/j.1365-2818.1993.tb03313.x
229
           https://imagej.net/media/manders.pdf
230
    .. [2] Dunn, K. W., Kamocka, M. M., & McDonald, J. H. (2011). A practical
231
           guide to evaluating colocalization in biological microscopy.
232
           American journal of physiology. Cell physiology, 300(4), C723–C742.
233
           https://doi.org/10.1152/ajpcell.00462.2010
234
    .. [3] Bolte, S. and Cordelières, F.P. (2006), A guided tour into
235
           subcellular colocalization analysis in light microscopy. Journal of
236
           Microscopy, 224: 213-232.
237
           https://doi.org/10.1111/j.1365-2818.2006.01
238
    .. [4] Adler J, Parmryd I. (2010), Quantifying colocalization by
239
           correlation: the Pearson correlation coefficient is
240
           superior to the Mander's overlap coefficient. Cytometry A.
241
           Aug;77(8):733-42.https://doi.org/10.1002/cyto.a.20896
242
    .. [5] Adler, J, Parmryd, I. Quantifying colocalization: The case for
243
           discarding the Manders overlap coefficient. Cytometry. 2021; 99:
244
           910– 920. https://doi.org/10.1002/cyto.a.24336
245

246
    """
247
    image0 = np.asarray(image0)
248
    image1 = np.asarray(image1)
249
    if mask is not None:
250
        mask = as_binary_ndarray(mask, variable_name="mask")
251
        check_shape_equality(image0, image1, mask)
252
        image0 = image0[mask]
253
        image1 = image1[mask]
254
    else:
255
        check_shape_equality(image0, image1)
256

257
    # check non-negative image
258
    if image0.min() < 0:
259
        raise ValueError("image0 contains negative values")
260
    if image1.min() < 0:
261
        raise ValueError("image1 contains negative values")
262

263
    denom = (np.sum(np.square(image0)) * (np.sum(np.square(image1)))) ** 0.5
264
    return np.sum(np.multiply(image0, image1)) / denom
265

266

267
def intersection_coeff(image0_mask, image1_mask, mask=None):
268
    r"""Fraction of a channel's segmented binary mask that overlaps with a
269
    second channel's segmented binary mask.
270

271
    Parameters
272
    ----------
273
    image0_mask : (M, N) ndarray of dtype bool
274
        Image mask of channel A.
275
    image1_mask : (M, N) ndarray of dtype bool
276
        Image mask of channel B.
277
        Must have same dimensions as `image0_mask`.
278
    mask : (M, N) ndarray of dtype bool, optional
279
        Only `image0_mask` and `image1_mask` pixels within this region of
280
        interest
281
        mask are included in the calculation.
282
        Must have same dimensions as `image0_mask`.
283

284
    Returns
285
    -------
286
    Intersection coefficient, float
287
        Fraction of `image0_mask` that overlaps with `image1_mask`.
288

289
    """
290
    image0_mask = as_binary_ndarray(image0_mask, variable_name="image0_mask")
291
    image1_mask = as_binary_ndarray(image1_mask, variable_name="image1_mask")
292
    if mask is not None:
293
        mask = as_binary_ndarray(mask, variable_name="mask")
294
        check_shape_equality(image0_mask, image1_mask, mask)
295
        image0_mask = image0_mask[mask]
296
        image1_mask = image1_mask[mask]
297
    else:
298
        check_shape_equality(image0_mask, image1_mask)
299

300
    nonzero_image0 = np.count_nonzero(image0_mask)
301
    if nonzero_image0 == 0:
302
        return 0
303
    nonzero_joint = np.count_nonzero(np.logical_and(image0_mask, image1_mask))
304
    return nonzero_joint / nonzero_image0
305

Использование cookies

Мы используем файлы cookie в соответствии с Политикой конфиденциальности и Политикой использования cookies.

Нажимая кнопку «Принимаю», Вы даете АО «СберТех» согласие на обработку Ваших персональных данных в целях совершенствования нашего веб-сайта и Сервиса GitVerse, а также повышения удобства их использования.

Запретить использование cookies Вы можете самостоятельно в настройках Вашего браузера.