Graphite  Version 3
An experimental 3D geometry processing program
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 
43 #include <geogram/basic/common.h>
44 #include <geogram/mesh/mesh.h>
45 #include <geogram/mesh/mesh_io.h>
47 #include <geogram/basic/process.h>
49 #include <geogram/basic/debug_stream.h>
50 #include <functional>
51 
58 namespace GEO {
59 
60  struct IsectInfo;
61 
62  /********************************************************************/
63 
68  class GEOGRAM_API MeshSurfaceIntersection {
69  public:
70 
71  typedef exact::vec3h ExactPoint;
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 
253  void set_build_skeleton(Mesh* skeleton, bool trim_fins=false) {
254  skeleton_ = skeleton;
255  skeleton_trim_fins_ = trim_fins;
256  }
257 
264  interpolate_attributes_ = x;
265  }
266 
267  protected:
277 
292 
311 
327  void intersect_epilogue(const vector<IsectInfo>& intersections);
328 
329 
338  void lock() {
340  }
341 
345  void unlock() {
347  }
348 
359 
371 
378  return mesh_;
379  }
380 
386  const Mesh& target_mesh() const {
387  return mesh_;
388  }
389 
402  const Mesh& readonly_mesh() const {
403  return mesh_copy_;
404  }
405 
406  class RadialSort;
407 
418 
422  void mark_external_shell(vector<index_t>& on_external_shell);
423 
424  protected:
425 
429  class GEOGRAM_API RadialSort {
430  public:
436  const MeshSurfaceIntersection& mesh
437  ) : mesh_(mesh),
438  h_ref_(index_t(-1)),
439  degenerate_(false)
440  {
441  }
442 
447  void init(index_t h_ref);
448 
455  bool operator()(index_t h1, index_t h2) const;
456 
463  bool degenerate() const {
464  return degenerate_;
465  }
466 
475  const ExactPoint& p1, const ExactPoint& p2
476  );
477 
486  const ExactPoint& p1, const ExactPoint& p2
487  );
488 
489  protected:
490 
499  Sign h_orient(index_t h1, index_t h2) const;
500 
508 
509  public:
510  void test(index_t h1, index_t h2) {
511  (*this)(h1,h2);
512  Sign o_ref1 = h_orient(h_ref_, h1);
513  Sign o_ref2 = h_orient(h_ref_, h2);
514  Sign oN_ref1 = h_refNorient(h1);
515  Sign oN_ref2 = h_refNorient(h2);
516  Sign o_12 = h_orient(h1,h2);
517  std::cerr
518  << " o_ref1=" << int(o_ref1) << " o_ref2=" << int(o_ref2)
519  << " oN_ref1=" << int(oN_ref1) << " oN_ref2=" << int(oN_ref2)
520  << " o_12=" << int(o_12)
521  << std::endl;
522  }
523 
524 
525  private:
526  const MeshSurfaceIntersection& mesh_;
527  index_t h_ref_; // ---reference halfedge
528  exact::vec3 U_ref_; // -.
529  exact::vec3 V_ref_; // +-reference basis
530  exact::vec3 N_ref_; // -'
531  vec3I U_ref_I_; // -.
532  vec3I V_ref_I_; // +-reference basis (interval arithmetics)
533  vec3I N_ref_I_; // _'
534  mutable vector< std::pair<index_t, Sign> > refNorient_cache_;
535  mutable bool degenerate_;
536  };
537 
538 
539  protected:
540  Process::spinlock lock_;
541  Mesh& mesh_;
542  Mesh mesh_copy_;
543  Attribute<const ExactPoint*> vertex_to_exact_point_;
544 
545 #ifdef GEOGRAM_USE_EXACT_NT
546  // Exact points are canonicalized
547  // (by Numeric::optimize_number_representation(vec3HEx)) so
548  // we can use this comparator that makes the global vertex map
549  // much much faster.
550  typedef vec3HExLexicoCompareCanonical ExactPointCompare;
551 #else
552  // Generic comparator for global vertex map.
553  typedef vec3HgLexicoCompare<exact::scalar> ExactPointCompare;
554 #endif
555  std::map<ExactPoint,index_t,ExactPointCompare> exact_point_to_vertex_;
556 
557  bool verbose_;
558  bool fine_verbose_;
559  bool delaunay_;
560  bool detect_intersecting_neighbors_;
561  bool use_radial_sort_;
562 
563  PCK::SOSMode SOS_bkp_;
564  bool rescale_;
565  bool normalize_;
566  vec3 normalize_center_;
567  double normalize_radius_;
568 
569  index_t monster_threshold_;
570  bool dry_run_;
571  friend class MeshInTriangle;
572  friend class CoplanarFacets;
573 
574  Mesh* skeleton_;
575  bool skeleton_trim_fins_;
576  bool interpolate_attributes_;
577 
578  /***************************************************/
579 
588  class Halfedges {
589  public:
590 
595  Halfedges(MeshSurfaceIntersection& I) : mesh_(I.mesh_) {
596  }
597 
602  // TODO: destroy alpha3 attribute (kept now for debugging
603  }
604 
609  void initialize() {
610  facet_corner_alpha3_.bind(
611  mesh_.facet_corners.attributes(), "alpha3"
612  );
613  }
614 
620  index_t nb() const {
621  return mesh_.facet_corners.nb();
622  }
623 
629  return index_as_iterator(0);
630  }
631 
637  return index_as_iterator(nb());
638  }
639 
645  index_t facet(index_t h) const {
646  return h/3;
647  }
648 
658  index_t alpha2(index_t h) const {
659  index_t t1 = facet(h);
660  index_t t2 = mesh_.facet_corners.adjacent_facet(h);
661  if(t2 == NO_INDEX) {
662  return NO_INDEX;
663  }
664  for(index_t h2: mesh_.facets.corners(t2)) {
665  if(mesh_.facet_corners.adjacent_facet(h2) == t1) {
666  return h2;
667  }
668  }
670  }
671 
681  index_t alpha3(index_t h) const {
682  return facet_corner_alpha3_[h];
683  }
684 
692  return alpha3(3*f)/3;
693  }
694 
705  index_t vertex(index_t h, index_t dlv) const {
706  index_t f = h/3;
707  index_t lv = (h+dlv)%3;
708  return mesh_.facets.vertex(f,lv);
709  }
710 
711 
720  void sew2(index_t h1, index_t h2) {
721  geo_debug_assert(vertex(h1,0) == vertex(h2,1));
722  geo_debug_assert(vertex(h2,0) == vertex(h1,1));
723  index_t t1 = h1/3;
724  index_t t2 = h2/3;
725  mesh_.facet_corners.set_adjacent_facet(h1,t2);
726  mesh_.facet_corners.set_adjacent_facet(h2,t1);
727  }
728 
737  void sew3(index_t h1, index_t h2) {
738  geo_debug_assert(vertex(h1,0) == vertex(h2,1));
739  geo_debug_assert(vertex(h2,0) == vertex(h1,1));
740  facet_corner_alpha3_[h1] = h2;
741  facet_corner_alpha3_[h2] = h1;
742  }
743 
744  private:
745  Mesh& mesh_;
746  Attribute<index_t> facet_corner_alpha3_;
747  } halfedges_;
748 
749  /***************************************************/
750 
757  public:
758 
763  RadialBundles(MeshSurfaceIntersection& I) : I_(I), mesh_(I.mesh_) {
764  }
765 
770  void initialize();
771 
775  index_t nb() const {
776  return bndl_start_.size() - 1;
777  }
778 
784  return index_as_iterator(0);
785  }
786 
792  return index_as_iterator(nb());
793  }
794 
801  geo_debug_assert(bndl < nb());
802  return bndl_start_[bndl+1] - bndl_start_[bndl];
803  }
804 
812  index_t halfedge(index_t bndl, index_t li) const {
813  geo_debug_assert(bndl < nb());
814  geo_debug_assert(li < nb_halfedges(bndl));
815  return H_[bndl_start_[bndl] + li];
816  }
817 
825  void set_halfedge(index_t bndl, index_t li, index_t h) {
826  geo_debug_assert(bndl < nb());
827  geo_debug_assert(li < nb_halfedges(bndl));
828  H_[bndl_start_[bndl] + li] = h;
829  }
830 
837  return index_ptr_range(
838  H_, bndl_start_[bndl], bndl_start_[bndl+1]
839  );
840  }
841 
848  return const_index_ptr_range(
849  H_, bndl_start_[bndl], bndl_start_[bndl+1]
850  );
851  }
852 
860  index_t vertex(index_t bndl, index_t lv) const {
861  geo_debug_assert(bndl_start_[bndl+1] - bndl_start_[bndl] > 0);
862  index_t h = H_[bndl_start_[bndl]];
863  return I_.halfedges_.vertex(h,lv);
864  }
865 
874  return v_first_bndl_[v];
875  }
876 
885  return bndl_next_around_v_[bndl];
886  }
887 
894  index_t result = 0;
895  for(
896  index_t bndl = vertex_first_bundle(v);
897  bndl != NO_INDEX;
898  bndl = next_around_vertex(bndl)
899  ) {
900  ++result;
901  }
902  return result;
903  }
904 
912  geo_debug_assert(bndl < nb());
913  return (bndl >= nb()/2) ? (bndl-nb()/2) : (bndl+nb()/2);
914  }
915 
922  index_t v = vertex(bndl,0);
923  if(nb_bundles_around_vertex(v) != 2) {
924  return NO_INDEX;
925  }
926  for(
927  index_t bndl2 = vertex_first_bundle(v);
928  bndl2 != NO_INDEX; bndl2 = next_around_vertex(bndl2)
929  ) {
930  if(bndl2 != bndl) {
931  return opposite(bndl2);
932  }
933  }
935  }
936 
943  index_t v = vertex(bndl,1);
944  if(nb_bundles_around_vertex(v) != 2) {
945  return NO_INDEX;
946  }
947  for(
948  index_t bndl2 = vertex_first_bundle(v);
949  bndl2 != NO_INDEX; bndl2 = next_around_vertex(bndl2)
950  ) {
951  if(opposite(bndl2) != bndl) {
952  return bndl2;
953  }
954  }
956  }
957 
966  bool radial_sort(index_t bndl, RadialSort& RS) {
967  geo_debug_assert(bndl < nb());
968  if(nb_halfedges(bndl) <= 2) {
969  return true;
970  }
971  auto b = H_.begin() + std::ptrdiff_t(bndl_start_[bndl]);
972  auto e = H_.begin() + std::ptrdiff_t(bndl_start_[bndl+1]);
973  RS.init(*b);
974  std::sort(
975  b, e,
976  [&](index_t h1, index_t h2)->bool {
977  return RS(h1,h2);
978  }
979  );
980  bool OK = !RS.degenerate();
981  bndl_is_sorted_[bndl] = OK;
982  return OK;
983  }
984 
994  index_t bndl, const vector<index_t>& halfedges
995  ) {
996  geo_debug_assert(halfedges.size() == nb_halfedges(bndl));
997  for(index_t i=0; i<halfedges.size(); ++i) {
998  set_halfedge(bndl, i, halfedges[i]);
999  }
1000  bndl_is_sorted_[bndl] = true;
1001  }
1002 
1008  typedef std::pair<index_t, index_t> ChartPos;
1009 
1018  index_t bndl, vector<ChartPos>& chart_pos
1019  );
1020 
1021  bool is_sorted(index_t bndl) const {
1022  geo_assert(bndl < nb());
1023  return bndl_is_sorted_[bndl];
1024  }
1025 
1026  // private:
1028  Mesh& mesh_;
1029  Attribute<index_t> facet_chart_;
1030  vector<index_t> H_;
1031  vector<index_t> bndl_start_;
1032  vector<index_t> v_first_bndl_;
1033  vector<index_t> bndl_next_around_v_;
1034  vector<bool> bndl_is_sorted_;
1035  } radial_bundles_;
1036 
1037  /***************************************************/
1038 
1040  public:
1045  RadialPolylines(MeshSurfaceIntersection& I) : I_(I), mesh_(I.mesh_) {
1046  }
1047 
1052  void initialize();
1053 
1059  void radial_sort();
1060 
1064  index_t nb() const {
1065  return polyline_start_.size() - 1;
1066  }
1067 
1073  return index_as_iterator(0);
1074  }
1075 
1081  return index_as_iterator(nb());
1082  }
1083 
1090  geo_debug_assert(polyline < nb());
1091  return const_index_ptr_range(
1092  B_, polyline_start_[polyline], polyline_start_[polyline+1]
1093  );
1094  }
1095 
1096  index_t nb_bundles(index_t polyline) const {
1097  geo_debug_assert(polyline < nb());
1098  return polyline_start_[polyline+1] - polyline_start_[polyline];
1099  }
1100 
1101  index_t bundle(index_t polyline, index_t li) const {
1102  geo_debug_assert(polyline < nb());
1103  geo_debug_assert(li < nb_bundles(polyline));
1104  return B_[polyline_start_[polyline] + li];
1105  }
1106 
1114  void get_skeleton(Mesh& to, bool trim_fins=false);
1115 
1116  private:
1118  Mesh& mesh_;
1119  vector<index_t> B_;
1120  vector<index_t> polyline_start_;
1121  } radial_polylines_;
1122  };
1123 
1124  /********************************************************************/
1125 
1136  void GEOGRAM_API mesh_boolean_operation(
1137  Mesh& result, Mesh& A, Mesh& B, const std::string& operation,
1138  bool verbose=false
1139  );
1140 
1150  void GEOGRAM_API mesh_union(
1151  Mesh& result, Mesh& A, Mesh& B, bool verbose=false
1152  );
1153 
1163  void GEOGRAM_API mesh_intersection(
1164  Mesh& result, Mesh& A, Mesh& B, bool verbose=false
1165  );
1166 
1176  void GEOGRAM_API mesh_difference(
1177  Mesh& result, Mesh& A, Mesh& B, bool verbose=false
1178  );
1179 
1186  void GEOGRAM_API mesh_remove_intersections(
1187  Mesh& M, index_t max_iter = 3, bool verbose=false
1188  );
1189 
1200  Mesh& M, index_t f1, index_t f2
1201  );
1202 
1203  /**************************************************************************/
1204 }
1205 
1206 #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.
Definition: attributes.h:1394
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:1412
Meshes a single triangle with the constraints that come from the intersections with the other triangl...
index_t nb() const
Gets the number of (sub-)elements.
Definition: mesh.h:89
AttributesManager & attributes() const
Gets the attributes manager.
Definition: mesh.h:100
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.
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.
Mesh & target_mesh()
Gets the target mesh.
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.
const Mesh & readonly_mesh() const
Gets a copy of the initial mesh passed to the constructor.
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.
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.
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.
const Mesh & target_mesh() const
Gets the target mesh.
void build_Weiler_model()
Builds the Weiler model.
void set_interpolate_attributes(bool x)
Specifies that attributes should be interpolated.
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:2693
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
Specialization of vector for elements of type bool.
Definition: memory.h:795
Vector with aligned memory allocation.
Definition: memory.h:635
index_t size() const
Gets the number of elements.
Definition: memory.h:674
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.
Definition: thread_sync.h:149
void release_spinlock(volatile spinlock &x)
Makes x available to other threads.
Definition: thread_sync.h:173
void acquire_spinlock(volatile spinlock &x)
Loops until x is available then reserves it.
Definition: thread_sync.h:155
Global Vorpaline namespace.
void mesh_boolean_operation(Mesh &result, Mesh &A, Mesh &B, const std::string &operation, bool verbose=false)
Computes a boolean operation with two surface meshes.
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_union(Mesh &result, Mesh &A, Mesh &B, bool verbose=false)
Computes the union of two surface meshes.
void mesh_intersection(Mesh &result, Mesh &A, Mesh &B, bool verbose=false)
Computes the intersection 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 sort(const ITERATOR &begin, const ITERATOR &end)
Sorts elements in parallel.
Definition: algorithm.h:90
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...
void mesh_difference(Mesh &result, Mesh &A, Mesh &B, bool verbose=false)
Computes the difference of two surface meshes.
Function and classes for process manipulation.