propagate_winding_numbers.cpp 9.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260
  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 "order_facets_around_edge.h"
  17. #include "outer_facet.h"
  18. #include "closest_facet.h"
  19. #include "assign_scalar.h"
  20. #include "extract_cells.h"
  21. #include <stdexcept>
  22. #include <limits>
  23. #include <vector>
  24. #include <tuple>
  25. #include <queue>
  26. namespace propagate_winding_numbers_helper {
  27. template<
  28. typename DerivedF,
  29. typename DeriveduE,
  30. typename uE2EType >
  31. bool is_orientable(
  32. const Eigen::PlainObjectBase<DerivedF>& F,
  33. const Eigen::PlainObjectBase<DeriveduE>& uE,
  34. const std::vector<std::vector<uE2EType> >& uE2E) {
  35. const size_t num_faces = F.rows();
  36. const size_t num_edges = uE.rows();
  37. auto edge_index_to_face_index = [&](size_t ei) {
  38. return ei % num_faces;
  39. };
  40. auto is_consistent = [&](size_t fid, size_t s, size_t d) {
  41. if ((size_t)F(fid, 0) == s && (size_t)F(fid, 1) == d) return true;
  42. if ((size_t)F(fid, 1) == s && (size_t)F(fid, 2) == d) return true;
  43. if ((size_t)F(fid, 2) == s && (size_t)F(fid, 0) == d) return true;
  44. if ((size_t)F(fid, 0) == d && (size_t)F(fid, 1) == s) return false;
  45. if ((size_t)F(fid, 1) == d && (size_t)F(fid, 2) == s) return false;
  46. if ((size_t)F(fid, 2) == d && (size_t)F(fid, 0) == s) return false;
  47. throw "Invalid face!!";
  48. };
  49. for (size_t i=0; i<num_edges; i++) {
  50. const size_t s = uE(i,0);
  51. const size_t d = uE(i,1);
  52. int count=0;
  53. for (const auto& ei : uE2E[i]) {
  54. const size_t fid = edge_index_to_face_index(ei);
  55. if (is_consistent(fid, s, d)) count++;
  56. else count--;
  57. }
  58. if (count != 0) {
  59. return false;
  60. }
  61. }
  62. return true;
  63. }
  64. }
  65. template<
  66. typename DerivedV,
  67. typename DerivedF,
  68. typename DerivedL,
  69. typename DerivedW>
  70. IGL_INLINE void igl::copyleft::cgal::propagate_winding_numbers(
  71. const Eigen::PlainObjectBase<DerivedV>& V,
  72. const Eigen::PlainObjectBase<DerivedF>& F,
  73. const Eigen::PlainObjectBase<DerivedL>& labels,
  74. Eigen::PlainObjectBase<DerivedW>& W) {
  75. const size_t num_faces = F.rows();
  76. //typedef typename DerivedF::Scalar Index;
  77. Eigen::MatrixXi E, uE;
  78. Eigen::VectorXi EMAP;
  79. std::vector<std::vector<size_t> > uE2E;
  80. igl::unique_edge_map(F, E, uE, EMAP, uE2E);
  81. assert(propagate_winding_numbers_helper::is_orientable(F, uE, uE2E));
  82. Eigen::VectorXi P;
  83. const size_t num_patches = igl::extract_manifold_patches(F, EMAP, uE2E, P);
  84. DerivedW per_patch_cells;
  85. const size_t num_cells =
  86. igl::copyleft::cgal::extract_cells(
  87. V, F, P, E, uE, uE2E, EMAP, per_patch_cells);
  88. typedef std::tuple<size_t, bool, size_t> CellConnection;
  89. std::vector<std::set<CellConnection> > cell_adjacency(num_cells);
  90. for (size_t i=0; i<num_patches; i++) {
  91. const int positive_cell = per_patch_cells(i,0);
  92. const int negative_cell = per_patch_cells(i,1);
  93. cell_adjacency[positive_cell].emplace(negative_cell, false, i);
  94. cell_adjacency[negative_cell].emplace(positive_cell, true, i);
  95. }
  96. auto save_cell = [&](const std::string& filename, size_t cell_id) {
  97. std::vector<size_t> faces;
  98. for (size_t i=0; i<num_patches; i++) {
  99. if ((per_patch_cells.row(i).array() == cell_id).any()) {
  100. for (size_t j=0; j<num_faces; j++) {
  101. if ((size_t)P[j] == i) {
  102. faces.push_back(j);
  103. }
  104. }
  105. }
  106. }
  107. Eigen::MatrixXi cell_faces(faces.size(), 3);
  108. for (size_t i=0; i<faces.size(); i++) {
  109. cell_faces.row(i) = F.row(faces[i]);
  110. }
  111. Eigen::MatrixXd vertices(V.rows(), 3);
  112. for (size_t i=0; i<(size_t)V.rows(); i++) {
  113. assign_scalar(V(i,0), vertices(i,0));
  114. assign_scalar(V(i,1), vertices(i,1));
  115. assign_scalar(V(i,2), vertices(i,2));
  116. }
  117. writePLY(filename, vertices, cell_faces);
  118. };
  119. #ifndef NDEBUG
  120. {
  121. // Check for odd cycle.
  122. Eigen::VectorXi cell_labels(num_cells);
  123. cell_labels.setZero();
  124. Eigen::VectorXi parents(num_cells);
  125. parents.setConstant(-1);
  126. auto trace_parents = [&](size_t idx) {
  127. std::list<size_t> path;
  128. path.push_back(idx);
  129. while ((size_t)parents[path.back()] != path.back()) {
  130. path.push_back(parents[path.back()]);
  131. }
  132. return path;
  133. };
  134. for (size_t i=0; i<num_cells; i++) {
  135. if (cell_labels[i] == 0) {
  136. cell_labels[i] = 1;
  137. std::queue<size_t> Q;
  138. Q.push(i);
  139. parents[i] = i;
  140. while (!Q.empty()) {
  141. size_t curr_idx = Q.front();
  142. Q.pop();
  143. int curr_label = cell_labels[curr_idx];
  144. for (const auto& neighbor : cell_adjacency[curr_idx]) {
  145. if (cell_labels[std::get<0>(neighbor)] == 0) {
  146. cell_labels[std::get<0>(neighbor)] = curr_label * -1;
  147. Q.push(std::get<0>(neighbor));
  148. parents[std::get<0>(neighbor)] = curr_idx;
  149. } else {
  150. if (cell_labels[std::get<0>(neighbor)] !=
  151. curr_label * -1) {
  152. std::cerr << "Odd cell cycle detected!" << std::endl;
  153. auto path = trace_parents(curr_idx);
  154. path.reverse();
  155. auto path2 = trace_parents(std::get<0>(neighbor));
  156. path.insert(path.end(),
  157. path2.begin(), path2.end());
  158. for (auto cell_id : path) {
  159. std::cout << cell_id << " ";
  160. std::stringstream filename;
  161. filename << "cell_" << cell_id << ".ply";
  162. save_cell(filename.str(), cell_id);
  163. }
  164. std::cout << std::endl;
  165. }
  166. assert(cell_labels[std::get<0>(neighbor)] == curr_label * -1);
  167. }
  168. }
  169. }
  170. }
  171. }
  172. }
  173. #endif
  174. size_t outer_facet;
  175. bool flipped;
  176. Eigen::VectorXi I;
  177. I.setLinSpaced(num_faces, 0, num_faces-1);
  178. igl::copyleft::cgal::outer_facet(V, F, I, outer_facet, flipped);
  179. const size_t outer_patch = P[outer_facet];
  180. const size_t infinity_cell = per_patch_cells(outer_patch, flipped?1:0);
  181. Eigen::VectorXi patch_labels(num_patches);
  182. const int INVALID = std::numeric_limits<int>::max();
  183. patch_labels.setConstant(INVALID);
  184. for (size_t i=0; i<num_faces; i++) {
  185. if (patch_labels[P[i]] == INVALID) {
  186. patch_labels[P[i]] = labels[i];
  187. } else {
  188. assert(patch_labels[P[i]] == labels[i]);
  189. }
  190. }
  191. assert((patch_labels.array() != INVALID).all());
  192. const size_t num_labels = patch_labels.maxCoeff()+1;
  193. Eigen::MatrixXi per_cell_W(num_cells, num_labels);
  194. per_cell_W.setConstant(INVALID);
  195. per_cell_W.row(infinity_cell).setZero();
  196. std::queue<size_t> Q;
  197. Q.push(infinity_cell);
  198. while (!Q.empty()) {
  199. size_t curr_cell = Q.front();
  200. Q.pop();
  201. for (const auto& neighbor : cell_adjacency[curr_cell]) {
  202. size_t neighbor_cell, patch_idx;
  203. bool direction;
  204. std::tie(neighbor_cell, direction, patch_idx) = neighbor;
  205. if ((per_cell_W.row(neighbor_cell).array() == INVALID).any()) {
  206. per_cell_W.row(neighbor_cell) = per_cell_W.row(curr_cell);
  207. for (size_t i=0; i<num_labels; i++) {
  208. int inc = (patch_labels[patch_idx] == (int)i) ?
  209. (direction ? -1:1) :0;
  210. per_cell_W(neighbor_cell, i) =
  211. per_cell_W(curr_cell, i) + inc;
  212. }
  213. Q.push(neighbor_cell);
  214. } else {
  215. #ifndef NDEBUG
  216. for (size_t i=0; i<num_labels; i++) {
  217. if ((int)i == patch_labels[patch_idx]) {
  218. int inc = direction ? -1:1;
  219. assert(per_cell_W(neighbor_cell, i) ==
  220. per_cell_W(curr_cell, i) + inc);
  221. } else {
  222. assert(per_cell_W(neighbor_cell, i) ==
  223. per_cell_W(curr_cell, i));
  224. }
  225. }
  226. #endif
  227. }
  228. }
  229. }
  230. W.resize(num_faces, num_labels*2);
  231. for (size_t i=0; i<num_faces; i++) {
  232. const size_t patch = P[i];
  233. const size_t positive_cell = per_patch_cells(patch, 0);
  234. const size_t negative_cell = per_patch_cells(patch, 1);
  235. for (size_t j=0; j<num_labels; j++) {
  236. W(i,j*2 ) = per_cell_W(positive_cell, j);
  237. W(i,j*2+1) = per_cell_W(negative_cell, j);
  238. }
  239. }
  240. }
  241. #ifdef IGL_STATIC_LIBRARY
  242. 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> >&);
  243. #endif