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