SelfIntersectMesh.h 32 KB

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