sparse_voxel_grid.cpp 6.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148
  1. #include "sparse_voxel_grid.h"
  2. #include <unordered_map>
  3. #include <array>
  4. #include <vector>
  5. template <typename DerivedP0, typename Func, typename DerivedS, typename DerivedV, typename DerivedI>
  6. IGL_INLINE void igl::sparse_voxel_grid(const Eigen::MatrixBase<DerivedP0>& p0,
  7. const Func& scalarFunc,
  8. const double eps,
  9. const int expected_number_of_cubes,
  10. Eigen::PlainObjectBase<DerivedS>& CS,
  11. Eigen::PlainObjectBase<DerivedV>& CV,
  12. Eigen::PlainObjectBase<DerivedI>& CI)
  13. {
  14. typedef typename DerivedV::Scalar ScalarV;
  15. typedef typename DerivedS::Scalar ScalarS;
  16. typedef typename DerivedI::Scalar ScalarI;
  17. typedef Eigen::Matrix<ScalarV, 1, 3> VertexRowVector;
  18. typedef Eigen::Matrix<ScalarI, 1, 8> IndexRowVector;
  19. struct IndexRowVectorHash {
  20. std::size_t operator()(const Eigen::RowVector3i& key) const {
  21. std::size_t seed = 0;
  22. std::hash<int> hasher;
  23. for (int i = 0; i < 3; i++) {
  24. seed ^= hasher(key[i]) + 0x9e3779b9 + (seed<<6) + (seed>>2); // Copied from boost::hash_combine
  25. }
  26. return seed;
  27. }
  28. };
  29. auto sgn = [](ScalarS val) -> int {
  30. return (ScalarS(0) < val) - (val < ScalarS(0));
  31. };
  32. ScalarV half_eps = 0.5 * eps;
  33. std::vector<IndexRowVector> CI_vector;
  34. std::vector<VertexRowVector> CV_vector;
  35. std::vector<ScalarS> CS_vector;
  36. CI_vector.reserve(expected_number_of_cubes);
  37. CV_vector.reserve(8 * expected_number_of_cubes);
  38. CS_vector.reserve(8 * expected_number_of_cubes);
  39. // Track visisted neighbors
  40. std::unordered_map<Eigen::RowVector3i, int, IndexRowVectorHash> visited;
  41. visited.reserve(6 * expected_number_of_cubes);
  42. visited.max_load_factor(0.5);
  43. // BFS Queue
  44. std::vector<Eigen::RowVector3i> queue;
  45. queue.reserve(expected_number_of_cubes * 8);
  46. queue.push_back(Eigen::RowVector3i(0, 0, 0));
  47. while (queue.size() > 0)
  48. {
  49. Eigen::RowVector3i pi = queue.back();
  50. queue.pop_back();
  51. VertexRowVector ctr = p0 + eps*pi.cast<ScalarV>(); // R^3 center of this cube
  52. // X, Y, Z basis vectors, and array of neighbor offsets used to construct cubes
  53. const Eigen::RowVector3i bx(1, 0, 0), by(0, 1, 0), bz(0, 0, -1);
  54. const std::array<Eigen::RowVector3i, 6> neighbors = {
  55. bx, -bx, by, -by, bz, -bz
  56. };
  57. // Compute the position of the cube corners and the scalar values at those corners
  58. std::array<VertexRowVector, 8> cubeCorners = {
  59. ctr+half_eps*(bx+by+bz).cast<ScalarV>(), ctr+half_eps*(bx+by-bz).cast<ScalarV>(), ctr+half_eps*(-bx+by-bz).cast<ScalarV>(), ctr+half_eps*(-bx+by+bz).cast<ScalarV>(),
  60. ctr+half_eps*(bx-by+bz).cast<ScalarV>(), ctr+half_eps*(bx-by-bz).cast<ScalarV>(), ctr+half_eps*(-bx-by-bz).cast<ScalarV>(), ctr+half_eps*(-bx-by+bz).cast<ScalarV>()
  61. };
  62. std::array<ScalarS, 8> cubeScalars;
  63. for (int i = 0; i < 8; i++) { cubeScalars[i] = scalarFunc(cubeCorners[i]); }
  64. // If this cube doesn't intersect the surface, disregard it
  65. bool validCube = false;
  66. int sign = sgn(cubeScalars[0]);
  67. for (int i = 1; i < 8; i++) {
  68. if (sign != sgn(cubeScalars[i])) {
  69. validCube = true;
  70. break;
  71. }
  72. }
  73. if (!validCube) {
  74. continue;
  75. }
  76. // Add the cube vertices and indices to the output arrays if they are not there already
  77. IndexRowVector cube;
  78. uint8_t vertexAlreadyAdded = 0; // This is a bimask. If a bit is 1, it has been visited already by the BFS
  79. constexpr std::array<uint8_t, 6> zv = {
  80. (1 << 0) | (1 << 1) | (1 << 4) | (1 << 5),
  81. (1 << 2) | (1 << 3) | (1 << 6) | (1 << 7),
  82. (1 << 0) | (1 << 1) | (1 << 2) | (1 << 3),
  83. (1 << 4) | (1 << 5) | (1 << 6) | (1 << 7),
  84. (1 << 0) | (1 << 3) | (1 << 4) | (1 << 7),
  85. (1 << 1) | (1 << 2) | (1 << 5) | (1 << 6), };
  86. constexpr std::array<std::array<int, 4>, 6> zvv {{
  87. {{0, 1, 4, 5}}, {{3, 2, 7, 6}}, {{0, 1, 2, 3}},
  88. {{4, 5, 6, 7}}, {{0, 3, 4, 7}}, {{1, 2, 5, 6}} }};
  89. for (int n = 0; n < 6; n++) { // For each neighbor, check the hash table to see if its been added before
  90. Eigen::RowVector3i nkey = pi + neighbors[n];
  91. auto nbr = visited.find(nkey);
  92. if (nbr != visited.end()) { // We've already visited this neighbor, use references to its vertices instead of duplicating them
  93. vertexAlreadyAdded |= zv[n];
  94. for (int i = 0; i < 4; i++) { cube[zvv[n][i]] = CI_vector[nbr->second][zvv[n % 2 == 0 ? n + 1 : n - 1][i]]; }
  95. } else {
  96. queue.push_back(nkey); // Otherwise, we have not visited the neighbor, put it in the BFS queue
  97. }
  98. }
  99. for (int i = 0; i < 8; i++) { // Add new, non-visited,2 vertices to the arrays
  100. if (0 == ((1 << i) & vertexAlreadyAdded)) {
  101. cube[i] = CS_vector.size();
  102. CV_vector.push_back(cubeCorners[i]);
  103. CS_vector.push_back(cubeScalars[i]);
  104. }
  105. }
  106. visited[pi] = CI_vector.size();
  107. CI_vector.push_back(cube);
  108. }
  109. CV.conservativeResize(CV_vector.size(), 3);
  110. CS.conservativeResize(CS_vector.size(), 1);
  111. CI.conservativeResize(CI_vector.size(), 8);
  112. // If you pass in column-major matrices, this is going to be slooooowwwww
  113. for (int i = 0; i < CV_vector.size(); i++) {
  114. CV.row(i) = CV_vector[i];
  115. }
  116. for (int i = 0; i < CS_vector.size(); i++) {
  117. CS[i] = CS_vector[i];
  118. }
  119. for (int i = 0; i < CI_vector.size(); i++) {
  120. CI.row(i) = CI_vector[i];
  121. }
  122. }
  123. #ifdef IGL_STATIC_LIBRARY
  124. template void igl::sparse_voxel_grid<class Eigen::Matrix<double, -1, -1, 0, -1, -1>, class std::function<double(class Eigen::Matrix<double, -1, -1, 0, -1, -1> const &)>, class Eigen::Matrix<double, -1, 1, 0, -1, 1>, class Eigen::Matrix<double, -1, -1, 0, -1, -1>, class Eigen::Matrix<int, -1, -1, 0, -1, -1> >(class Eigen::MatrixBase<class Eigen::Matrix<double, -1, -1, 0, -1, -1> > const &, class std::function<double(class Eigen::Matrix<double, -1, -1, 0, -1, -1> const &)> const &, double, int, class Eigen::PlainObjectBase<class Eigen::Matrix<double, -1, 1, 0, -1, 1> > &, class Eigen::PlainObjectBase<class Eigen::Matrix<double, -1, -1, 0, -1, -1> > &, class Eigen::PlainObjectBase<class Eigen::Matrix<int, -1, -1, 0, -1, -1> > &);
  125. template void igl::sparse_voxel_grid<Eigen::Matrix<double, 1, 3, 1, 1, 3>, std::function<double (Eigen::Matrix<double, 1, 3, 1, 1, 3> const&)>, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1> >(Eigen::MatrixBase<Eigen::Matrix<double, 1, 3, 1, 1, 3> > const&, std::function<double (Eigen::Matrix<double, 1, 3, 1, 1, 3> const&)> const&, double, int, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&);
  126. #endif