Graphite  Version 3
An experimental 3D geometry processing program
delaunay_3d.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_DELAUNAY_DELAUNAY_3D
41 #define GEOGRAM_DELAUNAY_DELAUNAY_3D
42 
43 #include <geogram/basic/common.h>
45 #include <geogram/delaunay/cavity.h>
47 #include <geogram/basic/geometry.h>
48 
49 #include <stack>
50 
56 namespace GEO {
57 
107  class GEOGRAM_API Delaunay3d : public Delaunay {
108  public:
122  Delaunay3d(coord_index_t dimension = 3);
123 
128  index_t nb_vertices, const double* vertices
129  ) override;
130 
134  index_t nearest_vertex(const double* p) const override;
135 
136 
137  protected:
138 
145  static const index_t NO_TETRAHEDRON = index_t(-1);
146 
160  index_t& iv0, index_t& iv1, index_t& iv2, index_t& iv3
161  );
162 
183  const double* p, index_t hint = NO_TETRAHEDRON,
184  bool thread_safe = false,
185  Sign* orient = nullptr
186  ) const;
187 
206  const double* p, index_t hint, index_t max_iter
207  ) const;
208 
217  index_t insert(index_t v, index_t hint = NO_TETRAHEDRON);
218 
243  index_t v,
244  index_t t, const Sign* orient,
245  index_t& t_bndry, index_t& f_bndry,
246  index_t& first, index_t& last
247  );
248 
265  const double* p, index_t t,
266  index_t& t_bndry, index_t& f_bndry,
267  index_t& first, index_t& last
268  );
269 
280 
281 
300  index_t v,
301  index_t t_bndry, index_t f_bndry,
302  index_t prev_f=index_t(-1)
303  );
304 
326  index_t t1,
327  index_t t1fborder,
328  index_t t1ft2,
329  index_t& t2,
330  index_t& t2fborder,
331  index_t& t2ft1
332  ) const {
333 
334  // Note: this function is a bit long for an inline function,
335  // but I observed a (modest) performance gain doing so.
336 
337  // Find two vertices that are both on facets new_f and f1
338  // (the edge around which we are turning)
339  // This uses duality as follows:
340  // Primal form (not used here):
341  // halfedge_facet_[v1][v2] returns a facet that is incident
342  // to both v1 and v2.
343  // Dual form (used here):
344  // halfedge_facet_[f1][f2] returns a vertex that both
345  // f1 and f2 are incident to.
346  signed_index_t ev1 =
347  tet_vertex(t1, index_t(halfedge_facet_[t1ft2][t1fborder]));
348  signed_index_t ev2 =
349  tet_vertex(t1, index_t(halfedge_facet_[t1fborder][t1ft2]));
350 
351  // Turn around edge [ev1,ev2] inside the conflict zone
352  // until we reach again the boundary of the conflict zone.
353  // Traversing inside the conflict zone is faster (as compared
354  // to outside) since it traverses a smaller number of tets.
355  index_t cur_t = t1;
356  index_t cur_f = t1ft2;
357  index_t next_t = index_t(tet_adjacent(cur_t,cur_f));
358  while(tet_is_in_list(next_t)) {
359  geo_debug_assert(next_t != t1);
360  cur_t = next_t;
361  cur_f = get_facet_by_halfedge(cur_t,ev1,ev2);
362  next_t = index_t(tet_adjacent(cur_t, cur_f));
363  }
364 
365  // At this point, cur_t is in conflict zone and
366  // next_t is outside the conflict zone.
367  index_t f12,f21;
368  get_facets_by_halfedge(next_t, ev1, ev2, f12, f21);
369  t2 = index_t(tet_adjacent(next_t,f21));
370  signed_index_t v_neigh_opposite = tet_vertex(next_t,f12);
371  t2ft1 = find_tet_vertex(t2, v_neigh_opposite);
372  t2fborder = cur_f;
373 
374  // Test whether the found neighboring tet was created
375  // (then return true) or is an old tet in conflict
376  // (then return false).
377  return(t2 != cur_t);
378  }
379 
380  // _________ Combinatorics - new and delete _________________________
381 
389  index_t max_t() const {
390  return cell_to_v_store_.size() / 4;
391  }
392 
393 
409  static const index_t NOT_IN_LIST = index_t(~0);
410 
426  static const index_t NOT_IN_LIST_BIT = index_t(1u << 31);
427 
433  static const index_t END_OF_LIST = ~(NOT_IN_LIST_BIT);
434 
435 
453  bool tet_is_in_list(index_t t) const {
454  geo_debug_assert(t < max_t());
455  return (cell_next_[t] & NOT_IN_LIST_BIT) == 0;
456  }
457 
469  geo_debug_assert(t < max_t());
470  geo_debug_assert(tet_is_in_list(t));
471  return cell_next_[t];
472  }
473 
484  void add_tet_to_list(index_t t, index_t& first, index_t& last) {
485  geo_debug_assert(t < max_t());
486  geo_debug_assert(!tet_is_in_list(t));
487  if(last == END_OF_LIST) {
488  geo_debug_assert(first == END_OF_LIST);
489  first = last = t;
490  cell_next_[t] = END_OF_LIST;
491  } else {
492  cell_next_[t] = first;
493  first = t;
494  }
495  }
496 
507  geo_debug_assert(t < max_t());
508  geo_debug_assert(tet_is_in_list(t));
509  cell_next_[t] = NOT_IN_LIST;
510  }
511 
518  static const signed_index_t VERTEX_AT_INFINITY = -1;
519 
530  bool tet_is_finite(index_t t) const {
531  return
532  cell_to_v_store_[4 * t] >= 0 &&
533  cell_to_v_store_[4 * t + 1] >= 0 &&
534  cell_to_v_store_[4 * t + 2] >= 0 &&
535  cell_to_v_store_[4 * t + 3] >= 0;
536  }
537 
549  bool tet_is_real(index_t t) const {
550  return !tet_is_free(t) && tet_is_finite(t);
551  }
552 
562  bool tet_is_virtual(index_t t) const {
563  return
564  !tet_is_free(t) && (
565  cell_to_v_store_[4 * t] == VERTEX_AT_INFINITY ||
566  cell_to_v_store_[4 * t + 1] == VERTEX_AT_INFINITY ||
567  cell_to_v_store_[4 * t + 2] == VERTEX_AT_INFINITY ||
568  cell_to_v_store_[4 * t + 3] == VERTEX_AT_INFINITY) ;
569  }
570 
581  bool tet_is_free(index_t t) const {
582  return tet_is_in_list(t);
583  }
584 
593  index_t result;
594  if(first_free_ == END_OF_LIST) {
595  cell_to_v_store_.resize(cell_to_v_store_.size() + 4, -1);
596  cell_to_cell_store_.resize(cell_to_cell_store_.size() + 4, -1);
597  // index_t(NOT_IN_LIST) is necessary, else with
598  // NOT_IN_LIST alone the compiler tries to generate a
599  // reference to NOT_IN_LIST resulting in a link error.
600  cell_next_.push_back(index_t(NOT_IN_LIST));
601  result = max_t() - 1;
602  } else {
603  result = first_free_;
604  first_free_ = tet_next(first_free_);
605  remove_tet_from_list(result);
606  }
607 
608  cell_to_cell_store_[4 * result] = -1;
609  cell_to_cell_store_[4 * result + 1] = -1;
610  cell_to_cell_store_[4 * result + 2] = -1;
611  cell_to_cell_store_[4 * result + 3] = -1;
612 
613  return result;
614  }
615 
631  ) {
632  index_t result = new_tetrahedron();
633  cell_to_v_store_[4 * result] = v1;
634  cell_to_v_store_[4 * result + 1] = v2;
635  cell_to_v_store_[4 * result + 2] = v3;
636  cell_to_v_store_[4 * result + 3] = v4;
637  return result;
638  }
639 
648  cur_stamp_ = (stamp | NOT_IN_LIST_BIT);
649  }
650 
666  bool tet_is_marked(index_t t) const {
667  return cell_next_[t] == cur_stamp_;
668  }
669 
683  void mark_tet(index_t t) {
684  cell_next_[t] = cur_stamp_;
685  }
686 
687  // _________ Combinatorics ___________________________________
688 
707  geo_debug_assert(f < 4);
708  geo_debug_assert(v < 3);
709  return index_t(tet_facet_vertex_[f][v]);
710  }
711 
720  geo_debug_assert(t < max_t());
721  geo_debug_assert(lv < 4);
722  return cell_to_v_store_[4 * t + lv];
723  }
724 
733  geo_debug_assert(t < max_t());
734  // Find local index of v in tetrahedron t vertices.
735  const signed_index_t* T = &(cell_to_v_store_[4 * t]);
736  return find_4(T,v);
737  }
738 
739 
748  geo_debug_assert(t < max_t());
749  geo_debug_assert(lv < 4);
750  geo_debug_assert(cell_to_v_store_[4 * t + lv] != -1);
751  return index_t(cell_to_v_store_[4 * t + lv]);
752  }
753 
761  geo_debug_assert(t < max_t());
762  geo_debug_assert(lv < 4);
763  cell_to_v_store_[4 * t + lv] = v;
764  }
765 
773  geo_debug_assert(t < max_t());
774  geo_debug_assert(lf < 4);
775  signed_index_t result = cell_to_cell_store_[4 * t + lf];
776  return result;
777  }
778 
787  geo_debug_assert(t1 < max_t());
788  geo_debug_assert(t2 < max_t());
789  geo_debug_assert(lf1 < 4);
790  cell_to_cell_store_[4 * t1 + lf1] = signed_index_t(t2);
791  }
792 
802  index_t t1, index_t t2_in
803  ) const {
804  geo_debug_assert(t1 < max_t());
805  geo_debug_assert(t2_in < max_t());
806  geo_debug_assert(t1 != t2_in);
807 
808  signed_index_t t2 = signed_index_t(t2_in);
809 
810  // Find local index of t2 in tetrahedron t1 adajcent tets.
811  const signed_index_t* T = &(cell_to_cell_store_[4 * t1]);
812  index_t result = find_4(T,t2);
813 
814  // Sanity check: make sure that t1 is adjacent to t2
815  // only once!
816  geo_debug_assert(tet_adjacent(t1,(result+1)%4) != t2);
817  geo_debug_assert(tet_adjacent(t1,(result+2)%4) != t2);
818  geo_debug_assert(tet_adjacent(t1,(result+3)%4) != t2);
819  return result;
820  }
821 
822 
823 
837  void set_tet(
838  index_t t,
841  index_t a0, index_t a1, index_t a2, index_t a3
842  ) {
843  geo_debug_assert(t < max_t());
844  cell_to_v_store_[4 * t] = v0;
845  cell_to_v_store_[4 * t + 1] = v1;
846  cell_to_v_store_[4 * t + 2] = v2;
847  cell_to_v_store_[4 * t + 3] = v3;
848  cell_to_cell_store_[4 * t] = signed_index_t(a0);
849  cell_to_cell_store_[4 * t + 1] = signed_index_t(a1);
850  cell_to_cell_store_[4 * t + 2] = signed_index_t(a2);
851  cell_to_cell_store_[4 * t + 3] = signed_index_t(a3);
852  }
853 
854  // _________ Combinatorics - traversals ______________________________
855 
867  ) const {
868  geo_debug_assert(t < max_t());
869  geo_debug_assert(v1 != v2);
870  // Find local index of v1 and v2 in tetrahedron t
871  const signed_index_t* T = &(cell_to_v_store_[4 * t]);
872  index_t lv1 = find_4(T,v1);
873  index_t lv2 = find_4(T,v2);
874  geo_debug_assert(lv1 != lv2);
875  return index_t(halfedge_facet_[lv1][lv2]);
876  }
877 
878 
892  index_t& f12, index_t& f21
893  ) const {
894  geo_debug_assert(t < max_t());
895  geo_debug_assert(v1 != v2);
896 
897  // Find local index of v1 and v2 in tetrahedron t
898  // The following expression is 10% faster than using
899  // if() statements (multiply by boolean result of test).
900  // Thank to Laurent Alonso for this idea.
901  const signed_index_t* T = &(cell_to_v_store_[4 * t]);
902 
903  signed_index_t lv1 =
904  (T[1] == v1) | ((T[2] == v1) * 2) | ((T[3] == v1) * 3);
905 
906  signed_index_t lv2 =
907  (T[1] == v2) | ((T[2] == v2) * 2) | ((T[3] == v2) * 3);
908 
909  geo_debug_assert(lv1 != 0 || T[0] == v1);
910  geo_debug_assert(lv2 != 0 || T[0] == v2);
911  geo_debug_assert(lv1 >= 0);
912  geo_debug_assert(lv2 >= 0);
913  geo_debug_assert(lv1 != lv2);
914 
915  f12 = index_t(halfedge_facet_[lv1][lv2]);
916  f21 = index_t(halfedge_facet_[lv2][lv1]);
917  }
918 
919 
931  ) const {
932  return (index_t)tet_adjacent(
933  t, get_facet_by_halfedge(t, v1, v2)
934  );
935  }
936 
937  // _________ Predicates _____________________________________________
938 
952  bool tet_is_conflict(index_t t, const double* p) const {
953 
954  // Lookup tetrahedron vertices
955  const double* pv[4];
956  for(index_t i=0; i<4; ++i) {
957  signed_index_t v = tet_vertex(t,i);
958  pv[i] = (v == -1) ? nullptr : vertex_ptr(index_t(v));
959  }
960 
961  // Check for virtual tetrahedra (then in_sphere()
962  // is replaced with orient3d())
963  for(index_t lf = 0; lf < 4; ++lf) {
964 
965  if(pv[lf] == nullptr) {
966 
967  // Facet of a virtual tetrahedron opposite to
968  // infinite vertex corresponds to
969  // the triangle on the convex hull of the points.
970  // Orientation is obtained by replacing vertex lf
971  // with p.
972  pv[lf] = p;
973  Sign sign = PCK::orient_3d(pv[0],pv[1],pv[2],pv[3]);
974 
975  if(sign > 0) {
976  return true;
977  }
978 
979  if(sign < 0) {
980  return false;
981  }
982 
983  // If sign is zero, we check the real tetrahedron
984  // adjacent to the facet on the convex hull.
985  geo_debug_assert(tet_adjacent(t, lf) >= 0);
986  index_t t2 = index_t(tet_adjacent(t, lf));
987  geo_debug_assert(!tet_is_virtual(t2));
988 
989  // If t2 is already chained in the conflict list,
990  // then it is conflict
991  if(tet_is_in_list(t2)) {
992  return true;
993  }
994 
995  // If t2 is marked, then it is not in conflict.
996  if(tet_is_marked(t2)) {
997  return false;
998  }
999 
1000  return tet_is_conflict(t2, p);
1001  }
1002  }
1003 
1004  // If the tetrahedron is a finite one, it is in conflict
1005  // if its circumscribed sphere contains the point (this is
1006  // the standard case).
1007 
1008  if(weighted_) {
1009  double h0 = heights_[finite_tet_vertex(t, 0)];
1010  double h1 = heights_[finite_tet_vertex(t, 1)];
1011  double h2 = heights_[finite_tet_vertex(t, 2)];
1012  double h3 = heights_[finite_tet_vertex(t, 3)];
1013  index_t pindex = index_t(
1014  (p - vertex_ptr(0)) / int(vertex_stride_)
1015  );
1016  double h = heights_[pindex];
1017  return (PCK::orient_3dlifted_SOS(
1018  pv[0],pv[1],pv[2],pv[3],p,h0,h1,h2,h3,h
1019  ) > 0) ;
1020  }
1021 
1022  return (PCK::in_sphere_3d_SOS(pv[0], pv[1], pv[2], pv[3], p) > 0);
1023  }
1024 
1025  protected:
1026 
1036  // The following expression is 10% faster than using
1037  // if() statements. This uses the C++ norm, that
1038  // ensures that the 'true' boolean value converted to
1039  // an int is always 1. With most compilers, this avoids
1040  // generating branching instructions.
1041  // Thank to Laurent Alonso for this idea.
1042  // Note: Laurent also has this version:
1043  // (T[0] != v)+(T[2]==v)+2*(T[3]==v)
1044  // that avoids a *3 multiply, but it is not faster in
1045  // practice.
1046  index_t result = index_t(
1047  (T[1] == v) | ((T[2] == v) * 2) | ((T[3] == v) * 3)
1048  );
1049  // Sanity check, important if it was T[0], not explicitly
1050  // tested (detects input that does not meet the precondition).
1051  geo_debug_assert(T[result] == v);
1052  return result;
1053  }
1054 
1058  ~Delaunay3d() override;
1059 
1064  void show_tet(index_t t) const;
1065 
1072  void show_tet_adjacent(index_t t, index_t lf) const;
1073 
1080  index_t first, const std::string& list_name
1081  ) const;
1082 
1086  void check_combinatorics(bool verbose = false) const;
1087 
1091  void check_geometry(bool verbose = false) const;
1092 
1093  private:
1094  vector<signed_index_t> cell_to_v_store_;
1095  vector<signed_index_t> cell_to_cell_store_;
1096  vector<index_t> cell_next_;
1097  vector<index_t> reorder_;
1098  index_t cur_stamp_; // used for marking
1099  index_t first_free_;
1100  bool weighted_;
1101  vector<double> heights_; // only used in weighted mode
1102 
1106  bool debug_mode_;
1107 
1111  bool verbose_debug_mode_;
1112 
1116  bool benchmark_mode_;
1117 
1126  static char tet_facet_vertex_[4][3];
1127 
1132  static char halfedge_facet_[4][4];
1133 
1138  std::stack<index_t> S_;
1139 
1145  class StellateConflictStack {
1146  public:
1147 
1157  void push(index_t t1, index_t t1fbord, index_t t1fprev) {
1158  store_.resize(store_.size()+1);
1159  top().t1 = t1;
1160  top().t1fbord = Numeric::uint8(t1fbord);
1161  top().t1fprev = Numeric::uint8(t1fprev);
1162  }
1163 
1170  void save_locals(index_t new_t, index_t t1ft2, index_t t2ft1) {
1171  geo_debug_assert(!empty());
1172  top().new_t = new_t;
1173  top().t1ft2 = Numeric::uint8(t1ft2);
1174  top().t2ft1 = Numeric::uint8(t2ft1);
1175  }
1176 
1186  void get_parameters(
1187  index_t& t1, index_t& t1fbord, index_t& t1fprev
1188  ) const {
1189  geo_debug_assert(!empty());
1190  t1 = top().t1;
1191  t1fbord = index_t(top().t1fbord);
1192  t1fprev = index_t(top().t1fprev);
1193  }
1194 
1195 
1202  void get_locals(
1203  index_t& new_t, index_t& t1ft2, index_t& t2ft1
1204  ) const {
1205  geo_debug_assert(!empty());
1206  new_t = top().new_t;
1207  t1ft2 = index_t(top().t1ft2);
1208  t2ft1 = index_t(top().t2ft1);
1209  }
1210 
1214  void pop() {
1215  geo_debug_assert(!empty());
1216  store_.pop_back();
1217  }
1218 
1224  bool empty() const {
1225  return store_.empty();
1226  }
1227 
1228  private:
1229 
1234  struct Frame {
1235  // Parameters
1236  index_t t1;
1237  index_t new_t;
1238  Numeric::uint8 t1fbord ;
1239 
1240  // Local variables
1241  Numeric::uint8 t1fprev ;
1242  Numeric::uint8 t1ft2 ;
1243  Numeric::uint8 t2ft1 ;
1244  };
1245 
1252  Frame& top() {
1253  geo_debug_assert(!empty());
1254  return *store_.rbegin();
1255  }
1256 
1263  const Frame& top() const {
1264  geo_debug_assert(!empty());
1265  return *store_.rbegin();
1266  }
1267 
1268  std::vector<Frame> store_;
1269  };
1270 
1275  StellateConflictStack S2_;
1276 
1277  Cavity cavity_;
1278  };
1279 
1280  /************************************************************************/
1281 
1293  class GEOGRAM_API RegularWeightedDelaunay3d : public Delaunay3d {
1294  public:
1305 
1306  protected:
1311  };
1312 }
1313 
1314 #endif
#define geo_debug_assert(x)
Verifies that a condition is met.
Definition: assert.h:196
Implementation of Delaunay in 3d.
Definition: delaunay_3d.h:107
index_t locate_inexact(const double *p, index_t hint, index_t max_iter) const
Finds the tetrahedron that (approximately) contains a point using inexact predicates.
void check_geometry(bool verbose=false) const
For debugging purposes, test some geometrical properties.
void show_tet_adjacent(index_t t, index_t lf) const
For debugging purposes, displays a tetrahedron adjacency.
signed_index_t tet_vertex(index_t t, index_t lv) const
Gets the index of a vertex of a tetrahedron.
Definition: delaunay_3d.h:719
void set_tet_adjacent(index_t t1, index_t lf1, index_t t2)
Sets a tetrahedron-to-tetrahedron adjacency.
Definition: delaunay_3d.h:786
index_t get_facet_by_halfedge(index_t t, signed_index_t v1, signed_index_t v2) const
Definition: delaunay_3d.h:865
bool get_neighbor_along_conflict_zone_border(index_t t1, index_t t1fborder, index_t t1ft2, index_t &t2, index_t &t2fborder, index_t &t2ft1) const
Finds the neighbor of a tetrahedron on the border of the conflict zone.
Definition: delaunay_3d.h:325
void set_tet(index_t t, signed_index_t v0, signed_index_t v1, signed_index_t v2, signed_index_t v3, index_t a0, index_t a1, index_t a2, index_t a3)
Sets the vertices and adjacent tetrahedra of a tetrahedron.
Definition: delaunay_3d.h:837
index_t next_around_halfedge(index_t &t, signed_index_t v1, signed_index_t v2) const
Gets the next tetrahedron around an oriented edge of a tetrahedron.
Definition: delaunay_3d.h:929
void mark_tet(index_t t)
Marks a tetrahedron.
Definition: delaunay_3d.h:683
index_t find_tet_adjacent(index_t t1, index_t t2_in) const
Finds the index of the facet accros which t1 is adjacent to t2_in.
Definition: delaunay_3d.h:801
static index_t find_4(const signed_index_t *T, signed_index_t v)
Finds the index of an integer in an array of four integers.
Definition: delaunay_3d.h:1035
index_t nearest_vertex(const double *p) const override
Computes the nearest vertex from a query point.
~Delaunay3d() override
Delaunay3d destructor.
bool tet_is_in_list(index_t t) const
Tests whether a tetrahedron belongs to a linked list.
Definition: delaunay_3d.h:453
index_t find_tet_vertex(index_t t, signed_index_t v) const
Finds the index of the vertex in a tetrahedron.
Definition: delaunay_3d.h:732
index_t tet_next(index_t t) const
Gets the index of a successor of a tetrahedron.
Definition: delaunay_3d.h:468
bool tet_is_finite(index_t t) const
Tests whether a given tetrahedron is a finite one.
Definition: delaunay_3d.h:530
void check_combinatorics(bool verbose=false) const
For debugging purposes, tests some combinatorial properties.
bool tet_is_marked(index_t t) const
Tests whether a tetrahedron is marked.
Definition: delaunay_3d.h:666
bool tet_is_free(index_t t) const
Tests whether a tetrahedron is in the free list.
Definition: delaunay_3d.h:581
bool tet_is_virtual(index_t t) const
Tests whether a tetrahedron is a virtual one.
Definition: delaunay_3d.h:562
index_t stellate_conflict_zone_iterative(index_t v, index_t t_bndry, index_t f_bndry, index_t prev_f=index_t(-1))
Creates a star of tetrahedra filling the conflict zone.
index_t insert(index_t v, index_t hint=NO_TETRAHEDRON)
Inserts a point in the triangulation.
void set_tet_mark_stamp(index_t stamp)
Generates a unique stamp for marking tets.
Definition: delaunay_3d.h:647
void set_tet_vertex(index_t t, index_t lv, signed_index_t v)
Sets a tetrahedron-to-vertex adjacency.
Definition: delaunay_3d.h:760
Delaunay3d(coord_index_t dimension=3)
Constructs a new Delaunay3d.
void set_vertices(index_t nb_vertices, const double *vertices) override
Sets the vertices of this Delaunay, and recomputes the cells.
void remove_tet_from_list(index_t t)
Removes a tetrahedron from the linked list it belongs to.
Definition: delaunay_3d.h:506
void add_tet_to_list(index_t t, index_t &first, index_t &last)
Adds a tetrahedron to a linked list.
Definition: delaunay_3d.h:484
signed_index_t tet_adjacent(index_t t, index_t lf) const
Gets the index of a tetrahedron adjacent to another one.
Definition: delaunay_3d.h:772
void show_tet(index_t t) const
For debugging purposes, displays a tetrahedron.
void find_conflict_zone_iterative(const double *p, index_t t, index_t &t_bndry, index_t &f_bndry, index_t &first, index_t &last)
This function is used to implement find_conflict_zone.
void get_facets_by_halfedge(index_t t, signed_index_t v1, signed_index_t v2, index_t &f12, index_t &f21) const
Definition: delaunay_3d.h:890
bool create_first_tetrahedron(index_t &iv0, index_t &iv1, index_t &iv2, index_t &iv3)
Finds in the pointset a set of four non-coplanar points.
void find_conflict_zone(index_t v, index_t t, const Sign *orient, index_t &t_bndry, index_t &f_bndry, index_t &first, index_t &last)
Determines the list of tetrahedra in conflict with a given point.
index_t new_tetrahedron()
Creates a new tetrahedron.
Definition: delaunay_3d.h:592
void show_list(index_t first, const std::string &list_name) const
For debugging purposes, displays a tetrahedron.
index_t new_tetrahedron(signed_index_t v1, signed_index_t v2, signed_index_t v3, signed_index_t v4)
Creates a new tetrahedron.
Definition: delaunay_3d.h:628
bool tet_is_conflict(index_t t, const double *p) const
Tests whether a given tetrahedron is in conflict with a given 3d point.
Definition: delaunay_3d.h:952
index_t finite_tet_vertex(index_t t, index_t lv) const
Gets the index of a vertex of a tetrahedron.
Definition: delaunay_3d.h:747
index_t stellate_cavity(index_t v)
Creates a star of tetrahedra filling the conflict zone.
bool tet_is_real(index_t t) const
Tests whether a tetrahedron is a real one.
Definition: delaunay_3d.h:549
index_t max_t() const
Maximum valid index for a tetrahedron.
Definition: delaunay_3d.h:389
static index_t tet_facet_vertex(index_t f, index_t v)
Returns the local index of a vertex by facet and by local vertex index in the facet.
Definition: delaunay_3d.h:706
index_t locate(const double *p, index_t hint=NO_TETRAHEDRON, bool thread_safe=false, Sign *orient=nullptr) const
Finds the tetrahedron that contains a point.
Abstract interface for Delaunay triangulation in Nd.
Definition: delaunay.h:71
Regular Delaunay triangulation of weighted points.
Definition: delaunay_3d.h:1293
RegularWeightedDelaunay3d(coord_index_t dimension=4)
Constructs a new Regular Delaunay3d triangulation.
~RegularWeightedDelaunay3d() override
RegularWeightedDelaunay3d destructor.
Abstract interface for Delaunay.
Common include file, providing basic definitions. Should be included before anything else by all head...
Geometric functions in 2d and 3d.
uint8_t uint8
Definition: numeric.h:135
Sign orient_3dlifted_SOS(const double *p0, const double *p1, const double *p2, const double *p3, const double *p4, double h0, double h1, double h2, double h3, double h4)
Computes the 4d orientation test with symbolic perturbation.
Sign in_sphere_3d_SOS(const double *p0, const double *p1, const double *p2, const double *p3, const double *p4)
Tests whether a 3d point is inside the circumscribed sphere of a 3d tetrahedron.
Sign orient_3d(const vec3HE &p0, const vec3HE &p1, const vec3HE &p2, const vec3HE &p3)
Computes the orientation predicate in 3d.
Global Vorpaline namespace.
geo_signed_index_t signed_index_t
The type for storing and manipulating indices differences.
Definition: numeric.h:343
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
geo_coord_index_t coord_index_t
The type for storing coordinate indices, and iterating on the coordinates of a point.
Definition: numeric.h:363
Filtered exact predicates for restricted Voronoi diagrams.