SelfIntersectMesh.h 33 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068106910701071107210731074107510761077107810791080108110821083108410851086108710881089109010911092109310941095109610971098109911001101110211031104110511061107110811091110111111121113111411151116111711181119112011211122112311241125112611271128112911301131113211331134113511361137113811391140114111421143114411451146114711481149115011511152115311541155115611571158115911601161116211631164116511661167116811691170117111721173117411751176117711781179118011811182118311841185118611871188118911901191
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2014 Alec Jacobson <alecjacobson@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. #ifndef IGL_SELFINTERSECTMESH_H
  9. #define IGL_SELFINTERSECTMESH_H
  10. #include "CGAL_includes.hpp"
  11. #include "remesh_self_intersections.h"
  12. #include <Eigen/Dense>
  13. #include <list>
  14. #include <map>
  15. #include <vector>
  16. #ifndef IGL_FIRST_HIT_EXCEPTION
  17. #define IGL_FIRST_HIT_EXCEPTION 10
  18. #endif
  19. // The easiest way to keep track of everything is to use a class
  20. namespace igl
  21. {
  22. // Kernel is a CGAL kernel like:
  23. // CGAL::Exact_predicates_inexact_constructions_kernel
  24. // or
  25. // CGAL::Exact_predicates_exact_constructions_kernel
  26. template <
  27. typename Kernel,
  28. typename DerivedV,
  29. typename DerivedF,
  30. typename DerivedVV,
  31. typename DerivedFF,
  32. typename DerivedIF,
  33. typename DerivedJ,
  34. typename DerivedIM>
  35. class SelfIntersectMesh
  36. {
  37. typedef
  38. SelfIntersectMesh<
  39. Kernel,
  40. DerivedV,
  41. DerivedF,
  42. DerivedVV,
  43. DerivedFF,
  44. DerivedIF,
  45. DerivedJ,
  46. DerivedIM> Self;
  47. public:
  48. // 3D Primitives
  49. typedef CGAL::Point_3<Kernel> Point_3;
  50. typedef CGAL::Segment_3<Kernel> Segment_3;
  51. typedef CGAL::Triangle_3<Kernel> Triangle_3;
  52. typedef CGAL::Plane_3<Kernel> Plane_3;
  53. typedef CGAL::Tetrahedron_3<Kernel> Tetrahedron_3;
  54. //typedef CGAL::Polyhedron_3<Kernel> Polyhedron_3;
  55. //typedef CGAL::Nef_polyhedron_3<Kernel> Nef_polyhedron_3;
  56. // 2D Primitives
  57. typedef CGAL::Point_2<Kernel> Point_2;
  58. typedef CGAL::Segment_2<Kernel> Segment_2;
  59. typedef CGAL::Triangle_2<Kernel> Triangle_2;
  60. // 2D Constrained Delaunay Triangulation types
  61. typedef CGAL::Triangulation_vertex_base_2<Kernel> TVB_2;
  62. typedef CGAL::Constrained_triangulation_face_base_2<Kernel> CTFB_2;
  63. typedef CGAL::Triangulation_data_structure_2<TVB_2,CTFB_2> TDS_2;
  64. typedef CGAL::Exact_intersections_tag Itag;
  65. typedef CGAL::Constrained_Delaunay_triangulation_2<Kernel,TDS_2,Itag>
  66. CDT_2;
  67. typedef CGAL::Constrained_triangulation_plus_2<CDT_2> CDT_plus_2;
  68. // Axis-align boxes for all-pairs self-intersection detection
  69. typedef std::vector<Triangle_3> Triangles;
  70. typedef typename Triangles::iterator TrianglesIterator;
  71. typedef typename Triangles::const_iterator TrianglesConstIterator;
  72. typedef
  73. CGAL::Box_intersection_d::Box_with_handle_d<double,3,TrianglesIterator>
  74. Box;
  75. // Input mesh
  76. const Eigen::PlainObjectBase<DerivedV> & V;
  77. const Eigen::PlainObjectBase<DerivedF> & F;
  78. // Number of self-intersecting triangle pairs
  79. typedef typename DerivedF::Index Index;
  80. Index count;
  81. typedef std::vector<CGAL::Object> ObjectList;
  82. std::vector<ObjectList > F_objects;
  83. Triangles T;
  84. typedef std::vector<Index> IndexList;
  85. IndexList lIF;
  86. std::vector<bool> offensive;
  87. std::vector<Index> offending_index;
  88. std::vector<Index> offending;
  89. // Make a short name for the edge map's key
  90. typedef std::pair<Index,Index> EMK;
  91. // Make a short name for the type stored at each edge, the edge map's
  92. // value
  93. typedef std::vector<Index> EMV;
  94. // Make a short name for the edge map
  95. typedef std::map<EMK,EMV> EdgeMap;
  96. EdgeMap edge2faces;
  97. public:
  98. RemeshSelfIntersectionsParam params;
  99. public:
  100. // Constructs (VV,FF) a new mesh with self-intersections of (V,F)
  101. // subdivided
  102. //
  103. // See also: remesh_self_intersections.h
  104. inline SelfIntersectMesh(
  105. const Eigen::PlainObjectBase<DerivedV> & V,
  106. const Eigen::PlainObjectBase<DerivedF> & F,
  107. const RemeshSelfIntersectionsParam & params,
  108. Eigen::PlainObjectBase<DerivedVV> & VV,
  109. Eigen::PlainObjectBase<DerivedFF> & FF,
  110. Eigen::PlainObjectBase<DerivedIF> & IF,
  111. Eigen::PlainObjectBase<DerivedJ> & J,
  112. Eigen::PlainObjectBase<DerivedIM> & IM);
  113. private:
  114. // Helper function to mark a face as offensive
  115. //
  116. // Inputs:
  117. // f index of face in F
  118. inline void mark_offensive(const Index f);
  119. // Helper function to count intersections between faces
  120. //
  121. // Input:
  122. // fa index of face A in F
  123. // fb index of face B in F
  124. inline void count_intersection( const Index fa, const Index fb);
  125. // Helper function for box_intersect. Intersect two triangles A and B,
  126. // append the intersection object (point,segment,triangle) to a running
  127. // list for A and B
  128. //
  129. // Inputs:
  130. // A triangle in 3D
  131. // B triangle in 3D
  132. // fa index of A in F (and F_objects)
  133. // fb index of A in F (and F_objects)
  134. // Returns true only if A intersects B
  135. //
  136. inline bool intersect(
  137. const Triangle_3 & A,
  138. const Triangle_3 & B,
  139. const Index fa,
  140. const Index fb);
  141. // Helper function for box_intersect. In the case where A and B have
  142. // already been identified to share a vertex, then we only want to add
  143. // possible segment intersections. Assumes truly duplicate triangles are
  144. // not given as input
  145. //
  146. // Inputs:
  147. // A triangle in 3D
  148. // B triangle in 3D
  149. // fa index of A in F (and F_objects)
  150. // fb index of B in F (and F_objects)
  151. // va index of shared vertex in A (and F_objects)
  152. // vb index of shared vertex in B (and F_objects)
  153. //// Returns object of intersection (should be Segment or point)
  154. // Returns true if intersection (besides shared point)
  155. //
  156. inline bool single_shared_vertex(
  157. const Triangle_3 & A,
  158. const Triangle_3 & B,
  159. const Index fa,
  160. const Index fb,
  161. const Index va,
  162. const Index vb);
  163. // Helper handling one direction
  164. inline bool single_shared_vertex(
  165. const Triangle_3 & A,
  166. const Triangle_3 & B,
  167. const Index fa,
  168. const Index fb,
  169. const Index va);
  170. // Helper function for box_intersect. In the case where A and B have
  171. // already been identified to share two vertices, then we only want to add
  172. // a possible coplanar (Triangle) intersection. Assumes truly degenerate
  173. // facets are not givin as input.
  174. inline bool double_shared_vertex(
  175. const Triangle_3 & A,
  176. const Triangle_3 & B,
  177. const Index fa,
  178. const Index fb);
  179. public:
  180. // Callback function called during box self intersections test. Means
  181. // boxes a and b intersect. This method then checks if the triangles in
  182. // each box intersect and if so, then processes the intersections
  183. //
  184. // Inputs:
  185. // a box containing a triangle
  186. // b box containing a triangle
  187. inline void box_intersect(const Box& a, const Box& b);
  188. private:
  189. // Compute 2D delaunay triangulation of a given 3d triangle and a list of
  190. // intersection objects (points,segments,triangles). CGAL uses an affine
  191. // projection rather than an isometric projection, so we're not
  192. // guaranteed that the 2D delaunay triangulation here will be a delaunay
  193. // triangulation in 3D.
  194. //
  195. // Inputs:
  196. // A triangle in 3D
  197. // A_objects_3 updated list of intersection objects for A
  198. // Outputs:
  199. // cdt Contrained delaunay triangulation in projected 2D plane
  200. inline void projected_delaunay(
  201. const Triangle_3 & A,
  202. const ObjectList & A_objects_3,
  203. CDT_plus_2 & cdt);
  204. // Getters:
  205. public:
  206. //const IndexList& get_lIF() const{ return lIF;}
  207. static inline void box_intersect_static(
  208. SelfIntersectMesh * SIM,
  209. const Box &a,
  210. const Box &b);
  211. };
  212. }
  213. // Implementation
  214. #include "mesh_to_cgal_triangle_list.h"
  215. #include <igl/REDRUM.h>
  216. #include <igl/get_seconds.h>
  217. #include <igl/C_STR.h>
  218. #include <functional>
  219. #include <algorithm>
  220. #include <exception>
  221. #include <cassert>
  222. #include <iostream>
  223. // References:
  224. // http://minregret.googlecode.com/svn/trunk/skyline/src/extern/CGAL-3.3.1/examples/Polyhedron/polyhedron_self_intersection.cpp
  225. // http://www.cgal.org/Manual/3.9/examples/Boolean_set_operations_2/do_intersect.cpp
  226. // Q: Should we be using CGAL::Polyhedron_3?
  227. // A: No! Input is just a list of unoriented triangles. Polyhedron_3 requires
  228. // a 2-manifold.
  229. // A: But! It seems we could use CGAL::Triangulation_3. Though it won't be easy
  230. // to take advantage of functions like insert_in_facet because we want to
  231. // constrain segments. Hmmm. Actualy Triangulation_3 doesn't look right...
  232. //static void box_intersect(SelfIntersectMesh * SIM,const Box & A, const Box & B)
  233. //{
  234. // return SIM->box_intersect(A,B);
  235. //}
  236. // CGAL's box_self_intersection_d uses C-style function callbacks without
  237. // userdata. This is a leapfrog method for calling a member function. It should
  238. // be bound as if the prototype was:
  239. // static void box_intersect(const Box &a, const Box &b)
  240. // using boost:
  241. // boost::function<void(const Box &a,const Box &b)> cb
  242. // = boost::bind(&::box_intersect, this, _1,_2);
  243. //
  244. template <
  245. typename Kernel,
  246. typename DerivedV,
  247. typename DerivedF,
  248. typename DerivedVV,
  249. typename DerivedFF,
  250. typename DerivedIF,
  251. typename DerivedJ,
  252. typename DerivedIM>
  253. inline void igl::SelfIntersectMesh<
  254. Kernel,
  255. DerivedV,
  256. DerivedF,
  257. DerivedVV,
  258. DerivedFF,
  259. DerivedIF,
  260. DerivedJ,
  261. DerivedIM>::box_intersect_static(
  262. Self * SIM,
  263. const typename Self::Box &a,
  264. const typename Self::Box &b)
  265. {
  266. SIM->box_intersect(a,b);
  267. }
  268. template <
  269. typename Kernel,
  270. typename DerivedV,
  271. typename DerivedF,
  272. typename DerivedVV,
  273. typename DerivedFF,
  274. typename DerivedIF,
  275. typename DerivedJ,
  276. typename DerivedIM>
  277. inline igl::SelfIntersectMesh<
  278. Kernel,
  279. DerivedV,
  280. DerivedF,
  281. DerivedVV,
  282. DerivedFF,
  283. DerivedIF,
  284. DerivedJ,
  285. DerivedIM>::SelfIntersectMesh(
  286. const Eigen::PlainObjectBase<DerivedV> & V,
  287. const Eigen::PlainObjectBase<DerivedF> & F,
  288. const RemeshSelfIntersectionsParam & params,
  289. Eigen::PlainObjectBase<DerivedVV> & VV,
  290. Eigen::PlainObjectBase<DerivedFF> & FF,
  291. Eigen::PlainObjectBase<DerivedIF> & IF,
  292. Eigen::PlainObjectBase<DerivedJ> & J,
  293. Eigen::PlainObjectBase<DerivedIM> & IM):
  294. V(V),
  295. F(F),
  296. count(0),
  297. F_objects(F.rows()),
  298. T(),
  299. lIF(),
  300. offensive(F.rows(),false),
  301. offending_index(F.rows(),-1),
  302. offending(),
  303. edge2faces(),
  304. params(params)
  305. {
  306. using namespace std;
  307. using namespace Eigen;
  308. //const auto & tictoc = []()
  309. //{
  310. // static double t_start = igl::get_seconds();
  311. // double diff = igl::get_seconds()-t_start;
  312. // t_start += diff;
  313. // return diff;
  314. //};
  315. //tictoc();
  316. // Compute and process self intersections
  317. mesh_to_cgal_triangle_list(V,F,T);
  318. //cout<<"mesh_to_cgal_triangle_list: "<<tictoc()<<endl;
  319. // http://www.cgal.org/Manual/latest/doc_html/cgal_manual/Box_intersection_d/Chapter_main.html#Section_63.5
  320. // Create the corresponding vector of bounding boxes
  321. std::vector<Box> boxes;
  322. boxes.reserve(T.size());
  323. for (
  324. TrianglesIterator tit = T.begin();
  325. tit != T.end();
  326. ++tit)
  327. {
  328. boxes.push_back(Box(tit->bbox(), tit));
  329. }
  330. // Leapfrog callback
  331. std::function<void(const Box &a,const Box &b)> cb =
  332. std::bind(&box_intersect_static, this,
  333. // Explicitly use std namespace to avoid confusion with boost (who puts
  334. // _1 etc. in global namespace)
  335. std::placeholders::_1,
  336. std::placeholders::_2);
  337. //cout<<"boxes and bind: "<<tictoc()<<endl;
  338. // Run the self intersection algorithm with all defaults
  339. try{
  340. CGAL::box_self_intersection_d(boxes.begin(), boxes.end(),cb);
  341. }catch(int e)
  342. {
  343. // Rethrow if not IGL_FIRST_HIT_EXCEPTION
  344. if(e != IGL_FIRST_HIT_EXCEPTION)
  345. {
  346. throw e;
  347. }
  348. // Otherwise just fall through
  349. }
  350. //cout<<"box_self_intersection_d: "<<tictoc()<<endl;
  351. // Convert lIF to Eigen matrix
  352. assert(lIF.size()%2 == 0);
  353. IF.resize(lIF.size()/2,2);
  354. {
  355. Index i=0;
  356. for(
  357. typename IndexList::const_iterator ifit = lIF.begin();
  358. ifit!=lIF.end();
  359. )
  360. {
  361. IF(i,0) = (*ifit);
  362. ifit++;
  363. IF(i,1) = (*ifit);
  364. ifit++;
  365. i++;
  366. }
  367. }
  368. //cout<<"IF: "<<tictoc()<<endl;
  369. if(params.detect_only)
  370. {
  371. return;
  372. }
  373. int NF_count = 0;
  374. // list of new faces, we'll fix F later
  375. vector<
  376. typename Eigen::Matrix<typename DerivedFF::Scalar,Dynamic,Dynamic>
  377. > NF(offending.size());
  378. // list of new vertices
  379. typedef vector<Point_3> Point_3List;
  380. Point_3List NV;
  381. Index NV_count = 0;
  382. vector<CDT_plus_2> cdt(offending.size());
  383. vector<Plane_3> P(offending.size());
  384. // Use map for *all* faces
  385. map<typename CDT_plus_2::Vertex_handle,Index> v2i;
  386. // Loop over offending triangles
  387. const size_t noff = offending.size();
  388. // Unfortunately it looks like CGAL has trouble allocating memory when
  389. // multiple openmp threads are running. Crashes durring CDT...
  390. //# pragma omp parallel for if (noff>1000)
  391. for(Index o = 0;o<(Index)noff;o++)
  392. {
  393. // index in F
  394. const Index f = offending[o];
  395. {
  396. projected_delaunay(T[f],F_objects[f],cdt[o]);
  397. }
  398. // Q: Is this also delaunay in 3D?
  399. // A: No, because the projection is affine and delaunay is not affine
  400. // invariant.
  401. // Q: Then, can't we first get the 2D delaunay triangulation, then lift it
  402. // to 3D and flip any offending edges?
  403. // Plane of projection (also used by projected_delaunay)
  404. P[o] = Plane_3(T[f].vertex(0),T[f].vertex(1),T[f].vertex(2));
  405. // Build index map
  406. {
  407. Index i=0;
  408. for(
  409. typename CDT_plus_2::Finite_vertices_iterator vit = cdt[o].finite_vertices_begin();
  410. vit != cdt[o].finite_vertices_end();
  411. ++vit)
  412. {
  413. if(i<3)
  414. {
  415. //cout<<T[f].vertex(i)<<
  416. // (T[f].vertex(i) == P[o].to_3d(vit->point())?" == ":" != ")<<
  417. // P[o].to_3d(vit->point())<<endl;
  418. #ifndef NDEBUG
  419. // I want to be sure that the original corners really show up as the
  420. // original corners of the CDT. I.e. I don't trust CGAL to maintain
  421. // the order
  422. assert(T[f].vertex(i) == P[o].to_3d(vit->point()));
  423. #endif
  424. // For first three, use original index in F
  425. //# pragma omp critical
  426. v2i[vit] = F(f,i);
  427. }else
  428. {
  429. const Point_3 vit_point_3 = P[o].to_3d(vit->point());
  430. // First look up each edge's neighbors to see if exact point has
  431. // already been added (This makes everything a bit quadratic)
  432. bool found = false;
  433. for(int e = 0; e<3 && !found;e++)
  434. {
  435. // Index of F's eth edge in V
  436. Index i = F(f,(e+1)%3);
  437. Index j = F(f,(e+2)%3);
  438. // Be sure that i<j
  439. if(i>j)
  440. {
  441. swap(i,j);
  442. }
  443. assert(edge2faces.count(EMK(i,j))==1);
  444. // loop over neighbors
  445. for(
  446. typename IndexList::const_iterator nit = edge2faces[EMK(i,j)].begin();
  447. nit != edge2faces[EMK(i,j)].end() && !found;
  448. nit++)
  449. {
  450. // don't consider self
  451. if(*nit == f)
  452. {
  453. continue;
  454. }
  455. // index of neighbor in offending (to find its cdt)
  456. Index no = offending_index[*nit];
  457. // Loop over vertices of that neighbor's cdt (might not have been
  458. // processed yet, but then it's OK because it'll just be empty)
  459. for(
  460. typename CDT_plus_2::Finite_vertices_iterator uit = cdt[no].finite_vertices_begin();
  461. uit != cdt[no].finite_vertices_end() && !found;
  462. ++uit)
  463. {
  464. if(vit_point_3 == P[no].to_3d(uit->point()))
  465. {
  466. assert(v2i.count(uit) == 1);
  467. //# pragma omp critical
  468. v2i[vit] = v2i[uit];
  469. found = true;
  470. }
  471. }
  472. }
  473. }
  474. if(!found)
  475. {
  476. //# pragma omp critical
  477. {
  478. v2i[vit] = V.rows()+NV_count;
  479. NV.push_back(vit_point_3);
  480. NV_count++;
  481. }
  482. }
  483. }
  484. i++;
  485. }
  486. }
  487. {
  488. Index i = 0;
  489. // Resize to fit new number of triangles
  490. NF[o].resize(cdt[o].number_of_faces(),3);
  491. //# pragma omp atomic
  492. NF_count+=NF[o].rows();
  493. // Append new faces to NF
  494. for(
  495. typename CDT_plus_2::Finite_faces_iterator fit = cdt[o].finite_faces_begin();
  496. fit != cdt[o].finite_faces_end();
  497. ++fit)
  498. {
  499. NF[o](i,0) = v2i[fit->vertex(0)];
  500. NF[o](i,1) = v2i[fit->vertex(1)];
  501. NF[o](i,2) = v2i[fit->vertex(2)];
  502. i++;
  503. }
  504. }
  505. }
  506. //cout<<"CDT: "<<tictoc()<<" "<<t_proj_del<<endl;
  507. assert(NV_count == (Index)NV.size());
  508. // Build output
  509. #ifndef NDEBUG
  510. {
  511. Index off_count = 0;
  512. for(Index f = 0;f<F.rows();f++)
  513. {
  514. off_count+= (offensive[f]?1:0);
  515. }
  516. assert(off_count==(Index)offending.size());
  517. }
  518. #endif
  519. // Append faces
  520. FF.resize(F.rows()-offending.size()+NF_count,3);
  521. J.resize(FF.rows());
  522. // First append non-offending original faces
  523. // There's an Eigen way to do this in one line but I forget
  524. Index off = 0;
  525. for(Index f = 0;f<F.rows();f++)
  526. {
  527. if(!offensive[f])
  528. {
  529. FF.row(off) = F.row(f);
  530. J(off) = f;
  531. off++;
  532. }
  533. }
  534. assert(off == (Index)(F.rows()-offending.size()));
  535. // Now append replacement faces for offending faces
  536. for(Index o = 0;o<(Index)offending.size();o++)
  537. {
  538. FF.block(off,0,NF[o].rows(),3) = NF[o];
  539. J.block(off,0,NF[o].rows(),1).setConstant(offending[o]);
  540. off += NF[o].rows();
  541. }
  542. // Append vertices
  543. VV.resize(V.rows()+NV_count,3);
  544. VV.block(0,0,V.rows(),3) = V;
  545. {
  546. Index i = 0;
  547. for(
  548. typename Point_3List::const_iterator nvit = NV.begin();
  549. nvit != NV.end();
  550. nvit++)
  551. {
  552. for(Index d = 0;d<3;d++)
  553. {
  554. const Point_3 & p = *nvit;
  555. VV(V.rows()+i,d) = CGAL::to_double(p[d]);
  556. // This distinction does not seem necessary:
  557. //#ifdef INEXACT_CONSTRUCTION
  558. // VV(V.rows()+i,d) = CGAL::to_double(p[d]);
  559. //#else
  560. // VV(V.rows()+i,d) = CGAL::to_double(p[d].exact());
  561. //#endif
  562. }
  563. i++;
  564. }
  565. }
  566. IM.resize(VV.rows(),1);
  567. map<Point_3,Index> vv2i;
  568. // Safe to check for duplicates using double for original vertices: if
  569. // incoming reps are different then the points are unique.
  570. for(Index v = 0;v<V.rows();v++)
  571. {
  572. const Point_3 p(V(v,0),V(v,1),V(v,2));
  573. if(vv2i.count(p)==0)
  574. {
  575. vv2i[p] = v;
  576. }
  577. assert(vv2i.count(p) == 1);
  578. IM(v) = vv2i[p];
  579. }
  580. // Must check for duplicates of new vertices using exact.
  581. {
  582. Index v = V.rows();
  583. for(
  584. typename Point_3List::const_iterator nvit = NV.begin();
  585. nvit != NV.end();
  586. nvit++)
  587. {
  588. const Point_3 & p = *nvit;
  589. if(vv2i.count(p)==0)
  590. {
  591. vv2i[p] = v;
  592. }
  593. assert(vv2i.count(p) == 1);
  594. IM(v) = vv2i[p];
  595. v++;
  596. }
  597. }
  598. //cout<<"Output + dupes: "<<tictoc()<<endl;
  599. // Q: Does this give the same result as TETGEN?
  600. // A: For the cow and beast, yes.
  601. // Q: Is tetgen faster than this CGAL implementation?
  602. // A: Well, yes. But Tetgen is only solving the detection (predicates)
  603. // problem. This is also appending the intersection objects (construction).
  604. // But maybe tetgen is still faster for the detection part. For example, this
  605. // CGAL implementation on the beast takes 98 seconds but tetgen detection
  606. // takes 14 seconds
  607. }
  608. template <
  609. typename Kernel,
  610. typename DerivedV,
  611. typename DerivedF,
  612. typename DerivedVV,
  613. typename DerivedFF,
  614. typename DerivedIF,
  615. typename DerivedJ,
  616. typename DerivedIM>
  617. inline void igl::SelfIntersectMesh<
  618. Kernel,
  619. DerivedV,
  620. DerivedF,
  621. DerivedVV,
  622. DerivedFF,
  623. DerivedIF,
  624. DerivedJ,
  625. DerivedIM>::mark_offensive(const Index f)
  626. {
  627. using namespace std;
  628. lIF.push_back(f);
  629. if(!offensive[f])
  630. {
  631. offensive[f]=true;
  632. offending_index[f]=offending.size();
  633. offending.push_back(f);
  634. // Add to edge map
  635. for(Index e = 0; e<3;e++)
  636. {
  637. // Index of F's eth edge in V
  638. Index i = F(f,(e+1)%3);
  639. Index j = F(f,(e+2)%3);
  640. // Be sure that i<j
  641. if(i>j)
  642. {
  643. swap(i,j);
  644. }
  645. // Create new list if there is no entry
  646. if(edge2faces.count(EMK(i,j))==0)
  647. {
  648. edge2faces[EMK(i,j)] = EMV();
  649. }
  650. // append to list
  651. edge2faces[EMK(i,j)].push_back(f);
  652. }
  653. }
  654. }
  655. template <
  656. typename Kernel,
  657. typename DerivedV,
  658. typename DerivedF,
  659. typename DerivedVV,
  660. typename DerivedFF,
  661. typename DerivedIF,
  662. typename DerivedJ,
  663. typename DerivedIM>
  664. inline void igl::SelfIntersectMesh<
  665. Kernel,
  666. DerivedV,
  667. DerivedF,
  668. DerivedVV,
  669. DerivedFF,
  670. DerivedIF,
  671. DerivedJ,
  672. DerivedIM>::count_intersection(
  673. const Index fa,
  674. const Index fb)
  675. {
  676. mark_offensive(fa);
  677. mark_offensive(fb);
  678. this->count++;
  679. // We found the first intersection
  680. if(params.first_only && this->count >= 1)
  681. {
  682. throw IGL_FIRST_HIT_EXCEPTION;
  683. }
  684. }
  685. template <
  686. typename Kernel,
  687. typename DerivedV,
  688. typename DerivedF,
  689. typename DerivedVV,
  690. typename DerivedFF,
  691. typename DerivedIF,
  692. typename DerivedJ,
  693. typename DerivedIM>
  694. inline bool igl::SelfIntersectMesh<
  695. Kernel,
  696. DerivedV,
  697. DerivedF,
  698. DerivedVV,
  699. DerivedFF,
  700. DerivedIF,
  701. DerivedJ,
  702. DerivedIM>::intersect(
  703. const Triangle_3 & A,
  704. const Triangle_3 & B,
  705. const Index fa,
  706. const Index fb)
  707. {
  708. // Determine whether there is an intersection
  709. if(!CGAL::do_intersect(A,B))
  710. {
  711. return false;
  712. }
  713. if(!params.detect_only)
  714. {
  715. // Construct intersection
  716. CGAL::Object result = CGAL::intersection(A,B);
  717. F_objects[fa].push_back(result);
  718. F_objects[fb].push_back(result);
  719. }
  720. count_intersection(fa,fb);
  721. return true;
  722. }
  723. template <
  724. typename Kernel,
  725. typename DerivedV,
  726. typename DerivedF,
  727. typename DerivedVV,
  728. typename DerivedFF,
  729. typename DerivedIF,
  730. typename DerivedJ,
  731. typename DerivedIM>
  732. inline bool igl::SelfIntersectMesh<
  733. Kernel,
  734. DerivedV,
  735. DerivedF,
  736. DerivedVV,
  737. DerivedFF,
  738. DerivedIF,
  739. DerivedJ,
  740. DerivedIM>::single_shared_vertex(
  741. const Triangle_3 & A,
  742. const Triangle_3 & B,
  743. const Index fa,
  744. const Index fb,
  745. const Index va,
  746. const Index vb)
  747. {
  748. ////using namespace std;
  749. //CGAL::Object result = CGAL::intersection(A,B);
  750. //if(CGAL::object_cast<Segment_3 >(&result))
  751. //{
  752. // // Append to each triangle's running list
  753. // F_objects[fa].push_back(result);
  754. // F_objects[fb].push_back(result);
  755. // count_intersection(fa,fb);
  756. //}else
  757. //{
  758. // // Then intersection must be at point
  759. // // And point must be at shared vertex
  760. // assert(CGAL::object_cast<Point_3>(&result));
  761. //}
  762. if(single_shared_vertex(A,B,fa,fb,va))
  763. {
  764. return true;
  765. }
  766. return single_shared_vertex(B,A,fb,fa,vb);
  767. }
  768. template <
  769. typename Kernel,
  770. typename DerivedV,
  771. typename DerivedF,
  772. typename DerivedVV,
  773. typename DerivedFF,
  774. typename DerivedIF,
  775. typename DerivedJ,
  776. typename DerivedIM>
  777. inline bool igl::SelfIntersectMesh<
  778. Kernel,
  779. DerivedV,
  780. DerivedF,
  781. DerivedVV,
  782. DerivedFF,
  783. DerivedIF,
  784. DerivedJ,
  785. DerivedIM>::single_shared_vertex(
  786. const Triangle_3 & A,
  787. const Triangle_3 & B,
  788. const Index fa,
  789. const Index fb,
  790. const Index va)
  791. {
  792. // This was not a good idea. It will not handle coplanar triangles well.
  793. using namespace std;
  794. Segment_3 sa(
  795. A.vertex((va+1)%3),
  796. A.vertex((va+2)%3));
  797. if(CGAL::do_intersect(sa,B))
  798. {
  799. CGAL::Object result = CGAL::intersection(sa,B);
  800. if(const Point_3 * p = CGAL::object_cast<Point_3 >(&result))
  801. {
  802. if(!params.detect_only)
  803. {
  804. // Single intersection --> segment from shared point to intersection
  805. CGAL::Object seg = CGAL::make_object(Segment_3(
  806. A.vertex(va),
  807. *p));
  808. F_objects[fa].push_back(seg);
  809. F_objects[fb].push_back(seg);
  810. }
  811. count_intersection(fa,fb);
  812. return true;
  813. }else if(CGAL::object_cast<Segment_3 >(&result))
  814. {
  815. //cerr<<REDRUM("Coplanar at: "<<fa<<" & "<<fb<<" (single shared).")<<endl;
  816. // Must be coplanar
  817. if(params.detect_only)
  818. {
  819. count_intersection(fa,fb);
  820. }else
  821. {
  822. // WRONG:
  823. //// Segment intersection --> triangle from shared point to intersection
  824. //CGAL::Object tri = CGAL::make_object(Triangle_3(
  825. // A.vertex(va),
  826. // s->vertex(0),
  827. // s->vertex(1)));
  828. //F_objects[fa].push_back(tri);
  829. //F_objects[fb].push_back(tri);
  830. //count_intersection(fa,fb);
  831. // Need to do full test. Intersection could be a general poly.
  832. bool test = intersect(A,B,fa,fb);
  833. ((void)test);
  834. assert(test);
  835. }
  836. return true;
  837. }else
  838. {
  839. cerr<<REDRUM("Segment ∩ triangle neither point nor segment?")<<endl;
  840. assert(false);
  841. }
  842. }
  843. return false;
  844. }
  845. template <
  846. typename Kernel,
  847. typename DerivedV,
  848. typename DerivedF,
  849. typename DerivedVV,
  850. typename DerivedFF,
  851. typename DerivedIF,
  852. typename DerivedJ,
  853. typename DerivedIM>
  854. inline bool igl::SelfIntersectMesh<
  855. Kernel,
  856. DerivedV,
  857. DerivedF,
  858. DerivedVV,
  859. DerivedFF,
  860. DerivedIF,
  861. DerivedJ,
  862. DerivedIM>::double_shared_vertex(
  863. const Triangle_3 & A,
  864. const Triangle_3 & B,
  865. const Index fa,
  866. const Index fb)
  867. {
  868. using namespace std;
  869. // Cheaper way to do this than calling do_intersect?
  870. if(
  871. // Can only have an intersection if co-planar
  872. (A.supporting_plane() == B.supporting_plane() ||
  873. A.supporting_plane() == B.supporting_plane().opposite()) &&
  874. CGAL::do_intersect(A,B))
  875. {
  876. // Construct intersection
  877. try
  878. {
  879. // This can fail for Epick but not Epeck
  880. CGAL::Object result = CGAL::intersection(A,B);
  881. if(!result.empty())
  882. {
  883. if(CGAL::object_cast<Segment_3 >(&result))
  884. {
  885. // not coplanar
  886. return false;
  887. } else if(CGAL::object_cast<Point_3 >(&result))
  888. {
  889. // this "shouldn't" happen but does for inexact
  890. return false;
  891. } else
  892. {
  893. if(!params.detect_only)
  894. {
  895. // Triangle object
  896. F_objects[fa].push_back(result);
  897. F_objects[fb].push_back(result);
  898. }
  899. count_intersection(fa,fb);
  900. //cerr<<REDRUM("Coplanar at: "<<fa<<" & "<<fb<<" (double shared).")<<endl;
  901. return true;
  902. }
  903. }else
  904. {
  905. // CGAL::intersection is disagreeing with do_intersect
  906. return false;
  907. }
  908. }catch(...)
  909. {
  910. // This catches some cgal assertion:
  911. // CGAL error: assertion violation!
  912. // Expression : is_finite(d)
  913. // File : /opt/local/include/CGAL/GMP/Gmpq_type.h
  914. // Line : 132
  915. // Explanation:
  916. // But only if NDEBUG is not defined, otherwise there's an uncaught
  917. // "Floating point exception: 8" SIGFPE
  918. return false;
  919. }
  920. }
  921. // No intersection.
  922. return false;
  923. }
  924. template <
  925. typename Kernel,
  926. typename DerivedV,
  927. typename DerivedF,
  928. typename DerivedVV,
  929. typename DerivedFF,
  930. typename DerivedIF,
  931. typename DerivedJ,
  932. typename DerivedIM>
  933. inline void igl::SelfIntersectMesh<
  934. Kernel,
  935. DerivedV,
  936. DerivedF,
  937. DerivedVV,
  938. DerivedFF,
  939. DerivedIF,
  940. DerivedJ,
  941. DerivedIM>::box_intersect(
  942. const Box& a,
  943. const Box& b)
  944. {
  945. using namespace std;
  946. // Could we write this as a static function of:
  947. //
  948. // F.row(fa)
  949. // F.row(fb)
  950. // A
  951. // B
  952. // index in F and T
  953. Index fa = a.handle()-T.begin();
  954. Index fb = b.handle()-T.begin();
  955. const Triangle_3 & A = *a.handle();
  956. const Triangle_3 & B = *b.handle();
  957. // I'm not going to deal with degenerate triangles, though at some point we
  958. // should
  959. assert(!a.handle()->is_degenerate());
  960. assert(!b.handle()->is_degenerate());
  961. // Number of combinatorially shared vertices
  962. Index comb_shared_vertices = 0;
  963. // Number of geometrically shared vertices (*not* including combinatorially
  964. // shared)
  965. Index geo_shared_vertices = 0;
  966. // Keep track of shared vertex indices (we only handles single shared
  967. // vertices as a special case, so just need last/first/only ones)
  968. Index va=-1,vb=-1;
  969. Index ea,eb;
  970. for(ea=0;ea<3;ea++)
  971. {
  972. for(eb=0;eb<3;eb++)
  973. {
  974. if(F(fa,ea) == F(fb,eb))
  975. {
  976. comb_shared_vertices++;
  977. va = ea;
  978. vb = eb;
  979. }else if(A.vertex(ea) == B.vertex(eb))
  980. {
  981. geo_shared_vertices++;
  982. va = ea;
  983. vb = eb;
  984. }
  985. }
  986. }
  987. const Index total_shared_vertices = comb_shared_vertices + geo_shared_vertices;
  988. if(comb_shared_vertices== 3)
  989. {
  990. //// Combinatorially duplicate face, these should be removed by preprocessing
  991. //cerr<<REDRUM("Facets "<<fa<<" and "<<fb<<" are combinatorial duplicates")<<endl;
  992. goto done;
  993. }
  994. if(total_shared_vertices== 3)
  995. {
  996. //// Geometrically duplicate face, these should be removed by preprocessing
  997. //cerr<<REDRUM("Facets "<<fa<<" and "<<fb<<" are geometrical duplicates")<<endl;
  998. goto done;
  999. }
  1000. //// SPECIAL CASES ARE BROKEN FOR COPLANAR TRIANGLES
  1001. //if(total_shared_vertices > 0)
  1002. //{
  1003. // bool coplanar =
  1004. // CGAL::coplanar(A.vertex(0),A.vertex(1),A.vertex(2),B.vertex(0)) &&
  1005. // CGAL::coplanar(A.vertex(0),A.vertex(1),A.vertex(2),B.vertex(1)) &&
  1006. // CGAL::coplanar(A.vertex(0),A.vertex(1),A.vertex(2),B.vertex(2));
  1007. // if(coplanar)
  1008. // {
  1009. // cerr<<MAGENTAGIN("Facets "<<fa<<" and "<<fb<<
  1010. // " are coplanar and share vertices")<<endl;
  1011. // goto full;
  1012. // }
  1013. //}
  1014. if(total_shared_vertices == 2)
  1015. {
  1016. // Q: What about coplanar?
  1017. //
  1018. // o o
  1019. // |\ /|
  1020. // | \/ |
  1021. // | /\ |
  1022. // |/ \|
  1023. // o----o
  1024. double_shared_vertex(A,B,fa,fb);
  1025. goto done;
  1026. }
  1027. assert(total_shared_vertices<=1);
  1028. if(total_shared_vertices==1)
  1029. {
  1030. assert(va>=0 && va<3);
  1031. assert(vb>=0 && vb<3);
  1032. //#ifndef NDEBUG
  1033. // CGAL::Object result =
  1034. //#endif
  1035. single_shared_vertex(A,B,fa,fb,va,vb);
  1036. //#ifndef NDEBUG
  1037. // if(!CGAL::object_cast<Segment_3 >(&result))
  1038. // {
  1039. // const Point_3 * p = CGAL::object_cast<Point_3 >(&result);
  1040. // assert(p);
  1041. // for(int ea=0;ea<3;ea++)
  1042. // {
  1043. // for(int eb=0;eb<3;eb++)
  1044. // {
  1045. // if(F(fa,ea) == F(fb,eb))
  1046. // {
  1047. // assert(*p==A.vertex(ea));
  1048. // assert(*p==B.vertex(eb));
  1049. // }
  1050. // }
  1051. // }
  1052. // }
  1053. //#endif
  1054. }else
  1055. {
  1056. //full:
  1057. // No geometrically shared vertices, do general intersect
  1058. intersect(*a.handle(),*b.handle(),fa,fb);
  1059. }
  1060. done:
  1061. return;
  1062. }
  1063. // Compute 2D delaunay triangulation of a given 3d triangle and a list of
  1064. // intersection objects (points,segments,triangles). CGAL uses an affine
  1065. // projection rather than an isometric projection, so we're not guaranteed
  1066. // that the 2D delaunay triangulation here will be a delaunay triangulation
  1067. // in 3D.
  1068. //
  1069. // Inputs:
  1070. // A triangle in 3D
  1071. // A_objects_3 updated list of intersection objects for A
  1072. // Outputs:
  1073. // cdt Contrained delaunay triangulation in projected 2D plane
  1074. template <
  1075. typename Kernel,
  1076. typename DerivedV,
  1077. typename DerivedF,
  1078. typename DerivedVV,
  1079. typename DerivedFF,
  1080. typename DerivedIF,
  1081. typename DerivedJ,
  1082. typename DerivedIM>
  1083. inline void igl::SelfIntersectMesh<
  1084. Kernel,
  1085. DerivedV,
  1086. DerivedF,
  1087. DerivedVV,
  1088. DerivedFF,
  1089. DerivedIF,
  1090. DerivedJ,
  1091. DerivedIM>::projected_delaunay(
  1092. const Triangle_3 & A,
  1093. const ObjectList & A_objects_3,
  1094. CDT_plus_2 & cdt)
  1095. {
  1096. using namespace std;
  1097. // http://www.cgal.org/Manual/3.2/doc_html/cgal_manual/Triangulation_2/Chapter_main.html#Section_2D_Triangulations_Constrained_Plus
  1098. // Plane of triangle A
  1099. Plane_3 P(A.vertex(0),A.vertex(1),A.vertex(2));
  1100. // Insert triangle into vertices
  1101. typename CDT_plus_2::Vertex_handle corners[3];
  1102. for(Index i = 0;i<3;i++)
  1103. {
  1104. corners[i] = cdt.insert(P.to_2d(A.vertex(i)));
  1105. }
  1106. // Insert triangle edges as constraints
  1107. for(Index i = 0;i<3;i++)
  1108. {
  1109. cdt.insert_constraint( corners[(i+1)%3], corners[(i+2)%3]);
  1110. }
  1111. // Insert constraints for intersection objects
  1112. for( const auto & obj : A_objects_3)
  1113. {
  1114. if(const Segment_3 *iseg = CGAL::object_cast<Segment_3 >(&obj))
  1115. {
  1116. // Add segment constraint
  1117. cdt.insert_constraint(P.to_2d(iseg->vertex(0)),P.to_2d(iseg->vertex(1)));
  1118. }else if(const Point_3 *ipoint = CGAL::object_cast<Point_3 >(&obj))
  1119. {
  1120. // Add point
  1121. cdt.insert(P.to_2d(*ipoint));
  1122. } else if(const Triangle_3 *itri = CGAL::object_cast<Triangle_3 >(&obj))
  1123. {
  1124. // Add 3 segment constraints
  1125. cdt.insert_constraint(P.to_2d(itri->vertex(0)),P.to_2d(itri->vertex(1)));
  1126. cdt.insert_constraint(P.to_2d(itri->vertex(1)),P.to_2d(itri->vertex(2)));
  1127. cdt.insert_constraint(P.to_2d(itri->vertex(2)),P.to_2d(itri->vertex(0)));
  1128. } else if(const std::vector<Point_3 > *polyp =
  1129. CGAL::object_cast< std::vector<Point_3 > >(&obj))
  1130. {
  1131. //cerr<<REDRUM("Poly...")<<endl;
  1132. const std::vector<Point_3 > & poly = *polyp;
  1133. const Index m = poly.size();
  1134. assert(m>=2);
  1135. for(Index p = 0;p<m;p++)
  1136. {
  1137. const Index np = (p+1)%m;
  1138. cdt.insert_constraint(P.to_2d(poly[p]),P.to_2d(poly[np]));
  1139. }
  1140. }else
  1141. {
  1142. cerr<<REDRUM("What is this object?!")<<endl;
  1143. assert(false);
  1144. }
  1145. }
  1146. }
  1147. #endif