SelfIntersectMesh.h 33 KB

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