blob.py 28 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723
  1. import math
  2. import numpy as np
  3. import scipy.ndimage as ndi
  4. from scipy import spatial
  5. from .._shared.filters import gaussian
  6. from .._shared.utils import _supported_float_type, check_nD
  7. from ..transform import integral_image
  8. from ..util import img_as_float
  9. from ._hessian_det_appx import _hessian_matrix_det
  10. from .peak import peak_local_max
  11. # This basic blob detection algorithm is based on:
  12. # http://www.cs.utah.edu/~jfishbau/advimproc/project1/ (04.04.2013)
  13. # Theory behind: https://en.wikipedia.org/wiki/Blob_detection (04.04.2013)
  14. def _compute_disk_overlap(d, r1, r2):
  15. """
  16. Compute fraction of surface overlap between two disks of radii
  17. ``r1`` and ``r2``, with centers separated by a distance ``d``.
  18. Parameters
  19. ----------
  20. d : float
  21. Distance between centers.
  22. r1 : float
  23. Radius of the first disk.
  24. r2 : float
  25. Radius of the second disk.
  26. Returns
  27. -------
  28. fraction: float
  29. Fraction of area of the overlap between the two disks.
  30. """
  31. ratio1 = (d**2 + r1**2 - r2**2) / (2 * d * r1)
  32. ratio1 = np.clip(ratio1, -1, 1)
  33. acos1 = math.acos(ratio1)
  34. ratio2 = (d**2 + r2**2 - r1**2) / (2 * d * r2)
  35. ratio2 = np.clip(ratio2, -1, 1)
  36. acos2 = math.acos(ratio2)
  37. a = -d + r2 + r1
  38. b = d - r2 + r1
  39. c = d + r2 - r1
  40. d = d + r2 + r1
  41. area = r1**2 * acos1 + r2**2 * acos2 - 0.5 * math.sqrt(abs(a * b * c * d))
  42. return area / (math.pi * (min(r1, r2) ** 2))
  43. def _compute_sphere_overlap(d, r1, r2):
  44. """
  45. Compute volume overlap fraction between two spheres of radii
  46. ``r1`` and ``r2``, with centers separated by a distance ``d``.
  47. Parameters
  48. ----------
  49. d : float
  50. Distance between centers.
  51. r1 : float
  52. Radius of the first sphere.
  53. r2 : float
  54. Radius of the second sphere.
  55. Returns
  56. -------
  57. fraction: float
  58. Fraction of volume of the overlap between the two spheres.
  59. Notes
  60. -----
  61. See for example http://mathworld.wolfram.com/Sphere-SphereIntersection.html
  62. for more details.
  63. """
  64. vol = (
  65. math.pi
  66. / (12 * d)
  67. * (r1 + r2 - d) ** 2
  68. * (d**2 + 2 * d * (r1 + r2) - 3 * (r1**2 + r2**2) + 6 * r1 * r2)
  69. )
  70. return vol / (4.0 / 3 * math.pi * min(r1, r2) ** 3)
  71. def _blob_overlap(blob1, blob2, *, sigma_dim=1):
  72. """Finds the overlapping area fraction between two blobs.
  73. Returns a float representing fraction of overlapped area. Note that 0.0
  74. is *always* returned for dimension greater than 3.
  75. Parameters
  76. ----------
  77. blob1 : sequence of arrays
  78. A sequence of ``(row, col, sigma)`` or ``(pln, row, col, sigma)``,
  79. where ``row, col`` (or ``(pln, row, col)``) are coordinates
  80. of blob and ``sigma`` is the standard deviation of the Gaussian kernel
  81. which detected the blob.
  82. blob2 : sequence of arrays
  83. A sequence of ``(row, col, sigma)`` or ``(pln, row, col, sigma)``,
  84. where ``row, col`` (or ``(pln, row, col)``) are coordinates
  85. of blob and ``sigma`` is the standard deviation of the Gaussian kernel
  86. which detected the blob.
  87. sigma_dim : int, optional
  88. The dimensionality of the sigma value. Can be 1 or the same as the
  89. dimensionality of the blob space (2 or 3).
  90. Returns
  91. -------
  92. f : float
  93. Fraction of overlapped area (or volume in 3D).
  94. """
  95. ndim = len(blob1) - sigma_dim
  96. if ndim > 3:
  97. return 0.0
  98. root_ndim = math.sqrt(ndim)
  99. # we divide coordinates by sigma * sqrt(ndim) to rescale space to isotropy,
  100. # giving spheres of radius = 1 or < 1.
  101. if blob1[-1] == blob2[-1] == 0:
  102. return 0.0
  103. elif blob1[-1] > blob2[-1]:
  104. max_sigma = blob1[-sigma_dim:]
  105. r1 = 1
  106. r2 = blob2[-1] / blob1[-1]
  107. else:
  108. max_sigma = blob2[-sigma_dim:]
  109. r2 = 1
  110. r1 = blob1[-1] / blob2[-1]
  111. pos1 = blob1[:ndim] / (max_sigma * root_ndim)
  112. pos2 = blob2[:ndim] / (max_sigma * root_ndim)
  113. d = np.sqrt(np.sum((pos2 - pos1) ** 2))
  114. if d > r1 + r2: # centers farther than sum of radii, so no overlap
  115. return 0.0
  116. # one blob is inside the other
  117. if d <= abs(r1 - r2):
  118. return 1.0
  119. if ndim == 2:
  120. return _compute_disk_overlap(d, r1, r2)
  121. else: # ndim=3 http://mathworld.wolfram.com/Sphere-SphereIntersection.html
  122. return _compute_sphere_overlap(d, r1, r2)
  123. def _prune_blobs(blobs_array, overlap, *, sigma_dim=1):
  124. """Eliminated blobs with area overlap.
  125. Parameters
  126. ----------
  127. blobs_array : ndarray
  128. A 2d array with each row representing 3 (or 4) values,
  129. ``(row, col, sigma)`` or ``(pln, row, col, sigma)`` in 3D,
  130. where ``(row, col)`` (``(pln, row, col)``) are coordinates of the blob
  131. and ``sigma`` is the standard deviation of the Gaussian kernel which
  132. detected the blob.
  133. This array must not have a dimension of size 0.
  134. overlap : float
  135. A value between 0 and 1. If the fraction of area overlapping for 2
  136. blobs is greater than `overlap` the smaller blob is eliminated.
  137. sigma_dim : int, optional
  138. The number of columns in ``blobs_array`` corresponding to sigmas rather
  139. than positions.
  140. Returns
  141. -------
  142. A : ndarray
  143. `array` with overlapping blobs removed.
  144. """
  145. sigma = blobs_array[:, -sigma_dim:].max()
  146. distance = 2 * sigma * math.sqrt(blobs_array.shape[1] - sigma_dim)
  147. tree = spatial.cKDTree(blobs_array[:, :-sigma_dim])
  148. pairs = np.array(list(tree.query_pairs(distance)))
  149. if len(pairs) == 0:
  150. return blobs_array
  151. else:
  152. for i, j in pairs:
  153. blob1, blob2 = blobs_array[i], blobs_array[j]
  154. if _blob_overlap(blob1, blob2, sigma_dim=sigma_dim) > overlap:
  155. # note: this test works even in the anisotropic case because
  156. # all sigmas increase together.
  157. if blob1[-1] > blob2[-1]:
  158. blob2[-1] = 0
  159. else:
  160. blob1[-1] = 0
  161. return np.stack([b for b in blobs_array if b[-1] > 0])
  162. def _format_exclude_border(img_ndim, exclude_border):
  163. """Format an ``exclude_border`` argument as a tuple of ints for calling
  164. ``peak_local_max``.
  165. """
  166. if isinstance(exclude_border, tuple):
  167. if len(exclude_border) != img_ndim:
  168. raise ValueError(
  169. "`exclude_border` should have the same length as the "
  170. "dimensionality of the image."
  171. )
  172. for exclude in exclude_border:
  173. if not isinstance(exclude, int):
  174. raise ValueError(
  175. "exclude border, when expressed as a tuple, must only "
  176. "contain ints."
  177. )
  178. return exclude_border + (0,)
  179. elif isinstance(exclude_border, int):
  180. return (exclude_border,) * img_ndim + (0,)
  181. elif exclude_border is True:
  182. raise ValueError("exclude_border cannot be True")
  183. elif exclude_border is False:
  184. return (0,) * (img_ndim + 1)
  185. else:
  186. raise ValueError(f'Unsupported value ({exclude_border}) for exclude_border')
  187. def blob_dog(
  188. image,
  189. min_sigma=1,
  190. max_sigma=50,
  191. sigma_ratio=1.6,
  192. threshold=0.5,
  193. overlap=0.5,
  194. *,
  195. threshold_rel=None,
  196. exclude_border=False,
  197. ):
  198. r"""Finds blobs in the given grayscale image.
  199. Blobs are found using the Difference of Gaussian (DoG) method [1]_, [2]_.
  200. For each blob found, the method returns its coordinates and the standard
  201. deviation of the Gaussian kernel that detected the blob.
  202. Parameters
  203. ----------
  204. image : ndarray
  205. Input grayscale image, blobs are assumed to be light on dark
  206. background (white on black).
  207. min_sigma : scalar or sequence of scalars, optional
  208. Minimum standard deviation for Gaussian kernel. Keep this value low to
  209. detect smaller blobs. The standard deviation of the Gaussian kernel
  210. is given either as a sequence for each axis, or as a single number, in
  211. which case it is equal for all axes.
  212. max_sigma : scalar or sequence of scalars, optional
  213. The maximum standard deviation for Gaussian kernel. Keep this high to
  214. detect larger blobs. The standard deviation of the Gaussian kernel
  215. is given either as a sequence for each axis, or as a single number, in
  216. which case it is equal for all axes.
  217. sigma_ratio : float, optional
  218. The ratio between the standard deviation of Gaussian Kernels used for
  219. computing the Difference of Gaussians
  220. threshold : float or None, optional
  221. The absolute lower bound for scale space maxima. Local maxima smaller
  222. than `threshold` are ignored. Reduce this to detect blobs with lower
  223. intensities. If `threshold_rel` is also specified, whichever threshold
  224. is larger will be used. If None, `threshold_rel` is used instead.
  225. overlap : float, optional
  226. A value between 0 and 1. If the area of two blobs overlaps by a
  227. fraction greater than `threshold`, the smaller blob is eliminated.
  228. threshold_rel : float or None, optional
  229. Minimum intensity of peaks, calculated as
  230. ``max(dog_space) * threshold_rel``, where ``dog_space`` refers to the
  231. stack of Difference-of-Gaussian (DoG) images computed internally. This
  232. should have a value between 0 and 1. If None, `threshold` is used
  233. instead.
  234. exclude_border : tuple of ints, int, or False, optional
  235. If tuple of ints, the length of the tuple must match the input array's
  236. dimensionality. Each element of the tuple will exclude peaks from
  237. within `exclude_border`-pixels of the border of the image along that
  238. dimension.
  239. If nonzero int, `exclude_border` excludes peaks from within
  240. `exclude_border`-pixels of the border of the image.
  241. If zero or False, peaks are identified regardless of their
  242. distance from the border.
  243. Returns
  244. -------
  245. A : (n, image.ndim + sigma) ndarray
  246. A 2d array with each row representing 2 coordinate values for a 2D
  247. image, or 3 coordinate values for a 3D image, plus the sigma(s) used.
  248. When a single sigma is passed, outputs are:
  249. ``(r, c, sigma)`` or ``(p, r, c, sigma)`` where ``(r, c)`` or
  250. ``(p, r, c)`` are coordinates of the blob and ``sigma`` is the standard
  251. deviation of the Gaussian kernel which detected the blob. When an
  252. anisotropic gaussian is used (sigmas per dimension), the detected sigma
  253. is returned for each dimension.
  254. See also
  255. --------
  256. skimage.filters.difference_of_gaussians
  257. References
  258. ----------
  259. .. [1] https://en.wikipedia.org/wiki/Blob_detection#The_difference_of_Gaussians_approach
  260. .. [2] Lowe, D. G. "Distinctive Image Features from Scale-Invariant
  261. Keypoints." International Journal of Computer Vision 60, 91–110 (2004).
  262. https://www.cs.ubc.ca/~lowe/papers/ijcv04.pdf
  263. :DOI:`10.1023/B:VISI.0000029664.99615.94`
  264. Examples
  265. --------
  266. >>> from skimage import data, feature
  267. >>> coins = data.coins()
  268. >>> feature.blob_dog(coins, threshold=.05, min_sigma=10, max_sigma=40)
  269. array([[128., 155., 10.],
  270. [198., 155., 10.],
  271. [124., 338., 10.],
  272. [127., 102., 10.],
  273. [193., 281., 10.],
  274. [126., 208., 10.],
  275. [267., 115., 10.],
  276. [197., 102., 10.],
  277. [198., 215., 10.],
  278. [123., 279., 10.],
  279. [126., 46., 10.],
  280. [259., 247., 10.],
  281. [196., 43., 10.],
  282. [ 54., 276., 10.],
  283. [267., 358., 10.],
  284. [ 58., 100., 10.],
  285. [259., 305., 10.],
  286. [185., 347., 16.],
  287. [261., 174., 16.],
  288. [ 46., 336., 16.],
  289. [ 54., 217., 10.],
  290. [ 55., 157., 10.],
  291. [ 57., 41., 10.],
  292. [260., 47., 16.]])
  293. Notes
  294. -----
  295. The radius of each blob is approximately :math:`\sqrt{2}\sigma` for
  296. a 2-D image and :math:`\sqrt{3}\sigma` for a 3-D image.
  297. """
  298. image = img_as_float(image)
  299. float_dtype = _supported_float_type(image.dtype)
  300. image = image.astype(float_dtype, copy=False)
  301. # if both min and max sigma are scalar, function returns only one sigma
  302. scalar_sigma = np.isscalar(max_sigma) and np.isscalar(min_sigma)
  303. # Gaussian filter requires that sequence-type sigmas have same
  304. # dimensionality as image. This broadcasts scalar kernels
  305. if np.isscalar(max_sigma):
  306. max_sigma = np.full(image.ndim, max_sigma, dtype=float_dtype)
  307. if np.isscalar(min_sigma):
  308. min_sigma = np.full(image.ndim, min_sigma, dtype=float_dtype)
  309. # Convert sequence types to array
  310. min_sigma = np.asarray(min_sigma, dtype=float_dtype)
  311. max_sigma = np.asarray(max_sigma, dtype=float_dtype)
  312. if sigma_ratio <= 1.0:
  313. raise ValueError('sigma_ratio must be > 1.0')
  314. # k such that min_sigma*(sigma_ratio**k) > max_sigma
  315. k = int(np.mean(np.log(max_sigma / min_sigma) / np.log(sigma_ratio) + 1))
  316. # a geometric progression of standard deviations for gaussian kernels
  317. sigma_list = np.array([min_sigma * (sigma_ratio**i) for i in range(k + 1)])
  318. # computing difference between two successive Gaussian blurred images
  319. # to obtain an approximation of the scale invariant Laplacian of the
  320. # Gaussian operator
  321. dog_image_cube = np.empty(image.shape + (k,), dtype=float_dtype)
  322. gaussian_previous = gaussian(image, sigma=sigma_list[0], mode='reflect')
  323. for i, s in enumerate(sigma_list[1:]):
  324. gaussian_current = gaussian(image, sigma=s, mode='reflect')
  325. dog_image_cube[..., i] = gaussian_previous - gaussian_current
  326. gaussian_previous = gaussian_current
  327. # normalization factor for consistency in DoG magnitude
  328. sf = 1 / (sigma_ratio - 1)
  329. dog_image_cube *= sf
  330. exclude_border = _format_exclude_border(image.ndim, exclude_border)
  331. local_maxima = peak_local_max(
  332. dog_image_cube,
  333. threshold_abs=threshold,
  334. threshold_rel=threshold_rel,
  335. exclude_border=exclude_border,
  336. footprint=np.ones((3,) * (image.ndim + 1)),
  337. )
  338. # Catch no peaks
  339. if local_maxima.size == 0:
  340. return np.empty((0, image.ndim + (1 if scalar_sigma else image.ndim)))
  341. # Convert local_maxima to float64
  342. lm = local_maxima.astype(float_dtype)
  343. # translate final column of lm, which contains the index of the
  344. # sigma that produced the maximum intensity value, into the sigma
  345. sigmas_of_peaks = sigma_list[local_maxima[:, -1]]
  346. if scalar_sigma:
  347. # select one sigma column, keeping dimension
  348. sigmas_of_peaks = sigmas_of_peaks[:, 0:1]
  349. # Remove sigma index and replace with sigmas
  350. lm = np.hstack([lm[:, :-1], sigmas_of_peaks])
  351. sigma_dim = sigmas_of_peaks.shape[1]
  352. return _prune_blobs(lm, overlap, sigma_dim=sigma_dim)
  353. def blob_log(
  354. image,
  355. min_sigma=1,
  356. max_sigma=50,
  357. num_sigma=10,
  358. threshold=0.2,
  359. overlap=0.5,
  360. log_scale=False,
  361. *,
  362. threshold_rel=None,
  363. exclude_border=False,
  364. ):
  365. r"""Finds blobs in the given grayscale image.
  366. Blobs are found using the Laplacian of Gaussian (LoG) method [1]_.
  367. For each blob found, the method returns its coordinates and the standard
  368. deviation of the Gaussian kernel that detected the blob.
  369. Parameters
  370. ----------
  371. image : ndarray
  372. Input grayscale image, blobs are assumed to be light on dark
  373. background (white on black).
  374. min_sigma : scalar or sequence of scalars, optional
  375. Minimum standard deviation for Gaussian kernel. Keep this value low to
  376. detect smaller blobs. The standard deviation of the Gaussian kernel
  377. is given either as a sequence for each axis, or as a single number, in
  378. which case it is equal for all axes.
  379. max_sigma : scalar or sequence of scalars, optional
  380. The maximum standard deviation for Gaussian kernel. Keep this high to
  381. detect larger blobs. The standard deviation of the Gaussian kernel
  382. is given either as a sequence for each axis, or as a single number, in
  383. which case it is equal for all axes.
  384. num_sigma : int, optional
  385. The number of evenly spaced values for standard deviation of the
  386. Gaussian kernel to consider on the closed interval
  387. ``[min_sigma, max_sigma]``.
  388. threshold : float or None, optional
  389. The absolute lower bound for scale space maxima. Local maxima smaller
  390. than `threshold` are ignored. Reduce this to detect blobs with lower
  391. intensities. If `threshold_rel` is also specified, whichever threshold
  392. is larger will be used. If None, `threshold_rel` is used instead.
  393. overlap : float, optional
  394. A value between 0 and 1. If the area of two blobs overlaps by a
  395. fraction greater than `threshold`, the smaller blob is eliminated.
  396. log_scale : bool, optional
  397. If set intermediate values of standard deviations are interpolated
  398. using a logarithmic scale to the base `10`. If not, linear
  399. interpolation is used.
  400. threshold_rel : float or None, optional
  401. Minimum intensity of peaks, calculated as
  402. ``max(log_space) * threshold_rel``, where ``log_space`` refers to the
  403. stack of Laplacian-of-Gaussian (LoG) images computed internally. This
  404. should have a value between 0 and 1. If None, `threshold` is used
  405. instead.
  406. exclude_border : tuple of ints, int, or False, optional
  407. If tuple of ints, the length of the tuple must match the input array's
  408. dimensionality. Each element of the tuple will exclude peaks from
  409. within `exclude_border`-pixels of the border of the image along that
  410. dimension.
  411. If nonzero int, `exclude_border` excludes peaks from within
  412. `exclude_border`-pixels of the border of the image.
  413. If zero or False, peaks are identified regardless of their
  414. distance from the border.
  415. Returns
  416. -------
  417. A : (n, image.ndim + sigma) ndarray
  418. A 2d array with each row representing 2 coordinate values for a 2D
  419. image, or 3 coordinate values for a 3D image, plus the sigma(s) used.
  420. When a single sigma is passed, outputs are:
  421. ``(r, c, sigma)`` or ``(p, r, c, sigma)`` where ``(r, c)`` or
  422. ``(p, r, c)`` are coordinates of the blob and ``sigma`` is the standard
  423. deviation of the Gaussian kernel which detected the blob. When an
  424. anisotropic gaussian is used (sigmas per dimension), the detected sigma
  425. is returned for each dimension.
  426. References
  427. ----------
  428. .. [1] https://en.wikipedia.org/wiki/Blob_detection#The_Laplacian_of_Gaussian
  429. Examples
  430. --------
  431. >>> from skimage import data, feature, exposure
  432. >>> img = data.coins()
  433. >>> img = exposure.equalize_hist(img) # improves detection
  434. >>> feature.blob_log(img, threshold = .3)
  435. array([[124. , 336. , 11.88888889],
  436. [198. , 155. , 11.88888889],
  437. [194. , 213. , 17.33333333],
  438. [121. , 272. , 17.33333333],
  439. [263. , 244. , 17.33333333],
  440. [194. , 276. , 17.33333333],
  441. [266. , 115. , 11.88888889],
  442. [128. , 154. , 11.88888889],
  443. [260. , 174. , 17.33333333],
  444. [198. , 103. , 11.88888889],
  445. [126. , 208. , 11.88888889],
  446. [127. , 102. , 11.88888889],
  447. [263. , 302. , 17.33333333],
  448. [197. , 44. , 11.88888889],
  449. [185. , 344. , 17.33333333],
  450. [126. , 46. , 11.88888889],
  451. [113. , 323. , 1. ]])
  452. Notes
  453. -----
  454. The radius of each blob is approximately :math:`\sqrt{2}\sigma` for
  455. a 2-D image and :math:`\sqrt{3}\sigma` for a 3-D image.
  456. """
  457. image = img_as_float(image)
  458. float_dtype = _supported_float_type(image.dtype)
  459. image = image.astype(float_dtype, copy=False)
  460. # if both min and max sigma are scalar, function returns only one sigma
  461. scalar_sigma = True if np.isscalar(max_sigma) and np.isscalar(min_sigma) else False
  462. # Gaussian filter requires that sequence-type sigmas have same
  463. # dimensionality as image. This broadcasts scalar kernels
  464. if np.isscalar(max_sigma):
  465. max_sigma = np.full(image.ndim, max_sigma, dtype=float_dtype)
  466. if np.isscalar(min_sigma):
  467. min_sigma = np.full(image.ndim, min_sigma, dtype=float_dtype)
  468. # Convert sequence types to array
  469. min_sigma = np.asarray(min_sigma, dtype=float_dtype)
  470. max_sigma = np.asarray(max_sigma, dtype=float_dtype)
  471. if log_scale:
  472. start = np.log10(min_sigma)
  473. stop = np.log10(max_sigma)
  474. sigma_list = np.logspace(start, stop, num_sigma)
  475. else:
  476. sigma_list = np.linspace(min_sigma, max_sigma, num_sigma)
  477. # computing gaussian laplace
  478. image_cube = np.empty(image.shape + (len(sigma_list),), dtype=float_dtype)
  479. for i, s in enumerate(sigma_list):
  480. # average s**2 provides scale invariance
  481. image_cube[..., i] = -ndi.gaussian_laplace(image, s) * np.mean(s) ** 2
  482. exclude_border = _format_exclude_border(image.ndim, exclude_border)
  483. local_maxima = peak_local_max(
  484. image_cube,
  485. threshold_abs=threshold,
  486. threshold_rel=threshold_rel,
  487. exclude_border=exclude_border,
  488. footprint=np.ones((3,) * (image.ndim + 1)),
  489. )
  490. # Catch no peaks
  491. if local_maxima.size == 0:
  492. return np.empty((0, image.ndim + (1 if scalar_sigma else image.ndim)))
  493. # Convert local_maxima to float64
  494. lm = local_maxima.astype(float_dtype)
  495. # translate final column of lm, which contains the index of the
  496. # sigma that produced the maximum intensity value, into the sigma
  497. sigmas_of_peaks = sigma_list[local_maxima[:, -1]]
  498. if scalar_sigma:
  499. # select one sigma column, keeping dimension
  500. sigmas_of_peaks = sigmas_of_peaks[:, 0:1]
  501. # Remove sigma index and replace with sigmas
  502. lm = np.hstack([lm[:, :-1], sigmas_of_peaks])
  503. sigma_dim = sigmas_of_peaks.shape[1]
  504. return _prune_blobs(lm, overlap, sigma_dim=sigma_dim)
  505. def blob_doh(
  506. image,
  507. min_sigma=1,
  508. max_sigma=30,
  509. num_sigma=10,
  510. threshold=0.01,
  511. overlap=0.5,
  512. log_scale=False,
  513. *,
  514. threshold_rel=None,
  515. ):
  516. """Finds blobs in the given grayscale image.
  517. Blobs are found using the Determinant of Hessian method [1]_. For each blob
  518. found, the method returns its coordinates and the standard deviation
  519. of the Gaussian Kernel used for the Hessian matrix whose determinant
  520. detected the blob. Determinant of Hessians is approximated using [2]_.
  521. Parameters
  522. ----------
  523. image : 2D ndarray
  524. Input grayscale image. Blobs can either be light on dark or vice versa.
  525. min_sigma : float, optional
  526. The minimum standard deviation for Gaussian Kernel used to compute
  527. Hessian matrix. Keep this value low to detect smaller blobs.
  528. The standard deviation of the Gaussian kernel is given either as a
  529. sequence for each axis, or as a single number, in which case it is
  530. equal for all axes.
  531. max_sigma : float, optional
  532. The maximum standard deviation for Gaussian Kernel used to compute
  533. Hessian matrix. Keep this value high to detect larger blobs.
  534. The standard deviation of the Gaussian kernel is given either as a
  535. sequence for each axis, or as a single number, in which case it is
  536. equal for all axes.
  537. num_sigma : int, optional
  538. The number of evenly spaced values for standard deviation of the
  539. Gaussian kernel to consider on the closed interval
  540. ``[min_sigma, max_sigma]``.
  541. threshold : float or None, optional
  542. The absolute lower bound for scale space maxima. Local maxima smaller
  543. than `threshold` are ignored. Reduce this to detect blobs with lower
  544. intensities. If `threshold_rel` is also specified, whichever threshold
  545. is larger will be used. If None, `threshold_rel` is used instead.
  546. overlap : float, optional
  547. A value between 0 and 1. If the area of two blobs overlaps by a
  548. fraction greater than `threshold`, the smaller blob is eliminated.
  549. log_scale : bool, optional
  550. If set intermediate values of standard deviations are interpolated
  551. using a logarithmic scale to the base `10`. If not, linear
  552. interpolation is used.
  553. threshold_rel : float or None, optional
  554. Minimum intensity of peaks, calculated as
  555. ``max(doh_space) * threshold_rel``, where ``doh_space`` refers to the
  556. stack of Determinant-of-Hessian (DoH) images computed internally. This
  557. should have a value between 0 and 1. If None, `threshold` is used
  558. instead.
  559. Returns
  560. -------
  561. A : (n, 3) ndarray
  562. A 2d array with each row representing 3 values, ``(y,x,sigma)``
  563. where ``(y,x)`` are coordinates of the blob and ``sigma`` is the
  564. standard deviation of the Gaussian kernel of the Hessian Matrix whose
  565. determinant detected the blob.
  566. References
  567. ----------
  568. .. [1] https://en.wikipedia.org/wiki/Blob_detection#The_determinant_of_the_Hessian
  569. .. [2] Herbert Bay, Andreas Ess, Tinne Tuytelaars, Luc Van Gool,
  570. "SURF: Speeded Up Robust Features"
  571. ftp://ftp.vision.ee.ethz.ch/publications/articles/eth_biwi_00517.pdf
  572. Examples
  573. --------
  574. >>> from skimage import data, feature
  575. >>> img = data.coins()
  576. >>> feature.blob_doh(img)
  577. array([[197. , 153. , 20.33333333],
  578. [124. , 336. , 20.33333333],
  579. [126. , 153. , 20.33333333],
  580. [195. , 100. , 23.55555556],
  581. [192. , 212. , 23.55555556],
  582. [121. , 271. , 30. ],
  583. [126. , 101. , 20.33333333],
  584. [193. , 275. , 23.55555556],
  585. [123. , 205. , 20.33333333],
  586. [270. , 363. , 30. ],
  587. [265. , 113. , 23.55555556],
  588. [262. , 243. , 23.55555556],
  589. [185. , 348. , 30. ],
  590. [156. , 302. , 30. ],
  591. [123. , 44. , 23.55555556],
  592. [260. , 173. , 30. ],
  593. [197. , 44. , 20.33333333]])
  594. Notes
  595. -----
  596. The radius of each blob is approximately `sigma`.
  597. Computation of Determinant of Hessians is independent of the standard
  598. deviation. Therefore detecting larger blobs won't take more time. In
  599. methods line :py:meth:`blob_dog` and :py:meth:`blob_log` the computation
  600. of Gaussians for larger `sigma` takes more time. The downside is that
  601. this method can't be used for detecting blobs of radius less than `3px`
  602. due to the box filters used in the approximation of Hessian Determinant.
  603. """
  604. check_nD(image, 2)
  605. image = img_as_float(image)
  606. float_dtype = _supported_float_type(image.dtype)
  607. image = image.astype(float_dtype, copy=False)
  608. image = integral_image(image)
  609. if log_scale:
  610. start, stop = math.log(min_sigma, 10), math.log(max_sigma, 10)
  611. sigma_list = np.logspace(start, stop, num_sigma)
  612. else:
  613. sigma_list = np.linspace(min_sigma, max_sigma, num_sigma)
  614. image_cube = np.empty(shape=image.shape + (len(sigma_list),), dtype=float_dtype)
  615. for j, s in enumerate(sigma_list):
  616. image_cube[..., j] = _hessian_matrix_det(image, s)
  617. local_maxima = peak_local_max(
  618. image_cube,
  619. threshold_abs=threshold,
  620. threshold_rel=threshold_rel,
  621. exclude_border=False,
  622. footprint=np.ones((3,) * image_cube.ndim),
  623. )
  624. # Catch no peaks
  625. if local_maxima.size == 0:
  626. return np.empty((0, 3))
  627. # Convert local_maxima to float64
  628. lm = local_maxima.astype(np.float64)
  629. # Convert the last index to its corresponding scale value
  630. lm[:, -1] = sigma_list[local_maxima[:, -1]]
  631. return _prune_blobs(lm, overlap)