8#ifndef IGL_COPYLEFT_CGAL_SELFINTERSECTMESH_H 
    9#define IGL_COPYLEFT_CGAL_SELFINTERSECTMESH_H 
   11#include "CGAL_includes.hpp" 
   13#include "../../unique.h" 
   14#include "../../default_num_threads.h" 
   25#ifndef IGL_FIRST_HIT_EXCEPTION 
   26#define IGL_FIRST_HIT_EXCEPTION 10 
   76          typedef CGAL::Exact_intersections_tag 
Itag;
 
   82            CGAL::Box_intersection_d::Box_with_handle_d<double,3,TrianglesIterator>
 
   86          const Eigen::MatrixBase<DerivedV> & 
V;
 
   87          const Eigen::MatrixBase<DerivedF> & 
F;
 
   89          typedef typename DerivedF::Index 
Index;
 
   91          typedef std::vector<std::pair<Index, CGAL::Object>> 
ObjectList;
 
  100          typedef std::pair<Index,Index> 
EMK;
 
  103          typedef std::vector<Index> 
EMV;
 
  107          std::vector<std::pair<TrianglesIterator, TrianglesIterator> >
 
  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);
 
  140          inline void mark_offensive(
const Index f);
 
  145          inline void count_intersection( 
const Index fa, 
const Index fb);
 
  156          inline bool intersect(
 
  174          inline bool single_shared_vertex(
 
  189          inline bool single_shared_vertex(
 
  206          inline bool double_shared_vertex(
 
  211              const std::vector<std::pair<Index,Index> > shared);
 
  236          std::mutex m_offending_lock;
 
 
  247#include "../../get_seconds.h" 
  248#include "../../C_STR.h" 
  293  DerivedIM>::box_intersect_static(
 
 
  318  DerivedIM>::SelfIntersectMesh(
 
  319  const Eigen::MatrixBase<DerivedV> & V,
 
  320  const Eigen::MatrixBase<DerivedF> & F,
 
  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):
 
  336  using namespace Eigen;
 
  338#ifdef IGL_SELFINTERSECTMESH_TIMING 
  339  const auto & tictoc = []() -> 
double 
  346  const auto log_time = [&](
const std::string& label) -> 
void{
 
  347    printf(
"%50s: %0.5lf\n",
 
  348      C_STR(
"SelfIntersectMesh." << label),tictoc());
 
  355#ifdef IGL_SELFINTERSECTMESH_TIMING 
  356  log_time(
"convert_to_triangle_list");
 
  360  std::vector<Box> boxes;
 
  361  boxes.reserve(
T.size());
 
  367    if (!tit->is_degenerate())
 
  369      boxes.push_back(
Box(tit->bbox(), tit));
 
  373  std::function<void(
const Box &a,
const Box &b)> cb =
 
  377      std::placeholders::_1,
 
  378      std::placeholders::_2);
 
  379#ifdef IGL_SELFINTERSECTMESH_TIMING 
  380  log_time(
"box_and_bind");
 
  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");
 
  398#ifdef IGL_SELFINTERSECTMESH_TIMING 
  399  log_time(
"resolve_intersection");
 
  403  assert(
lIF.size()%2 == 0);
 
  404  IF.resize(
lIF.size()/2,2);
 
  408      typename IndexList::const_iterator ifit = 
lIF.begin();
 
  419#ifdef IGL_SELFINTERSECTMESH_TIMING 
  420  log_time(
"store_intersecting_face_pairs");
 
  432#ifdef IGL_SELFINTERSECTMESH_TIMING 
  433  log_time(
"remesh_intersection");
 
 
  455  DerivedIM>::mark_offensive(
const Index f)
 
  459  if(offending.count(f) == 0)
 
  483  DerivedIM>::count_intersection(
 
  487  std::lock_guard<std::mutex> guard(m_offending_lock);
 
  492  if(params.first_only && this->count >= 1)
 
  517  const Triangle_3 & A,
 
  518  const Triangle_3 & B,
 
  523  if(!CGAL::do_intersect(A,B))
 
  527  count_intersection(fa,fb);
 
  528  if(!params.detect_only)
 
  531    CGAL::Object result = CGAL::intersection(A,B);
 
  534    std::lock_guard<std::mutex> guard(m_offending_lock);
 
  535    offending[fa].push_back({fb, result});
 
  536    offending[fb].push_back({fa, result});
 
  558  DerivedIM>::single_shared_vertex(
 
  559  const Triangle_3 & A,
 
  560  const Triangle_3 & B,
 
  566  if(single_shared_vertex(A,B,fa,fb,va))
 
  570  return single_shared_vertex(B,A,fb,fa,vb);
 
  590  DerivedIM>::single_shared_vertex(
 
  591  const Triangle_3 & A,
 
  592  const Triangle_3 & B,
 
  603  if(CGAL::do_intersect(sa,B))
 
  607    if(params.detect_only)
 
  609      count_intersection(fa,fb);
 
  612    CGAL::Object result = CGAL::intersection(sa,B);
 
  613    if(
const Point_3 * p = CGAL::object_cast<Point_3 >(&result))
 
  616      CGAL::Object seg = CGAL::make_object(Segment_3(
 
  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});
 
  624    }
else if(CGAL::object_cast<Segment_3 >(&result))
 
  629      assert(test && 
"intersect should agree with do_intersect");
 
  633      cerr<<
"Segment ∩ triangle neither point nor segment?"<<endl;
 
  659  DerivedIM>::double_shared_vertex(
 
  660  const Triangle_3 & A,
 
  661  const Triangle_3 & B,
 
  664  const std::vector<std::pair<Index,Index> > shared)
 
  668  auto opposite_vertex = [](
const Index a0, 
const Index a1) {
 
  681  Index a2 = opposite_vertex(shared[0].first, shared[1].first);
 
  682  if (! B.supporting_plane().has_on(A.vertex(a2)))
 
  685  Index b2 = opposite_vertex(shared[0].second, shared[1].second);
 
  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)
 
  699  const auto & opposite_point_inside = [](
 
  700    const Triangle_3 & A, 
const Index a2, 
const Triangle_3 & B) 
 
  703    return CGAL::do_intersect(A.vertex(a2),B);
 
  708  const auto & opposite_edges_intersect = [](
 
  709    const Triangle_3 & A, 
const Index va,
 
  710    const Triangle_3 & B, 
const Index vb) -> 
bool 
  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);
 
  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))
 
  728  count_intersection(fa,fb);
 
  729  if(params.detect_only)
 
  737    CGAL::Object result = CGAL::intersection(A,B);
 
  740      if(CGAL::object_cast<Segment_3 >(&result))
 
  744          "Co-planar non-degenerate triangles should intersect over triangle");
 
  746      } 
else if(CGAL::object_cast<Point_3 >(&result))
 
  750          "Co-planar non-degenerate triangles should intersect over triangle");
 
  755        std::lock_guard<std::mutex> guard(m_offending_lock);
 
  756        offending[fa].push_back({fb, result});
 
  757        offending[fb].push_back({fa, result});
 
  763      assert(
false && 
"CGAL::intersection should agree with predicate tests");
 
  799  DerivedIM>::box_intersect(
 
  803  candidate_triangle_pairs.push_back({a.handle(), b.handle()});
 
 
  823  DerivedIM>::process_intersecting_boxes()
 
  825  std::mutex exception_mutex;
 
  826  bool exception_fired = 
false;
 
  833  auto process_chunk = [&]( 
const size_t first, 
const size_t last) -> 
void 
  837      assert(last >= first);
 
  839      for (
size_t i=first; i<last; i++)
 
  841        if(exception_fired) 
return;
 
  842        Index fa=T.size(), fb=T.size();
 
  844          const auto& tri_pair = candidate_triangle_pairs[i];
 
  845          fa = tri_pair.first - T.begin();
 
  846          fb = tri_pair.second - T.begin();
 
  848        assert(fa < T.size());
 
  849        assert(fb < T.size());
 
  851        if(exception_fired) 
return;
 
  857        Index comb_shared_vertices = 0;
 
  860        Index geo_shared_vertices = 0;
 
  862        std::vector<std::pair<Index,Index> > shared;
 
  868            if(F(fa,ea) == F(fb,eb))
 
  870              comb_shared_vertices++;
 
  871              shared.emplace_back(ea,eb);
 
  872            }
else if(A.vertex(ea) == B.vertex(eb))
 
  874              geo_shared_vertices++;
 
  875              shared.emplace_back(ea,eb);
 
  879        const Index total_shared_vertices =
 
  880          comb_shared_vertices + geo_shared_vertices;
 
  881        if(exception_fired) 
return;
 
  883        if(comb_shared_vertices== 3)
 
  885          assert(shared.size() == 3);
 
  890        if(total_shared_vertices== 3)
 
  892          assert(shared.size() == 3);
 
  897        if(total_shared_vertices == 2)
 
  899          assert(shared.size() == 2);
 
  908          double_shared_vertex(A,B,fa,fb,shared);
 
  911        assert(total_shared_vertices<=1);
 
  912        if(total_shared_vertices==1)
 
  914          single_shared_vertex(A,B,fa,fb,shared[0].first,shared[0].second);
 
  922      std::lock_guard<std::mutex> exception_lock(exception_mutex);
 
  923      exception_fired = 
true;
 
  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++)
 
  934    threads.emplace_back(process_chunk, i*chunk_size, (i+1)*chunk_size);
 
  937  process_chunk((num_threads-1)*chunk_size, num_pairs);
 
  938  for (
auto& t : threads)
 
  940    if (t.joinable()) t.join();
 
  942  if(exception_fired) 
throw exception;
 
 
#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.
 
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