main.cpp 3.4 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126
  1. #include <igl/marching_cubes.h>
  2. #include <igl/signed_distance.h>
  3. #include <igl/read_triangle_mesh.h>
  4. #include <igl/voxel_grid.h>
  5. #include <igl/opengl/glfw/Viewer.h>
  6. #include <Eigen/Core>
  7. #include <iostream>
  8. int main(int argc, char * argv[])
  9. {
  10. using namespace Eigen;
  11. using namespace std;
  12. using namespace igl;
  13. MatrixXi F;
  14. MatrixXd V;
  15. // Read in inputs as double precision floating point meshes
  16. read_triangle_mesh(
  17. TUTORIAL_SHARED_PATH "/armadillo.obj",V,F);
  18. cout<<"Creating grid..."<<endl;
  19. // number of vertices on the largest side
  20. const int s = 100;
  21. // create grid
  22. MatrixXd GV;
  23. Eigen::RowVector3i res;
  24. igl::voxel_grid(V,0,s,1,GV,res);
  25. // compute values
  26. cout<<"Computing distances..."<<endl;
  27. // Batch based function for evaluating implicit
  28. const auto batch_implicit = [&V,&F](const MatrixXd & Q)->VectorXd
  29. {
  30. VectorXd S;
  31. {
  32. VectorXi I;
  33. MatrixXd C,N;
  34. signed_distance(Q,V,F,SIGNED_DISTANCE_TYPE_FAST_WINDING_NUMBER,S,I,C,N);
  35. // Extremely flatten out near zero
  36. S = S.array().sign() * S.array().abs().exp();
  37. }
  38. return S;
  39. };
  40. VectorXd S = batch_implicit(GV);
  41. cout<<"Marching cubes..."<<endl;
  42. MatrixXd SV;
  43. MatrixXi SF;
  44. igl::marching_cubes(S,GV,res(0),res(1),res(2),0,SV,SF);
  45. std::unordered_map<std::int64_t,int> E2V;
  46. igl::marching_cubes(S,GV,res(0),res(1),res(2),0,SV,SF,E2V);
  47. // Initialize min and max for root finding bisection
  48. assert(E2V.size() == SV.rows());
  49. Eigen::MatrixXd T(SV.rows(),2);
  50. Eigen::VectorXi I(SV.rows());
  51. // Precompute slices for end points
  52. Eigen::MatrixXd GVi(SV.rows(),3);
  53. Eigen::MatrixXd GVj(SV.rows(),3);
  54. // This is only used for the assertion below
  55. const auto ij2key = [](std::int32_t i,std::int32_t j)
  56. {
  57. if(i>j){ std::swap(i,j); }
  58. std::int64_t ret = 0;
  59. ret |= i;
  60. ret |= static_cast<std::int64_t>(j) << 32;
  61. return ret;
  62. };
  63. const auto key2ij = [](const std::int64_t & key, std::int32_t & i, std::int32_t & j)
  64. {
  65. i = key & 0xFFFFFFFF;
  66. j = key >> 32;
  67. };
  68. for (const auto& e2v: E2V)
  69. {
  70. const std::int64_t key = e2v.first;
  71. const std::int32_t v = e2v.second;
  72. std::int32_t i,j;
  73. key2ij(key,i,j);
  74. const std::int64_t key0 = ij2key(i,j);
  75. assert(key0 != key);
  76. T.row(v) << 0,1;
  77. // (i,j) is ordered so that i<j, but let's order so that S(i)<S(j)
  78. if(S(i)>S(j)) { std::swap(i,j); }
  79. GVi.row(v) = GV.row(i);
  80. GVj.row(v) = GV.row(j);
  81. }
  82. const auto root_find_iteration = [&SV,&I,&T,&GVi,&GVj,&batch_implicit]()
  83. {
  84. Eigen::VectorXd T_mid = (T.col(0)+T.col(1))/2;
  85. // Use this midpoint guess for the current visualization
  86. SV = GVi.array().colwise() * (1.0 - T_mid.array()) + GVj.array().colwise() * T_mid.array();
  87. // Compute values at midpoints
  88. VectorXd S_mid = batch_implicit(SV);
  89. // Update bounds
  90. T.col(1) = (S_mid.array() > 0).select(T_mid, T.col(1));
  91. T.col(0) = (S_mid.array() <= 0).select(T_mid, T.col(0));
  92. };
  93. cout<<R"(Usage:
  94. ' ' Conduct a bisection iteration
  95. )";
  96. igl::opengl::glfw::Viewer viewer;
  97. viewer.data().set_mesh(SV,SF);
  98. viewer.data().show_lines = false;
  99. viewer.callback_key_down =
  100. [&](igl::opengl::glfw::Viewer & viewer, unsigned char key, int mod)->bool
  101. {
  102. switch(key)
  103. {
  104. default:
  105. return false;
  106. case ' ':
  107. root_find_iteration();
  108. viewer.data().set_vertices(SV);
  109. viewer.data().compute_normals();
  110. break;
  111. }
  112. return true;
  113. };
  114. viewer.launch();
  115. }