Geogram Version 1.9.9
A programming library of geometric algorithms
Loading...
Searching...
No Matches
mesh_surface_intersection.h
Go to the documentation of this file.
1/*
2 * Copyright (c) 2000-2022 Inria
3 * All rights reserved.
4 *
5 * Redistribution and use in source and binary forms, with or without
6 * modification, are permitted provided that the following conditions are met:
7 *
8 * * Redistributions of source code must retain the above copyright notice,
9 * this list of conditions and the following disclaimer.
10 * * Redistributions in binary form must reproduce the above copyright notice,
11 * this list of conditions and the following disclaimer in the documentation
12 * and/or other materials provided with the distribution.
13 * * Neither the name of the ALICE Project-Team nor the names of its
14 * contributors may be used to endorse or promote products derived from this
15 * software without specific prior written permission.
16 *
17 * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
18 * AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
19 * IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
20 * ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
21 * LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
22 * CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
23 * SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
24 * INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
25 * CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
26 * ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
27 * POSSIBILITY OF SUCH DAMAGE.
28 *
29 * Contact: Bruno Levy
30 *
31 * https://www.inria.fr/fr/bruno-levy
32 *
33 * Inria,
34 * Domaine de Voluceau,
35 * 78150 Le Chesnay - Rocquencourt
36 * FRANCE
37 *
38 */
39
40#ifndef GEOGRAM_MESH_MESH_SURFACE_INTERSECTION
41#define GEOGRAM_MESH_MESH_SURFACE_INTERSECTION
42
44#include <geogram/mesh/mesh.h>
49#include <geogram/basic/debug_stream.h>
50#include <functional>
51#include <tuple>
52
59namespace GEO {
60
61 struct IsectInfo;
62
63 /********************************************************************/
64
69 class GEOGRAM_API MeshSurfaceIntersection {
70 public:
71
73
76
84 void intersect();
85
91
97
98
99 void remove_fins();
100
124 void classify(const std::string& expr);
125
126
135
149 index_t component, index_t v
150 );
151
163 index_t component, index_t v
164 );
165
166
173 void simplify_coplanar_facets(double angle_tolerance = 0.0);
174
179 void set_verbose(bool x) {
180 verbose_ = x;
181 if(!verbose_ && fine_verbose_) {
182 fine_verbose_ = false;
183 }
184 }
185
190 void set_fine_verbose(bool x) {
191 fine_verbose_ = x;
192 if(fine_verbose_ && !verbose_) {
193 verbose_ = true;
194 }
195 }
196
205 monster_threshold_ = nb;
206 }
207
213 void set_dry_run(bool x) {
214 dry_run_ = x;
215 }
216
223 void set_delaunay(bool x) {
224 delaunay_ = x;
225 }
226
233 detect_intersecting_neighbors_ = x;
234 }
235
243 void set_radial_sort(bool x) {
244 use_radial_sort_ = x;
245 }
246
254 void set_normalize(bool x) {
255 normalize_ = x;
256 }
257
258
268 void set_build_skeleton(Mesh* skeleton, bool trim_fins=false) {
269 skeleton_ = skeleton;
270 skeleton_trim_fins_ = trim_fins;
271 }
272
279 interpolate_attributes_ = x;
280 }
281
282 protected:
292
307
326
342 void intersect_epilogue(const vector<IsectInfo>& intersections);
343
344
353 void lock() {
354 Process::acquire_spinlock(lock_);
355 }
356
360 void unlock() {
361 Process::release_spinlock(lock_);
362 }
363
374
383 return (vertex_to_exact_point_[v] == nullptr);
384 }
385
397
404 return mesh_;
405 }
406
412 const Mesh& target_mesh() const {
413 return mesh_;
414 }
415
428 const Mesh& readonly_mesh() const {
429 return mesh_copy_;
430 }
431
432 class RadialSort;
433
444
448 void mark_external_shell(vector<index_t>& on_external_shell);
449
457 return original_facet_id_[f];
458 }
459
470 return f_is_flipped_[f];
471 }
472
483 std::tuple<vec3, vec3, vec3> get_initial_facet_vertices(index_t f) const {
484 // All facets are duplicated before radial sort. If f is one of
485 // the duplicated facets then its orientation is flipped as
486 // compared to initial facet.
487 index_t orig_f = original_facet_id_[f];
488 vec3 p1 = mesh_copy_.facets.point(orig_f,0);
489 vec3 p2 = mesh_copy_.facets.point(orig_f,1);
490 vec3 p3 = mesh_copy_.facets.point(orig_f,2);
491 if(f_is_flipped_[f]) {
492 std::swap(p1,p3);
493 }
494 return std::make_tuple(p1,p2,p3);
495 }
496
497 protected:
498
502 class GEOGRAM_API RadialSort {
503 public:
509 I_(I),
510 mesh_(I_.target_mesh()),
511 h_ref_(NO_INDEX),
512 degenerate_(false)
513 {
514 }
515
520 void init(index_t h_ref);
521
528 bool operator()(index_t h1, index_t h2) const;
529
536 bool degenerate() const {
537 return degenerate_;
538 }
539
540 protected:
550
558
565
574 void report_problem(const char* message) const;
575
576 private:
577 const MeshSurfaceIntersection& I_;
578 const Mesh& mesh_;
579 index_t h_ref_; // reference halfedge
580 exact::vec3 N_ref_; // normal to reference triangle (exact)
581 mutable bool degenerate_;
582 };
583
584 protected:
585 Process::spinlock lock_;
586 Mesh& mesh_;
587 Mesh mesh_copy_;
588 Attribute<const ExactPoint*> vertex_to_exact_point_;
589 Attribute<index_t> original_facet_id_; // mesh_ facet to mesh_copy_ facet
590 Attribute<bool> f_is_flipped_; // mesh_ facet is flipped wrt mesh_copy_
591
592#if defined(GEOGRAM_USE_EXACT_NT) && defined(GEOGRAM_EXACT_NT_IS_MPF_NT)
593 // Exact points are canonicalized
594 // (by Numeric::optimize_number_representation(vec3HEx)) so
595 // we can use this comparator that makes the global vertex map
596 // much much faster.
597 typedef vec3HExLexicoCompareCanonical ExactPointCompare;
598#else
599 // Generic comparator for global vertex map.
601#endif
602 std::map<ExactPoint,index_t,ExactPointCompare> exact_point_to_vertex_;
603
604 bool verbose_;
605 bool fine_verbose_;
606 bool delaunay_;
607 bool detect_intersecting_neighbors_;
608 bool use_radial_sort_;
609
610 PCK::SOSMode SOS_bkp_;
611 bool rescale_;
612 bool normalize_;
613 vec3 normalize_center_;
614 double normalize_radius_;
615
616 index_t monster_threshold_;
617 bool dry_run_;
618 friend class MeshInTriangle;
619 friend class CoplanarFacets;
620
621 Mesh* skeleton_;
622 bool skeleton_trim_fins_;
623 bool interpolate_attributes_;
624
625 /***************************************************/
626
635 class Halfedges {
636 public:
637
642 Halfedges(MeshSurfaceIntersection& I) : mesh_(I.mesh_) {
643 }
644
649 // TODO: destroy alpha3 attribute (kept now for debugging
650 }
651
656 void initialize() {
657 facet_corner_alpha3_.bind(
658 mesh_.facet_corners.attributes(), "alpha3"
659 );
660 }
661
667 index_t nb() const {
668 return mesh_.facet_corners.nb();
669 }
670
676 return index_as_iterator(0);
677 }
678
684 return index_as_iterator(nb());
685 }
686
693 return h/3;
694 }
695
706 index_t t1 = facet(h);
707 index_t t2 = mesh_.facet_corners.adjacent_facet(h);
708 if(t2 == NO_INDEX) {
709 return NO_INDEX;
710 }
711 for(index_t h2: mesh_.facets.corners(t2)) {
712 if(mesh_.facet_corners.adjacent_facet(h2) == t1) {
713 return h2;
714 }
715 }
717 }
718
729 return facet_corner_alpha3_[h];
730 }
731
739 return alpha3(3*f)/3;
740 }
741
753 index_t f = h/3;
754 index_t lv = (h+dlv)%3;
755 return mesh_.facets.vertex(f,lv);
756 }
757
758
767 void sew2(index_t h1, index_t h2) {
768 geo_debug_assert(vertex(h1,0) == vertex(h2,1));
769 geo_debug_assert(vertex(h2,0) == vertex(h1,1));
770 index_t t1 = h1/3;
771 index_t t2 = h2/3;
772 mesh_.facet_corners.set_adjacent_facet(h1,t2);
773 mesh_.facet_corners.set_adjacent_facet(h2,t1);
774 }
775
784 void sew3(index_t h1, index_t h2) {
785 geo_debug_assert(vertex(h1,0) == vertex(h2,1));
786 geo_debug_assert(vertex(h2,0) == vertex(h1,1));
787 facet_corner_alpha3_[h1] = h2;
788 facet_corner_alpha3_[h2] = h1;
789 }
790
791 private:
792 Mesh& mesh_;
793 Attribute<index_t> facet_corner_alpha3_;
794 } halfedges_;
795
796 /***************************************************/
797
804 public:
805
810 RadialBundles(MeshSurfaceIntersection& I) : I_(I), mesh_(I.mesh_) {
811 }
812
818
822 index_t nb() const {
823 return bndl_start_.size() - 1;
824 }
825
831 return index_as_iterator(0);
832 }
833
839 return index_as_iterator(nb());
840 }
841
848 geo_debug_assert(bndl < nb());
849 return bndl_start_[bndl+1] - bndl_start_[bndl];
850 }
851
859 index_t halfedge(index_t bndl, index_t li) const {
860 geo_debug_assert(bndl < nb());
861 geo_debug_assert(li < nb_halfedges(bndl));
862 return H_[bndl_start_[bndl] + li];
863 }
864
873 geo_debug_assert(bndl < nb());
874 geo_debug_assert(li < nb_halfedges(bndl));
875 H_[bndl_start_[bndl] + li] = h;
876 }
877
884 return index_ptr_range(
885 H_, bndl_start_[bndl], bndl_start_[bndl+1]
886 );
887 }
888
896 H_, bndl_start_[bndl], bndl_start_[bndl+1]
897 );
898 }
899
907 index_t vertex(index_t bndl, index_t lv) const {
908 geo_debug_assert(bndl_start_[bndl+1] - bndl_start_[bndl] > 0);
909 index_t h = H_[bndl_start_[bndl]];
910 return I_.halfedges_.vertex(h,lv);
911 }
912
921 return v_first_bndl_[v];
922 }
923
932 return bndl_next_around_v_[bndl];
933 }
934
941 index_t result = 0;
942 for(
943 index_t bndl = vertex_first_bundle(v);
944 bndl != NO_INDEX;
945 bndl = next_around_vertex(bndl)
946 ) {
947 ++result;
948 }
949 return result;
950 }
951
959 geo_debug_assert(bndl < nb());
960 return (bndl >= nb()/2) ? (bndl-nb()/2) : (bndl+nb()/2);
961 }
962
969 index_t v = vertex(bndl,0);
970 if(nb_bundles_around_vertex(v) != 2) {
971 return NO_INDEX;
972 }
973 for(
974 index_t bndl2 = vertex_first_bundle(v);
975 bndl2 != NO_INDEX; bndl2 = next_around_vertex(bndl2)
976 ) {
977 if(bndl2 != bndl) {
978 return opposite(bndl2);
979 }
980 }
982 }
983
990 index_t v = vertex(bndl,1);
991 if(nb_bundles_around_vertex(v) != 2) {
992 return NO_INDEX;
993 }
994 for(
995 index_t bndl2 = vertex_first_bundle(v);
996 bndl2 != NO_INDEX; bndl2 = next_around_vertex(bndl2)
997 ) {
998 if(opposite(bndl2) != bndl) {
999 return bndl2;
1000 }
1001 }
1003 }
1004
1014 geo_debug_assert(bndl < nb());
1015 if(nb_halfedges(bndl) <= 2) {
1016 return true;
1017 }
1018 auto b = H_.begin() + std::ptrdiff_t(bndl_start_[bndl]);
1019 auto e = H_.begin() + std::ptrdiff_t(bndl_start_[bndl+1]);
1020 RS.init(*b);
1021 std::sort(
1022 b, e, [&RS](index_t h1, index_t h2) {
1023 return RS(h1,h2);
1024 }
1025 );
1026 bool OK = !RS.degenerate();
1027 bndl_is_sorted_[bndl] = OK;
1028 return OK;
1029 }
1030
1040 index_t bndl, const vector<index_t>& halfedges
1041 ) {
1042 geo_debug_assert(halfedges.size() == nb_halfedges(bndl));
1043 for(index_t i=0; i<halfedges.size(); ++i) {
1044 set_halfedge(bndl, i, halfedges[i]);
1045 }
1046 bndl_is_sorted_[bndl] = true;
1047 }
1048
1054 typedef std::pair<index_t, index_t> ChartPos;
1055
1064 index_t bndl, vector<ChartPos>& chart_pos
1065 );
1066
1067 bool is_sorted(index_t bndl) const {
1068 geo_debug_assert(bndl < nb());
1069 return bndl_is_sorted_[bndl];
1070 }
1071
1072 // private:
1074 Mesh& mesh_;
1075 Attribute<index_t> facet_chart_;
1076 vector<index_t> H_;
1077 vector<index_t> bndl_start_;
1078 vector<index_t> v_first_bndl_;
1079 vector<index_t> bndl_next_around_v_;
1080 vector<char> bndl_is_sorted_; // not vector<bool>, multithread! (#308)
1081 } radial_bundles_;
1082
1083 /***************************************************/
1084
1086 public:
1091 RadialPolylines(MeshSurfaceIntersection& I) : I_(I), mesh_(I.mesh_) {
1092 }
1093
1099
1106
1110 index_t nb() const {
1111 return polyline_start_.size() - 1;
1112 }
1113
1119 return index_as_iterator(0);
1120 }
1121
1127 return index_as_iterator(nb());
1128 }
1129
1136 geo_debug_assert(polyline < nb());
1137 return const_index_ptr_range(
1138 B_, polyline_start_[polyline], polyline_start_[polyline+1]
1139 );
1140 }
1141
1142 index_t nb_bundles(index_t polyline) const {
1143 geo_debug_assert(polyline < nb());
1144 return polyline_start_[polyline+1] - polyline_start_[polyline];
1145 }
1146
1147 index_t bundle(index_t polyline, index_t li) const {
1148 geo_debug_assert(polyline < nb());
1149 geo_debug_assert(li < nb_bundles(polyline));
1150 return B_[polyline_start_[polyline] + li];
1151 }
1152
1160 void get_skeleton(Mesh& to, bool trim_fins=false);
1161
1162 private:
1164 Mesh& mesh_;
1165 vector<index_t> B_;
1166 vector<index_t> polyline_start_;
1167 } radial_polylines_;
1168 };
1169
1170 /********************************************************************/
1171
1172 enum MeshBooleanOperationFlags {
1173 MESH_BOOL_OPS_DEFAULT = 0,
1174 MESH_BOOL_OPS_VERBOSE = 1,
1175 MESH_BOOL_OPS_ATTRIBS = 2,
1176 MESH_BOOL_OPS_NO_SIMPLIFY = 4
1177 };
1178
1192 void GEOGRAM_API mesh_boolean_operation(
1193 Mesh& result, const Mesh& A, const Mesh& B, const std::string& operation,
1194 MeshBooleanOperationFlags flags = MESH_BOOL_OPS_DEFAULT
1195 );
1196
1208 Mesh& result, const Mesh& A, const Mesh& B, const std::string& operation,
1209 bool verbose
1210 ) {
1212 result, A, B, operation,
1213 verbose ? MESH_BOOL_OPS_VERBOSE : MESH_BOOL_OPS_DEFAULT
1214 );
1215 }
1216
1229 inline void mesh_union(
1230 Mesh& result, const Mesh& A, const Mesh& B,
1231 MeshBooleanOperationFlags flags=MESH_BOOL_OPS_DEFAULT
1232 ) {
1233 mesh_boolean_operation(result, A, B, "A+B", flags);
1234 }
1235
1245 inline void mesh_union(
1246 Mesh& result, const Mesh& A, const Mesh& B, bool verbose
1247 ) {
1248 mesh_boolean_operation(result, A, B, "A+B", verbose);
1249 }
1250
1251
1265 Mesh& result, const Mesh& A, const Mesh& B,
1266 MeshBooleanOperationFlags flags=MESH_BOOL_OPS_DEFAULT
1267 ) {
1268 mesh_boolean_operation(result, A, B, "A*B", flags);
1269 }
1270
1281 Mesh& result, const Mesh& A, const Mesh& B, bool verbose
1282 ) {
1283 mesh_boolean_operation(result, A, B, "A*B", verbose);
1284 }
1285
1298 inline void mesh_difference(
1299 Mesh& result, const Mesh& A, const Mesh& B,
1300 MeshBooleanOperationFlags flags=MESH_BOOL_OPS_DEFAULT
1301 ) {
1302 mesh_boolean_operation(result, A, B, "A-B", flags);
1303 }
1304
1314 inline void mesh_difference(
1315 Mesh& result, const Mesh& A, const Mesh& B, bool verbose
1316 ) {
1317 mesh_boolean_operation(result, A, B, "A-B", verbose);
1318 }
1319
1327 Mesh& M, index_t max_iter = 3, bool verbose=false
1328 );
1329
1340 Mesh& M, index_t f1, index_t f2
1341 );
1342
1343 /**************************************************************************/
1344}
1345
1346#endif
#define geo_assert_not_reached
Sets a non reachable point in the program.
Definition assert.h:177
#define geo_debug_assert(x)
Verifies that a condition is met.
Definition assert.h:196
Generic mechanism for attributes.
Manages an attribute attached to a set of object.
Detects and retriangulates a set of coplanar facets for MeshSurfaceIntersection.
index_t adjacent_facet(index_t c) const
Gets the facet that a corner is adjacent to.
Definition mesh.h:970
void set_adjacent_facet(index_t c, index_t f)
Sets the facet that a corner is adjacent to.
Definition mesh.h:1033
index_t vertex(index_t f, index_t lv) const
Gets a vertex by facet and local vertex index.
Definition mesh.h:1162
index_range corners(index_t f) const
Gets the corners of a facet.
Definition mesh.h:1557
Meshes a single triangle with the constraints that come from the intersections with the other triangl...
AttributesManager & attributes() const
Gets the attributes manager.
Definition mesh.h:109
index_t nb() const
Gets the number of (sub-)elements.
Definition mesh.h:98
Halfedfge-like API wrappers on top of a triangulated mesh.
index_t facet(index_t h) const
Gets the facet associated to a halfedge.
void sew2(index_t h1, index_t h2)
Creates a surfacic link between two halfedges.
void sew3(index_t h1, index_t h2)
Creates a volumetric link between two halfedges.
index_t alpha3(index_t h) const
gets the volumetric neighbor of a halfedge
index_t vertex(index_t h, index_t dlv) const
gets a vertex of an halfedge
index_t alpha2(index_t h) const
gets the surfacic neighbor of a halfedge
index_t facet_alpha3(index_t f) const
gets the volumetric neighbor of a facet
Halfedges(MeshSurfaceIntersection &I)
Halfedges constructor.
index_as_iterator begin() const
used by range-based for
index_t nb() const
Gets the number of halfedegs in the map.
index_as_iterator end() const
used by range-based for
Represents the set of radial halfedge bundles.
index_as_iterator begin() const
used by range-based for
index_t nb() const
Gets the number of bundles.
index_t vertex(index_t bndl, index_t lv) const
gets one of the vertices at the two extremities of a bundle
index_t vertex_first_bundle(index_t v) const
gets the first bundle starting from a vertex
void get_sorted_incident_charts(index_t bndl, vector< ChartPos > &chart_pos)
Gets the sorted list of charts around bundle.
index_as_iterator end() const
used by range-based for
const_index_ptr_range halfedges(index_t bndl) const
gets the halfedges in a bundle
index_t next_around_vertex(index_t bndl) const
gets the next bundle around a vertex
bool radial_sort(index_t bndl, RadialSort &RS)
Sorts the halfedges of the bundle in-place.
index_ptr_range halfedges(index_t bndl)
gets the halfedges in a bundle
index_t prev_along_polyline(index_t bndl)
gets the predecessor of a bundle along its polyline
index_t nb_halfedges(index_t bndl) const
Gets the number of halfedges in a bundle.
void initialize()
Initializes the structure.
index_t next_along_polyline(index_t bndl)
gets the successor of a bundle along its polyline
void set_sorted_halfedges(index_t bndl, const vector< index_t > &halfedges)
Sets the halfedges of a bundle.
std::pair< index_t, index_t > ChartPos
Indicates where to find a chart in a bundle.
index_t nb_bundles_around_vertex(index_t v) const
gets the bumber of bundles around a vertex
RadialBundles(MeshSurfaceIntersection &I)
RadialBundles constructor.
index_t opposite(index_t bndl)
gets the opposite bundle
index_t halfedge(index_t bndl, index_t li) const
Gets a halfedge in a bundle from local index.
void set_halfedge(index_t bndl, index_t li, index_t h)
Sets a halfedge in a bundle.
void radial_sort()
Sorts all the bundles of all polylines.
index_as_iterator begin() const
used by range-based for
const_index_ptr_range bundles(index_t polyline) const
gets the bundles in a polyline
void get_skeleton(Mesh &to, bool trim_fins=false)
Copies the set of polylines to a mesh.
void initialize()
Initializes the structure.
index_t nb() const
Gets the number of polylines.
index_as_iterator end() const
used by range-based for
RadialPolylines(MeshSurfaceIntersection &I)
RadialPolylines constructor.
Sign h_orient(index_t h1, index_t h2) const
Computes the relative orientations of two halfedges.
Sign h_refNorient(index_t h2) const
Computes the normal orientation of a halfedge relative to h_ref.
void init(index_t h_ref)
Initializes radial sorting around a given halfedge.
bool degenerate() const
Tests if a degeneracy was encountered.
RadialSort(const MeshSurfaceIntersection &I)
RadialSort constructor.
void report_problem(const char *message) const
This function is called whenever radial sort encounters a configuration not supposed to happen....
exact::vec3 normal(index_t h) const
Computes the normal to a facet with exact coordinates.
bool operator()(index_t h1, index_t h2) const
Compares two halfedges.
Computes surface intersections.
void set_verbose(bool x)
Display information while computing the intersection. Default is unset.
index_t tentatively_classify_component_vertex_fast(index_t component, index_t v)
Classifies a vertex of the computed intersection.
void set_monster_threshold(index_t nb)
Sets the threshold from which triangle is considered to be a monster.
void simplify_coplanar_facets(double angle_tolerance=0.0)
Merge coplanar facets and retriangulate them using a Constrained Delaunay triangulation.
void remove_internal_shells()
Removes all the facets that are not on the outer boundary.
index_t classify_component(index_t component, index_t v)
Classifies a connected component.
void unlock()
Releases the lock associated with this mesh.
void set_detect_intersecting_neighbors(bool x)
detect and compute intersections between facets that share a facet or an edge. Set to false if input ...
void set_dry_run(bool x)
In dry run mode, the computed local triangulations are not inserted in the global mesh....
index_t find_or_create_exact_vertex(const ExactPoint &p)
Finds or creates a vertex in the mesh, by exact coordinates.
void set_fine_verbose(bool x)
Display detailed information while computing the intersection. Default is unset.
void set_normalize(bool x)
Specifies whether coordinates should be normalized during computation. If set, original coordinates a...
void intersect_prologue()
substep of intersect(), prepares the mesh
void remove_external_shell()
Removes all the facets that are on the outer boundary.
bool initial_facet_is_flipped(index_t f) const
Indicates whether the initial facet (in mesh_copy_) that supports a facet in the intersection mesh (i...
const Mesh & readonly_mesh() const
Gets a copy of the initial mesh passed to the constructor.
void set_build_skeleton(Mesh *skeleton, bool trim_fins=false)
Optionally save the skeleton (that is, the collection of non-manifold edges) to a given mesh....
const Mesh & target_mesh() const
Gets the target mesh.
ExactPoint exact_vertex(index_t v) const
Gets the exact point associated with a vertex.
void set_radial_sort(bool x)
Specifies whether surfaces should be duplicated and radial edges sorted in order to create the volume...
void intersect_epilogue(const vector< IsectInfo > &intersections)
subset of intersect(), cleans the resulting mesh and undoes optional geometric normalization.
void intersect_get_intersections(vector< IsectInfo > &intersections)
substep of intersect(), finds all the intersection points and segments.
void classify(const std::string &expr)
Classifies the facets and keep only the ones on the boundary of a combination of regions defined by a...
index_t tentatively_classify_component_vertex(index_t component, index_t v)
Classifies a vertex of the computed intersection.
void lock()
Acquires a lock on this mesh.
bool is_original_vertex(index_t v) const
Tests whether a given vertex is an original mesh vertex or an intersection.
void intersect_remesh_intersections(vector< IsectInfo > &intersections)
substep of intersect(), inserts the intersection points and segments into the triangles.
index_t get_initial_facet(index_t f) const
Gets the vertices of the initial facet (in mesh_copy_) that supports a facet in the intersection mesh...
void mark_external_shell(vector< index_t > &on_external_shell)
Marks all the facets that are on the external shell.
void build_Weiler_model()
Builds the Weiler model.
void set_interpolate_attributes(bool x)
Specifies that attributes should be interpolated.
Mesh & target_mesh()
Gets the target mesh.
void set_delaunay(bool x)
If set, compute constrained Delaunay triangulation in the intersected triangles. If there are interse...
std::tuple< vec3, vec3, vec3 > get_initial_facet_vertices(index_t f) const
Gets the vertices of the initial facet (in mesh_copy_) that supports a facet in the intersection mesh...
Represents a mesh.
Definition mesh.h:3132
Wraps an integer to be used with the range-based for construct.
Definition range.h:66
Comparator class for vec3Hg \detail Used to create maps indexed by vec3Hg or SOS symbolic perturbatio...
Definition vechg.h:293
3d vector with homogeneous coordinates
Definition vechg.h:185
Vector with aligned memory allocation.
Definition memory.h:710
index_t size() const
Gets the number of elements.
Definition memory.h:756
Exact predicates and constructs.
Common include file, providing basic definitions. Should be included before anything else by all head...
The class that represents a mesh.
Functions to load and save meshes.
SOSMode
Mode for symbolic perturbations.
Definition predicates.h:67
std::atomic_flag spinlock
A lightweight synchronization structure.
Global Vorpaline namespace.
Definition basic.h:55
bool mesh_facets_have_intersection(Mesh &M, index_t f1, index_t f2)
Tests whether two mesh facets have a non-degenerate intersection.
void mesh_intersection(Mesh &result, const Mesh &A, const Mesh &B, MeshBooleanOperationFlags flags=MESH_BOOL_OPS_DEFAULT)
Computes the intersection of two surface meshes.
void mesh_difference(Mesh &result, const Mesh &A, const Mesh &B, MeshBooleanOperationFlags flags=MESH_BOOL_OPS_DEFAULT)
Computes the difference of two surface meshes.
void mesh_boolean_operation(Mesh &result, const Mesh &A, const Mesh &B, const std::string &operation, MeshBooleanOperationFlags flags=MESH_BOOL_OPS_DEFAULT)
Computes a boolean operation with two surface meshes.
void mesh_union(Mesh &result, const Mesh &A, const Mesh &B, MeshBooleanOperationFlags flags=MESH_BOOL_OPS_DEFAULT)
Computes the union of two surface meshes.
geo_index_t index_t
The type for storing and manipulating indices.
Definition numeric.h:330
Sign
Integer constants that represent the sign of a value.
Definition numeric.h:69
void mesh_remove_intersections(Mesh &M, index_t max_iter=3, bool verbose=false)
Attempts to make a surface mesh conformal by removing intersecting facets and re-triangulating the ho...
Function and classes for process manipulation.