Graphite Version 3
An experimental 3D geometry processing program
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
58namespace GEO {
59
60 struct IsectInfo;
61
62 /********************************************************************/
63
68 class GEOGRAM_API MeshSurfaceIntersection {
69 public:
70
72
75
83 void intersect();
84
90
96
97
98 void remove_fins();
99
123 void classify(const std::string& expr);
124
125
134
148 index_t component, index_t v
149 );
150
162 index_t component, index_t v
163 );
164
165
172 void simplify_coplanar_facets(double angle_tolerance = 0.0);
173
178 void set_verbose(bool x) {
179 verbose_ = x;
180 fine_verbose_ = x;
181 }
182
191 monster_threshold_ = nb;
192 }
193
199 void set_dry_run(bool x) {
200 dry_run_ = x;
201 }
202
209 void set_delaunay(bool x) {
210 delaunay_ = x;
211 }
212
219 detect_intersecting_neighbors_ = x;
220 }
221
229 void set_radial_sort(bool x) {
230 use_radial_sort_ = x;
231 }
232
240 void set_normalize(bool x) {
241 normalize_ = x;
242 }
243
244
254 void set_build_skeleton(Mesh* skeleton, bool trim_fins=false) {
255 skeleton_ = skeleton;
256 skeleton_trim_fins_ = trim_fins;
257 }
258
265 interpolate_attributes_ = x;
266 }
267
268 protected:
278
293
312
328 void intersect_epilogue(const vector<IsectInfo>& intersections);
329
330
339 void lock() {
340 Process::acquire_spinlock(lock_);
341 }
342
346 void unlock() {
347 Process::release_spinlock(lock_);
348 }
349
360
372
379 return mesh_;
380 }
381
387 const Mesh& target_mesh() const {
388 return mesh_;
389 }
390
403 const Mesh& readonly_mesh() const {
404 return mesh_copy_;
405 }
406
407 class RadialSort;
408
419
423 void mark_external_shell(vector<index_t>& on_external_shell);
424
425 protected:
426
430 class GEOGRAM_API RadialSort {
431 public:
437 const MeshSurfaceIntersection& mesh
438 ) : mesh_(mesh),
439 h_ref_(index_t(-1)),
440 degenerate_(false)
441 {
442 }
443
448 void init(index_t h_ref);
449
456 bool operator()(index_t h1, index_t h2) const;
457
464 bool degenerate() const {
465 return degenerate_;
466 }
467
476 const ExactPoint& p1, const ExactPoint& p2
477 );
478
487 const ExactPoint& p1, const ExactPoint& p2
488 );
489
490 protected:
491
501
509
510 public:
511 void test(index_t h1, index_t h2) {
512 (*this)(h1,h2);
513 Sign o_ref1 = h_orient(h_ref_, h1);
514 Sign o_ref2 = h_orient(h_ref_, h2);
515 Sign oN_ref1 = h_refNorient(h1);
516 Sign oN_ref2 = h_refNorient(h2);
517 Sign o_12 = h_orient(h1,h2);
518 std::cerr
519 << " o_ref1=" << int(o_ref1) << " o_ref2=" << int(o_ref2)
520 << " oN_ref1=" << int(oN_ref1) << " oN_ref2=" << int(oN_ref2)
521 << " o_12=" << int(o_12)
522 << std::endl;
523 }
524
525
526 private:
527 const MeshSurfaceIntersection& mesh_;
528 index_t h_ref_; // ---reference halfedge
529 exact::vec3 U_ref_; // -.
530 exact::vec3 V_ref_; // +-reference basis
531 exact::vec3 N_ref_; // -'
532 vec3I U_ref_I_; // -.
533 vec3I V_ref_I_; // +-reference basis (interval arithmetics)
534 vec3I N_ref_I_; // _'
535 mutable vector< std::pair<index_t, Sign> > refNorient_cache_;
536 mutable bool degenerate_;
537 };
538
539
540 protected:
541 Process::spinlock lock_;
542 Mesh& mesh_;
543 Mesh mesh_copy_;
544 Attribute<const ExactPoint*> vertex_to_exact_point_;
545
546#ifdef GEOGRAM_USE_EXACT_NT
547 // Exact points are canonicalized
548 // (by Numeric::optimize_number_representation(vec3HEx)) so
549 // we can use this comparator that makes the global vertex map
550 // much much faster.
551 typedef vec3HExLexicoCompareCanonical ExactPointCompare;
552#else
553 // Generic comparator for global vertex map.
554 typedef vec3HgLexicoCompare<exact::scalar> ExactPointCompare;
555#endif
556 std::map<ExactPoint,index_t,ExactPointCompare> exact_point_to_vertex_;
557
558 bool verbose_;
559 bool fine_verbose_;
560 bool delaunay_;
561 bool detect_intersecting_neighbors_;
562 bool use_radial_sort_;
563
564 PCK::SOSMode SOS_bkp_;
565 bool rescale_;
566 bool normalize_;
567 vec3 normalize_center_;
568 double normalize_radius_;
569
570 index_t monster_threshold_;
571 bool dry_run_;
572 friend class MeshInTriangle;
573 friend class CoplanarFacets;
574
575 Mesh* skeleton_;
576 bool skeleton_trim_fins_;
577 bool interpolate_attributes_;
578
579 /***************************************************/
580
589 class Halfedges {
590 public:
591
596 Halfedges(MeshSurfaceIntersection& I) : mesh_(I.mesh_) {
597 }
598
603 // TODO: destroy alpha3 attribute (kept now for debugging
604 }
605
610 void initialize() {
611 facet_corner_alpha3_.bind(
612 mesh_.facet_corners.attributes(), "alpha3"
613 );
614 }
615
621 index_t nb() const {
622 return mesh_.facet_corners.nb();
623 }
624
630 return index_as_iterator(0);
631 }
632
638 return index_as_iterator(nb());
639 }
640
647 return h/3;
648 }
649
660 index_t t1 = facet(h);
661 index_t t2 = mesh_.facet_corners.adjacent_facet(h);
662 if(t2 == NO_INDEX) {
663 return NO_INDEX;
664 }
665 for(index_t h2: mesh_.facets.corners(t2)) {
666 if(mesh_.facet_corners.adjacent_facet(h2) == t1) {
667 return h2;
668 }
669 }
671 }
672
683 return facet_corner_alpha3_[h];
684 }
685
693 return alpha3(3*f)/3;
694 }
695
707 index_t f = h/3;
708 index_t lv = (h+dlv)%3;
709 return mesh_.facets.vertex(f,lv);
710 }
711
712
721 void sew2(index_t h1, index_t h2) {
722 geo_debug_assert(vertex(h1,0) == vertex(h2,1));
723 geo_debug_assert(vertex(h2,0) == vertex(h1,1));
724 index_t t1 = h1/3;
725 index_t t2 = h2/3;
726 mesh_.facet_corners.set_adjacent_facet(h1,t2);
727 mesh_.facet_corners.set_adjacent_facet(h2,t1);
728 }
729
738 void sew3(index_t h1, index_t h2) {
739 geo_debug_assert(vertex(h1,0) == vertex(h2,1));
740 geo_debug_assert(vertex(h2,0) == vertex(h1,1));
741 facet_corner_alpha3_[h1] = h2;
742 facet_corner_alpha3_[h2] = h1;
743 }
744
745 private:
746 Mesh& mesh_;
747 Attribute<index_t> facet_corner_alpha3_;
748 } halfedges_;
749
750 /***************************************************/
751
758 public:
759
764 RadialBundles(MeshSurfaceIntersection& I) : I_(I), mesh_(I.mesh_) {
765 }
766
772
776 index_t nb() const {
777 return bndl_start_.size() - 1;
778 }
779
785 return index_as_iterator(0);
786 }
787
793 return index_as_iterator(nb());
794 }
795
802 geo_debug_assert(bndl < nb());
803 return bndl_start_[bndl+1] - bndl_start_[bndl];
804 }
805
813 index_t halfedge(index_t bndl, index_t li) const {
814 geo_debug_assert(bndl < nb());
815 geo_debug_assert(li < nb_halfedges(bndl));
816 return H_[bndl_start_[bndl] + li];
817 }
818
827 geo_debug_assert(bndl < nb());
828 geo_debug_assert(li < nb_halfedges(bndl));
829 H_[bndl_start_[bndl] + li] = h;
830 }
831
838 return index_ptr_range(
839 H_, bndl_start_[bndl], bndl_start_[bndl+1]
840 );
841 }
842
850 H_, bndl_start_[bndl], bndl_start_[bndl+1]
851 );
852 }
853
861 index_t vertex(index_t bndl, index_t lv) const {
862 geo_debug_assert(bndl_start_[bndl+1] - bndl_start_[bndl] > 0);
863 index_t h = H_[bndl_start_[bndl]];
864 return I_.halfedges_.vertex(h,lv);
865 }
866
875 return v_first_bndl_[v];
876 }
877
886 return bndl_next_around_v_[bndl];
887 }
888
895 index_t result = 0;
896 for(
897 index_t bndl = vertex_first_bundle(v);
898 bndl != NO_INDEX;
899 bndl = next_around_vertex(bndl)
900 ) {
901 ++result;
902 }
903 return result;
904 }
905
913 geo_debug_assert(bndl < nb());
914 return (bndl >= nb()/2) ? (bndl-nb()/2) : (bndl+nb()/2);
915 }
916
923 index_t v = vertex(bndl,0);
924 if(nb_bundles_around_vertex(v) != 2) {
925 return NO_INDEX;
926 }
927 for(
928 index_t bndl2 = vertex_first_bundle(v);
929 bndl2 != NO_INDEX; bndl2 = next_around_vertex(bndl2)
930 ) {
931 if(bndl2 != bndl) {
932 return opposite(bndl2);
933 }
934 }
936 }
937
944 index_t v = vertex(bndl,1);
945 if(nb_bundles_around_vertex(v) != 2) {
946 return NO_INDEX;
947 }
948 for(
949 index_t bndl2 = vertex_first_bundle(v);
950 bndl2 != NO_INDEX; bndl2 = next_around_vertex(bndl2)
951 ) {
952 if(opposite(bndl2) != bndl) {
953 return bndl2;
954 }
955 }
957 }
958
968 geo_debug_assert(bndl < nb());
969 if(nb_halfedges(bndl) <= 2) {
970 return true;
971 }
972 auto b = H_.begin() + std::ptrdiff_t(bndl_start_[bndl]);
973 auto e = H_.begin() + std::ptrdiff_t(bndl_start_[bndl+1]);
974 RS.init(*b);
975 std::sort(
976 b, e,
977 [&](index_t h1, index_t h2)->bool {
978 return RS(h1,h2);
979 }
980 );
981 bool OK = !RS.degenerate();
982 bndl_is_sorted_[bndl] = OK;
983 return OK;
984 }
985
995 index_t bndl, const vector<index_t>& halfedges
996 ) {
997 geo_debug_assert(halfedges.size() == nb_halfedges(bndl));
998 for(index_t i=0; i<halfedges.size(); ++i) {
999 set_halfedge(bndl, i, halfedges[i]);
1000 }
1001 bndl_is_sorted_[bndl] = true;
1002 }
1003
1009 typedef std::pair<index_t, index_t> ChartPos;
1010
1019 index_t bndl, vector<ChartPos>& chart_pos
1020 );
1021
1022 bool is_sorted(index_t bndl) const {
1023 geo_assert(bndl < nb());
1024 return bndl_is_sorted_[bndl];
1025 }
1026
1027 // private:
1029 Mesh& mesh_;
1030 Attribute<index_t> facet_chart_;
1031 vector<index_t> H_;
1032 vector<index_t> bndl_start_;
1033 vector<index_t> v_first_bndl_;
1034 vector<index_t> bndl_next_around_v_;
1035 vector<bool> bndl_is_sorted_;
1036 } radial_bundles_;
1037
1038 /***************************************************/
1039
1041 public:
1046 RadialPolylines(MeshSurfaceIntersection& I) : I_(I), mesh_(I.mesh_) {
1047 }
1048
1054
1061
1065 index_t nb() const {
1066 return polyline_start_.size() - 1;
1067 }
1068
1074 return index_as_iterator(0);
1075 }
1076
1082 return index_as_iterator(nb());
1083 }
1084
1091 geo_debug_assert(polyline < nb());
1092 return const_index_ptr_range(
1093 B_, polyline_start_[polyline], polyline_start_[polyline+1]
1094 );
1095 }
1096
1097 index_t nb_bundles(index_t polyline) const {
1098 geo_debug_assert(polyline < nb());
1099 return polyline_start_[polyline+1] - polyline_start_[polyline];
1100 }
1101
1102 index_t bundle(index_t polyline, index_t li) const {
1103 geo_debug_assert(polyline < nb());
1104 geo_debug_assert(li < nb_bundles(polyline));
1105 return B_[polyline_start_[polyline] + li];
1106 }
1107
1115 void get_skeleton(Mesh& to, bool trim_fins=false);
1116
1117 private:
1119 Mesh& mesh_;
1120 vector<index_t> B_;
1121 vector<index_t> polyline_start_;
1122 } radial_polylines_;
1123 };
1124
1125 /********************************************************************/
1126
1127 enum MeshBooleanOperationFlags {
1128 MESH_BOOL_OPS_DEFAULT = 0,
1129 MESH_BOOL_OPS_VERBOSE = 1,
1130 MESH_BOOL_OPS_ATTRIBS = 2,
1131 MESH_BOOL_OPS_NO_SIMPLIFY = 4
1132 };
1133
1147 void GEOGRAM_API mesh_boolean_operation(
1148 Mesh& result, const Mesh& A, const Mesh& B, const std::string& operation,
1149 MeshBooleanOperationFlags flags = MESH_BOOL_OPS_DEFAULT
1150 );
1151
1163 Mesh& result, const Mesh& A, const Mesh& B, const std::string& operation,
1164 bool verbose
1165 ) {
1167 result, A, B, operation,
1168 verbose ? MESH_BOOL_OPS_VERBOSE : MESH_BOOL_OPS_DEFAULT
1169 );
1170 }
1171
1184 inline void mesh_union(
1185 Mesh& result, const Mesh& A, const Mesh& B,
1186 MeshBooleanOperationFlags flags=MESH_BOOL_OPS_DEFAULT
1187 ) {
1188 mesh_boolean_operation(result, A, B, "A+B", flags);
1189 }
1190
1200 inline void mesh_union(
1201 Mesh& result, const Mesh& A, const Mesh& B, bool verbose
1202 ) {
1203 mesh_boolean_operation(result, A, B, "A+B", verbose);
1204 }
1205
1206
1220 Mesh& result, const Mesh& A, const Mesh& B,
1221 MeshBooleanOperationFlags flags=MESH_BOOL_OPS_DEFAULT
1222 ) {
1223 mesh_boolean_operation(result, A, B, "A*B", flags);
1224 }
1225
1236 Mesh& result, const Mesh& A, const Mesh& B, bool verbose
1237 ) {
1238 mesh_boolean_operation(result, A, B, "A*B", verbose);
1239 }
1240
1253 inline void mesh_difference(
1254 Mesh& result, const Mesh& A, const Mesh& B,
1255 MeshBooleanOperationFlags flags=MESH_BOOL_OPS_DEFAULT
1256 ) {
1257 mesh_boolean_operation(result, A, B, "A-B", flags);
1258 }
1259
1269 inline void mesh_difference(
1270 Mesh& result, const Mesh& A, const Mesh& B, bool verbose
1271 ) {
1272 mesh_boolean_operation(result, A, B, "A-B", verbose);
1273 }
1274
1282 Mesh& M, index_t max_iter = 3, bool verbose=false
1283 );
1284
1295 Mesh& M, index_t f1, index_t f2
1296 );
1297
1298 /**************************************************************************/
1299}
1300
1301#endif
#define geo_assert_not_reached
Sets a non reachable point in the program.
Definition assert.h:177
#define geo_assert(x)
Verifies that a condition is met.
Definition assert.h:149
#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:887
void set_adjacent_facet(index_t c, index_t f)
Sets the facet that a corner is adjacent to.
Definition mesh.h:950
index_t vertex(index_t f, index_t lv) const
Gets a vertex by facet and local vertex index.
Definition mesh.h:1054
index_range corners(index_t f) const
Gets the corners of a facet.
Definition mesh.h:1420
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:100
index_t nb() const
Gets the number of (sub-)elements.
Definition mesh.h:89
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.
static exact::vec3 exact_direction(const ExactPoint &p1, const ExactPoint &p2)
Computes a vector of arbitrary length with its direction given by two points.
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 &mesh)
RadialSort constructor.
bool operator()(index_t h1, index_t h2) const
Compares two halfedges.
static vec3I exact_direction_I(const ExactPoint &p1, const ExactPoint &p2)
Computes an interval vector of arbitrary length with its direction given by two points.
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_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.
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.
void intersect_remesh_intersections(vector< IsectInfo > &intersections)
substep of intersect(), inserts the intersection points and segments into the triangles.
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...
Represents a mesh.
Definition mesh.h:2701
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:660
index_t size() const
Gets the number of elements.
Definition memory.h:706
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.
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:329
Sign
Integer constants that represent the sign of a value.
Definition numeric.h:68
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.