_marching_cubes_lewiner.py 13 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352
  1. import base64
  2. import numpy as np
  3. from . import _marching_cubes_lewiner_luts as mcluts
  4. from . import _marching_cubes_lewiner_cy
  5. def marching_cubes(
  6. volume,
  7. level=None,
  8. *,
  9. spacing=(1.0, 1.0, 1.0),
  10. gradient_direction='descent',
  11. step_size=1,
  12. allow_degenerate=True,
  13. method='lewiner',
  14. mask=None,
  15. ):
  16. """Marching cubes algorithm to find surfaces in 3d volumetric data.
  17. In contrast with Lorensen et al. approach [2]_, Lewiner et
  18. al. algorithm is faster, resolves ambiguities, and guarantees
  19. topologically correct results. Therefore, this algorithm generally
  20. a better choice.
  21. Parameters
  22. ----------
  23. volume : (M, N, P) ndarray
  24. Input data volume to find isosurfaces. Will internally be
  25. converted to float32 if necessary.
  26. level : float, optional
  27. Contour value to search for isosurfaces in `volume`. If not
  28. given or None, the average of the min and max of vol is used.
  29. spacing : length-3 tuple of floats, optional
  30. Voxel spacing in spatial dimensions corresponding to numpy array
  31. indexing dimensions (M, N, P) as in `volume`.
  32. gradient_direction : {'descent', 'ascent'}, optional
  33. Controls if the mesh was generated from an isosurface with gradient
  34. descent toward objects of interest (the default), or the opposite,
  35. considering the *left-hand* rule.
  36. The two options are:
  37. * descent : Object was greater than exterior
  38. * ascent : Exterior was greater than object
  39. step_size : int, optional
  40. Step size in voxels. Default 1. Larger steps yield faster but
  41. coarser results. The result will always be topologically correct
  42. though.
  43. allow_degenerate : bool, optional
  44. Whether to allow degenerate (i.e. zero-area) triangles in the
  45. end-result. Default True. If False, degenerate triangles are
  46. removed, at the cost of making the algorithm slower.
  47. method : {'lewiner', 'lorensen'}, optional
  48. Whether the method of Lewiner et al. or Lorensen et al. will be used.
  49. mask : (M, N, P) array, optional
  50. Boolean array. The marching cube algorithm will be computed only on
  51. True elements. This will save computational time when interfaces
  52. are located within certain region of the volume M, N, P-e.g. the top
  53. half of the cube-and also allow to compute finite surfaces-i.e. open
  54. surfaces that do not end at the border of the cube.
  55. Returns
  56. -------
  57. verts : (V, 3) array
  58. Spatial coordinates for V unique mesh vertices. Coordinate order
  59. matches input `volume` (M, N, P). If ``allow_degenerate`` is set to
  60. True, then the presence of degenerate triangles in the mesh can make
  61. this array have duplicate vertices.
  62. faces : (F, 3) array
  63. Define triangular faces via referencing vertex indices from ``verts``.
  64. This algorithm specifically outputs triangles, so each face has
  65. exactly three indices.
  66. normals : (V, 3) array
  67. The normal direction at each vertex, as calculated from the
  68. data.
  69. values : (V,) array
  70. Gives a measure for the maximum value of the data in the local region
  71. near each vertex. This can be used by visualization tools to apply
  72. a colormap to the mesh.
  73. See Also
  74. --------
  75. skimage.measure.mesh_surface_area
  76. skimage.measure.find_contours
  77. Notes
  78. -----
  79. The algorithm [1]_ is an improved version of Chernyaev's Marching
  80. Cubes 33 algorithm. It is an efficient algorithm that relies on
  81. heavy use of lookup tables to handle the many different cases,
  82. keeping the algorithm relatively easy. This implementation is
  83. written in Cython, ported from Lewiner's C++ implementation.
  84. To quantify the area of an isosurface generated by this algorithm, pass
  85. verts and faces to `skimage.measure.mesh_surface_area`.
  86. Regarding visualization of algorithm output, to contour a volume
  87. named `myvolume` about the level 0.0, using the ``mayavi`` package::
  88. >>>
  89. >> from mayavi import mlab
  90. >> verts, faces, _, _ = marching_cubes(myvolume, 0.0)
  91. >> mlab.triangular_mesh([vert[0] for vert in verts],
  92. [vert[1] for vert in verts],
  93. [vert[2] for vert in verts],
  94. faces)
  95. >> mlab.show()
  96. Similarly using the ``visvis`` package::
  97. >>>
  98. >> import visvis as vv
  99. >> verts, faces, normals, values = marching_cubes(myvolume, 0.0)
  100. >> vv.mesh(np.fliplr(verts), faces, normals, values)
  101. >> vv.use().Run()
  102. To reduce the number of triangles in the mesh for better performance,
  103. see this `example
  104. <https://docs.enthought.com/mayavi/mayavi/auto/example_julia_set_decimation.html#example-julia-set-decimation>`_
  105. using the ``mayavi`` package.
  106. References
  107. ----------
  108. .. [1] Thomas Lewiner, Helio Lopes, Antonio Wilson Vieira and Geovan
  109. Tavares. Efficient implementation of Marching Cubes' cases with
  110. topological guarantees. Journal of Graphics Tools 8(2)
  111. pp. 1-15 (december 2003).
  112. :DOI:`10.1080/10867651.2003.10487582`
  113. .. [2] Lorensen, William and Harvey E. Cline. Marching Cubes: A High
  114. Resolution 3D Surface Construction Algorithm. Computer Graphics
  115. (SIGGRAPH 87 Proceedings) 21(4) July 1987, p. 163-170).
  116. :DOI:`10.1145/37401.37422`
  117. """
  118. use_classic = False
  119. if method == 'lorensen':
  120. use_classic = True
  121. elif method != 'lewiner':
  122. raise ValueError("method should be either 'lewiner' or 'lorensen'")
  123. return _marching_cubes_lewiner(
  124. volume,
  125. level,
  126. spacing,
  127. gradient_direction,
  128. step_size,
  129. allow_degenerate,
  130. use_classic=use_classic,
  131. mask=mask,
  132. )
  133. def _marching_cubes_lewiner(
  134. volume,
  135. level,
  136. spacing,
  137. gradient_direction,
  138. step_size,
  139. allow_degenerate,
  140. use_classic,
  141. mask,
  142. ):
  143. """Lewiner et al. algorithm for marching cubes. See
  144. marching_cubes_lewiner for documentation.
  145. """
  146. # Check volume and ensure its in the format that the alg needs
  147. if not isinstance(volume, np.ndarray) or (volume.ndim != 3):
  148. raise ValueError('Input volume should be a 3D numpy array.')
  149. if volume.shape[0] < 2 or volume.shape[1] < 2 or volume.shape[2] < 2:
  150. raise ValueError("Input array must be at least 2x2x2.")
  151. volume = np.ascontiguousarray(volume, np.float32) # no copy if not necessary
  152. # Check/convert other inputs:
  153. # level
  154. if level is None:
  155. level = 0.5 * (volume.min() + volume.max())
  156. else:
  157. level = float(level)
  158. if level < volume.min() or level > volume.max():
  159. raise ValueError("Surface level must be within volume data range.")
  160. # spacing
  161. if len(spacing) != 3:
  162. raise ValueError("`spacing` must consist of three floats.")
  163. # step_size
  164. step_size = int(step_size)
  165. if step_size < 1:
  166. raise ValueError('step_size must be at least one.')
  167. # use_classic
  168. use_classic = bool(use_classic)
  169. # Get LutProvider class (reuse if possible)
  170. L = _get_mc_luts()
  171. # Check if a mask array is passed
  172. if mask is not None:
  173. if not mask.shape == volume.shape:
  174. raise ValueError('volume and mask must have the same shape.')
  175. # Apply algorithm
  176. func = _marching_cubes_lewiner_cy.marching_cubes
  177. vertices, faces, normals, values = func(
  178. volume, level, L, step_size, use_classic, mask
  179. )
  180. if not len(vertices):
  181. raise RuntimeError('No surface found at the given iso value.')
  182. # Output in z-y-x order, as is common in skimage
  183. vertices = np.fliplr(vertices)
  184. normals = np.fliplr(normals)
  185. # Finishing touches to output
  186. faces.shape = -1, 3
  187. if gradient_direction == 'descent':
  188. # MC implementation is right-handed, but gradient_direction is
  189. # left-handed
  190. faces = np.fliplr(faces)
  191. elif not gradient_direction == 'ascent':
  192. raise ValueError(
  193. f"Incorrect input {gradient_direction} in `gradient_direction`, "
  194. "see docstring."
  195. )
  196. if not np.array_equal(spacing, (1, 1, 1)):
  197. vertices = vertices * np.r_[spacing]
  198. if allow_degenerate:
  199. return vertices, faces, normals, values
  200. else:
  201. fun = _marching_cubes_lewiner_cy.remove_degenerate_faces
  202. return fun(vertices.astype(np.float32), faces, normals, values)
  203. def _to_array(args):
  204. shape, text = args
  205. byts = base64.decodebytes(text.encode('utf-8'))
  206. ar = np.frombuffer(byts, dtype='int8')
  207. ar.shape = shape
  208. return ar
  209. # Map an edge-index to two relative pixel positions. The edge index
  210. # represents a point that lies somewhere in between these pixels.
  211. # Linear interpolation should be used to determine where it is exactly.
  212. # 0
  213. # 3 1 -> 0x
  214. # 2 xx
  215. # fmt: off
  216. EDGETORELATIVEPOSX = np.array([ [0,1],[1,1],[1,0],[0,0], [0,1],[1,1],[1,0],[0,0], [0,0],[1,1],[1,1],[0,0] ], 'int8')
  217. EDGETORELATIVEPOSY = np.array([ [0,0],[0,1],[1,1],[1,0], [0,0],[0,1],[1,1],[1,0], [0,0],[0,0],[1,1],[1,1] ], 'int8')
  218. EDGETORELATIVEPOSZ = np.array([ [0,0],[0,0],[0,0],[0,0], [1,1],[1,1],[1,1],[1,1], [0,1],[0,1],[0,1],[0,1] ], 'int8')
  219. # fmt: on
  220. def _get_mc_luts():
  221. """Kind of lazy obtaining of the luts."""
  222. if not hasattr(mcluts, 'THE_LUTS'):
  223. mcluts.THE_LUTS = _marching_cubes_lewiner_cy.LutProvider(
  224. EDGETORELATIVEPOSX,
  225. EDGETORELATIVEPOSY,
  226. EDGETORELATIVEPOSZ,
  227. _to_array(mcluts.CASESCLASSIC),
  228. _to_array(mcluts.CASES),
  229. _to_array(mcluts.TILING1),
  230. _to_array(mcluts.TILING2),
  231. _to_array(mcluts.TILING3_1),
  232. _to_array(mcluts.TILING3_2),
  233. _to_array(mcluts.TILING4_1),
  234. _to_array(mcluts.TILING4_2),
  235. _to_array(mcluts.TILING5),
  236. _to_array(mcluts.TILING6_1_1),
  237. _to_array(mcluts.TILING6_1_2),
  238. _to_array(mcluts.TILING6_2),
  239. _to_array(mcluts.TILING7_1),
  240. _to_array(mcluts.TILING7_2),
  241. _to_array(mcluts.TILING7_3),
  242. _to_array(mcluts.TILING7_4_1),
  243. _to_array(mcluts.TILING7_4_2),
  244. _to_array(mcluts.TILING8),
  245. _to_array(mcluts.TILING9),
  246. _to_array(mcluts.TILING10_1_1),
  247. _to_array(mcluts.TILING10_1_1_),
  248. _to_array(mcluts.TILING10_1_2),
  249. _to_array(mcluts.TILING10_2),
  250. _to_array(mcluts.TILING10_2_),
  251. _to_array(mcluts.TILING11),
  252. _to_array(mcluts.TILING12_1_1),
  253. _to_array(mcluts.TILING12_1_1_),
  254. _to_array(mcluts.TILING12_1_2),
  255. _to_array(mcluts.TILING12_2),
  256. _to_array(mcluts.TILING12_2_),
  257. _to_array(mcluts.TILING13_1),
  258. _to_array(mcluts.TILING13_1_),
  259. _to_array(mcluts.TILING13_2),
  260. _to_array(mcluts.TILING13_2_),
  261. _to_array(mcluts.TILING13_3),
  262. _to_array(mcluts.TILING13_3_),
  263. _to_array(mcluts.TILING13_4),
  264. _to_array(mcluts.TILING13_5_1),
  265. _to_array(mcluts.TILING13_5_2),
  266. _to_array(mcluts.TILING14),
  267. _to_array(mcluts.TEST3),
  268. _to_array(mcluts.TEST4),
  269. _to_array(mcluts.TEST6),
  270. _to_array(mcluts.TEST7),
  271. _to_array(mcluts.TEST10),
  272. _to_array(mcluts.TEST12),
  273. _to_array(mcluts.TEST13),
  274. _to_array(mcluts.SUBCONFIG13),
  275. )
  276. return mcluts.THE_LUTS
  277. def mesh_surface_area(verts, faces):
  278. """Compute surface area, given vertices and triangular faces.
  279. Parameters
  280. ----------
  281. verts : (V, 3) array of floats
  282. Array containing coordinates for V unique mesh vertices.
  283. faces : (F, 3) array of ints
  284. List of length-3 lists of integers, referencing vertex coordinates as
  285. provided in `verts`.
  286. Returns
  287. -------
  288. area : float
  289. Surface area of mesh. Units now [coordinate units] ** 2.
  290. Notes
  291. -----
  292. The arguments expected by this function are the first two outputs from
  293. `skimage.measure.marching_cubes`. For unit correct output, ensure correct
  294. `spacing` was passed to `skimage.measure.marching_cubes`.
  295. This algorithm works properly only if the ``faces`` provided are all
  296. triangles.
  297. See Also
  298. --------
  299. skimage.measure.marching_cubes
  300. """
  301. # Fancy indexing to define two vector arrays from triangle vertices
  302. actual_verts = verts[faces]
  303. a = actual_verts[:, 0, :] - actual_verts[:, 1, :]
  304. b = actual_verts[:, 0, :] - actual_verts[:, 2, :]
  305. del actual_verts
  306. # Area of triangle in 3D = 1/2 * Euclidean norm of cross product
  307. return ((np.cross(a, b) ** 2).sum(axis=1) ** 0.5).sum() / 2.0