propagate_winding_numbers.cpp 9.5 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2015 Qingnan Zhou <qnzhou@gmail.com>
  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. //
  9. #include "propagate_winding_numbers.h"
  10. #include "../../extract_manifold_patches.h"
  11. #include "../../extract_non_manifold_edge_curves.h"
  12. #include "../../facet_components.h"
  13. #include "../../unique_edge_map.h"
  14. #include "../../piecewise_constant_winding_number.h"
  15. #include "../../writeOBJ.h"
  16. #include "../../writePLY.h"
  17. #include "../../get_seconds.h"
  18. #include "order_facets_around_edge.h"
  19. #include "outer_facet.h"
  20. #include "closest_facet.h"
  21. #include "assign_scalar.h"
  22. #include "extract_cells.h"
  23. #include <stdexcept>
  24. #include <limits>
  25. #include <vector>
  26. #include <tuple>
  27. #include <queue>
  28. //#define PROPAGATE_WINDING_NUMBER_TIMING
  29. template<
  30. typename DerivedV,
  31. typename DerivedF,
  32. typename DerivedL,
  33. typename DerivedW>
  34. IGL_INLINE void igl::copyleft::cgal::propagate_winding_numbers(
  35. const Eigen::PlainObjectBase<DerivedV>& V,
  36. const Eigen::PlainObjectBase<DerivedF>& F,
  37. const Eigen::PlainObjectBase<DerivedL>& labels,
  38. Eigen::PlainObjectBase<DerivedW>& W) {
  39. #ifdef PROPAGATE_WINDING_NUMBER_TIMING
  40. const auto & tictoc = []()
  41. {
  42. static double t_start = igl::get_seconds();
  43. double diff = igl::get_seconds()-t_start;
  44. t_start += diff;
  45. return diff;
  46. };
  47. tictoc();
  48. #endif
  49. const size_t num_faces = F.rows();
  50. //typedef typename DerivedF::Scalar Index;
  51. Eigen::MatrixXi E, uE;
  52. Eigen::VectorXi EMAP;
  53. std::vector<std::vector<size_t> > uE2E;
  54. igl::unique_edge_map(F, E, uE, EMAP, uE2E);
  55. if (!piecewise_constant_winding_number(F, uE, uE2E))
  56. {
  57. std::cerr << "Input mesh is not orientable!" << std::endl;
  58. }
  59. Eigen::VectorXi P;
  60. const size_t num_patches = igl::extract_manifold_patches(F, EMAP, uE2E, P);
  61. #ifdef PROPAGATE_WINDING_NUMBER_TIMING
  62. std::cout << "extract manifold patches: " << tictoc() << std::endl;
  63. #endif
  64. DerivedW per_patch_cells;
  65. const size_t num_cells =
  66. igl::copyleft::cgal::extract_cells(
  67. V, F, P, E, uE, uE2E, EMAP, per_patch_cells);
  68. #ifdef PROPAGATE_WINDING_NUMBER_TIMING
  69. std::cout << "extract cells: " << tictoc() << std::endl;;
  70. #endif
  71. typedef std::tuple<size_t, bool, size_t> CellConnection;
  72. std::vector<std::set<CellConnection> > cell_adjacency(num_cells);
  73. for (size_t i=0; i<num_patches; i++) {
  74. const int positive_cell = per_patch_cells(i,0);
  75. const int negative_cell = per_patch_cells(i,1);
  76. cell_adjacency[positive_cell].emplace(negative_cell, false, i);
  77. cell_adjacency[negative_cell].emplace(positive_cell, true, i);
  78. }
  79. #ifdef PROPAGATE_WINDING_NUMBER_TIMING
  80. std::cout << "cell connection: " << tictoc() << std::endl;
  81. #endif
  82. auto save_cell = [&](const std::string& filename, size_t cell_id) {
  83. std::vector<size_t> faces;
  84. for (size_t i=0; i<num_patches; i++) {
  85. if ((per_patch_cells.row(i).array() == cell_id).any()) {
  86. for (size_t j=0; j<num_faces; j++) {
  87. if ((size_t)P[j] == i) {
  88. faces.push_back(j);
  89. }
  90. }
  91. }
  92. }
  93. Eigen::MatrixXi cell_faces(faces.size(), 3);
  94. for (size_t i=0; i<faces.size(); i++) {
  95. cell_faces.row(i) = F.row(faces[i]);
  96. }
  97. Eigen::MatrixXd vertices(V.rows(), 3);
  98. for (size_t i=0; i<(size_t)V.rows(); i++) {
  99. assign_scalar(V(i,0), vertices(i,0));
  100. assign_scalar(V(i,1), vertices(i,1));
  101. assign_scalar(V(i,2), vertices(i,2));
  102. }
  103. writePLY(filename, vertices, cell_faces);
  104. };
  105. #ifndef NDEBUG
  106. {
  107. // Check for odd cycle.
  108. Eigen::VectorXi cell_labels(num_cells);
  109. cell_labels.setZero();
  110. Eigen::VectorXi parents(num_cells);
  111. parents.setConstant(-1);
  112. auto trace_parents = [&](size_t idx) {
  113. std::list<size_t> path;
  114. path.push_back(idx);
  115. while ((size_t)parents[path.back()] != path.back()) {
  116. path.push_back(parents[path.back()]);
  117. }
  118. return path;
  119. };
  120. for (size_t i=0; i<num_cells; i++) {
  121. if (cell_labels[i] == 0) {
  122. cell_labels[i] = 1;
  123. std::queue<size_t> Q;
  124. Q.push(i);
  125. parents[i] = i;
  126. while (!Q.empty()) {
  127. size_t curr_idx = Q.front();
  128. Q.pop();
  129. int curr_label = cell_labels[curr_idx];
  130. for (const auto& neighbor : cell_adjacency[curr_idx]) {
  131. if (cell_labels[std::get<0>(neighbor)] == 0) {
  132. cell_labels[std::get<0>(neighbor)] = curr_label * -1;
  133. Q.push(std::get<0>(neighbor));
  134. parents[std::get<0>(neighbor)] = curr_idx;
  135. } else {
  136. if (cell_labels[std::get<0>(neighbor)] !=
  137. curr_label * -1) {
  138. std::cerr << "Odd cell cycle detected!" << std::endl;
  139. auto path = trace_parents(curr_idx);
  140. path.reverse();
  141. auto path2 = trace_parents(std::get<0>(neighbor));
  142. path.insert(path.end(),
  143. path2.begin(), path2.end());
  144. for (auto cell_id : path) {
  145. std::cout << cell_id << " ";
  146. std::stringstream filename;
  147. filename << "cell_" << cell_id << ".ply";
  148. save_cell(filename.str(), cell_id);
  149. }
  150. std::cout << std::endl;
  151. }
  152. // Do not fail when odd cycle is detected because the resulting
  153. // integer winding number field, although inconsistent, may still
  154. // be used if the problem region is local and embedded within a
  155. // valid volume.
  156. //assert(cell_labels[std::get<0>(neighbor)] == curr_label * -1);
  157. }
  158. }
  159. }
  160. }
  161. }
  162. #ifdef PROPAGATE_WINDING_NUMBER_TIMING
  163. std::cout << "check for odd cycle: " << tictoc() << std::endl;
  164. #endif
  165. }
  166. #endif
  167. size_t outer_facet;
  168. bool flipped;
  169. Eigen::VectorXi I;
  170. I.setLinSpaced(num_faces, 0, num_faces-1);
  171. igl::copyleft::cgal::outer_facet(V, F, I, outer_facet, flipped);
  172. #ifdef PROPAGATE_WINDING_NUMBER_TIMING
  173. std::cout << "outer facet: " << tictoc() << std::endl;
  174. #endif
  175. const size_t outer_patch = P[outer_facet];
  176. const size_t infinity_cell = per_patch_cells(outer_patch, flipped?1:0);
  177. Eigen::VectorXi patch_labels(num_patches);
  178. const int INVALID = std::numeric_limits<int>::max();
  179. patch_labels.setConstant(INVALID);
  180. for (size_t i=0; i<num_faces; i++) {
  181. if (patch_labels[P[i]] == INVALID) {
  182. patch_labels[P[i]] = labels[i];
  183. } else {
  184. assert(patch_labels[P[i]] == labels[i]);
  185. }
  186. }
  187. assert((patch_labels.array() != INVALID).all());
  188. const size_t num_labels = patch_labels.maxCoeff()+1;
  189. Eigen::MatrixXi per_cell_W(num_cells, num_labels);
  190. per_cell_W.setConstant(INVALID);
  191. per_cell_W.row(infinity_cell).setZero();
  192. std::queue<size_t> Q;
  193. Q.push(infinity_cell);
  194. while (!Q.empty()) {
  195. size_t curr_cell = Q.front();
  196. Q.pop();
  197. for (const auto& neighbor : cell_adjacency[curr_cell]) {
  198. size_t neighbor_cell, patch_idx;
  199. bool direction;
  200. std::tie(neighbor_cell, direction, patch_idx) = neighbor;
  201. if ((per_cell_W.row(neighbor_cell).array() == INVALID).any()) {
  202. per_cell_W.row(neighbor_cell) = per_cell_W.row(curr_cell);
  203. for (size_t i=0; i<num_labels; i++) {
  204. int inc = (patch_labels[patch_idx] == (int)i) ?
  205. (direction ? -1:1) :0;
  206. per_cell_W(neighbor_cell, i) =
  207. per_cell_W(curr_cell, i) + inc;
  208. }
  209. Q.push(neighbor_cell);
  210. } else {
  211. #ifndef NDEBUG
  212. // Checking for winding number consistency.
  213. // This check would inevitably fail for meshes that contain open
  214. // boundary or non-orientable. However, the inconsistent winding number
  215. // field would still be useful in some cases such as when problem region
  216. // is local and embedded within the volume. This, unfortunately, is the
  217. // best we can do because the problem of computing integer winding
  218. // number is ill-defined for open and non-orientable surfaces.
  219. for (size_t i=0; i<num_labels; i++) {
  220. if ((int)i == patch_labels[patch_idx]) {
  221. int inc = direction ? -1:1;
  222. //assert(per_cell_W(neighbor_cell, i) ==
  223. // per_cell_W(curr_cell, i) + inc);
  224. } else {
  225. //assert(per_cell_W(neighbor_cell, i) ==
  226. // per_cell_W(curr_cell, i));
  227. }
  228. }
  229. #endif
  230. }
  231. }
  232. }
  233. #ifdef PROPAGATE_WINDING_NUMBER_TIMING
  234. std::cout << "propagate winding number: " << tictoc() << std::endl;
  235. #endif
  236. W.resize(num_faces, num_labels*2);
  237. for (size_t i=0; i<num_faces; i++) {
  238. const size_t patch = P[i];
  239. const size_t positive_cell = per_patch_cells(patch, 0);
  240. const size_t negative_cell = per_patch_cells(patch, 1);
  241. for (size_t j=0; j<num_labels; j++) {
  242. W(i,j*2 ) = per_cell_W(positive_cell, j);
  243. W(i,j*2+1) = per_cell_W(negative_cell, j);
  244. }
  245. }
  246. #ifdef PROPAGATE_WINDING_NUMBER_TIMING
  247. std::cout << "save result: " << tictoc() << std::endl;
  248. #endif
  249. }
  250. #ifdef IGL_STATIC_LIBRARY
  251. template void igl::copyleft::cgal::propagate_winding_numbers<Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, 1, 0, -1, 1>, Eigen::Matrix<int, -1, -1, 0, -1, -1> >(Eigen::PlainObjectBase<Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&);
  252. #endif