SelfIntersectMesh.h 27 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2014 Alec Jacobson <[email protected]>
  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_COPYLEFT_CGAL_SELFINTERSECTMESH_H
  9. #define IGL_COPYLEFT_CGAL_SELFINTERSECTMESH_H
  10. #include "CGAL_includes.hpp"
  11. #include "RemeshSelfIntersectionsParam.h"
  12. #include "../../unique.h"
  13. #include "../../default_num_threads.h"
  14. #include <Eigen/Dense>
  15. #include <list>
  16. #include <map>
  17. #include <vector>
  18. #include <thread>
  19. #include <mutex>
  20. #include <cstdio>
  21. //#define IGL_SELFINTERSECTMESH_TIMING
  22. #ifndef IGL_FIRST_HIT_EXCEPTION
  23. #define IGL_FIRST_HIT_EXCEPTION 10
  24. #endif
  25. // The easiest way to keep track of everything is to use a class
  26. namespace igl
  27. {
  28. namespace copyleft
  29. {
  30. namespace cgal
  31. {
  32. /// Class for computing the self-intersections of a mesh
  33. ///
  34. /// @tparam Kernel is a CGAL kernel like:
  35. /// CGAL::Exact_predicates_inexact_constructions_kernel
  36. /// or
  37. /// CGAL::Exact_predicates_exact_constructions_kernel
  38. template <
  39. typename Kernel,
  40. typename DerivedV,
  41. typename DerivedF,
  42. typename DerivedVV,
  43. typename DerivedFF,
  44. typename DerivedIF,
  45. typename DerivedJ,
  46. typename DerivedIM>
  47. class SelfIntersectMesh
  48. {
  49. typedef
  50. SelfIntersectMesh<
  51. Kernel,
  52. DerivedV,
  53. DerivedF,
  54. DerivedVV,
  55. DerivedFF,
  56. DerivedIF,
  57. DerivedJ,
  58. DerivedIM> Self;
  59. public:
  60. // 3D Primitives
  61. typedef CGAL::Point_3<Kernel> Point_3;
  62. typedef CGAL::Segment_3<Kernel> Segment_3;
  63. typedef CGAL::Triangle_3<Kernel> Triangle_3;
  64. typedef CGAL::Plane_3<Kernel> Plane_3;
  65. typedef CGAL::Tetrahedron_3<Kernel> Tetrahedron_3;
  66. // 2D Primitives
  67. typedef CGAL::Point_2<Kernel> Point_2;
  68. typedef CGAL::Segment_2<Kernel> Segment_2;
  69. typedef CGAL::Triangle_2<Kernel> Triangle_2;
  70. // 2D Constrained Delaunay Triangulation types
  71. typedef CGAL::Exact_intersections_tag Itag;
  72. // Axis-align boxes for all-pairs self-intersection detection
  73. typedef std::vector<Triangle_3> Triangles;
  74. typedef typename Triangles::iterator TrianglesIterator;
  75. typedef typename Triangles::const_iterator TrianglesConstIterator;
  76. typedef
  77. CGAL::Box_intersection_d::Box_with_handle_d<double,3,TrianglesIterator>
  78. Box;
  79. // Input mesh
  80. const Eigen::MatrixBase<DerivedV> & V;
  81. const Eigen::MatrixBase<DerivedF> & F;
  82. // Number of self-intersecting triangle pairs
  83. typedef typename DerivedF::Index Index;
  84. Index count;
  85. typedef std::vector<std::pair<Index, CGAL::Object>> ObjectList;
  86. // Using a vector here makes this **not** output sensitive
  87. Triangles T;
  88. typedef std::vector<Index> IndexList;
  89. IndexList lIF;
  90. // #F-long list of faces with intersections mapping to the order in
  91. // which they were first found
  92. std::map<Index,ObjectList> offending;
  93. // Make a short name for the edge map's key
  94. typedef std::pair<Index,Index> EMK;
  95. // Make a short name for the type stored at each edge, the edge map's
  96. // value
  97. typedef std::vector<Index> EMV;
  98. // Make a short name for the edge map
  99. typedef std::map<EMK,EMV> EdgeMap;
  100. // Maps edges of offending faces to all incident offending faces
  101. std::vector<std::pair<TrianglesIterator, TrianglesIterator> >
  102. candidate_triangle_pairs;
  103. public:
  104. RemeshSelfIntersectionsParam params;
  105. public:
  106. /// Constructs (VV,FF) a new mesh with self-intersections of (V,F)
  107. /// subdivided
  108. ///
  109. /// @param[in] V #V by 3 list of vertex positions
  110. /// @param[in] F #F by 3 list of triangle indices into V
  111. /// @param[in] params parameters
  112. /// @param[out] VV #VV by 3 list of vertex positions
  113. /// @param[out] FF #FF by 3 list of triangle indices into VV
  114. /// @param[out] IF #IF by 2 list of edge indices into VV
  115. /// @param[out] J #F list of indices into FF of birth parents
  116. /// @param[out] IM #VV list of indices into V of birth parents
  117. ///
  118. ///
  119. /// \see remesh_self_intersections.h
  120. inline SelfIntersectMesh(
  121. const Eigen::MatrixBase<DerivedV> & V,
  122. const Eigen::MatrixBase<DerivedF> & F,
  123. const RemeshSelfIntersectionsParam & params,
  124. Eigen::PlainObjectBase<DerivedVV> & VV,
  125. Eigen::PlainObjectBase<DerivedFF> & FF,
  126. Eigen::PlainObjectBase<DerivedIF> & IF,
  127. Eigen::PlainObjectBase<DerivedJ> & J,
  128. Eigen::PlainObjectBase<DerivedIM> & IM);
  129. private:
  130. /// Helper function to mark a face as offensive
  131. ///
  132. /// @param[in] f index of face in F
  133. inline void mark_offensive(const Index f);
  134. /// Helper function to count intersections between faces
  135. ///
  136. /// @param[in] fa index of face A in F
  137. /// @param[in] fb index of face B in F
  138. inline void count_intersection( const Index fa, const Index fb);
  139. /// Helper function for box_intersect. Intersect two triangles A and B,
  140. /// append the intersection object (point,segment,triangle) to a running
  141. /// list for A and B
  142. ///
  143. /// @param[in] A triangle in 3D
  144. /// @param[in] B triangle in 3D
  145. /// @param[in] fa index of A in F (and key into offending)
  146. /// @param[in] fb index of B in F (and key into offending)
  147. /// @return true only if A intersects B
  148. ///
  149. inline bool intersect(
  150. const Triangle_3 & A,
  151. const Triangle_3 & B,
  152. const Index fa,
  153. const Index fb);
  154. /// Helper function for box_intersect. In the case where A and B have
  155. /// already been identified to share a vertex, then we only want to
  156. /// add possible segment intersections. Assumes truly duplicate
  157. /// triangles are not given as input
  158. ///
  159. /// @param[in] A triangle in 3D
  160. /// @param[in] B triangle in 3D
  161. /// @param[in] fa index of A in F (and key into offending)
  162. /// @param[in] fb index of B in F (and key into offending)
  163. /// @param[in] va index of shared vertex in A (and key into offending)
  164. /// @param[in] vb index of shared vertex in B (and key into offending)
  165. /// @return true if intersection (besides shared point)
  166. ///
  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. const Index vb);
  174. //// Helper handling one direction
  175. ///
  176. /// @param[in] A triangle in 3D
  177. /// @param[in] B triangle in 3D
  178. /// @param[in] fa index of A in F (and key into offending)
  179. /// @param[in] fb index of B in F (and key into offending)
  180. /// @param[in] va index of shared vertex in A (and key into offending)
  181. /// @return true if intersection (besides shared point)
  182. inline bool single_shared_vertex(
  183. const Triangle_3 & A,
  184. const Triangle_3 & B,
  185. const Index fa,
  186. const Index fb,
  187. const Index va);
  188. /// Helper function for box_intersect. In the case where A and B have
  189. /// already been identified to share two vertices, then we only want
  190. /// to add a possible coplanar (Triangle) intersection. Assumes truly
  191. /// degenerate facets are not givin as input.
  192. ///
  193. /// @param[in] A triangle in 3D
  194. /// @param[in] B triangle in 3D
  195. /// @param[in] fa index of A in F (and key into offending)
  196. /// @param[in] fb index of B in F (and key into offending)
  197. /// @param[in] shared list of pairs of indices of shared vertices
  198. /// @return true if intersection (besides shared point)
  199. inline bool double_shared_vertex(
  200. const Triangle_3 & A,
  201. const Triangle_3 & B,
  202. const Index fa,
  203. const Index fb,
  204. const std::vector<std::pair<Index,Index> > shared);
  205. public:
  206. /// Callback function called during box self intersections test. Means
  207. /// boxes a and b intersect. This method then checks if the triangles
  208. /// in each box intersect and if so, then processes the intersections
  209. ///
  210. /// @param[in] a box containing a triangle
  211. /// @param[in] b box containing a triangle
  212. inline void box_intersect(const Box& a, const Box& b);
  213. /// Process all of the intersecting boxes
  214. inline void process_intersecting_boxes();
  215. public:
  216. // Getters:
  217. //const IndexList& get_lIF() const{ return lIF;}
  218. /// Static function that captures a SelfIntersectMesh instance to pass
  219. /// to cgal.
  220. /// @param[in] SIM pointer to SelfIntersectMesh instance
  221. /// @param[in] a box containing a triangle
  222. /// @param[in] b box containing a triangle
  223. static inline void box_intersect_static(
  224. SelfIntersectMesh * SIM,
  225. const Box &a,
  226. const Box &b);
  227. private:
  228. std::mutex m_offending_lock;
  229. };
  230. }
  231. }
  232. }
  233. // Implementation
  234. #include "mesh_to_cgal_triangle_list.h"
  235. #include "remesh_intersections.h"
  236. #include "../../get_seconds.h"
  237. #include "../../C_STR.h"
  238. #include <functional>
  239. #include <algorithm>
  240. #include <exception>
  241. #include <cassert>
  242. // References:
  243. // http://minregret.googlecode.com/svn/trunk/skyline/src/extern/CGAL-3.3.1/examples/Polyhedron/polyhedron_self_intersection.cpp
  244. // http://www.cgal.org/Manual/3.9/examples/Boolean_set_operations_2/do_intersect.cpp
  245. // Q: Should we be using CGAL::Polyhedron_3?
  246. // A: No! Input is just a list of unoriented triangles. Polyhedron_3 requires
  247. // a 2-manifold.
  248. // A: But! It seems we could use CGAL::Triangulation_3. Though it won't be easy
  249. // to take advantage of functions like insert_in_facet because we want to
  250. // constrain segments. Hmmm. Actually Triangulation_3 doesn't look right...
  251. // CGAL's box_self_intersection_d uses C-style function callbacks without
  252. // userdata. This is a leapfrog method for calling a member function. It should
  253. // be bound as if the prototype was:
  254. // static void box_intersect(const Box &a, const Box &b)
  255. // using boost:
  256. // boost::function<void(const Box &a,const Box &b)> cb
  257. // = boost::bind(&::box_intersect, this, _1,_2);
  258. //
  259. template <
  260. typename Kernel,
  261. typename DerivedV,
  262. typename DerivedF,
  263. typename DerivedVV,
  264. typename DerivedFF,
  265. typename DerivedIF,
  266. typename DerivedJ,
  267. typename DerivedIM>
  268. inline void igl::copyleft::cgal::SelfIntersectMesh<
  269. Kernel,
  270. DerivedV,
  271. DerivedF,
  272. DerivedVV,
  273. DerivedFF,
  274. DerivedIF,
  275. DerivedJ,
  276. DerivedIM>::box_intersect_static(
  277. Self * SIM,
  278. const typename Self::Box &a,
  279. const typename Self::Box &b)
  280. {
  281. SIM->box_intersect(a,b);
  282. }
  283. template <
  284. typename Kernel,
  285. typename DerivedV,
  286. typename DerivedF,
  287. typename DerivedVV,
  288. typename DerivedFF,
  289. typename DerivedIF,
  290. typename DerivedJ,
  291. typename DerivedIM>
  292. inline igl::copyleft::cgal::SelfIntersectMesh<
  293. Kernel,
  294. DerivedV,
  295. DerivedF,
  296. DerivedVV,
  297. DerivedFF,
  298. DerivedIF,
  299. DerivedJ,
  300. DerivedIM>::SelfIntersectMesh(
  301. const Eigen::MatrixBase<DerivedV> & V,
  302. const Eigen::MatrixBase<DerivedF> & F,
  303. const RemeshSelfIntersectionsParam & params,
  304. Eigen::PlainObjectBase<DerivedVV> & VV,
  305. Eigen::PlainObjectBase<DerivedFF> & FF,
  306. Eigen::PlainObjectBase<DerivedIF> & IF,
  307. Eigen::PlainObjectBase<DerivedJ> & J,
  308. Eigen::PlainObjectBase<DerivedIM> & IM):
  309. V(V),
  310. F(F),
  311. count(0),
  312. T(),
  313. lIF(),
  314. offending(),
  315. params(params)
  316. {
  317. #ifdef IGL_SELFINTERSECTMESH_TIMING
  318. const auto & tictoc = []() -> double
  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. const auto log_time = [&](const std::string& label) -> void{
  326. std::printf("%50s: %0.5lf\n",
  327. C_STR("SelfIntersectMesh." << label),tictoc());
  328. };
  329. tictoc();
  330. #endif
  331. // Compute and process self intersections
  332. mesh_to_cgal_triangle_list(V,F,T);
  333. #ifdef IGL_SELFINTERSECTMESH_TIMING
  334. log_time("convert_to_triangle_list");
  335. #endif
  336. // http://www.cgal.org/Manual/latest/doc_html/cgal_manual/Box_intersection_d/Chapter_main.html#Section_63.5
  337. // Create the corresponding vector of bounding boxes
  338. std::vector<Box> boxes;
  339. boxes.reserve(T.size());
  340. for (
  341. TrianglesIterator tit = T.begin();
  342. tit != T.end();
  343. ++tit)
  344. {
  345. if (!tit->is_degenerate())
  346. {
  347. boxes.push_back(Box(tit->bbox(), tit));
  348. }
  349. }
  350. // Leapfrog callback
  351. std::function<void(const Box &a,const Box &b)> cb =
  352. std::bind(&box_intersect_static, this,
  353. // Explicitly use std namespace to avoid confusion with boost (who puts
  354. // _1 etc. in global namespace)
  355. std::placeholders::_1,
  356. std::placeholders::_2);
  357. #ifdef IGL_SELFINTERSECTMESH_TIMING
  358. log_time("box_and_bind");
  359. #endif
  360. // Run the self intersection algorithm with given cutoff size
  361. CGAL::box_self_intersection_d(boxes.begin(), boxes.end(),cb,std::ptrdiff_t(params.cutoff));
  362. #ifdef IGL_SELFINTERSECTMESH_TIMING
  363. log_time("box_intersection_d");
  364. #endif
  365. try{
  366. process_intersecting_boxes();
  367. }catch(int e)
  368. {
  369. // Rethrow if not IGL_FIRST_HIT_EXCEPTION
  370. if(e != IGL_FIRST_HIT_EXCEPTION)
  371. {
  372. throw e;
  373. }
  374. // Otherwise just fall through
  375. }
  376. #ifdef IGL_SELFINTERSECTMESH_TIMING
  377. log_time("resolve_intersection");
  378. #endif
  379. // Convert lIF to Eigen matrix
  380. assert(lIF.size()%2 == 0);
  381. IF.resize(lIF.size()/2,2);
  382. {
  383. Index i=0;
  384. for(
  385. typename IndexList::const_iterator ifit = lIF.begin();
  386. ifit!=lIF.end();
  387. )
  388. {
  389. IF(i,0) = (*ifit);
  390. ifit++;
  391. IF(i,1) = (*ifit);
  392. ifit++;
  393. i++;
  394. }
  395. }
  396. #ifdef IGL_SELFINTERSECTMESH_TIMING
  397. log_time("store_intersecting_face_pairs");
  398. #endif
  399. if(params.detect_only)
  400. {
  401. return;
  402. }
  403. remesh_intersections(
  404. V,F,T,offending,
  405. params.stitch_all,params.slow_and_more_precise_rounding,VV,FF,J,IM);
  406. #ifdef IGL_SELFINTERSECTMESH_TIMING
  407. log_time("remesh_intersection");
  408. #endif
  409. }
  410. template <
  411. typename Kernel,
  412. typename DerivedV,
  413. typename DerivedF,
  414. typename DerivedVV,
  415. typename DerivedFF,
  416. typename DerivedIF,
  417. typename DerivedJ,
  418. typename DerivedIM>
  419. inline void igl::copyleft::cgal::SelfIntersectMesh<
  420. Kernel,
  421. DerivedV,
  422. DerivedF,
  423. DerivedVV,
  424. DerivedFF,
  425. DerivedIF,
  426. DerivedJ,
  427. DerivedIM>::mark_offensive(const Index f)
  428. {
  429. lIF.push_back(f);
  430. if(offending.count(f) == 0)
  431. {
  432. // first time marking, initialize with new id and empty list
  433. offending[f] = {};
  434. }
  435. }
  436. template <
  437. typename Kernel,
  438. typename DerivedV,
  439. typename DerivedF,
  440. typename DerivedVV,
  441. typename DerivedFF,
  442. typename DerivedIF,
  443. typename DerivedJ,
  444. typename DerivedIM>
  445. inline void igl::copyleft::cgal::SelfIntersectMesh<
  446. Kernel,
  447. DerivedV,
  448. DerivedF,
  449. DerivedVV,
  450. DerivedFF,
  451. DerivedIF,
  452. DerivedJ,
  453. DerivedIM>::count_intersection(
  454. const Index fa,
  455. const Index fb)
  456. {
  457. std::lock_guard<std::mutex> guard(m_offending_lock);
  458. mark_offensive(fa);
  459. mark_offensive(fb);
  460. this->count++;
  461. // We found the first intersection
  462. if(params.first_only && this->count >= 1)
  463. {
  464. throw IGL_FIRST_HIT_EXCEPTION;
  465. }
  466. }
  467. template <
  468. typename Kernel,
  469. typename DerivedV,
  470. typename DerivedF,
  471. typename DerivedVV,
  472. typename DerivedFF,
  473. typename DerivedIF,
  474. typename DerivedJ,
  475. typename DerivedIM>
  476. inline bool igl::copyleft::cgal::SelfIntersectMesh<
  477. Kernel,
  478. DerivedV,
  479. DerivedF,
  480. DerivedVV,
  481. DerivedFF,
  482. DerivedIF,
  483. DerivedJ,
  484. DerivedIM>::intersect(
  485. const Triangle_3 & A,
  486. const Triangle_3 & B,
  487. const Index fa,
  488. const Index fb)
  489. {
  490. // Determine whether there is an intersection
  491. if(!CGAL::do_intersect(A,B))
  492. {
  493. return false;
  494. }
  495. count_intersection(fa,fb);
  496. if(!params.detect_only)
  497. {
  498. // Construct intersection
  499. CGAL::Object result = CGAL::intersection(A,B);
  500. // Could avoid this mutex if `offending` was per-thread and passed as input
  501. // reference.
  502. std::lock_guard<std::mutex> guard(m_offending_lock);
  503. offending[fa].push_back({fb, result});
  504. offending[fb].push_back({fa, result});
  505. }
  506. return true;
  507. }
  508. template <
  509. typename Kernel,
  510. typename DerivedV,
  511. typename DerivedF,
  512. typename DerivedVV,
  513. typename DerivedFF,
  514. typename DerivedIF,
  515. typename DerivedJ,
  516. typename DerivedIM>
  517. inline bool igl::copyleft::cgal::SelfIntersectMesh<
  518. Kernel,
  519. DerivedV,
  520. DerivedF,
  521. DerivedVV,
  522. DerivedFF,
  523. DerivedIF,
  524. DerivedJ,
  525. DerivedIM>::single_shared_vertex(
  526. const Triangle_3 & A,
  527. const Triangle_3 & B,
  528. const Index fa,
  529. const Index fb,
  530. const Index va,
  531. const Index vb)
  532. {
  533. if(single_shared_vertex(A,B,fa,fb,va))
  534. {
  535. return true;
  536. }
  537. return single_shared_vertex(B,A,fb,fa,vb);
  538. }
  539. template <
  540. typename Kernel,
  541. typename DerivedV,
  542. typename DerivedF,
  543. typename DerivedVV,
  544. typename DerivedFF,
  545. typename DerivedIF,
  546. typename DerivedJ,
  547. typename DerivedIM>
  548. inline bool igl::copyleft::cgal::SelfIntersectMesh<
  549. Kernel,
  550. DerivedV,
  551. DerivedF,
  552. DerivedVV,
  553. DerivedFF,
  554. DerivedIF,
  555. DerivedJ,
  556. DerivedIM>::single_shared_vertex(
  557. const Triangle_3 & A,
  558. const Triangle_3 & B,
  559. const Index fa,
  560. const Index fb,
  561. const Index va)
  562. {
  563. // This was not a good idea. It will not handle coplanar triangles well.
  564. Segment_3 sa(
  565. A.vertex((va+1)%3),
  566. A.vertex((va+2)%3));
  567. if(CGAL::do_intersect(sa,B))
  568. {
  569. // can't put count_intersection(fa,fb) here since we use intersect below
  570. // and then it will be counted twice.
  571. if(params.detect_only)
  572. {
  573. count_intersection(fa,fb);
  574. return true;
  575. }
  576. CGAL::Object result = CGAL::intersection(sa,B);
  577. if(const Point_3 * p = CGAL::object_cast<Point_3 >(&result))
  578. {
  579. // Single intersection --> segment from shared point to intersection
  580. CGAL::Object seg = CGAL::make_object(Segment_3(
  581. A.vertex(va),
  582. *p));
  583. count_intersection(fa,fb);
  584. std::lock_guard<std::mutex> guard(m_offending_lock);
  585. offending[fa].push_back({fb, seg});
  586. offending[fb].push_back({fa, seg});
  587. return true;
  588. }else if(CGAL::object_cast<Segment_3 >(&result))
  589. {
  590. // Need to do full test. Intersection could be a general poly.
  591. bool test = intersect(A,B,fa,fb);
  592. ((void)test);
  593. assert(test && "intersect should agree with do_intersect");
  594. return true;
  595. }else
  596. {
  597. // Should never happen.
  598. assert(false && "Segment ∩ triangle neither point nor segment?");
  599. }
  600. }
  601. return false;
  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 bool igl::copyleft::cgal::SelfIntersectMesh<
  613. Kernel,
  614. DerivedV,
  615. DerivedF,
  616. DerivedVV,
  617. DerivedFF,
  618. DerivedIF,
  619. DerivedJ,
  620. DerivedIM>::double_shared_vertex(
  621. const Triangle_3 & A,
  622. const Triangle_3 & B,
  623. const Index fa,
  624. const Index fb,
  625. const std::vector<std::pair<Index,Index> > shared)
  626. {
  627. auto opposite_vertex = [](const Index a0, const Index a1) {
  628. // get opposite index of A
  629. int a2=-1;
  630. for(int c=0;c<3;++c)
  631. if(c!=a0 && c!=a1) {
  632. a2 = c;
  633. break;
  634. }
  635. assert(a2 != -1);
  636. return a2;
  637. };
  638. // must be co-planar
  639. Index a2 = opposite_vertex(shared[0].first, shared[1].first);
  640. if (! B.supporting_plane().has_on(A.vertex(a2)))
  641. return false;
  642. Index b2 = opposite_vertex(shared[0].second, shared[1].second);
  643. if (int(CGAL::coplanar_orientation(A.vertex(shared[0].first), A.vertex(shared[1].first), A.vertex(a2))) *
  644. int(CGAL::coplanar_orientation(B.vertex(shared[0].second), B.vertex(shared[1].second), B.vertex(b2))) < 0)
  645. // There is certainly no self intersection as the non-shared triangle vertices lie on opposite sides of the shared edge.
  646. return false;
  647. // Since A and B are non-degenerate the intersection must be a polygon
  648. // (triangle). Either
  649. // - the vertex of A (B) opposite the shared edge of lies on B (A), or
  650. // - an edge of A intersects and edge of B without sharing a vertex
  651. //
  652. // Determine if the vertex opposite edge (a0,a1) in triangle A lies in
  653. // (intersects) triangle B
  654. const auto & opposite_point_inside = [](
  655. const Triangle_3 & A, const Index a2, const Triangle_3 & B)
  656. -> bool
  657. {
  658. return CGAL::do_intersect(A.vertex(a2),B);
  659. };
  660. // Determine if edge opposite vertex va in triangle A intersects edge
  661. // opposite vertex vb in triangle B.
  662. const auto & opposite_edges_intersect = [](
  663. const Triangle_3 & A, const Index va,
  664. const Triangle_3 & B, const Index vb) -> bool
  665. {
  666. Segment_3 sa( A.vertex((va+1)%3), A.vertex((va+2)%3));
  667. Segment_3 sb( B.vertex((vb+1)%3), B.vertex((vb+2)%3));
  668. bool ret = CGAL::do_intersect(sa,sb);
  669. return ret;
  670. };
  671. if(
  672. !opposite_point_inside(A,a2,B) &&
  673. !opposite_point_inside(B,b2,A) &&
  674. !opposite_edges_intersect(A,shared[0].first,B,shared[1].second) &&
  675. !opposite_edges_intersect(A,shared[1].first,B,shared[0].second))
  676. {
  677. return false;
  678. }
  679. // there is an intersection indeed
  680. count_intersection(fa,fb);
  681. if(params.detect_only)
  682. {
  683. return true;
  684. }
  685. // Construct intersection
  686. try
  687. {
  688. // This can fail for Epick but not Epeck
  689. CGAL::Object result = CGAL::intersection(A,B);
  690. if(!result.empty())
  691. {
  692. if(CGAL::object_cast<Segment_3 >(&result))
  693. {
  694. // not coplanar
  695. assert(false &&
  696. "Co-planar non-degenerate triangles should intersect over triangle");
  697. return false;
  698. } else if(CGAL::object_cast<Point_3 >(&result))
  699. {
  700. // this "shouldn't" happen but does for inexact
  701. assert(false &&
  702. "Co-planar non-degenerate triangles should intersect over triangle");
  703. return false;
  704. } else
  705. {
  706. // Triangle object
  707. std::lock_guard<std::mutex> guard(m_offending_lock);
  708. offending[fa].push_back({fb, result});
  709. offending[fb].push_back({fa, result});
  710. return true;
  711. }
  712. }else
  713. {
  714. // CGAL::intersection is disagreeing with do_intersect
  715. assert(false && "CGAL::intersection should agree with predicate tests");
  716. return false;
  717. }
  718. }catch(...)
  719. {
  720. // This catches some cgal assertion:
  721. // CGAL error: assertion violation!
  722. // Expression : is_finite(d)
  723. // File : /opt/local/include/CGAL/GMP/Gmpq_type.h
  724. // Line : 132
  725. // Explanation:
  726. // But only if NDEBUG is not defined, otherwise there's an uncaught
  727. // "Floating point exception: 8" SIGFPE
  728. return false;
  729. }
  730. // No intersection.
  731. return false;
  732. }
  733. template <
  734. typename Kernel,
  735. typename DerivedV,
  736. typename DerivedF,
  737. typename DerivedVV,
  738. typename DerivedFF,
  739. typename DerivedIF,
  740. typename DerivedJ,
  741. typename DerivedIM>
  742. inline void igl::copyleft::cgal::SelfIntersectMesh<
  743. Kernel,
  744. DerivedV,
  745. DerivedF,
  746. DerivedVV,
  747. DerivedFF,
  748. DerivedIF,
  749. DerivedJ,
  750. DerivedIM>::box_intersect(
  751. const Box& a,
  752. const Box& b)
  753. {
  754. candidate_triangle_pairs.push_back({a.handle(), b.handle()});
  755. }
  756. template <
  757. typename Kernel,
  758. typename DerivedV,
  759. typename DerivedF,
  760. typename DerivedVV,
  761. typename DerivedFF,
  762. typename DerivedIF,
  763. typename DerivedJ,
  764. typename DerivedIM>
  765. inline void igl::copyleft::cgal::SelfIntersectMesh<
  766. Kernel,
  767. DerivedV,
  768. DerivedF,
  769. DerivedVV,
  770. DerivedFF,
  771. DerivedIF,
  772. DerivedJ,
  773. DerivedIM>::process_intersecting_boxes()
  774. {
  775. std::mutex exception_mutex;
  776. bool exception_fired = false;
  777. int exception = -1;
  778. // Eventually switching to igl::parallel_for would be good, but currently
  779. // igl::parallel_for does not provide a way to catch exceptions fired on a
  780. // spawned thread _outside_ of its loop-chunk which is the mechanism used here
  781. // to bail out early when `first_only=true` to avoid
  782. // O(#candidate_triangle_pairs) behavior.
  783. auto process_chunk = [&]( const size_t first, const size_t last) -> void
  784. {
  785. try
  786. {
  787. assert(last >= first);
  788. for (size_t i=first; i<last; i++)
  789. {
  790. if(exception_fired) return;
  791. Index fa=T.size(), fb=T.size();
  792. {
  793. const auto& tri_pair = candidate_triangle_pairs[i];
  794. fa = tri_pair.first - T.begin();
  795. fb = tri_pair.second - T.begin();
  796. }
  797. assert(fa < T.size());
  798. assert(fb < T.size());
  799. if(exception_fired) return;
  800. const Triangle_3& A = T[fa];
  801. const Triangle_3& B = T[fb];
  802. // Number of combinatorially shared vertices
  803. Index comb_shared_vertices = 0;
  804. // Number of geometrically shared vertices (*not* including
  805. // combinatorially shared)
  806. Index geo_shared_vertices = 0;
  807. // Keep track of shared vertex indices
  808. std::vector<std::pair<Index,Index> > shared;
  809. Index ea,eb;
  810. for(ea=0;ea<3;ea++)
  811. {
  812. for(eb=0;eb<3;eb++)
  813. {
  814. if(F(fa,ea) == F(fb,eb))
  815. {
  816. comb_shared_vertices++;
  817. shared.emplace_back(ea,eb);
  818. }else if(A.vertex(ea) == B.vertex(eb))
  819. {
  820. geo_shared_vertices++;
  821. shared.emplace_back(ea,eb);
  822. }
  823. }
  824. }
  825. const Index total_shared_vertices =
  826. comb_shared_vertices + geo_shared_vertices;
  827. if(exception_fired) return;
  828. if(comb_shared_vertices== 3)
  829. {
  830. assert(shared.size() == 3);
  831. // Combinatorially duplicate face, these should be removed by
  832. // preprocessing
  833. continue;
  834. }
  835. if(total_shared_vertices== 3)
  836. {
  837. assert(shared.size() == 3);
  838. // Geometrically duplicate face, these should be removed by
  839. // preprocessing
  840. continue;
  841. }
  842. if(total_shared_vertices == 2)
  843. {
  844. assert(shared.size() == 2);
  845. // Q: What about coplanar?
  846. //
  847. // o o
  848. // |\ /|
  849. // | \/ |
  850. // | /\ |
  851. // |/ \|
  852. // o----o
  853. double_shared_vertex(A,B,fa,fb,shared);
  854. continue;
  855. }
  856. assert(total_shared_vertices<=1);
  857. if(total_shared_vertices==1)
  858. {
  859. single_shared_vertex(A,B,fa,fb,shared[0].first,shared[0].second);
  860. }else
  861. {
  862. intersect(A,B,fa,fb);
  863. }
  864. }
  865. }catch(int e)
  866. {
  867. std::lock_guard<std::mutex> exception_lock(exception_mutex);
  868. exception_fired = true;
  869. exception = e;
  870. }
  871. };
  872. const size_t num_threads = default_num_threads();
  873. assert(num_threads > 0);
  874. const size_t num_pairs = candidate_triangle_pairs.size();
  875. const size_t chunk_size = num_pairs / num_threads;
  876. std::vector<std::thread> threads;
  877. for (size_t i=0; i<num_threads-1; i++)
  878. {
  879. threads.emplace_back(process_chunk, i*chunk_size, (i+1)*chunk_size);
  880. }
  881. // Do some work in the master thread.
  882. process_chunk((num_threads-1)*chunk_size, num_pairs);
  883. for (auto& t : threads)
  884. {
  885. if (t.joinable()) t.join();
  886. }
  887. if(exception_fired) throw exception;
  888. //process_chunk(0, candidate_triangle_pairs.size());
  889. }
  890. #endif