active_contour_model.py 7.7 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250
  1. import numpy as np
  2. from scipy.interpolate import RectBivariateSpline
  3. from .._shared.utils import _supported_float_type
  4. from ..util import img_as_float
  5. from ..filters import sobel
  6. def active_contour(
  7. image,
  8. snake,
  9. alpha=0.01,
  10. beta=0.1,
  11. w_line=0.0,
  12. w_edge=1,
  13. gamma=0.01,
  14. max_px_move=1.0,
  15. max_num_iter=2500,
  16. convergence=0.1,
  17. *,
  18. boundary_condition='periodic',
  19. ):
  20. """Active contour model.
  21. Active contours by fitting snakes to features of images. Supports single
  22. and multichannel 2D images. Snakes can be periodic (for segmentation) or
  23. have fixed and/or free ends.
  24. The output snake has the same length as the input boundary.
  25. As the number of points is constant, make sure that the initial snake
  26. has enough points to capture the details of the final contour.
  27. Parameters
  28. ----------
  29. image : (M, N) or (M, N, 3) ndarray
  30. Input image.
  31. snake : (K, 2) ndarray
  32. Initial snake coordinates. For periodic boundary conditions, endpoints
  33. must not be duplicated.
  34. alpha : float, optional
  35. Snake length shape parameter. Higher values makes snake contract
  36. faster.
  37. beta : float, optional
  38. Snake smoothness shape parameter. Higher values makes snake smoother.
  39. w_line : float, optional
  40. Controls attraction to brightness. Use negative values to attract
  41. toward dark regions.
  42. w_edge : float, optional
  43. Controls attraction to edges. Use negative values to repel snake from
  44. edges.
  45. gamma : float, optional
  46. Explicit time stepping parameter.
  47. max_px_move : float, optional
  48. Maximum pixel distance to move per iteration.
  49. max_num_iter : int, optional
  50. Maximum iterations to optimize snake shape.
  51. convergence : float, optional
  52. Convergence criteria.
  53. boundary_condition : str, optional
  54. Boundary conditions for the contour. Can be one of 'periodic',
  55. 'free', 'fixed', 'free-fixed', or 'fixed-free'. 'periodic' attaches
  56. the two ends of the snake, 'fixed' holds the end-points in place,
  57. and 'free' allows free movement of the ends. 'fixed' and 'free' can
  58. be combined by parsing 'fixed-free', 'free-fixed'. Parsing
  59. 'fixed-fixed' or 'free-free' yields same behaviour as 'fixed' and
  60. 'free', respectively.
  61. Returns
  62. -------
  63. snake : (K, 2) ndarray
  64. Optimised snake, same shape as input parameter.
  65. References
  66. ----------
  67. .. [1] Kass, M.; Witkin, A.; Terzopoulos, D. "Snakes: Active contour
  68. models". International Journal of Computer Vision 1 (4): 321
  69. (1988). :DOI:`10.1007/BF00133570`
  70. Examples
  71. --------
  72. >>> from skimage.draw import circle_perimeter
  73. >>> from skimage.filters import gaussian
  74. Create and smooth image:
  75. >>> img = np.zeros((100, 100))
  76. >>> rr, cc = circle_perimeter(35, 45, 25)
  77. >>> img[rr, cc] = 1
  78. >>> img = gaussian(img, sigma=2, preserve_range=False)
  79. Initialize spline:
  80. >>> s = np.linspace(0, 2*np.pi, 100)
  81. >>> init = 50 * np.array([np.sin(s), np.cos(s)]).T + 50
  82. Fit spline to image:
  83. >>> snake = active_contour(img, init, w_edge=0, w_line=1) # doctest: +SKIP
  84. >>> dist = np.sqrt((45-snake[:, 0])**2 + (35-snake[:, 1])**2) # doctest: +SKIP
  85. >>> int(np.mean(dist)) # doctest: +SKIP
  86. 25
  87. """
  88. max_num_iter = int(max_num_iter)
  89. if max_num_iter <= 0:
  90. raise ValueError("max_num_iter should be >0.")
  91. convergence_order = 10
  92. valid_bcs = [
  93. 'periodic',
  94. 'free',
  95. 'fixed',
  96. 'free-fixed',
  97. 'fixed-free',
  98. 'fixed-fixed',
  99. 'free-free',
  100. ]
  101. if boundary_condition not in valid_bcs:
  102. raise ValueError(
  103. "Invalid boundary condition.\n"
  104. + "Should be one of: "
  105. + ", ".join(valid_bcs)
  106. + '.'
  107. )
  108. img = img_as_float(image)
  109. float_dtype = _supported_float_type(image.dtype)
  110. img = img.astype(float_dtype, copy=False)
  111. RGB = img.ndim == 3
  112. # Find edges using sobel:
  113. if w_edge != 0:
  114. if RGB:
  115. edge = [sobel(img[:, :, 0]), sobel(img[:, :, 1]), sobel(img[:, :, 2])]
  116. else:
  117. edge = [sobel(img)]
  118. else:
  119. edge = [0]
  120. # Superimpose intensity and edge images:
  121. if RGB:
  122. img = w_line * np.sum(img, axis=2) + w_edge * sum(edge)
  123. else:
  124. img = w_line * img + w_edge * edge[0]
  125. # Interpolate for smoothness:
  126. intp = RectBivariateSpline(
  127. np.arange(img.shape[1]), np.arange(img.shape[0]), img.T, kx=2, ky=2, s=0
  128. )
  129. snake_xy = snake[:, ::-1]
  130. x = snake_xy[:, 0].astype(float_dtype)
  131. y = snake_xy[:, 1].astype(float_dtype)
  132. n = len(x)
  133. xsave = np.empty((convergence_order, n), dtype=float_dtype)
  134. ysave = np.empty((convergence_order, n), dtype=float_dtype)
  135. # Build snake shape matrix for Euler equation in double precision
  136. eye_n = np.eye(n, dtype=float)
  137. a = (
  138. np.roll(eye_n, -1, axis=0) + np.roll(eye_n, -1, axis=1) - 2 * eye_n
  139. ) # second order derivative, central difference
  140. b = (
  141. np.roll(eye_n, -2, axis=0)
  142. + np.roll(eye_n, -2, axis=1)
  143. - 4 * np.roll(eye_n, -1, axis=0)
  144. - 4 * np.roll(eye_n, -1, axis=1)
  145. + 6 * eye_n
  146. ) # fourth order derivative, central difference
  147. A = -alpha * a + beta * b
  148. # Impose boundary conditions different from periodic:
  149. sfixed = False
  150. if boundary_condition.startswith('fixed'):
  151. A[0, :] = 0
  152. A[1, :] = 0
  153. A[1, :3] = [1, -2, 1]
  154. sfixed = True
  155. efixed = False
  156. if boundary_condition.endswith('fixed'):
  157. A[-1, :] = 0
  158. A[-2, :] = 0
  159. A[-2, -3:] = [1, -2, 1]
  160. efixed = True
  161. sfree = False
  162. if boundary_condition.startswith('free'):
  163. A[0, :] = 0
  164. A[0, :3] = [1, -2, 1]
  165. A[1, :] = 0
  166. A[1, :4] = [-1, 3, -3, 1]
  167. sfree = True
  168. efree = False
  169. if boundary_condition.endswith('free'):
  170. A[-1, :] = 0
  171. A[-1, -3:] = [1, -2, 1]
  172. A[-2, :] = 0
  173. A[-2, -4:] = [-1, 3, -3, 1]
  174. efree = True
  175. # Only one inversion is needed for implicit spline energy minimization:
  176. inv = np.linalg.inv(A + gamma * eye_n)
  177. # can use float_dtype once we have computed the inverse in double precision
  178. inv = inv.astype(float_dtype, copy=False)
  179. # Explicit time stepping for image energy minimization:
  180. for i in range(max_num_iter):
  181. # RectBivariateSpline always returns float64, so call astype here
  182. fx = intp(x, y, dx=1, grid=False).astype(float_dtype, copy=False)
  183. fy = intp(x, y, dy=1, grid=False).astype(float_dtype, copy=False)
  184. if sfixed:
  185. fx[0] = 0
  186. fy[0] = 0
  187. if efixed:
  188. fx[-1] = 0
  189. fy[-1] = 0
  190. if sfree:
  191. fx[0] *= 2
  192. fy[0] *= 2
  193. if efree:
  194. fx[-1] *= 2
  195. fy[-1] *= 2
  196. xn = inv @ (gamma * x + fx)
  197. yn = inv @ (gamma * y + fy)
  198. # Movements are capped to max_px_move per iteration:
  199. dx = max_px_move * np.tanh(xn - x)
  200. dy = max_px_move * np.tanh(yn - y)
  201. if sfixed:
  202. dx[0] = 0
  203. dy[0] = 0
  204. if efixed:
  205. dx[-1] = 0
  206. dy[-1] = 0
  207. x += dx
  208. y += dy
  209. # Convergence criteria needs to compare to a number of previous
  210. # configurations since oscillations can occur.
  211. j = i % (convergence_order + 1)
  212. if j < convergence_order:
  213. xsave[j, :] = x
  214. ysave[j, :] = y
  215. else:
  216. dist = np.min(
  217. np.max(np.abs(xsave - x[None, :]) + np.abs(ysave - y[None, :]), 1)
  218. )
  219. if dist < convergence:
  220. break
  221. return np.stack([y, x], axis=1)