Easy3D 2.5.3
self_intersection.h
1/********************************************************************
2 * Copyright (C) 2018-2021 by Liangliang Nan <liangliang.nan@gmail.com>
3 * Copyright (C) 2015 Alec Jacobson <alecjacobson@gmail.com>
4 *
5 * The code in this file is partly from libigl (version 3 Mar 2016) with
6 * modifications and enhancement:
7 * https://libigl.github.io/
8 * The original code was distributed under the Mozilla Public License v2.
9 ********************************************************************/
10
11#ifndef EASY3D_ALGO_EXT_SELF_INTERSECTION_H
12#define EASY3D_ALGO_EXT_SELF_INTERSECTION_H
13
14#include <vector>
15#include <unordered_map>
16
17// Use this instead to mute errors resulting from bad CGAL assertions
18//#define CGAL_KERNEL_NO_ASSERTIONS
19
20#include <CGAL/Exact_predicates_exact_constructions_kernel.h>
21#include <CGAL/intersections.h> // for triangle-triangle intersection
22#include <CGAL/box_intersection_d.h>
23#include <CGAL/Constrained_Delaunay_triangulation_2.h>
24#include <CGAL/Constrained_triangulation_plus_2.h>
25
26#include <easy3d/core/surface_mesh.h>
27
28namespace easy3d {
29
48
50 public:
53
60 std::vector<std::pair<SurfaceMesh::Face, SurfaceMesh::Face> >
61 detect(SurfaceMesh *mesh, bool construct = false);
62
70 bool remesh(SurfaceMesh *mesh, bool stitch);
71
72 private:
73
74 typedef CGAL::Exact_predicates_exact_constructions_kernel Kernel;
75
76 // 3D Primitives
77 typedef CGAL::Point_3<Kernel> Point_3;
78 typedef CGAL::Vector_3<Kernel> Vector_3;
79 typedef CGAL::Triangle_3<Kernel> Triangle_3;
80 typedef CGAL::Segment_3<Kernel> Segment_3;
81 typedef CGAL::Plane_3<Kernel> Plane_3;
82
83 // Constrained Delaunay Triangulation
84 typedef CGAL::Triangulation_vertex_base_2<Kernel> TVB_2;
85 typedef CGAL::Constrained_triangulation_face_base_2<Kernel> CTFB_2;
86 typedef CGAL::Triangulation_data_structure_2<TVB_2, CTFB_2> TDS_2;
87 typedef CGAL::Exact_intersections_tag Itag;
88 typedef CGAL::Constrained_Delaunay_triangulation_2<Kernel, TDS_2, Itag> CDT_2;
89 typedef CGAL::Constrained_triangulation_plus_2<CDT_2> CDT_plus_2;
90
91 typedef std::vector<std::pair<int, CGAL::Object> > ObjectList;
92
93 struct Triangle {
94 public:
95 Triangle(const Point_3 &a, const Point_3 &b, const Point_3 &c, SurfaceMesh::Face f) : triangle(a, b, c),
96 face(f) {}
97
98 Triangle_3 triangle;
100 int index; // face index
101 std::vector<SurfaceMesh::Vertex> vertices;
102 };
103
104 // Axis-align boxes for all-pairs self-intersection detection
105 typedef std::vector<Triangle> Triangles;
106 typedef typename Triangles::iterator TrianglesIterator;
107 typedef CGAL::Box_intersection_d::Box_with_handle_d<double, 3, TrianglesIterator> Box;
108
109 private:
110
111 void mesh_to_cgal_triangle_list(SurfaceMesh *mesh);
112
113 // test if two triangles intersect
114 bool do_intersect(const Triangle &A, const Triangle &B);
115
116 // Given a list of objects (e.g., resulting from intersecting a triangle
117 // with many other triangles), construct a constrained Delaunay
118 // triangulation on a given plane (P), by inserting constraints for each
119 // object projected onto that plane.
120 // Inputs:
121 // objects list of objects. This should lie on the given plane (P),
122 // otherwise they are added to the cdt _after_ their non-trivial
123 // projection
124 // P plane upon which all objects lie and upon which the CDT is
125 // conducted
126 // Outputs:
127 // vertices list of vertices of the CDT mesh _back on the 3D plane_
128 // faces list of triangle indices into vertices
129 //
130 void projected_cdt(
131 const std::vector<CGAL::Object> &objects,
132 const Plane_3 &P,
133 std::vector<Point_3> &vertices,
134 std::vector<std::vector<int> > &faces
135 );
136
137 // Given a current 2D constrained Delaunay triangulation (cdt), insert a
138 // 3D "object" (e.g., resulting from intersecting two triangles) into the
139 // cdt, by projecting it via the given plane (P) and adding appropriate
140 // constraints.
141 // Inputs:
142 // obj CGAL::Object representing a vertex, segment, or (convex)
143 // polygon. All vertices should lie on the plane P. If not, then this
144 // adds the _projection_ of this object to the cdt and that might not
145 // be what you wanted to do.
146 // P plane obj lies on and upon which the cdt is being performed
147 // cdt current CDT, see output
148 // Outputs:
149 // cdt CDT updated to contain constraints for the given object
150 //
151 void insert_into_cdt(const CGAL::Object &obj, const Plane_3 &P, CDT_plus_2 &cdt);
152
153 private:
154 // Helper function to mark a face as offensive
155 // Inputs:
156 // f index of face in F
157 inline void mark_offensive(int f);
158
159 // Helper function to count intersections between faces
160 inline void count_intersection(int fa, int fb);
161
162 // Helper function for box_intersect. Intersect two triangles A and B,
163 // append the intersection object (point,segment,triangle) to a running
164 // list for A and B
165 // Inputs:
166 // A triangle in 3D
167 // B triangle in 3D
168 // Returns true only if A intersects B
169 inline bool intersect(const Triangle &A, const Triangle &B);
170
171 // Helper function for box_intersect. In the case where A and B have
172 // already been identified to share a vertex, then we only want to
173 // add possible segment intersections. Assumes truly duplicate
174 // triangles are not given as input
175 // Inputs:
176 // A triangle in 3D
177 // B triangle in 3D
178 // va shared vertex in A (and key into offending)
179 // vb shared vertex in B (and key into offending)
180 // Returns true if intersection (besides shared point)
181 inline bool single_shared_vertex(const Triangle &A, const Triangle &B, int va, int vb);
182
183 // Helper handling one direction
184 inline bool single_shared_vertex(const Triangle &A, const Triangle &B, int va);
185
186 // Helper function for box_intersect. In the case where A and B have
187 // already been identified to share two vertices, then we only want
188 // to add a possible coplanar (Triangle) intersection. Assumes truly
189 // degenerate facets are not given as input.
190 inline bool double_shared_vertex(
191 const Triangle &A,
192 const Triangle &B,
193 const std::vector<std::pair<int, int> > &shared
194 );
195
196 private:
197 SurfaceMesh* mesh_;
198
199 bool construct_intersection_;
200
201 Triangles triangle_faces_;
202
203 // index in 'triangle_faces_' (degenerate faces removed) -> original face
204 std::unordered_map<int, SurfaceMesh::Face> original_face;
205
206 // list of faces with intersections mapping to the order in
207 // which they were first found
208 std::unordered_map<int, ObjectList> offending_;
209
210 // statistics on duplicate faces (to offer the user some feedback)
211 int total_comb_duplicate_face_; // Combinatorial duplicate faces
212 int total_geom_duplicate_face_; // Geometrically duplicate faces
213 };
214
215} // namespace easy3d
216
217#endif // EASY3D_ALGO_EXT_SELF_INTERSECTION_H
Detects and resolves self-intersection for surface mesh.
Definition: self_intersection.h:49
std::vector< std::pair< SurfaceMesh::Face, SurfaceMesh::Face > > detect(SurfaceMesh *mesh, bool construct=false)
Detect intersecting face pairs.
Definition: self_intersection.cpp:40
bool remesh(SurfaceMesh *mesh, bool stitch)
Detect and remesh the intersecting faces.
Definition: self_intersection.cpp:427
A halfedge data structure for polygonal meshes of 2-manifold.
Definition: surface_mesh.h:52
Definition: collider.cpp:182
Definition: surface_mesh.h:134