remesh_intersections_beta.cpp 13 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338
  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 "remesh_intersections_beta.h"
  10. #include <vector>
  11. #include <map>
  12. #include <unordered_map>
  13. #include <iostream>
  14. template <
  15. typename DerivedV,
  16. typename DerivedF,
  17. typename Kernel,
  18. typename DerivedVV,
  19. typename DerivedFF,
  20. typename DerivedJ,
  21. typename DerivedIM>
  22. IGL_INLINE void igl::cgal::remesh_intersections_beta(
  23. const Eigen::PlainObjectBase<DerivedV> & V,
  24. const Eigen::PlainObjectBase<DerivedF> & F,
  25. const std::vector<CGAL::Triangle_3<Kernel> > & T,
  26. const std::map<
  27. typename DerivedF::Index,
  28. std::pair<typename DerivedF::Index,
  29. std::vector<CGAL::Object> > > & offending,
  30. const std::map<
  31. std::pair<typename DerivedF::Index,typename DerivedF::Index>,
  32. std::vector<typename DerivedF::Index> > & edge2faces,
  33. Eigen::PlainObjectBase<DerivedVV> & VV,
  34. Eigen::PlainObjectBase<DerivedFF> & FF,
  35. Eigen::PlainObjectBase<DerivedJ> & J,
  36. Eigen::PlainObjectBase<DerivedIM> & IM) {
  37. typedef CGAL::Point_3<Kernel> Point_3;
  38. typedef CGAL::Segment_3<Kernel> Segment_3;
  39. typedef CGAL::Triangle_3<Kernel> Triangle_3;
  40. typedef CGAL::Plane_3<Kernel> Plane_3;
  41. typedef CGAL::Point_2<Kernel> Point_2;
  42. typedef CGAL::Segment_2<Kernel> Segment_2;
  43. typedef CGAL::Triangle_2<Kernel> Triangle_2;
  44. typedef CGAL::Triangulation_vertex_base_2<Kernel> TVB_2;
  45. typedef CGAL::Constrained_triangulation_face_base_2<Kernel> CTFB_2;
  46. typedef CGAL::Triangulation_data_structure_2<TVB_2,CTFB_2> TDS_2;
  47. typedef CGAL::Exact_intersections_tag Itag;
  48. typedef CGAL::Constrained_Delaunay_triangulation_2<Kernel,TDS_2,Itag>
  49. CDT_2;
  50. typedef CGAL::Constrained_triangulation_plus_2<CDT_2> CDT_plus_2;
  51. typedef typename DerivedF::Index Index;
  52. typedef std::pair<Index, Index> Edge;
  53. struct EdgeHash {
  54. size_t operator()(const Edge& e) const {
  55. return (e.first * 805306457) ^ (e.second * 201326611);
  56. }
  57. };
  58. typedef std::unordered_map<Edge, std::vector<Index>, EdgeHash > EdgeMap;
  59. auto normalize_plane_coeff = [](const Plane_3& P) {
  60. std::vector<typename Kernel::FT> coeffs = {
  61. P.a(), P.b(), P.c(), P.d()
  62. };
  63. const auto max_itr = std::max_element(coeffs.begin(), coeffs.end());
  64. const auto min_itr = std::min_element(coeffs.begin(), coeffs.end());
  65. typename Kernel::FT max_coeff;
  66. if (*max_itr < -1 * *min_itr) {
  67. max_coeff = *min_itr;
  68. } else {
  69. max_coeff = *max_itr;
  70. }
  71. std::transform(coeffs.begin(), coeffs.end(), coeffs.begin(),
  72. [&](const typename Kernel::FT& val)
  73. {return val / max_coeff; } );
  74. return coeffs;
  75. };
  76. auto plane_comp = [&](const Plane_3& p1, const Plane_3& p2) {
  77. const auto p1_coeffs = normalize_plane_coeff(p1);
  78. const auto p2_coeffs = normalize_plane_coeff(p2);
  79. if (p1_coeffs[0] != p2_coeffs[0])
  80. return p1_coeffs[0] < p2_coeffs[0];
  81. if (p1_coeffs[1] != p2_coeffs[1])
  82. return p1_coeffs[1] < p2_coeffs[1];
  83. if (p1_coeffs[2] != p2_coeffs[2])
  84. return p1_coeffs[2] < p2_coeffs[2];
  85. if (p1_coeffs[3] != p2_coeffs[3])
  86. return p1_coeffs[3] < p2_coeffs[3];
  87. return false;
  88. };
  89. std::map<Plane_3, std::vector<Index>, decltype(plane_comp)>
  90. unique_planes(plane_comp);
  91. const size_t num_faces = F.rows();
  92. const size_t num_base_vertices = V.rows();
  93. assert(num_faces == T.size());
  94. std::vector<bool> is_offending(num_faces, false);
  95. for (const auto itr : offending) {
  96. const auto& fid = itr.first;
  97. is_offending[fid] = true;
  98. Plane_3 key = T[fid].supporting_plane();
  99. assert(!key.is_degenerate());
  100. const auto jtr = unique_planes.find(key);
  101. if (jtr == unique_planes.end()) {
  102. unique_planes.insert({key, {fid}});
  103. } else {
  104. jtr->second.push_back(fid);
  105. }
  106. }
  107. const size_t INVALID = std::numeric_limits<size_t>::max();
  108. std::vector<std::vector<Index> > resolved_faces;
  109. std::vector<Index> source_faces;
  110. std::vector<Point_3> new_vertices;
  111. EdgeMap edge_vertices;
  112. /**
  113. * Run constraint Delaunay triangulation on the plane.
  114. */
  115. auto run_delaunay_triangulation = [&](const Plane_3& P,
  116. const std::vector<Index>& involved_faces,
  117. std::vector<Point_3>& vertices,
  118. std::vector<std::vector<Index> >& faces) {
  119. CDT_plus_2 cdt;
  120. for (const auto& fid : involved_faces) {
  121. const auto itr = offending.find(fid);
  122. const auto& triangle = T[fid];
  123. cdt.insert_constraint(P.to_2d(triangle[0]), P.to_2d(triangle[1]));
  124. cdt.insert_constraint(P.to_2d(triangle[1]), P.to_2d(triangle[2]));
  125. cdt.insert_constraint(P.to_2d(triangle[2]), P.to_2d(triangle[0]));
  126. if (itr == offending.end()) continue;
  127. for (const auto& obj : itr->second.second) {
  128. if(const Segment_3 *iseg = CGAL::object_cast<Segment_3 >(&obj)) {
  129. // Add segment constraint
  130. cdt.insert_constraint(
  131. P.to_2d(iseg->vertex(0)),P.to_2d(iseg->vertex(1)));
  132. }else if(const Point_3 *ipoint = CGAL::object_cast<Point_3 >(&obj)) {
  133. // Add point
  134. cdt.insert(P.to_2d(*ipoint));
  135. } else if(const Triangle_3 *itri = CGAL::object_cast<Triangle_3 >(&obj)) {
  136. // Add 3 segment constraints
  137. cdt.insert_constraint(
  138. P.to_2d(itri->vertex(0)),P.to_2d(itri->vertex(1)));
  139. cdt.insert_constraint(
  140. P.to_2d(itri->vertex(1)),P.to_2d(itri->vertex(2)));
  141. cdt.insert_constraint(
  142. P.to_2d(itri->vertex(2)),P.to_2d(itri->vertex(0)));
  143. } else if(const std::vector<Point_3 > *polyp =
  144. CGAL::object_cast< std::vector<Point_3 > >(&obj)) {
  145. //cerr<<REDRUM("Poly...")<<endl;
  146. const std::vector<Point_3 > & poly = *polyp;
  147. const Index m = poly.size();
  148. assert(m>=2);
  149. for(Index p = 0;p<m;p++)
  150. {
  151. const Index np = (p+1)%m;
  152. cdt.insert_constraint(P.to_2d(poly[p]),P.to_2d(poly[np]));
  153. }
  154. }else {
  155. throw std::runtime_error("Unknown intersection object!");
  156. }
  157. }
  158. }
  159. const size_t num_vertices = cdt.number_of_vertices();
  160. const size_t num_faces = cdt.number_of_faces();
  161. std::map<typename CDT_plus_2::Vertex_handle,Index> v2i;
  162. size_t count=0;
  163. for (auto itr = cdt.finite_vertices_begin();
  164. itr != cdt.finite_vertices_end(); itr++) {
  165. vertices.push_back(P.to_3d(itr->point()));
  166. v2i[itr] = count;
  167. count++;
  168. }
  169. for (auto itr = cdt.finite_faces_begin();
  170. itr != cdt.finite_faces_end(); itr++) {
  171. faces.push_back( {
  172. v2i[itr->vertex(0)],
  173. v2i[itr->vertex(1)],
  174. v2i[itr->vertex(2)] });
  175. }
  176. };
  177. /**
  178. * Given p on triangle indexed by ori_f, determine the index of p.
  179. */
  180. auto determine_point_index = [&](
  181. const Point_3& p, const size_t ori_f) -> Index {
  182. const auto& triangle = T[ori_f];
  183. const auto& f = F.row(ori_f).eval();
  184. // Check if p is one of the triangle corners.
  185. for (size_t i=0; i<3; i++) {
  186. if (p == triangle[i]) return f[i];
  187. }
  188. // Check if p is on one of the edges.
  189. for (size_t i=0; i<3; i++) {
  190. const Point_3 curr_corner = triangle[i];
  191. const Point_3 next_corner = triangle[(i+1)%3];
  192. const Segment_3 edge(curr_corner, next_corner);
  193. if (edge.has_on(p)) {
  194. const Index curr = f[i];
  195. const Index next = f[(i+1)%3];
  196. Edge key;
  197. key.first = curr<next?curr:next;
  198. key.second = curr<next?next:curr;
  199. auto itr = edge_vertices.find(key);
  200. if (itr == edge_vertices.end()) {
  201. const Index index =
  202. num_base_vertices + new_vertices.size();
  203. edge_vertices.insert({key, {index}});
  204. new_vertices.push_back(p);
  205. return index;
  206. } else {
  207. for (const auto vid : itr->second) {
  208. if (p == new_vertices[vid - num_base_vertices]) {
  209. return vid;
  210. }
  211. }
  212. const size_t index = num_base_vertices + new_vertices.size();
  213. new_vertices.push_back(p);
  214. itr->second.push_back(index);
  215. return index;
  216. }
  217. }
  218. }
  219. // p must be in the middle of the triangle.
  220. const size_t index = num_base_vertices + new_vertices.size();
  221. new_vertices.push_back(p);
  222. return index;
  223. };
  224. /**
  225. * Determine the vertex indices for each corner of each output triangle.
  226. */
  227. auto post_triangulation_process = [&](
  228. const std::vector<Point_3>& vertices,
  229. const std::vector<std::vector<Index> >& faces,
  230. const std::vector<Index>& involved_faces) {
  231. for (const auto& f : faces) {
  232. const Point_3& v0 = vertices[f[0]];
  233. const Point_3& v1 = vertices[f[1]];
  234. const Point_3& v2 = vertices[f[2]];
  235. Point_3 center(
  236. (v0[0] + v1[0] + v2[0]) / 3.0,
  237. (v0[1] + v1[1] + v2[1]) / 3.0,
  238. (v0[2] + v1[2] + v2[2]) / 3.0);
  239. for (const auto& ori_f : involved_faces) {
  240. const auto& triangle = T[ori_f];
  241. const Plane_3 P = triangle.supporting_plane();
  242. if (triangle.has_on(center)) {
  243. std::vector<Index> corners(3);
  244. corners[0] = determine_point_index(v0, ori_f);
  245. corners[1] = determine_point_index(v1, ori_f);
  246. corners[2] = determine_point_index(v2, ori_f);
  247. if (CGAL::orientation(
  248. P.to_2d(v0), P.to_2d(v1), P.to_2d(v2))
  249. == CGAL::RIGHT_TURN) {
  250. std::swap(corners[0], corners[1]);
  251. }
  252. resolved_faces.emplace_back(corners);
  253. source_faces.push_back(ori_f);
  254. }
  255. }
  256. }
  257. };
  258. // Process un-touched faces.
  259. for (size_t i=0; i<num_faces; i++) {
  260. if (!is_offending[i]) {
  261. resolved_faces.push_back(
  262. { F(i,0), F(i,1), F(i,2) } );
  263. source_faces.push_back(i);
  264. }
  265. }
  266. // Process self-intersecting faces.
  267. for (const auto itr : unique_planes) {
  268. Plane_3 P = itr.first;
  269. const auto& involved_faces = itr.second;
  270. std::vector<Point_3> vertices;
  271. std::vector<std::vector<Index> > faces;
  272. run_delaunay_triangulation(P, involved_faces, vertices, faces);
  273. post_triangulation_process(vertices, faces, involved_faces);
  274. }
  275. // Output resolved mesh.
  276. const size_t num_out_vertices = new_vertices.size() + num_base_vertices;
  277. VV.resize(num_out_vertices, 3);
  278. for (size_t i=0; i<num_base_vertices; i++) {
  279. assign_scalar(V(i,0), VV(i,0));
  280. assign_scalar(V(i,1), VV(i,1));
  281. assign_scalar(V(i,2), VV(i,2));
  282. }
  283. for (size_t i=num_base_vertices; i<num_out_vertices; i++) {
  284. assign_scalar(new_vertices[i-num_base_vertices][0], VV(i,0));
  285. assign_scalar(new_vertices[i-num_base_vertices][1], VV(i,1));
  286. assign_scalar(new_vertices[i-num_base_vertices][2], VV(i,2));
  287. }
  288. const size_t num_out_faces = resolved_faces.size();
  289. FF.resize(num_out_faces, 3);
  290. for (size_t i=0; i<num_out_faces; i++) {
  291. FF(i,0) = resolved_faces[i][0];
  292. FF(i,1) = resolved_faces[i][1];
  293. FF(i,2) = resolved_faces[i][2];
  294. }
  295. J.resize(num_out_faces);
  296. std::copy(source_faces.begin(), source_faces.end(), J.data());
  297. // Extract unique vertex indices.
  298. IM.resize(VV.rows(),1);
  299. std::map<Point_3,Index> vv2i;
  300. // Safe to check for duplicates using double for original vertices: if
  301. // incoming reps are different then the points are unique.
  302. for(Index v = 0;v<VV.rows();v++) {
  303. typename Kernel::FT p0,p1,p2;
  304. assign_scalar(VV(v,0),p0);
  305. assign_scalar(VV(v,1),p1);
  306. assign_scalar(VV(v,2),p2);
  307. const Point_3 p(p0,p1,p2);
  308. if(vv2i.count(p)==0) {
  309. vv2i[p] = v;
  310. }
  311. assert(vv2i.count(p) == 1);
  312. IM(v) = vv2i[p];
  313. }
  314. }