contours.cpp 6.0 KB

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