morphsnakes.py 15 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449
  1. from itertools import cycle
  2. import numpy as np
  3. from scipy import ndimage as ndi
  4. from .._shared.utils import check_nD
  5. __all__ = [
  6. 'morphological_chan_vese',
  7. 'morphological_geodesic_active_contour',
  8. 'inverse_gaussian_gradient',
  9. 'disk_level_set',
  10. 'checkerboard_level_set',
  11. ]
  12. class _fcycle:
  13. def __init__(self, iterable):
  14. """Call functions from the iterable each time it is called."""
  15. self.funcs = cycle(iterable)
  16. def __call__(self, *args, **kwargs):
  17. f = next(self.funcs)
  18. return f(*args, **kwargs)
  19. # SI and IS operators for 2D and 3D.
  20. _P2 = [
  21. np.eye(3),
  22. np.array([[0, 1, 0]] * 3),
  23. np.flipud(np.eye(3)),
  24. np.rot90([[0, 1, 0]] * 3),
  25. ]
  26. _P3 = [np.zeros((3, 3, 3)) for i in range(9)]
  27. _P3[0][:, :, 1] = 1
  28. _P3[1][:, 1, :] = 1
  29. _P3[2][1, :, :] = 1
  30. _P3[3][:, [0, 1, 2], [0, 1, 2]] = 1
  31. _P3[4][:, [0, 1, 2], [2, 1, 0]] = 1
  32. _P3[5][[0, 1, 2], :, [0, 1, 2]] = 1
  33. _P3[6][[0, 1, 2], :, [2, 1, 0]] = 1
  34. _P3[7][[0, 1, 2], [0, 1, 2], :] = 1
  35. _P3[8][[0, 1, 2], [2, 1, 0], :] = 1
  36. def sup_inf(u):
  37. """SI operator."""
  38. if np.ndim(u) == 2:
  39. P = _P2
  40. elif np.ndim(u) == 3:
  41. P = _P3
  42. else:
  43. raise ValueError("u has an invalid number of dimensions " "(should be 2 or 3)")
  44. erosions = []
  45. for P_i in P:
  46. erosions.append(ndi.binary_erosion(u, P_i).astype(np.int8))
  47. return np.stack(erosions, axis=0).max(0)
  48. def inf_sup(u):
  49. """IS operator."""
  50. if np.ndim(u) == 2:
  51. P = _P2
  52. elif np.ndim(u) == 3:
  53. P = _P3
  54. else:
  55. raise ValueError("u has an invalid number of dimensions " "(should be 2 or 3)")
  56. dilations = []
  57. for P_i in P:
  58. dilations.append(ndi.binary_dilation(u, P_i).astype(np.int8))
  59. return np.stack(dilations, axis=0).min(0)
  60. _curvop = _fcycle(
  61. [lambda u: sup_inf(inf_sup(u)), lambda u: inf_sup(sup_inf(u))] # SIoIS
  62. ) # ISoSI
  63. def _check_input(image, init_level_set):
  64. """Check that shapes of `image` and `init_level_set` match."""
  65. check_nD(image, [2, 3])
  66. if len(image.shape) != len(init_level_set.shape):
  67. raise ValueError(
  68. "The dimensions of the initial level set do not "
  69. "match the dimensions of the image."
  70. )
  71. def _init_level_set(init_level_set, image_shape):
  72. """Auxiliary function for initializing level sets with a string.
  73. If `init_level_set` is not a string, it is returned as is.
  74. """
  75. if isinstance(init_level_set, str):
  76. if init_level_set == 'checkerboard':
  77. res = checkerboard_level_set(image_shape)
  78. elif init_level_set == 'disk':
  79. res = disk_level_set(image_shape)
  80. else:
  81. raise ValueError("`init_level_set` not in " "['checkerboard', 'disk']")
  82. else:
  83. res = init_level_set
  84. return res
  85. def disk_level_set(image_shape, *, center=None, radius=None):
  86. """Create a disk level set with binary values.
  87. Parameters
  88. ----------
  89. image_shape : tuple of positive integers
  90. Shape of the image
  91. center : tuple of positive integers, optional
  92. Coordinates of the center of the disk given in (row, column). If not
  93. given, it defaults to the center of the image.
  94. radius : float, optional
  95. Radius of the disk. If not given, it is set to the 75% of the
  96. smallest image dimension.
  97. Returns
  98. -------
  99. out : array with shape `image_shape`
  100. Binary level set of the disk with the given `radius` and `center`.
  101. See Also
  102. --------
  103. checkerboard_level_set
  104. """
  105. if center is None:
  106. center = tuple(i // 2 for i in image_shape)
  107. if radius is None:
  108. radius = min(image_shape) * 3.0 / 8.0
  109. grid = np.mgrid[[slice(i) for i in image_shape]]
  110. grid = (grid.T - center).T
  111. phi = radius - np.sqrt(np.sum((grid) ** 2, 0))
  112. res = np.int8(phi > 0)
  113. return res
  114. def checkerboard_level_set(image_shape, square_size=5):
  115. """Create a checkerboard level set with binary values.
  116. Parameters
  117. ----------
  118. image_shape : tuple of positive integers
  119. Shape of the image.
  120. square_size : int, optional
  121. Size of the squares of the checkerboard. It defaults to 5.
  122. Returns
  123. -------
  124. out : array with shape `image_shape`
  125. Binary level set of the checkerboard.
  126. See Also
  127. --------
  128. disk_level_set
  129. """
  130. grid = np.mgrid[[slice(i) for i in image_shape]]
  131. grid = grid // square_size
  132. # Alternate 0/1 for even/odd numbers.
  133. grid = grid & 1
  134. checkerboard = np.bitwise_xor.reduce(grid, axis=0)
  135. res = np.int8(checkerboard)
  136. return res
  137. def inverse_gaussian_gradient(image, alpha=100.0, sigma=5.0):
  138. """Inverse of gradient magnitude.
  139. Compute the magnitude of the gradients in the image and then inverts the
  140. result in the range [0, 1]. Flat areas are assigned values close to 1,
  141. while areas close to borders are assigned values close to 0.
  142. This function or a similar one defined by the user should be applied over
  143. the image as a preprocessing step before calling
  144. `morphological_geodesic_active_contour`.
  145. Parameters
  146. ----------
  147. image : (M, N) or (L, M, N) array
  148. Grayscale image or volume.
  149. alpha : float, optional
  150. Controls the steepness of the inversion. A larger value will make the
  151. transition between the flat areas and border areas steeper in the
  152. resulting array.
  153. sigma : float, optional
  154. Standard deviation of the Gaussian filter applied over the image.
  155. Returns
  156. -------
  157. gimage : (M, N) or (L, M, N) array
  158. Preprocessed image (or volume) suitable for
  159. `morphological_geodesic_active_contour`.
  160. """
  161. gradnorm = ndi.gaussian_gradient_magnitude(image, sigma, mode='nearest')
  162. return 1.0 / np.sqrt(1.0 + alpha * gradnorm)
  163. def morphological_chan_vese(
  164. image,
  165. num_iter,
  166. init_level_set='checkerboard',
  167. smoothing=1,
  168. lambda1=1,
  169. lambda2=1,
  170. iter_callback=lambda x: None,
  171. ):
  172. """Morphological Active Contours without Edges (MorphACWE)
  173. Active contours without edges implemented with morphological operators. It
  174. can be used to segment objects in images and volumes without well defined
  175. borders. It is required that the inside of the object looks different on
  176. average than the outside (i.e., the inner area of the object should be
  177. darker or lighter than the outer area on average).
  178. Parameters
  179. ----------
  180. image : (M, N) or (L, M, N) array
  181. Grayscale image or volume to be segmented.
  182. num_iter : uint
  183. Number of num_iter to run
  184. init_level_set : str, (M, N) array, or (L, M, N) array
  185. Initial level set. If an array is given, it will be binarized and used
  186. as the initial level set. If a string is given, it defines the method
  187. to generate a reasonable initial level set with the shape of the
  188. `image`. Accepted values are 'checkerboard' and 'disk'. See the
  189. documentation of `checkerboard_level_set` and `disk_level_set`
  190. respectively for details about how these level sets are created.
  191. smoothing : uint, optional
  192. Number of times the smoothing operator is applied per iteration.
  193. Reasonable values are around 1-4. Larger values lead to smoother
  194. segmentations.
  195. lambda1 : float, optional
  196. Weight parameter for the outer region. If `lambda1` is larger than
  197. `lambda2`, the outer region will contain a larger range of values than
  198. the inner region.
  199. lambda2 : float, optional
  200. Weight parameter for the inner region. If `lambda2` is larger than
  201. `lambda1`, the inner region will contain a larger range of values than
  202. the outer region.
  203. iter_callback : function, optional
  204. If given, this function is called once per iteration with the current
  205. level set as the only argument. This is useful for debugging or for
  206. plotting intermediate results during the evolution.
  207. Returns
  208. -------
  209. out : (M, N) or (L, M, N) array
  210. Final segmentation (i.e., the final level set)
  211. See Also
  212. --------
  213. disk_level_set, checkerboard_level_set
  214. Notes
  215. -----
  216. This is a version of the Chan-Vese algorithm that uses morphological
  217. operators instead of solving a partial differential equation (PDE) for the
  218. evolution of the contour. The set of morphological operators used in this
  219. algorithm are proved to be infinitesimally equivalent to the Chan-Vese PDE
  220. (see [1]_). However, morphological operators are do not suffer from the
  221. numerical stability issues typically found in PDEs (it is not necessary to
  222. find the right time step for the evolution), and are computationally
  223. faster.
  224. The algorithm and its theoretical derivation are described in [1]_.
  225. References
  226. ----------
  227. .. [1] A Morphological Approach to Curvature-based Evolution of Curves and
  228. Surfaces, Pablo Márquez-Neila, Luis Baumela, Luis Álvarez. In IEEE
  229. Transactions on Pattern Analysis and Machine Intelligence (PAMI),
  230. 2014, :DOI:`10.1109/TPAMI.2013.106`
  231. """
  232. init_level_set = _init_level_set(init_level_set, image.shape)
  233. _check_input(image, init_level_set)
  234. u = np.int8(init_level_set > 0)
  235. iter_callback(u)
  236. for _ in range(num_iter):
  237. # inside = u > 0
  238. # outside = u <= 0
  239. c0 = (image * (1 - u)).sum() / float((1 - u).sum() + 1e-8)
  240. c1 = (image * u).sum() / float(u.sum() + 1e-8)
  241. # Image attachment
  242. du = np.gradient(u)
  243. abs_du = np.abs(du).sum(0)
  244. aux = abs_du * (lambda1 * (image - c1) ** 2 - lambda2 * (image - c0) ** 2)
  245. u[aux < 0] = 1
  246. u[aux > 0] = 0
  247. # Smoothing
  248. for _ in range(smoothing):
  249. u = _curvop(u)
  250. iter_callback(u)
  251. return u
  252. def morphological_geodesic_active_contour(
  253. gimage,
  254. num_iter,
  255. init_level_set='disk',
  256. smoothing=1,
  257. threshold='auto',
  258. balloon=0,
  259. iter_callback=lambda x: None,
  260. ):
  261. """Morphological Geodesic Active Contours (MorphGAC).
  262. Geodesic active contours implemented with morphological operators. It can
  263. be used to segment objects with visible but noisy, cluttered, broken
  264. borders.
  265. Parameters
  266. ----------
  267. gimage : (M, N) or (L, M, N) array
  268. Preprocessed image or volume to be segmented. This is very rarely the
  269. original image. Instead, this is usually a preprocessed version of the
  270. original image that enhances and highlights the borders (or other
  271. structures) of the object to segment.
  272. :func:`morphological_geodesic_active_contour` will try to stop the contour
  273. evolution in areas where `gimage` is small. See
  274. :func:`inverse_gaussian_gradient` as an example function to
  275. perform this preprocessing. Note that the quality of
  276. :func:`morphological_geodesic_active_contour` might greatly depend on this
  277. preprocessing.
  278. num_iter : uint
  279. Number of num_iter to run.
  280. init_level_set : str, (M, N) array, or (L, M, N) array
  281. Initial level set. If an array is given, it will be binarized and used
  282. as the initial level set. If a string is given, it defines the method
  283. to generate a reasonable initial level set with the shape of the
  284. `image`. Accepted values are 'checkerboard' and 'disk'. See the
  285. documentation of `checkerboard_level_set` and `disk_level_set`
  286. respectively for details about how these level sets are created.
  287. smoothing : uint, optional
  288. Number of times the smoothing operator is applied per iteration.
  289. Reasonable values are around 1-4. Larger values lead to smoother
  290. segmentations.
  291. threshold : float, optional
  292. Areas of the image with a value smaller than this threshold will be
  293. considered borders. The evolution of the contour will stop in these
  294. areas.
  295. balloon : float, optional
  296. Balloon force to guide the contour in non-informative areas of the
  297. image, i.e., areas where the gradient of the image is too small to push
  298. the contour towards a border. A negative value will shrink the contour,
  299. while a positive value will expand the contour in these areas. Setting
  300. this to zero will disable the balloon force.
  301. iter_callback : function, optional
  302. If given, this function is called once per iteration with the current
  303. level set as the only argument. This is useful for debugging or for
  304. plotting intermediate results during the evolution.
  305. Returns
  306. -------
  307. out : (M, N) or (L, M, N) array
  308. Final segmentation (i.e., the final level set)
  309. See Also
  310. --------
  311. inverse_gaussian_gradient, disk_level_set, checkerboard_level_set
  312. Notes
  313. -----
  314. This is a version of the Geodesic Active Contours (GAC) algorithm that uses
  315. morphological operators instead of solving partial differential equations
  316. (PDEs) for the evolution of the contour. The set of morphological operators
  317. used in this algorithm are proved to be infinitesimally equivalent to the
  318. GAC PDEs (see [1]_). However, morphological operators are do not suffer
  319. from the numerical stability issues typically found in PDEs (e.g., it is
  320. not necessary to find the right time step for the evolution), and are
  321. computationally faster.
  322. The algorithm and its theoretical derivation are described in [1]_.
  323. References
  324. ----------
  325. .. [1] A Morphological Approach to Curvature-based Evolution of Curves and
  326. Surfaces, Pablo Márquez-Neila, Luis Baumela, Luis Álvarez. In IEEE
  327. Transactions on Pattern Analysis and Machine Intelligence (PAMI),
  328. 2014, :DOI:`10.1109/TPAMI.2013.106`
  329. """
  330. image = gimage
  331. init_level_set = _init_level_set(init_level_set, image.shape)
  332. _check_input(image, init_level_set)
  333. if threshold == 'auto':
  334. threshold = np.percentile(image, 40)
  335. structure = np.ones((3,) * len(image.shape), dtype=np.int8)
  336. dimage = np.gradient(image)
  337. # threshold_mask = image > threshold
  338. if balloon != 0:
  339. threshold_mask_balloon = image > threshold / np.abs(balloon)
  340. u = np.int8(init_level_set > 0)
  341. iter_callback(u)
  342. for _ in range(num_iter):
  343. # Balloon
  344. if balloon > 0:
  345. aux = ndi.binary_dilation(u, structure)
  346. elif balloon < 0:
  347. aux = ndi.binary_erosion(u, structure)
  348. if balloon != 0:
  349. u[threshold_mask_balloon] = aux[threshold_mask_balloon]
  350. # Image attachment
  351. aux = np.zeros_like(image)
  352. du = np.gradient(u)
  353. for el1, el2 in zip(dimage, du):
  354. aux += el1 * el2
  355. u[aux > 0] = 1
  356. u[aux < 0] = 0
  357. # Smoothing
  358. for _ in range(smoothing):
  359. u = _curvop(u)
  360. iter_callback(u)
  361. return u