SelfIntersectMesh.h 35 KB

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