40#ifndef GEOGRAM_MESH_MESH_SURFACE_INTERSECTION
41#define GEOGRAM_MESH_MESH_SURFACE_INTERSECTION
49#include <geogram/basic/debug_stream.h>
181 if(!verbose_ && fine_verbose_) {
182 fine_verbose_ =
false;
192 if(fine_verbose_ && !verbose_) {
205 monster_threshold_ = nb;
233 detect_intersecting_neighbors_ = x;
244 use_radial_sort_ = x;
269 skeleton_ = skeleton;
270 skeleton_trim_fins_ = trim_fins;
279 interpolate_attributes_ = x;
354 Process::acquire_spinlock(lock_);
361 Process::release_spinlock(lock_);
383 return (vertex_to_exact_point_[v] ==
nullptr);
457 return original_facet_id_[f];
470 return f_is_flipped_[f];
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]) {
494 return std::make_tuple(p1,p2,p3);
510 mesh_(I_.target_mesh()),
581 mutable bool degenerate_;
592#if defined(GEOGRAM_USE_EXACT_NT) && defined(GEOGRAM_EXACT_NT_IS_MPF_NT)
602 std::map<ExactPoint,index_t,ExactPointCompare> exact_point_to_vertex_;
607 bool detect_intersecting_neighbors_;
608 bool use_radial_sort_;
613 vec3 normalize_center_;
614 double normalize_radius_;
622 bool skeleton_trim_fins_;
623 bool interpolate_attributes_;
657 facet_corner_alpha3_.bind(
668 return mesh_.facet_corners.
nb();
729 return facet_corner_alpha3_[h];
739 return alpha3(3*f)/3;
755 return mesh_.facets.
vertex(f,lv);
787 facet_corner_alpha3_[h1] = h2;
788 facet_corner_alpha3_[h2] = h1;
823 return bndl_start_.size() - 1;
849 return bndl_start_[bndl+1] - bndl_start_[bndl];
862 return H_[bndl_start_[bndl] + li];
875 H_[bndl_start_[bndl] + li] = h;
885 H_, bndl_start_[bndl], bndl_start_[bndl+1]
896 H_, bndl_start_[bndl], bndl_start_[bndl+1]
909 index_t h = H_[bndl_start_[bndl]];
910 return I_.halfedges_.vertex(h,lv);
921 return v_first_bndl_[v];
932 return bndl_next_around_v_[bndl];
943 index_t bndl = vertex_first_bundle(v);
945 bndl = next_around_vertex(bndl)
960 return (bndl >= nb()/2) ? (bndl-nb()/2) : (bndl+nb()/2);
970 if(nb_bundles_around_vertex(v) != 2) {
974 index_t bndl2 = vertex_first_bundle(v);
975 bndl2 != NO_INDEX; bndl2 = next_around_vertex(bndl2)
978 return opposite(bndl2);
991 if(nb_bundles_around_vertex(v) != 2) {
995 index_t bndl2 = vertex_first_bundle(v);
996 bndl2 != NO_INDEX; bndl2 = next_around_vertex(bndl2)
998 if(opposite(bndl2) != bndl) {
1015 if(nb_halfedges(bndl) <= 2) {
1018 auto b = H_.begin() + std::ptrdiff_t(bndl_start_[bndl]);
1019 auto e = H_.begin() + std::ptrdiff_t(bndl_start_[bndl+1]);
1027 bndl_is_sorted_[bndl] = OK;
1044 set_halfedge(bndl, i, halfedges[i]);
1046 bndl_is_sorted_[bndl] =
true;
1067 bool is_sorted(
index_t bndl)
const {
1069 return bndl_is_sorted_[bndl];
1111 return polyline_start_.size() - 1;
1138 B_, polyline_start_[polyline], polyline_start_[polyline+1]
1144 return polyline_start_[polyline+1] - polyline_start_[polyline];
1150 return B_[polyline_start_[polyline] + li];
1167 } radial_polylines_;
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
1193 Mesh& result,
const Mesh& A,
const Mesh& B,
const std::string& operation,
1194 MeshBooleanOperationFlags flags = MESH_BOOL_OPS_DEFAULT
1208 Mesh& result,
const Mesh& A,
const Mesh& B,
const std::string& operation,
1212 result, A, B, operation,
1213 verbose ? MESH_BOOL_OPS_VERBOSE : MESH_BOOL_OPS_DEFAULT
1231 MeshBooleanOperationFlags flags=MESH_BOOL_OPS_DEFAULT
1246 Mesh& result,
const Mesh& A,
const Mesh& B,
bool verbose
1266 MeshBooleanOperationFlags flags=MESH_BOOL_OPS_DEFAULT
1281 Mesh& result,
const Mesh& A,
const Mesh& B,
bool verbose
1300 MeshBooleanOperationFlags flags=MESH_BOOL_OPS_DEFAULT
1315 Mesh& result,
const Mesh& A,
const Mesh& B,
bool verbose
#define geo_assert_not_reached
Sets a non reachable point in the program.
#define geo_debug_assert(x)
Verifies that a condition is met.
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 vertex(index_t f, index_t lv) const
Gets a vertex by facet and local vertex index.
index_range corners(index_t f) const
Gets the corners of a facet.
Meshes a single triangle with the constraints that come from the intersections with the other triangl...
AttributesManager & attributes() const
Gets the attributes manager.
index_t nb() const
Gets the number of (sub-)elements.
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.
void initialize()
Initializes the structure.
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.
~Halfedges()
Halfedges destructor.
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...
Wraps an integer to be used with the range-based for construct.
Comparator class for vec3Hg \detail Used to create maps indexed by vec3Hg or SOS symbolic perturbatio...
3d vector with homogeneous coordinates
Vector with aligned memory allocation.
index_t size() const
Gets the number of elements.
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.
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.
Sign
Integer constants that represent the sign of a value.
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.