marching_cubes.cpp 12 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247
  1. /*===========================================================================*\
  2. * *
  3. * IsoEx *
  4. * Copyright (C) 2002 by Computer Graphics Group, RWTH Aachen *
  5. * www.rwth-graphics.de *
  6. * *
  7. *---------------------------------------------------------------------------*
  8. * *
  9. * License *
  10. * *
  11. * This library is free software; you can redistribute it and/or modify it *
  12. * under the terms of the GNU Library General Public License as published *
  13. * by the Free Software Foundation, version 2. *
  14. * *
  15. * This library is distributed in the hope that it will be useful, but *
  16. * WITHOUT ANY WARRANTY; without even the implied warranty of *
  17. * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU *
  18. * Library General Public License for more details. *
  19. * *
  20. * You should have received a copy of the GNU Library General Public *
  21. * License along with this library; if not, write to the Free Software *
  22. * Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA. *
  23. * *
  24. \*===========================================================================*/
  25. #include "marching_cubes.h"
  26. #include "marching_cubes_tables.h"
  27. #include <map>
  28. extern const int edgeTable[256];
  29. extern const int triTable[256][2][17];
  30. extern const int polyTable[8][16];
  31. class EdgeKey
  32. {
  33. public:
  34. EdgeKey(unsigned i0, unsigned i1) {
  35. if (i0 < i1) { i0_ = i0; i1_ = i1; }
  36. else { i0_ = i1; i1_ = i0; }
  37. }
  38. bool operator<(const EdgeKey& _rhs) const
  39. {
  40. if (i0_ != _rhs.i0_)
  41. return (i0_ < _rhs.i0_);
  42. else
  43. return (i1_ < _rhs.i1_);
  44. }
  45. private:
  46. unsigned i0_, i1_;
  47. };
  48. template <typename DerivedV, typename DerivedF>
  49. class MarchingCubes
  50. {
  51. typedef Eigen::PlainObjectBase<DerivedV> PointMatrixType;
  52. typedef Eigen::PlainObjectBase<DerivedF> FaceMatrixType;
  53. typedef std::map<EdgeKey, unsigned> MyMap;
  54. typedef typename MyMap::const_iterator MyMapIterator;
  55. public:
  56. MarchingCubes(const Eigen::Matrix<typename DerivedV::Scalar, Eigen::Dynamic, 1> &values,
  57. const Eigen::Matrix<typename DerivedV::Scalar, Eigen::Dynamic, 3> &points,
  58. const unsigned x_res,
  59. const unsigned y_res,
  60. const unsigned z_res,
  61. PointMatrixType &vertices,
  62. FaceMatrixType &faces)
  63. {
  64. if(x_res <2 || y_res<2 ||z_res<2)
  65. return;
  66. faces.resize(10000,3);
  67. int num_faces = 0;
  68. vertices.resize(10000,3);
  69. int num_vertices = 0;
  70. unsigned n_cubes = (x_res-1) * (y_res-1) * (z_res-1);
  71. assert(unsigned(points.rows()) == x_res * y_res * z_res);
  72. unsigned int offsets_[8];
  73. offsets_[0] = 0;
  74. offsets_[1] = 1;
  75. offsets_[2] = 1 + x_res;
  76. offsets_[3] = x_res;
  77. offsets_[4] = x_res*y_res;
  78. offsets_[5] = 1 + x_res*y_res;
  79. offsets_[6] = 1 + x_res + x_res*y_res;
  80. offsets_[7] = x_res + x_res*y_res;
  81. for (unsigned cube_it =0 ; cube_it < n_cubes; ++cube_it)
  82. {
  83. unsigned corner[8];
  84. typename DerivedF::Scalar samples[12];
  85. unsigned char cubetype(0);
  86. unsigned int i;
  87. // get point indices of corner vertices
  88. for (i=0; i<8; ++i)
  89. {
  90. // get cube coordinates
  91. unsigned int _idx = cube_it;
  92. unsigned int X(x_res-1), Y(y_res-1);
  93. unsigned int x = _idx % X; _idx /= X;
  94. unsigned int y = _idx % Y; _idx /= Y;
  95. unsigned int z = _idx;
  96. // transform to point coordinates
  97. _idx = x + y*x_res + z*x_res*y_res;
  98. // add offset
  99. corner[i] = _idx + offsets_[i];
  100. }
  101. // determine cube type
  102. for (i=0; i<8; ++i)
  103. if (values[corner[i]] > 0.0)
  104. cubetype |= (1<<i);
  105. // trivial reject ?
  106. if (cubetype == 0 || cubetype == 255)
  107. continue;
  108. // compute samples on cube's edges
  109. if (edgeTable[cubetype]&1)
  110. samples[0] = add_vertex(values, points, corner[0], corner[1], vertices, num_vertices, edge2vertex);
  111. if (edgeTable[cubetype]&2)
  112. samples[1] = add_vertex(values, points, corner[1], corner[2], vertices, num_vertices, edge2vertex);
  113. if (edgeTable[cubetype]&4)
  114. samples[2] = add_vertex(values, points, corner[3], corner[2], vertices, num_vertices, edge2vertex);
  115. if (edgeTable[cubetype]&8)
  116. samples[3] = add_vertex(values, points, corner[0], corner[3], vertices, num_vertices, edge2vertex);
  117. if (edgeTable[cubetype]&16)
  118. samples[4] = add_vertex(values, points, corner[4], corner[5], vertices, num_vertices, edge2vertex);
  119. if (edgeTable[cubetype]&32)
  120. samples[5] = add_vertex(values, points, corner[5], corner[6], vertices, num_vertices, edge2vertex);
  121. if (edgeTable[cubetype]&64)
  122. samples[6] = add_vertex(values, points, corner[7], corner[6], vertices, num_vertices, edge2vertex);
  123. if (edgeTable[cubetype]&128)
  124. samples[7] = add_vertex(values, points, corner[4], corner[7], vertices, num_vertices, edge2vertex);
  125. if (edgeTable[cubetype]&256)
  126. samples[8] = add_vertex(values, points, corner[0], corner[4], vertices, num_vertices, edge2vertex);
  127. if (edgeTable[cubetype]&512)
  128. samples[9] = add_vertex(values, points, corner[1], corner[5], vertices, num_vertices, edge2vertex);
  129. if (edgeTable[cubetype]&1024)
  130. samples[10] = add_vertex(values, points, corner[2], corner[6], vertices, num_vertices, edge2vertex);
  131. if (edgeTable[cubetype]&2048)
  132. samples[11] = add_vertex(values, points, corner[3], corner[7], vertices, num_vertices, edge2vertex);
  133. // connect samples by triangles
  134. for (i=0; triTable[cubetype][0][i] != -1; i+=3 )
  135. {
  136. num_faces++;
  137. if (num_faces > faces.rows())
  138. faces.conservativeResize(faces.rows()+10000, Eigen::NoChange);
  139. faces.row(num_faces-1) <<
  140. samples[triTable[cubetype][0][i ]],
  141. samples[triTable[cubetype][0][i+1]],
  142. samples[triTable[cubetype][0][i+2]];
  143. }
  144. }
  145. vertices.conservativeResize(num_vertices, Eigen::NoChange);
  146. faces.conservativeResize(num_faces, Eigen::NoChange);
  147. };
  148. static typename DerivedF::Scalar add_vertex(const Eigen::Matrix<typename DerivedV::Scalar, Eigen::Dynamic, 1> &values,
  149. const Eigen::Matrix<typename DerivedV::Scalar, Eigen::Dynamic, 3> &points,
  150. unsigned int i0,
  151. unsigned int i1,
  152. PointMatrixType &vertices,
  153. int &num_vertices,
  154. MyMap &edge2vertex)
  155. {
  156. // find vertex if it has been computed already
  157. MyMapIterator it = edge2vertex.find(EdgeKey(i0, i1));
  158. if (it != edge2vertex.end())
  159. return it->second;
  160. ;
  161. // generate new vertex
  162. const Eigen::Matrix<typename DerivedV::Scalar, 1, 3> & p0 = points.row(i0);
  163. const Eigen::Matrix<typename DerivedV::Scalar, 1, 3> & p1 = points.row(i1);
  164. typename DerivedV::Scalar s0 = fabs(values[i0]);
  165. typename DerivedV::Scalar s1 = fabs(values[i1]);
  166. typename DerivedV::Scalar t = s0 / (s0+s1);
  167. num_vertices++;
  168. if (num_vertices > vertices.rows())
  169. vertices.conservativeResize(vertices.rows()+10000, Eigen::NoChange);
  170. vertices.row(num_vertices-1) = (1.0f-t)*p0 + t*p1;
  171. edge2vertex[EdgeKey(i0, i1)] = num_vertices-1;
  172. return num_vertices-1;
  173. }
  174. ;
  175. // maps an edge to the sample vertex generated on it
  176. MyMap edge2vertex;
  177. };
  178. template <typename DerivedV, typename DerivedF>
  179. IGL_INLINE void igl::copyleft::marching_cubes(
  180. const Eigen::Matrix<typename DerivedV::Scalar, Eigen::Dynamic, 1> &values,
  181. const Eigen::Matrix<typename DerivedV::Scalar, Eigen::Dynamic, 3> &points,
  182. const unsigned x_res,
  183. const unsigned y_res,
  184. const unsigned z_res,
  185. Eigen::PlainObjectBase<DerivedV> &vertices,
  186. Eigen::PlainObjectBase<DerivedF> &faces)
  187. {
  188. MarchingCubes<DerivedV, DerivedF> mc(values,
  189. points,
  190. x_res,
  191. y_res,
  192. z_res,
  193. vertices,
  194. faces);
  195. }
  196. #ifdef IGL_STATIC_LIBRARY
  197. // Explicit template specialization
  198. template void igl::copyleft::marching_cubes<Eigen::Matrix<float, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3> >(Eigen::Matrix<Eigen::Matrix<float, -1, 3, 0, -1, 3>::Scalar, -1, 1, 0, -1, 1> const&, Eigen::Matrix<Eigen::Matrix<float, -1, 3, 0, -1, 3>::Scalar, -1, 3, 0, -1, 3> const&, unsigned int, unsigned int, unsigned int, Eigen::PlainObjectBase<Eigen::Matrix<float, -1, 3, 0, -1, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> >&);
  199. template void igl::copyleft::marching_cubes<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3> >(Eigen::Matrix<Eigen::Matrix<double, -1, 3, 0, -1, 3>::Scalar, -1, 1, 0, -1, 1> const&, Eigen::Matrix<Eigen::Matrix<double, -1, 3, 0, -1, 3>::Scalar, -1, 3, 0, -1, 3> const&, unsigned int, unsigned int, unsigned int, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> >&);
  200. template void igl::copyleft::marching_cubes<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<unsigned int, -1, 3, 0, -1, 3> >(Eigen::Matrix<Eigen::Matrix<double, -1, 3, 0, -1, 3>::Scalar, -1, 1, 0, -1, 1> const&, Eigen::Matrix<Eigen::Matrix<double, -1, 3, 0, -1, 3>::Scalar, -1, 3, 0, -1, 3> const&, unsigned int, unsigned int, unsigned int, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<unsigned int, -1, 3, 0, -1, 3> >&);
  201. template void igl::copyleft::marching_cubes<Eigen::Matrix<float, -1, 3, 0, -1, 3>, Eigen::Matrix<unsigned int, -1, 3, 0, -1, 3> >(Eigen::Matrix<Eigen::Matrix<float, -1, 3, 0, -1, 3>::Scalar, -1, 1, 0, -1, 1> const&, Eigen::Matrix<Eigen::Matrix<float, -1, 3, 0, -1, 3>::Scalar, -1, 3, 0, -1, 3> const&, unsigned int, unsigned int, unsigned int, Eigen::PlainObjectBase<Eigen::Matrix<float, -1, 3, 0, -1, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<unsigned int, -1, 3, 0, -1, 3> >&);
  202. template void igl::copyleft::marching_cubes<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1> >(Eigen::Matrix<Eigen::Matrix<double, -1, -1, 0, -1, -1>::Scalar, -1, 1, 0, -1, 1> const&, Eigen::Matrix<Eigen::Matrix<double, -1, -1, 0, -1, -1>::Scalar, -1, 3, 0, -1, 3> const&, unsigned int, unsigned int, unsigned int, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&);
  203. #endif