marching_tets.cpp 6.4 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2018 Francis Williams <francis@fwilliams.info>
  4. //
  5. // This Source Code Form is subject to the terms of the Mozilla Public License
  6. // v. 2.0. If a copy of the MPL was not distributed with this file, You can
  7. // obtain one at http://mozilla.org/MPL/2.0/.
  8. #include "marching_tets.h"
  9. #include <unordered_map>
  10. #include <vector>
  11. #include <utility>
  12. #include <cstdint>
  13. #include <iostream>
  14. template <typename DerivedTV,
  15. typename DerivedTT,
  16. typename DerivedS,
  17. typename DerivedSV,
  18. typename DerivedSF,
  19. typename DerivedJ,
  20. typename BCType>
  21. void igl::marching_tets(
  22. const Eigen::PlainObjectBase<DerivedTV>& TV,
  23. const Eigen::PlainObjectBase<DerivedTT>& TT,
  24. const Eigen::PlainObjectBase<DerivedS>& isovals,
  25. double isovalue,
  26. Eigen::PlainObjectBase<DerivedSV>& outV,
  27. Eigen::PlainObjectBase<DerivedSF>& outF,
  28. Eigen::PlainObjectBase<DerivedJ>& J,
  29. Eigen::SparseMatrix<BCType>& BC)
  30. {
  31. using namespace std;
  32. // We're hashing edges to deduplicate using 64 bit ints. The upper and lower
  33. // 32 bits of a key are the indices of vertices in the mesh. The implication is
  34. // that you can only have 2^32 vertices which I have deemed sufficient for
  35. // anything reasonable.
  36. const auto make_edge_key = [](const pair<int32_t, int32_t>& p) -> int64_t
  37. {
  38. std::int64_t ret = 0;
  39. ret |= p.first;
  40. ret |= static_cast<std::int64_t>(p.second) << 32;
  41. return ret;
  42. };
  43. const int mt_cell_lookup[16][4] =
  44. {
  45. { -1, -1, -1, -1 },
  46. { 0, 2, 1, -1 },
  47. { 0, 3, 4, -1 },
  48. { 2, 1, 3, 4 },
  49. { 5, 3, 1, -1 },
  50. { 0, 2, 5, 3 },
  51. { 0, 1, 5, 4 },
  52. { 2, 5, 4, -1 },
  53. { 4, 5, 2, -1 },
  54. { 0, 4, 5, 1 },
  55. { 0, 3, 5, 2 },
  56. { 1, 3, 5, -1 },
  57. { 4, 3, 1, 2 },
  58. { 0, 4, 3, -1 },
  59. { 0, 1, 2, -1 },
  60. { -1, -1, -1, -1 },
  61. };
  62. const int mt_edge_lookup[6][2] =
  63. {
  64. {0, 1},
  65. {0, 2},
  66. {0, 3},
  67. {1, 2},
  68. {1, 3},
  69. {2, 3},
  70. };
  71. // Store the faces and the tet they are in
  72. vector<pair<Eigen::RowVector3i, int>> faces;
  73. // Store the edges in the tet mesh which we add vertices on
  74. // so we can deduplicate
  75. vector<pair<int, int>> edge_table;
  76. assert(TT.cols() == 4 && TT.rows() >= 1);
  77. assert(TV.cols() == 3 && TV.rows() >= 4);
  78. assert(isovals.cols() == 1);
  79. // For each tet
  80. for (int i = 0; i < TT.rows(); i++)
  81. {
  82. uint8_t key = 0;
  83. for (int v = 0; v < 4; v++)
  84. {
  85. const int vid = TT(i, v);
  86. const uint8_t flag = isovals[vid] > isovalue;
  87. key |= flag << v;
  88. }
  89. // This will contain the index in TV of each vertex in the tet
  90. int v_ids[4] = {-1, -1, -1, -1};
  91. // Insert any vertices if the tet intersects the level surface
  92. for (int e = 0; e < 4 && mt_cell_lookup[key][e] != -1; e++)
  93. {
  94. const int tv1_idx = TT(i, mt_edge_lookup[mt_cell_lookup[key][e]][0]);
  95. const int tv2_idx = TT(i, mt_edge_lookup[mt_cell_lookup[key][e]][1]);
  96. const int vertex_id = edge_table.size();
  97. edge_table.push_back(make_pair(min(tv1_idx, tv2_idx), max(tv1_idx, tv2_idx)));
  98. v_ids[e] = vertex_id;
  99. }
  100. // Insert the corresponding faces
  101. if (v_ids[0] != -1)
  102. {
  103. bool is_quad = mt_cell_lookup[key][3] != -1;
  104. if (is_quad)
  105. {
  106. const Eigen::RowVector3i f1(v_ids[0], v_ids[1], v_ids[3]);
  107. const Eigen::RowVector3i f2(v_ids[1], v_ids[2], v_ids[3]);
  108. faces.push_back(make_pair(f1, i));
  109. faces.push_back(make_pair(f2, i));
  110. }
  111. else
  112. {
  113. const Eigen::RowVector3i f(v_ids[0], v_ids[1], v_ids[2]);
  114. faces.push_back(make_pair(f, i));
  115. }
  116. }
  117. }
  118. int num_unique = 0;
  119. outV.resize(edge_table.size(), 3);
  120. outF.resize(faces.size(), 3);
  121. J.resize(faces.size());
  122. // Sparse matrix triplets for BC
  123. vector<Eigen::Triplet<BCType>> bc_triplets;
  124. bc_triplets.reserve(edge_table.size());
  125. // Deduplicate vertices
  126. unordered_map<int64_t, int> emap;
  127. emap.max_load_factor(0.5);
  128. emap.reserve(edge_table.size());
  129. for (int f = 0; f < faces.size(); f++)
  130. {
  131. for (int v = 0; v < 3; v++)
  132. {
  133. const int vi = faces[f].first[v];
  134. const int ti = faces[f].second;
  135. const pair<int32_t, int32_t> edge = edge_table[vi];
  136. const int64_t key = make_edge_key(edge);
  137. auto it = emap.find(key);
  138. if (it == emap.end()) // New unique vertex, insert it
  139. {
  140. // Typedef to make sure we handle floats properly
  141. typedef Eigen::Matrix<typename DerivedTV::Scalar, 1, 3, Eigen::RowMajor, 1, 3> RowVector;
  142. const RowVector v1 = TV.row(edge.first);
  143. const RowVector v2 = TV.row(edge.second);
  144. const double a = fabs(isovals[edge.first] - isovalue);
  145. const double b = fabs(isovals[edge.second] - isovalue);
  146. const double w = a / (a+b);
  147. // Create a casted copy in case BCType is a float and we need to downcast
  148. const BCType bc_w = static_cast<BCType>(w);
  149. bc_triplets.push_back(Eigen::Triplet<BCType>(num_unique, edge.first, 1-bc_w));
  150. bc_triplets.push_back(Eigen::Triplet<BCType>(num_unique, edge.second, bc_w));
  151. // Create a casted copy in case DerivedTV::Scalar is a float and we need to downcast
  152. const typename DerivedTV::Scalar v_w = static_cast<typename DerivedTV::Scalar>(w);
  153. outV.row(num_unique) = (1-v_w)*v1 + v_w*v2;
  154. outF(f, v) = num_unique;
  155. J[f] = ti;
  156. emap.emplace(key, num_unique);
  157. num_unique += 1;
  158. } else {
  159. outF(f, v) = it->second;
  160. }
  161. }
  162. }
  163. outV.conservativeResize(num_unique, 3);
  164. J.conservativeResize(num_unique, 1);
  165. BC.resize(num_unique, TV.rows());
  166. BC.setFromTriplets(bc_triplets.begin(), bc_triplets.end());
  167. }
  168. #ifdef IGL_STATIC_LIBRARY
  169. template void igl::marching_tets<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, 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::Matrix<int, -1, 1, 0, -1, 1>, double>(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, double, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&, Eigen::SparseMatrix<double, 0, int>&);
  170. #endif // IGL_STATIC_LIBRARY