_colocalization.py 12 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305
  1. import numpy as np
  2. from scipy.stats import pearsonr
  3. from .._shared.utils import check_shape_equality, as_binary_ndarray
  4. __all__ = [
  5. 'pearson_corr_coeff',
  6. 'manders_coloc_coeff',
  7. 'manders_overlap_coeff',
  8. 'intersection_coeff',
  9. ]
  10. def pearson_corr_coeff(image0, image1, mask=None):
  11. r"""Calculate Pearson's Correlation Coefficient between pixel intensities
  12. in channels.
  13. Parameters
  14. ----------
  15. image0 : (M, N) ndarray
  16. Image of channel A.
  17. image1 : (M, N) ndarray
  18. Image of channel 2 to be correlated with channel B.
  19. Must have same dimensions as `image0`.
  20. mask : (M, N) ndarray of dtype bool, optional
  21. Only `image0` and `image1` pixels within this region of interest mask
  22. are included in the calculation. Must have same dimensions as `image0`.
  23. Returns
  24. -------
  25. pcc : float
  26. Pearson's correlation coefficient of the pixel intensities between
  27. the two images, within the mask if provided.
  28. p-value : float
  29. Two-tailed p-value.
  30. Notes
  31. -----
  32. Pearson's Correlation Coefficient (PCC) measures the linear correlation
  33. between the pixel intensities of the two images. Its value ranges from -1
  34. for perfect linear anti-correlation to +1 for perfect linear correlation.
  35. The calculation of the p-value assumes that the intensities of pixels in
  36. each input image are normally distributed.
  37. Scipy's implementation of Pearson's correlation coefficient is used. Please
  38. refer to it for further information and caveats [1]_.
  39. .. math::
  40. r = \frac{\sum (A_i - m_A_i) (B_i - m_B_i)}
  41. {\sqrt{\sum (A_i - m_A_i)^2 \sum (B_i - m_B_i)^2}}
  42. where
  43. :math:`A_i` is the value of the :math:`i^{th}` pixel in `image0`
  44. :math:`B_i` is the value of the :math:`i^{th}` pixel in `image1`,
  45. :math:`m_A_i` is the mean of the pixel values in `image0`
  46. :math:`m_B_i` is the mean of the pixel values in `image1`
  47. A low PCC value does not necessarily mean that there is no correlation
  48. between the two channel intensities, just that there is no linear
  49. correlation. You may wish to plot the pixel intensities of each of the two
  50. channels in a 2D scatterplot and use Spearman's rank correlation if a
  51. non-linear correlation is visually identified [2]_. Also consider if you
  52. are interested in correlation or co-occurence, in which case a method
  53. involving segmentation masks (e.g. MCC or intersection coefficient) may be
  54. more suitable [3]_ [4]_.
  55. Providing the mask of only relevant sections of the image (e.g., cells, or
  56. particular cellular compartments) and removing noise is important as the
  57. PCC is sensitive to these measures [3]_ [4]_.
  58. References
  59. ----------
  60. .. [1] https://docs.scipy.org/doc/scipy/reference/generated/scipy.stats.pearsonr.html
  61. .. [2] https://docs.scipy.org/doc/scipy/reference/generated/scipy.stats.spearmanr.html
  62. .. [3] Dunn, K. W., Kamocka, M. M., & McDonald, J. H. (2011). A practical
  63. guide to evaluating colocalization in biological microscopy.
  64. American journal of physiology. Cell physiology, 300(4), C723–C742.
  65. https://doi.org/10.1152/ajpcell.00462.2010
  66. .. [4] Bolte, S. and Cordelières, F.P. (2006), A guided tour into
  67. subcellular colocalization analysis in light microscopy. Journal of
  68. Microscopy, 224: 213-232.
  69. https://doi.org/10.1111/j.1365-2818.2006.01706.x
  70. """
  71. image0 = np.asarray(image0)
  72. image1 = np.asarray(image1)
  73. if mask is not None:
  74. mask = as_binary_ndarray(mask, variable_name="mask")
  75. check_shape_equality(image0, image1, mask)
  76. image0 = image0[mask]
  77. image1 = image1[mask]
  78. else:
  79. check_shape_equality(image0, image1)
  80. # scipy pearsonr function only takes flattened arrays
  81. image0 = image0.reshape(-1)
  82. image1 = image1.reshape(-1)
  83. return tuple(float(v) for v in pearsonr(image0, image1))
  84. def manders_coloc_coeff(image0, image1_mask, mask=None):
  85. r"""Manders' colocalization coefficient between two image channels.
  86. Parameters
  87. ----------
  88. image0 : (M, N) ndarray
  89. Input image (first channel). All pixel values should be non-negative.
  90. image1_mask : (M, N) ndarray of dtype bool
  91. Binary image giving the regions of interest in the second channel.
  92. Must have same shape as `image0`.
  93. mask : (M, N) ndarray of dtype bool, optional
  94. Only `image0` pixel values within `mask` are included in the calculation.
  95. Must have same shape as `image0`.
  96. Returns
  97. -------
  98. mcc : float
  99. Manders' colocalization coefficient.
  100. Notes
  101. -----
  102. Manders' colocalization coefficient (MCC) was developed in the context of
  103. confocal biological microscopy, to measure the fraction of colocalizing
  104. objects in each component of a dual-channel image. Out of the total
  105. intensity of, say, channel A, how much is found within the features
  106. (objects) of, say, channel B [1]_? The measure thus ranges from 0 for no
  107. colocalization to 1 for complete colocalization.
  108. MCC is commonly used to measure the colocalization of a particular protein
  109. in a subcelullar compartment. Typically, the mask for channel B is
  110. obtained by thresholding, to segment the features from the background.
  111. In this implementation, channel B is passed directly as a mask
  112. (`image1_mask`), leaving the segmentation step to the user (upstream).
  113. The implemented equation is:
  114. .. math::
  115. mcc = \frac{\sum_i A_{i,coloc}}{\sum_i A_i}
  116. where
  117. - :math:`A_i` is the value of the :math:`i^{th}` pixel in `image0`, and
  118. - :math:`A_{i, coloc} = A_i B_i`, considering that :math:`B_i` is the
  119. (``True`` or ``False``) value of the :math:`i^{th}` pixel in
  120. `image1_mask` cast into int or float (``1`` or ``0``, respectively).
  121. MCC is sensitive to noise, with diffuse signal in the first channel
  122. inflating its value. Therefore, images should be processed beforehand to
  123. remove out-of-focus and background light [2]_.
  124. References
  125. ----------
  126. .. [1] Manders, E.M.M., Verbeek, F.J. and Aten, J.A. (1993), Measurement of
  127. co-localization of objects in dual-colour confocal images. Journal
  128. of Microscopy, 169: 375-382.
  129. https://doi.org/10.1111/j.1365-2818.1993.tb03313.x
  130. https://imagej.net/media/manders.pdf
  131. .. [2] Dunn, K. W., Kamocka, M. M., & McDonald, J. H. (2011). A practical
  132. guide to evaluating colocalization in biological microscopy.
  133. American journal of physiology. Cell physiology, 300(4), C723–C742.
  134. https://doi.org/10.1152/ajpcell.00462.2010
  135. """
  136. image0 = np.asarray(image0)
  137. image1_mask = as_binary_ndarray(image1_mask, variable_name="image1_mask")
  138. if mask is not None:
  139. mask = as_binary_ndarray(mask, variable_name="mask")
  140. check_shape_equality(image0, image1_mask, mask)
  141. image0 = image0[mask]
  142. image1_mask = image1_mask[mask]
  143. else:
  144. check_shape_equality(image0, image1_mask)
  145. # check non-negative image
  146. if image0.min() < 0:
  147. raise ValueError("image contains negative values")
  148. sum = np.sum(image0)
  149. if sum == 0:
  150. return 0
  151. return np.sum(image0 * image1_mask) / sum
  152. def manders_overlap_coeff(image0, image1, mask=None):
  153. r"""Manders' overlap coefficient
  154. Parameters
  155. ----------
  156. image0 : (M, N) ndarray
  157. Image of channel A. All pixel values should be non-negative.
  158. image1 : (M, N) ndarray
  159. Image of channel B. All pixel values should be non-negative.
  160. Must have same dimensions as `image0`
  161. mask : (M, N) ndarray of dtype bool, optional
  162. Only `image0` and `image1` pixel values within this region of interest
  163. mask are included in the calculation.
  164. Must have ♣same dimensions as `image0`.
  165. Returns
  166. -------
  167. moc: float
  168. Manders' Overlap Coefficient of pixel intensities between the two
  169. images.
  170. Notes
  171. -----
  172. Manders' Overlap Coefficient (MOC) is given by the equation [1]_:
  173. .. math::
  174. r = \frac{\sum A_i B_i}{\sqrt{\sum A_i^2 \sum B_i^2}}
  175. where
  176. :math:`A_i` is the value of the :math:`i^{th}` pixel in `image0`
  177. :math:`B_i` is the value of the :math:`i^{th}` pixel in `image1`
  178. It ranges between 0 for no colocalization and 1 for complete colocalization
  179. of all pixels.
  180. MOC does not take into account pixel intensities, just the fraction of
  181. pixels that have positive values for both channels[2]_ [3]_. Its usefulness
  182. has been criticized as it changes in response to differences in both
  183. co-occurence and correlation and so a particular MOC value could indicate
  184. a wide range of colocalization patterns [4]_ [5]_.
  185. References
  186. ----------
  187. .. [1] Manders, E.M.M., Verbeek, F.J. and Aten, J.A. (1993), Measurement of
  188. co-localization of objects in dual-colour confocal images. Journal
  189. of Microscopy, 169: 375-382.
  190. https://doi.org/10.1111/j.1365-2818.1993.tb03313.x
  191. https://imagej.net/media/manders.pdf
  192. .. [2] Dunn, K. W., Kamocka, M. M., & McDonald, J. H. (2011). A practical
  193. guide to evaluating colocalization in biological microscopy.
  194. American journal of physiology. Cell physiology, 300(4), C723–C742.
  195. https://doi.org/10.1152/ajpcell.00462.2010
  196. .. [3] Bolte, S. and Cordelières, F.P. (2006), A guided tour into
  197. subcellular colocalization analysis in light microscopy. Journal of
  198. Microscopy, 224: 213-232.
  199. https://doi.org/10.1111/j.1365-2818.2006.01
  200. .. [4] Adler J, Parmryd I. (2010), Quantifying colocalization by
  201. correlation: the Pearson correlation coefficient is
  202. superior to the Mander's overlap coefficient. Cytometry A.
  203. Aug;77(8):733-42.https://doi.org/10.1002/cyto.a.20896
  204. .. [5] Adler, J, Parmryd, I. Quantifying colocalization: The case for
  205. discarding the Manders overlap coefficient. Cytometry. 2021; 99:
  206. 910– 920. https://doi.org/10.1002/cyto.a.24336
  207. """
  208. image0 = np.asarray(image0)
  209. image1 = np.asarray(image1)
  210. if mask is not None:
  211. mask = as_binary_ndarray(mask, variable_name="mask")
  212. check_shape_equality(image0, image1, mask)
  213. image0 = image0[mask]
  214. image1 = image1[mask]
  215. else:
  216. check_shape_equality(image0, image1)
  217. # check non-negative image
  218. if image0.min() < 0:
  219. raise ValueError("image0 contains negative values")
  220. if image1.min() < 0:
  221. raise ValueError("image1 contains negative values")
  222. denom = (np.sum(np.square(image0)) * (np.sum(np.square(image1)))) ** 0.5
  223. return np.sum(np.multiply(image0, image1)) / denom
  224. def intersection_coeff(image0_mask, image1_mask, mask=None):
  225. r"""Fraction of a channel's segmented binary mask that overlaps with a
  226. second channel's segmented binary mask.
  227. Parameters
  228. ----------
  229. image0_mask : (M, N) ndarray of dtype bool
  230. Image mask of channel A.
  231. image1_mask : (M, N) ndarray of dtype bool
  232. Image mask of channel B.
  233. Must have same dimensions as `image0_mask`.
  234. mask : (M, N) ndarray of dtype bool, optional
  235. Only `image0_mask` and `image1_mask` pixels within this region of
  236. interest
  237. mask are included in the calculation.
  238. Must have same dimensions as `image0_mask`.
  239. Returns
  240. -------
  241. Intersection coefficient, float
  242. Fraction of `image0_mask` that overlaps with `image1_mask`.
  243. """
  244. image0_mask = as_binary_ndarray(image0_mask, variable_name="image0_mask")
  245. image1_mask = as_binary_ndarray(image1_mask, variable_name="image1_mask")
  246. if mask is not None:
  247. mask = as_binary_ndarray(mask, variable_name="mask")
  248. check_shape_equality(image0_mask, image1_mask, mask)
  249. image0_mask = image0_mask[mask]
  250. image1_mask = image1_mask[mask]
  251. else:
  252. check_shape_equality(image0_mask, image1_mask)
  253. nonzero_image0 = np.count_nonzero(image0_mask)
  254. if nonzero_image0 == 0:
  255. return 0
  256. nonzero_joint = np.count_nonzero(np.logical_and(image0_mask, image1_mask))
  257. return nonzero_joint / nonzero_image0