propagate_winding_numbers.cpp 11 KB

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