main.cpp 7.1 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242
  1. #include <igl/copyleft/marching_cubes.h>
  2. #include <igl/unique_simplices.h>
  3. #include <igl/dual_contouring.h>
  4. #include <igl/get_seconds.h>
  5. #include <igl/grid.h>
  6. #include <igl/marching_cubes.h>
  7. #include <igl/per_corner_normals.h>
  8. #include <igl/parallel_for.h>
  9. #include <igl/per_corner_normals.h>
  10. #include <igl/per_face_normals.h>
  11. #include <igl/polygon_corners.h>
  12. #include <igl/slice.h>
  13. #include <igl/sparse_voxel_grid.h>
  14. #include <igl/opengl/glfw/Viewer.h>
  15. int main(int argc, char * argv[])
  16. {
  17. const auto & tictoc = []()
  18. {
  19. static double t_start = igl::get_seconds();
  20. double diff = igl::get_seconds()-t_start;
  21. t_start += diff;
  22. return diff;
  23. };
  24. // Create an interesting shape with sharp features using SDF CSG with spheres.
  25. const auto & sphere = [](
  26. const Eigen::RowVector3d & c,
  27. const double r,
  28. const Eigen::RowVector3d & x)->double
  29. {
  30. return (x-c).norm() - r;
  31. };
  32. const std::function<double(const Eigen::RowVector3d & x)> f =
  33. [&](const Eigen::RowVector3d & x)->double
  34. {
  35. return
  36. std::min(
  37. std::min(std::max(
  38. sphere(Eigen::RowVector3d(-0.2,0,-0.2),0.5,x),
  39. -sphere(Eigen::RowVector3d(+0.2,0,0.2),0.5,x)),
  40. sphere(Eigen::RowVector3d(-0.15,0,-0.15),0.3,x)
  41. ),
  42. std::max(
  43. std::max(
  44. sphere(Eigen::RowVector3d(-0.2,-0.5,-0.2),0.6,x),x(1)+0.45),-0.6-x(1))
  45. );
  46. };
  47. Eigen::RowVector3d p0(-0.2,0.5,-0.2);
  48. assert(abs(f(p0)) < 1e-10 && "p0 should be on zero level-set");
  49. // Simple finite difference gradients
  50. const auto & fd = [](
  51. const std::function<double(const Eigen::RowVector3d&)> & f,
  52. const Eigen::RowVector3d & x)
  53. {
  54. const double eps = 1e-10;
  55. Eigen::RowVector3d g;
  56. for(int c = 0;c<3;c++)
  57. {
  58. const Eigen::RowVector3d xp = x+eps*Eigen::RowVector3d(c==0,c==1,c==2);
  59. const double fp = f(xp);
  60. const Eigen::RowVector3d xn = x-eps*Eigen::RowVector3d(c==0,c==1,c==2);
  61. const double fn = f(xn);
  62. g(c) = (fp-fn)/(2*eps);
  63. }
  64. return g;
  65. };
  66. const auto & f_grad = [&fd,&f](const Eigen::RowVector3d & x)
  67. {
  68. return fd(f,x).normalized();
  69. };
  70. Eigen::MatrixXd V;
  71. Eigen::MatrixXi Q,F;
  72. Eigen::MatrixXd mcV,mcN;
  73. Eigen::MatrixXi mcF;
  74. // Grid parameters
  75. const Eigen::RowVector3d min_corner(-2,-2,-2);
  76. const Eigen::RowVector3d max_corner(+2,+2,+2);
  77. const int s = 256;
  78. int nx = s+1;
  79. int ny = s+1;
  80. int nz = s+1;
  81. const Eigen::RowVector3d step =
  82. (max_corner-min_corner).array()/(Eigen::RowVector3d(nx,ny,nz).array()-1);
  83. // Sparse grid below assumes regular grid
  84. assert((step(0) == step(1))&&(step(0) == step(2)));
  85. // Dual contouring parameters
  86. bool constrained = false;
  87. bool triangles = false;
  88. bool root_finding = true;
  89. for(int pass = 0;pass<2;pass++)
  90. {
  91. const bool sparse = pass == 1;
  92. printf("Using %s grid..\n",sparse?"sparse":"dense");
  93. if(sparse)
  94. {
  95. // igl::sparse_voxel_grid assumes (0,0,0) lies on the grid. But dense igl::grid
  96. // below won't necessarily do that depending on nx,ny,nz.
  97. tictoc();
  98. Eigen::MatrixXd GV;
  99. Eigen::VectorXd Gf;
  100. Eigen::Matrix<int,Eigen::Dynamic,8> GI;
  101. igl::sparse_voxel_grid(p0,f,step(0),16.*pow(step(0),-2.),Gf,GV,GI);
  102. const auto t_Gf = tictoc();
  103. printf(" %5f secs to populate sparse grid of %ld cells\n",t_Gf+tictoc(),GI.rows());
  104. // Dual contouring requires list of sparse edges (not cells)
  105. // extract _all_ edges from sparse_voxel_grid (conservative)
  106. Eigen::Matrix<int,Eigen::Dynamic,2> GI2;
  107. {
  108. Eigen::Matrix<int,Eigen::Dynamic,2> all_GI2(GI.rows()*12,2);
  109. all_GI2 <<
  110. // front
  111. GI.col(0),GI.col(1),
  112. GI.col(1),GI.col(2),
  113. GI.col(2),GI.col(3),
  114. GI.col(3),GI.col(0),
  115. // back
  116. GI.col(4+0),GI.col(4+1),
  117. GI.col(4+1),GI.col(4+2),
  118. GI.col(4+2),GI.col(4+3),
  119. GI.col(4+3),GI.col(4+0),
  120. // sides
  121. GI.col(0),GI.col(4+0),
  122. GI.col(1),GI.col(4+1),
  123. GI.col(2),GI.col(4+2),
  124. GI.col(3),GI.col(4+3);
  125. Eigen::VectorXi _1,_2;
  126. igl::unique_simplices(all_GI2,GI2,_1,_2);
  127. }
  128. tictoc();
  129. Eigen::RowVector3d step =
  130. (max_corner-min_corner).array()/(Eigen::RowVector3d(nx,ny,nz).array()-1);
  131. igl::dual_contouring(
  132. f,f_grad,step,Gf,GV,GI2,constrained,triangles,root_finding,V,Q);
  133. printf(" %5f secs dual contouring\n",t_Gf+tictoc());
  134. // Could use igl::marching_cubes once
  135. // https://github.com/libigl/libigl/pull/1687 is merged
  136. tictoc();
  137. igl::copyleft::marching_cubes(Gf,GV,GI,mcV, mcF);
  138. printf(" %5f secs marching cubes\n",t_Gf+tictoc());
  139. }else
  140. {
  141. tictoc();
  142. igl::dual_contouring(
  143. f,f_grad,min_corner,max_corner,nx,ny,nz,constrained,triangles,root_finding,V,Q);
  144. printf(" %5f secs dual contouring\n",tictoc());
  145. // build and sample grid
  146. tictoc();
  147. Eigen::MatrixXd GV;
  148. igl::grid(Eigen::RowVector3i(nx,ny,nz),GV);
  149. Eigen::VectorXd Gf(GV.rows());
  150. igl::parallel_for(GV.rows(),[&](const int i)
  151. {
  152. GV.row(i).array() *= (max_corner-min_corner).array();
  153. GV.row(i) += min_corner;
  154. Gf(i) = f(GV.row(i));
  155. },1000ul);
  156. const auto t_grid = tictoc();
  157. igl::marching_cubes(Gf,GV,nx,ny,nz,0,mcV,mcF);
  158. const auto t_mc = tictoc();
  159. printf(" %5f secs (%5f + %5f) marching cubes\n",t_grid+t_mc,t_grid,t_mc);
  160. }
  161. }
  162. // Crisp (as possible) rendering of resulting MC triangle mesh
  163. igl::per_corner_normals(mcV,mcF,20,mcN);
  164. // Crisp rendering of resulting DC quad mesh with edges
  165. Eigen::MatrixXi E;
  166. Eigen::MatrixXd VV,N,NN;
  167. Eigen::VectorXi J;
  168. Eigen::MatrixXi FF;
  169. if(triangles)
  170. {
  171. VV = V;
  172. FF = Q;
  173. E.resize(Q.rows()*3,2);
  174. E<<
  175. Q.col(0), Q.col(1),
  176. Q.col(1), Q.col(2),
  177. Q.col(2), Q.col(0);
  178. }else
  179. {
  180. Eigen::VectorXi I,C;
  181. igl::polygon_corners(Q,I,C);
  182. E.resize(Q.rows()*4,2);
  183. E<<
  184. Q.col(0), Q.col(1),
  185. Q.col(1), Q.col(2),
  186. Q.col(2), Q.col(3),
  187. Q.col(3), Q.col(0);
  188. igl::per_face_normals(V,I,C,N,VV,FF,J);
  189. igl::slice(N,J,1,NN);
  190. igl::per_corner_normals(V,I,C,20,N,VV,FF,J,NN);
  191. }
  192. igl::opengl::glfw::Viewer vr;
  193. bool show_edges = true;
  194. bool use_dc = true;
  195. const auto update = [&]()
  196. {
  197. const bool was_face_based = vr.data().face_based ;
  198. vr.data().clear();
  199. if(use_dc)
  200. {
  201. vr.data().set_mesh(VV,FF);
  202. vr.data().show_lines = false;
  203. vr.data().set_normals(NN);
  204. if(show_edges)
  205. {
  206. vr.data().clear_edges();
  207. vr.data().set_edges(V,E,Eigen::RowVector3d(0,0,0));
  208. }
  209. }else
  210. {
  211. vr.data().set_mesh(mcV,mcF);
  212. vr.data().set_normals(mcN);
  213. vr.data().show_lines = show_edges;
  214. }
  215. vr.data().face_based = was_face_based;
  216. };
  217. update();
  218. vr.data().face_based = true;
  219. vr.callback_key_pressed = [&](decltype(vr) &,unsigned int key, int mod)
  220. {
  221. switch(key)
  222. {
  223. case ' ': use_dc=!use_dc; update();return true;
  224. case 'L': case 'l': show_edges=!show_edges; update();return true;
  225. }
  226. return false;
  227. };
  228. std::cout<<R"(
  229. [space] Toggle between dual contouring and marching cubes
  230. )";
  231. vr.launch();
  232. }