Prusa Slicer 2.6.0
Loading...
Searching...
No Matches
igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM > Class Template Reference

#include <src/libigl/igl/copyleft/cgal/SelfIntersectMesh.h>

+ Collaboration diagram for igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >:

Public Types

typedef CGAL::Point_3< KernelPoint_3
 
typedef CGAL::Segment_3< KernelSegment_3
 
typedef CGAL::Triangle_3< KernelTriangle_3
 
typedef CGAL::Plane_3< KernelPlane_3
 
typedef CGAL::Tetrahedron_3< KernelTetrahedron_3
 
typedef CGAL::Point_2< KernelPoint_2
 
typedef CGAL::Segment_2< KernelSegment_2
 
typedef CGAL::Triangle_2< KernelTriangle_2
 
typedef CGAL::Exact_intersections_tag Itag
 
typedef std::vector< Triangle_3Triangles
 
typedef Triangles::iterator TrianglesIterator
 
typedef Triangles::const_iterator TrianglesConstIterator
 
typedef CGAL::Box_intersection_d::Box_with_handle_d< double, 3, TrianglesIteratorBox
 
typedef DerivedF::Index Index
 
typedef std::vector< std::pair< Index, CGAL::Object > > ObjectList
 
typedef std::vector< IndexIndexList
 
typedef std::pair< Index, IndexEMK
 
typedef std::vector< IndexEMV
 
typedef std::map< EMK, EMVEdgeMap
 

Public Member Functions

 SelfIntersectMesh (const Eigen::MatrixBase< DerivedV > &V, const Eigen::MatrixBase< DerivedF > &F, const RemeshSelfIntersectionsParam &params, Eigen::PlainObjectBase< DerivedVV > &VV, Eigen::PlainObjectBase< DerivedFF > &FF, Eigen::PlainObjectBase< DerivedIF > &IF, Eigen::PlainObjectBase< DerivedJ > &J, Eigen::PlainObjectBase< DerivedIM > &IM)
 
void box_intersect (const Box &a, const Box &b)
 
void process_intersecting_boxes ()
 

Static Public Member Functions

static void box_intersect_static (SelfIntersectMesh *SIM, const Box &a, const Box &b)
 

Public Attributes

const Eigen::MatrixBase< DerivedV > & V
 
const Eigen::MatrixBase< DerivedF > & F
 
Index count
 
Triangles T
 
IndexList lIF
 
std::map< Index, ObjectListoffending
 
std::vector< std::pair< TrianglesIterator, TrianglesIterator > > candidate_triangle_pairs
 
RemeshSelfIntersectionsParam params
 

Private Types

typedef SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM > Self
 

Private Member Functions

void mark_offensive (const Index f)
 
void count_intersection (const Index fa, const Index fb)
 
bool intersect (const Triangle_3 &A, const Triangle_3 &B, const Index fa, const Index fb)
 
bool single_shared_vertex (const Triangle_3 &A, const Triangle_3 &B, const Index fa, const Index fb, const Index va, const Index vb)
 
bool single_shared_vertex (const Triangle_3 &A, const Triangle_3 &B, const Index fa, const Index fb, const Index va)
 
bool double_shared_vertex (const Triangle_3 &A, const Triangle_3 &B, const Index fa, const Index fb, const std::vector< std::pair< Index, Index > > shared)
 

Private Attributes

std::mutex m_offending_lock
 

Detailed Description

template<typename Kernel, typename DerivedV, typename DerivedF, typename DerivedVV, typename DerivedFF, typename DerivedIF, typename DerivedJ, typename DerivedIM>
class igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >

Member Typedef Documentation

◆ Box

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
typedef CGAL::Box_intersection_d::Box_with_handle_d<double,3,TrianglesIterator> igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::Box

◆ EdgeMap

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
typedef std::map<EMK,EMV> igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::EdgeMap

◆ EMK

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
typedef std::pair<Index,Index> igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::EMK

◆ EMV

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
typedef std::vector<Index> igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::EMV

◆ Index

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
typedef DerivedF::Index igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::Index

◆ IndexList

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
typedef std::vector<Index> igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::IndexList

◆ Itag

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
typedef CGAL::Exact_intersections_tag igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::Itag

◆ ObjectList

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
typedef std::vector<std::pair<Index, CGAL::Object> > igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::ObjectList

◆ Plane_3

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
typedef CGAL::Plane_3<Kernel> igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::Plane_3

◆ Point_2

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
typedef CGAL::Point_2<Kernel> igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::Point_2

◆ Point_3

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
typedef CGAL::Point_3<Kernel> igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::Point_3

◆ Segment_2

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
typedef CGAL::Segment_2<Kernel> igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::Segment_2

◆ Segment_3

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
typedef CGAL::Segment_3<Kernel> igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::Segment_3

◆ Self

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
typedef SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM> igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::Self
private

◆ Tetrahedron_3

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
typedef CGAL::Tetrahedron_3<Kernel> igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::Tetrahedron_3

◆ Triangle_2

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
typedef CGAL::Triangle_2<Kernel> igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::Triangle_2

◆ Triangle_3

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
typedef CGAL::Triangle_3<Kernel> igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::Triangle_3

◆ Triangles

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
typedef std::vector<Triangle_3> igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::Triangles

◆ TrianglesConstIterator

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
typedef Triangles::const_iterator igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::TrianglesConstIterator

◆ TrianglesIterator

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
typedef Triangles::iterator igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::TrianglesIterator

Constructor & Destructor Documentation

◆ SelfIntersectMesh()

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::SelfIntersectMesh ( const Eigen::MatrixBase< DerivedV > &  V,
const Eigen::MatrixBase< DerivedF > &  F,
const RemeshSelfIntersectionsParam params,
Eigen::PlainObjectBase< DerivedVV > &  VV,
Eigen::PlainObjectBase< DerivedFF > &  FF,
Eigen::PlainObjectBase< DerivedIF > &  IF,
Eigen::PlainObjectBase< DerivedJ > &  J,
Eigen::PlainObjectBase< DerivedIM > &  IM 
)
inline
299 :
300 V(V),
301 F(F),
302 count(0),
303 T(),
304 lIF(),
305 offending(),
307{
308 using namespace std;
309 using namespace Eigen;
310
311#ifdef IGL_SELFINTERSECTMESH_DEBUG
312 const auto & tictoc = []() -> double
313 {
314 static double t_start = igl::get_seconds();
315 double diff = igl::get_seconds()-t_start;
316 t_start += diff;
317 return diff;
318 };
319 const auto log_time = [&](const std::string& label) -> void{
320 std::cout << "SelfIntersectMesh." << label << ": "
321 << tictoc() << std::endl;
322 };
323 tictoc();
324#endif
325
326 // Compute and process self intersections
328#ifdef IGL_SELFINTERSECTMESH_DEBUG
329 log_time("convert_to_triangle_list");
330#endif
331 // http://www.cgal.org/Manual/latest/doc_html/cgal_manual/Box_intersection_d/Chapter_main.html#Section_63.5
332 // Create the corresponding vector of bounding boxes
333 std::vector<Box> boxes;
334 boxes.reserve(T.size());
335 for (
336 TrianglesIterator tit = T.begin();
337 tit != T.end();
338 ++tit)
339 {
340 if (!tit->is_degenerate())
341 {
342 boxes.push_back(Box(tit->bbox(), tit));
343 }
344 }
345 // Leapfrog callback
346 std::function<void(const Box &a,const Box &b)> cb =
347 std::bind(&box_intersect_static, this,
348 // Explicitly use std namespace to avoid confusion with boost (who puts
349 // _1 etc. in global namespace)
350 std::placeholders::_1,
351 std::placeholders::_2);
352#ifdef IGL_SELFINTERSECTMESH_DEBUG
353 log_time("box_and_bind");
354#endif
355 // Run the self intersection algorithm with all defaults
356 CGAL::box_self_intersection_d(boxes.begin(), boxes.end(),cb);
357#ifdef IGL_SELFINTERSECTMESH_DEBUG
358 log_time("box_intersection_d");
359#endif
360 try{
362 }catch(int e)
363 {
364 // Rethrow if not IGL_FIRST_HIT_EXCEPTION
366 {
367 throw e;
368 }
369 // Otherwise just fall through
370 }
371#ifdef IGL_SELFINTERSECTMESH_DEBUG
372 log_time("resolve_intersection");
373#endif
374
375 // Convert lIF to Eigen matrix
376 assert(lIF.size()%2 == 0);
377 IF.resize(lIF.size()/2,2);
378 {
379 Index i=0;
380 for(
381 typename IndexList::const_iterator ifit = lIF.begin();
382 ifit!=lIF.end();
383 )
384 {
385 IF(i,0) = (*ifit);
386 ifit++;
387 IF(i,1) = (*ifit);
388 ifit++;
389 i++;
390 }
391 }
392#ifdef IGL_SELFINTERSECTMESH_DEBUG
393 log_time("store_intersecting_face_pairs");
394#endif
395
397 {
398 return;
399 }
400
402 V,F,T,offending,params.stitch_all,VV,FF,J,IM);
403
404#ifdef IGL_SELFINTERSECTMESH_DEBUG
405 log_time("remesh_intersection");
406#endif
407}
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
Definition PlainObjectBase.h:279
IndexList lIF
Definition SelfIntersectMesh.h:92
Index count
Definition SelfIntersectMesh.h:87
void process_intersecting_boxes()
Definition SelfIntersectMesh.h:793
CGAL::Box_intersection_d::Box_with_handle_d< double, 3, TrianglesIterator > Box
Definition SelfIntersectMesh.h:80
const Eigen::MatrixBase< DerivedF > & F
Definition SelfIntersectMesh.h:84
RemeshSelfIntersectionsParam params
Definition SelfIntersectMesh.h:108
Triangles T
Definition SelfIntersectMesh.h:90
const Eigen::MatrixBase< DerivedV > & V
Definition SelfIntersectMesh.h:83
static void box_intersect_static(SelfIntersectMesh *SIM, const Box &a, const Box &b)
Definition SelfIntersectMesh.h:266
Triangles::iterator TrianglesIterator
Definition SelfIntersectMesh.h:76
std::map< Index, ObjectList > offending
Definition SelfIntersectMesh.h:95
typedef void(GLAPIENTRYP _GLUfuncptr)(void)
#define IGL_FIRST_HIT_EXCEPTION
Definition intersect_other.cpp:16
Definition LDLT.h:16
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:33
Slic3r::Polygons diff(const Slic3r::Polygon &subject, const Slic3r::Polygon &clip, ApplySafetyOffset do_safety_offset)
Definition ClipperUtils.cpp:672
IGL_INLINE 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, Eigen::PlainObjectBase< DerivedVV > &VV, Eigen::PlainObjectBase< DerivedFF > &FF, Eigen::PlainObjectBase< DerivedJ > &J, Eigen::PlainObjectBase< DerivedIM > &IM)
Definition remesh_intersections.cpp:58
IGL_INLINE void mesh_to_cgal_triangle_list(const Eigen::MatrixBase< DerivedV > &V, const Eigen::MatrixBase< DerivedF > &F, std::vector< CGAL::Triangle_3< Kernel > > &T)
Definition mesh_to_cgal_triangle_list.cpp:17
IGL_INLINE double get_seconds()
Definition get_seconds.cpp:10
STL namespace.
bool detect_only
Definition RemeshSelfIntersectionsParam.h:21
bool stitch_all
Definition RemeshSelfIntersectionsParam.h:23

References igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::box_intersect_static(), igl::count(), igl::copyleft::cgal::RemeshSelfIntersectionsParam::detect_only, igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::F, igl::get_seconds(), IGL_FIRST_HIT_EXCEPTION, igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::lIF, igl::copyleft::cgal::mesh_to_cgal_triangle_list(), igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::offending, igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::params, igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::process_intersecting_boxes(), igl::copyleft::cgal::remesh_intersections(), Eigen::PlainObjectBase< Derived >::resize(), igl::copyleft::cgal::RemeshSelfIntersectionsParam::stitch_all, igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::T, igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::V, and void().

+ Here is the call graph for this function:

Member Function Documentation

◆ box_intersect()

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
void igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::box_intersect ( const Box a,
const Box b 
)
inline
772{
773 candidate_triangle_pairs.push_back({a.handle(), b.handle()});
774}
std::vector< std::pair< TrianglesIterator, TrianglesIterator > > candidate_triangle_pairs
Definition SelfIntersectMesh.h:105

Referenced by igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::box_intersect_static().

+ Here is the caller graph for this function:

◆ box_intersect_static()

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
void igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::box_intersect_static ( SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM > *  SIM,
const Box a,
const Box b 
)
inlinestatic
270{
271 SIM->box_intersect(a,b);
272}

References igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::box_intersect().

Referenced by igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::SelfIntersectMesh().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ count_intersection()

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
void igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::count_intersection ( const Index  fa,
const Index  fb 
)
inlineprivate
458{
459 std::lock_guard<std::mutex> guard(m_offending_lock);
460 mark_offensive(fa);
461 mark_offensive(fb);
462 this->count++;
463 // We found the first intersection
464 if(params.first_only && this->count >= 1)
465 {
467 }
468
469}
std::mutex m_offending_lock
Definition SelfIntersectMesh.h:208
void mark_offensive(const Index f)
Definition SelfIntersectMesh.h:427
bool first_only
Definition RemeshSelfIntersectionsParam.h:22

References igl::count(), and IGL_FIRST_HIT_EXCEPTION.

+ Here is the call graph for this function:

◆ double_shared_vertex()

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
bool igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::double_shared_vertex ( const Triangle_3 A,
const Triangle_3 B,
const Index  fa,
const Index  fb,
const std::vector< std::pair< Index, Index > >  shared 
)
inlineprivate
635{
636 using namespace std;
637
638 auto opposite_vertex = [](const Index a0, const Index a1) {
639 // get opposite index of A
640 int a2=-1;
641 for(int c=0;c<3;++c)
642 if(c!=a0 && c!=a1) {
643 a2 = c;
644 break;
645 }
646 assert(a2 != -1);
647 return a2;
648 };
649
650 // must be co-planar
651 Index a2 = opposite_vertex(shared[0].first, shared[1].first);
652 if (! B.supporting_plane().has_on(A.vertex(a2)))
653 return false;
654
655 Index b2 = opposite_vertex(shared[0].second, shared[1].second);
656
657 if (int(CGAL::coplanar_orientation(A.vertex(shared[0].first), A.vertex(shared[1].first), A.vertex(a2))) *
658 int(CGAL::coplanar_orientation(B.vertex(shared[0].second), B.vertex(shared[1].second), B.vertex(b2))) < 0)
659 // There is certainly no self intersection as the non-shared triangle vertices lie on opposite sides of the shared edge.
660 return false;
661
662 // Since A and B are non-degenerate the intersection must be a polygon
663 // (triangle). Either
664 // - the vertex of A (B) opposite the shared edge of lies on B (A), or
665 // - an edge of A intersects and edge of B without sharing a vertex
666 //
667 // Determine if the vertex opposite edge (a0,a1) in triangle A lies in
668 // (intersects) triangle B
669 const auto & opposite_point_inside = [](
670 const Triangle_3 & A, const Index a2, const Triangle_3 & B)
671 -> bool
672 {
673 return CGAL::do_intersect(A.vertex(a2),B);
674 };
675
676 // Determine if edge opposite vertex va in triangle A intersects edge
677 // opposite vertex vb in triangle B.
678 const auto & opposite_edges_intersect = [](
679 const Triangle_3 & A, const Index va,
680 const Triangle_3 & B, const Index vb) -> bool
681 {
682 Segment_3 sa( A.vertex((va+1)%3), A.vertex((va+2)%3));
683 Segment_3 sb( B.vertex((vb+1)%3), B.vertex((vb+2)%3));
684 bool ret = CGAL::do_intersect(sa,sb);
685 return ret;
686 };
687
688 if(
689 !opposite_point_inside(A,a2,B) &&
690 !opposite_point_inside(B,b2,A) &&
691 !opposite_edges_intersect(A,shared[0].first,B,shared[1].second) &&
692 !opposite_edges_intersect(A,shared[1].first,B,shared[0].second))
693 {
694 return false;
695 }
696
697 // there is an intersection indeed
698 count_intersection(fa,fb);
700 {
701 return true;
702 }
703 // Construct intersection
704 try
705 {
706 // This can fail for Epick but not Epeck
707 CGAL::Object result = CGAL::intersection(A,B);
708 if(!result.empty())
709 {
710 if(CGAL::object_cast<Segment_3 >(&result))
711 {
712 // not coplanar
713 assert(false &&
714 "Co-planar non-degenerate triangles should intersect over triangle");
715 return false;
716 } else if(CGAL::object_cast<Point_3 >(&result))
717 {
718 // this "shouldn't" happen but does for inexact
719 assert(false &&
720 "Co-planar non-degenerate triangles should intersect over triangle");
721 return false;
722 } else
723 {
724 // Triangle object
725 std::lock_guard<std::mutex> guard(m_offending_lock);
726 offending[fa].push_back({fb, result});
727 offending[fb].push_back({fa, result});
728 return true;
729 }
730 }else
731 {
732 // CGAL::intersection is disagreeing with do_intersect
733 assert(false && "CGAL::intersection should agree with predicate tests");
734 return false;
735 }
736 }catch(...)
737 {
738 // This catches some cgal assertion:
739 // CGAL error: assertion violation!
740 // Expression : is_finite(d)
741 // File : /opt/local/include/CGAL/GMP/Gmpq_type.h
742 // Line : 132
743 // Explanation:
744 // But only if NDEBUG is not defined, otherwise there's an uncaught
745 // "Floating point exception: 8" SIGFPE
746 return false;
747 }
748 // No intersection.
749 return false;
750}
CGAL::Triangle_3< Kernel > Triangle_3
Definition SelfIntersectMesh.h:65
CGAL::Segment_3< Kernel > Segment_3
Definition SelfIntersectMesh.h:64
DerivedF::Index Index
Definition SelfIntersectMesh.h:86
void count_intersection(const Index fa, const Index fb)
Definition SelfIntersectMesh.h:455

◆ intersect()

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
bool igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::intersect ( const Triangle_3 A,
const Triangle_3 B,
const Index  fa,
const Index  fb 
)
inlineprivate
493{
494 // Determine whether there is an intersection
495 if(!CGAL::do_intersect(A,B))
496 {
497 return false;
498 }
499 count_intersection(fa,fb);
501 {
502 // Construct intersection
503 CGAL::Object result = CGAL::intersection(A,B);
504 std::lock_guard<std::mutex> guard(m_offending_lock);
505 offending[fa].push_back({fb, result});
506 offending[fb].push_back({fa, result});
507 }
508 return true;
509}

◆ mark_offensive()

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
void igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::mark_offensive ( const Index  f)
inlineprivate
428{
429 using namespace std;
430 lIF.push_back(f);
431 if(offending.count(f) == 0)
432 {
433 // first time marking, initialize with new id and empty list
434 offending[f] = {};
435 }
436}

◆ process_intersecting_boxes()

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
void igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::process_intersecting_boxes
inline
794{
795 std::vector<std::mutex> triangle_locks(T.size());
796 std::vector<std::mutex> vertex_locks(V.rows());
797 std::mutex index_lock;
798 std::mutex exception_mutex;
799 bool exception_fired = false;
800 int exception = -1;
801 auto process_chunk =
802 [&](
803 const size_t first,
804 const size_t last) -> void
805 {
806 try
807 {
808 assert(last >= first);
809
810 for (size_t i=first; i<last; i++)
811 {
812 if(exception_fired) return;
813 Index fa=T.size(), fb=T.size();
814 {
815 // Before knowing which triangles are involved, we need to lock
816 // everything to prevent race condition in updating reference
817 // counters.
818 std::lock_guard<std::mutex> guard(index_lock);
819 const auto& tri_pair = candidate_triangle_pairs[i];
820 fa = tri_pair.first - T.begin();
821 fb = tri_pair.second - T.begin();
822 }
823 assert(fa < T.size());
824 assert(fb < T.size());
825
826 // Lock triangles
827 std::lock_guard<std::mutex> guard_A(triangle_locks[fa]);
828 std::lock_guard<std::mutex> guard_B(triangle_locks[fb]);
829
830 // Lock vertices
831 std::list<std::lock_guard<std::mutex> > guard_vertices;
832 {
833 std::vector<typename DerivedF::Scalar> unique_vertices;
834 std::vector<size_t> tmp1, tmp2;
835 igl::unique({F(fa,0), F(fa,1), F(fa,2), F(fb,0), F(fb,1), F(fb,2)},
836 unique_vertices, tmp1, tmp2);
837 std::for_each(unique_vertices.begin(), unique_vertices.end(),
838 [&](const typename DerivedF::Scalar& vi) {
839 guard_vertices.emplace_back(vertex_locks[vi]);
840 });
841 }
842 if(exception_fired) return;
843
844 const Triangle_3& A = T[fa];
845 const Triangle_3& B = T[fb];
846
847 // Number of combinatorially shared vertices
848 Index comb_shared_vertices = 0;
849 // Number of geometrically shared vertices (*not* including
850 // combinatorially shared)
851 Index geo_shared_vertices = 0;
852 // Keep track of shared vertex indices
853 std::vector<std::pair<Index,Index> > shared;
854 Index ea,eb;
855 for(ea=0;ea<3;ea++)
856 {
857 for(eb=0;eb<3;eb++)
858 {
859 if(F(fa,ea) == F(fb,eb))
860 {
861 comb_shared_vertices++;
862 shared.emplace_back(ea,eb);
863 }else if(A.vertex(ea) == B.vertex(eb))
864 {
865 geo_shared_vertices++;
866 shared.emplace_back(ea,eb);
867 }
868 }
869 }
870 const Index total_shared_vertices =
871 comb_shared_vertices + geo_shared_vertices;
872 if(exception_fired) return;
873
874 if(comb_shared_vertices== 3)
875 {
876 assert(shared.size() == 3);
877 // Combinatorially duplicate face, these should be removed by
878 // preprocessing
879 continue;
880 }
881 if(total_shared_vertices== 3)
882 {
883 assert(shared.size() == 3);
884 // Geometrically duplicate face, these should be removed by
885 // preprocessing
886 continue;
887 }
888 if(total_shared_vertices == 2)
889 {
890 assert(shared.size() == 2);
891 // Q: What about coplanar?
892 //
893 // o o
894 // |\ /|
895 // | \/ |
896 // | /\ |
897 // |/ \|
898 // o----o
899 double_shared_vertex(A,B,fa,fb,shared);
900 continue;
901 }
902 assert(total_shared_vertices<=1);
903 if(total_shared_vertices==1)
904 {
905 single_shared_vertex(A,B,fa,fb,shared[0].first,shared[0].second);
906 }else
907 {
908 intersect(A,B,fa,fb);
909 }
910 }
911 }catch(int e)
912 {
913 std::lock_guard<std::mutex> exception_lock(exception_mutex);
914 exception_fired = true;
915 exception = e;
916 }
917 };
918 size_t num_threads=0;
919 const size_t hardware_limit = std::thread::hardware_concurrency();
920 if (const char* igl_num_threads = std::getenv("LIBIGL_NUM_THREADS")) {
921 num_threads = atoi(igl_num_threads);
922 }
923 if (num_threads == 0 || num_threads > hardware_limit) {
924 num_threads = hardware_limit;
925 }
926 assert(num_threads > 0);
927 const size_t num_pairs = candidate_triangle_pairs.size();
928 const size_t chunk_size = num_pairs / num_threads;
929 std::vector<std::thread> threads;
930 for (size_t i=0; i<num_threads-1; i++)
931 {
932 threads.emplace_back(process_chunk, i*chunk_size, (i+1)*chunk_size);
933 }
934 // Do some work in the master thread.
935 process_chunk((num_threads-1)*chunk_size, num_pairs);
936 for (auto& t : threads)
937 {
938 if (t.joinable()) t.join();
939 }
940 if(exception_fired) throw exception;
941 //process_chunk(0, candidate_triangle_pairs.size());
942}
bool single_shared_vertex(const Triangle_3 &A, const Triangle_3 &B, const Index fa, const Index fb, const Index va, const Index vb)
Definition SelfIntersectMesh.h:528
bool double_shared_vertex(const Triangle_3 &A, const Triangle_3 &B, const Index fa, const Index fb, const std::vector< std::pair< Index, Index > > shared)
Definition SelfIntersectMesh.h:629
bool intersect(const Triangle_3 &A, const Triangle_3 &B, const Index fa, const Index fb)
Definition SelfIntersectMesh.h:488
IGL_INLINE void unique(const std::vector< T > &A, std::vector< T > &C, std::vector< size_t > &IA, std::vector< size_t > &IC)
Definition unique.cpp:21

References igl::intersect(), and igl::unique().

Referenced by igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::SelfIntersectMesh().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ single_shared_vertex() [1/2]

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
bool igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::single_shared_vertex ( const Triangle_3 A,
const Triangle_3 B,
const Index  fa,
const Index  fb,
const Index  va 
)
inlineprivate
566{
567 // This was not a good idea. It will not handle coplanar triangles well.
568 using namespace std;
569 Segment_3 sa(
570 A.vertex((va+1)%3),
571 A.vertex((va+2)%3));
572
573 if(CGAL::do_intersect(sa,B))
574 {
575 // can't put count_intersection(fa,fb) here since we use intersect below
576 // and then it will be counted twice.
578 {
579 count_intersection(fa,fb);
580 return true;
581 }
582 CGAL::Object result = CGAL::intersection(sa,B);
583 if(const Point_3 * p = CGAL::object_cast<Point_3 >(&result))
584 {
585 // Single intersection --> segment from shared point to intersection
586 CGAL::Object seg = CGAL::make_object(Segment_3(
587 A.vertex(va),
588 *p));
589 count_intersection(fa,fb);
590 std::lock_guard<std::mutex> guard(m_offending_lock);
591 offending[fa].push_back({fb, seg});
592 offending[fb].push_back({fa, seg});
593 return true;
594 }else if(CGAL::object_cast<Segment_3 >(&result))
595 {
596 // Need to do full test. Intersection could be a general poly.
597 bool test = intersect(A,B,fa,fb);
598 ((void)test);
599 assert(test && "intersect should agree with do_intersect");
600 return true;
601 }else
602 {
603 cerr<<REDRUM("Segment ∩ triangle neither point nor segment?")<<endl;
604 assert(false);
605 }
606 }
607
608 return false;
609}
#define REDRUM(X)
Definition REDRUM.h:40
CGAL::Point_3< Kernel > Point_3
Definition SelfIntersectMesh.h:63

References igl::intersect(), REDRUM, and void().

+ Here is the call graph for this function:

◆ single_shared_vertex() [2/2]

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
bool igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::single_shared_vertex ( const Triangle_3 A,
const Triangle_3 B,
const Index  fa,
const Index  fb,
const Index  va,
const Index  vb 
)
inlineprivate
535{
536 if(single_shared_vertex(A,B,fa,fb,va))
537 {
538 return true;
539 }
540 return single_shared_vertex(B,A,fb,fa,vb);
541}

Member Data Documentation

◆ candidate_triangle_pairs

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
std::vector<std::pair<TrianglesIterator, TrianglesIterator> > igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::candidate_triangle_pairs

◆ count

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
Index igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::count

◆ F

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
const Eigen::MatrixBase<DerivedF>& igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::F

◆ lIF

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
IndexList igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::lIF

◆ m_offending_lock

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
std::mutex igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::m_offending_lock
private

◆ offending

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
std::map<Index,ObjectList> igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::offending

◆ params

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
RemeshSelfIntersectionsParam igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::params

◆ T

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
Triangles igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::T

◆ V

template<typename Kernel , typename DerivedV , typename DerivedF , typename DerivedVV , typename DerivedFF , typename DerivedIF , typename DerivedJ , typename DerivedIM >
const Eigen::MatrixBase<DerivedV>& igl::copyleft::cgal::SelfIntersectMesh< Kernel, DerivedV, DerivedF, DerivedVV, DerivedFF, DerivedIF, DerivedJ, DerivedIM >::V

The documentation for this class was generated from the following file: