Geogram Version 1.9.6-rc
A programming library of geometric algorithms
Loading...
Searching...
No Matches
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
45#include <geogram/delaunay/cavity.h>
48
49#include <stack>
50
56namespace GEO {
57
107 class GEOGRAM_API Delaunay3d : public Delaunay {
108 public:
122 Delaunay3d(coord_index_t dimension = 3);
123
127 void set_vertices(index_t nb_vertices, const double* vertices) override;
128
132 index_t nearest_vertex(const double* p) const override;
133
134 protected:
135
142 static constexpr index_t NO_TETRAHEDRON = NO_INDEX;
143
157 index_t& iv0, index_t& iv1, index_t& iv2, index_t& iv3
158 );
159
180 const double* p, index_t hint = NO_TETRAHEDRON,
181 bool thread_safe = false,
182 Sign* orient = nullptr
183 ) const;
184
203 const double* p, index_t hint, index_t max_iter
204 ) const;
205
214 index_t insert(index_t v, index_t hint = NO_TETRAHEDRON);
215
240 index_t v,
241 index_t t, const Sign* orient,
242 index_t& t_bndry, index_t& f_bndry,
243 index_t& first, index_t& last
244 );
245
262 const double* p, index_t t,
263 index_t& t_bndry, index_t& f_bndry,
264 index_t& first, index_t& last
265 );
266
277
278
297 index_t v,
298 index_t t_bndry, index_t f_bndry,
299 index_t prev_f=NO_INDEX
300 );
301
323 index_t t1,
324 index_t t1fborder,
325 index_t t1ft2,
326 index_t& t2,
327 index_t& t2fborder,
328 index_t& t2ft1
329 ) const {
330
331 // Note: this function is a bit long for an inline function,
332 // but I observed a (modest) performance gain doing so.
333
334 // Find two vertices that are both on facets new_f and f1
335 // (the edge around which we are turning)
336 // This uses duality as follows:
337 // Primal form (not used here):
338 // halfedge_facet_[v1][v2] returns a facet that is incident
339 // to both v1 and v2.
340 // Dual form (used here):
341 // halfedge_facet_[f1][f2] returns a vertex that both
342 // f1 and f2 are incident to.
343 index_t ev1 =
344 tet_vertex(t1, index_t(halfedge_facet_[t1ft2][t1fborder]));
345 index_t ev2 =
346 tet_vertex(t1, index_t(halfedge_facet_[t1fborder][t1ft2]));
347
348 // Turn around edge [ev1,ev2] inside the conflict zone
349 // until we reach again the boundary of the conflict zone.
350 // Traversing inside the conflict zone is faster (as compared
351 // to outside) since it traverses a smaller number of tets.
352 index_t cur_t = t1;
353 index_t cur_f = t1ft2;
354 index_t next_t = tet_adjacent(cur_t,cur_f);
355 while(tet_is_in_list(next_t)) {
356 geo_debug_assert(next_t != t1);
357 cur_t = next_t;
358 cur_f = get_facet_by_halfedge(cur_t,ev1,ev2);
359 next_t = tet_adjacent(cur_t, cur_f);
360 }
361
362 // At this point, cur_t is in conflict zone and
363 // next_t is outside the conflict zone.
364 index_t f12,f21;
365 get_facets_by_halfedge(next_t, ev1, ev2, f12, f21);
366 t2 = tet_adjacent(next_t,f21);
367 index_t v_neigh_opposite = tet_vertex(next_t,f12);
368 t2ft1 = find_tet_vertex(t2, v_neigh_opposite);
369 t2fborder = cur_f;
370
371 // Test whether the found neighboring tet was created
372 // (then return true) or is an old tet in conflict
373 // (then return false).
374 return(t2 != cur_t);
375 }
376
377 /****** Combinatorics - new and delete ***************************/
378
386 index_t max_t() const {
387 return cell_to_v_store_.size() / 4;
388 }
389
405 static constexpr index_t NOT_IN_LIST = ~index_t(0);
406
422 static constexpr index_t NOT_IN_LIST_BIT =
423 index_t(1) << (sizeof(index_t)*8-1) ;
424
430 static constexpr index_t END_OF_LIST = ~NOT_IN_LIST_BIT;
431
432
450 bool tet_is_in_list(index_t t) const {
451 geo_debug_assert(t < max_t());
452 return (cell_next_[t] & NOT_IN_LIST_BIT) == 0;
453 }
454
466 geo_debug_assert(t < max_t());
467 geo_debug_assert(tet_is_in_list(t));
468 return cell_next_[t];
469 }
470
481 void add_tet_to_list(index_t t, index_t& first, index_t& last) {
482 geo_debug_assert(t < max_t());
483 geo_debug_assert(!tet_is_in_list(t));
484 if(last == END_OF_LIST) {
485 geo_debug_assert(first == END_OF_LIST);
486 first = last = t;
487 cell_next_[t] = END_OF_LIST;
488 } else {
489 cell_next_[t] = first;
490 first = t;
491 }
492 }
493
504 geo_debug_assert(t < max_t());
505 geo_debug_assert(tet_is_in_list(t));
506 cell_next_[t] = NOT_IN_LIST;
507 }
508
515 static constexpr index_t VERTEX_AT_INFINITY = NO_INDEX;
516
527 bool tet_is_finite(index_t t) const {
528 return
529 cell_to_v_store_[4 * t] != NO_INDEX &&
530 cell_to_v_store_[4 * t + 1] != NO_INDEX &&
531 cell_to_v_store_[4 * t + 2] != NO_INDEX &&
532 cell_to_v_store_[4 * t + 3] != NO_INDEX;
533 }
534
546 bool tet_is_real(index_t t) const {
547 return !tet_is_free(t) && tet_is_finite(t);
548 }
549
559 bool tet_is_virtual(index_t t) const {
560 return
561 !tet_is_free(t) && (
562 cell_to_v_store_[4 * t] == VERTEX_AT_INFINITY ||
563 cell_to_v_store_[4 * t + 1] == VERTEX_AT_INFINITY ||
564 cell_to_v_store_[4 * t + 2] == VERTEX_AT_INFINITY ||
565 cell_to_v_store_[4 * t + 3] == VERTEX_AT_INFINITY) ;
566 }
567
578 bool tet_is_free(index_t t) const {
579 return tet_is_in_list(t);
580 }
581
590 index_t result;
591 if(first_free_ == END_OF_LIST) {
592 cell_to_v_store_.resize(
593 cell_to_v_store_.size() + 4, NO_INDEX
594 );
595 cell_to_cell_store_.resize(
596 cell_to_cell_store_.size() + 4, NO_INDEX
597 );
598 // index_t(NOT_IN_LIST) is necessary, else with
599 // NOT_IN_LIST alone the compiler tries to generate a
600 // reference to NOT_IN_LIST resulting in a link error.
601 cell_next_.push_back(index_t(NOT_IN_LIST));
602 result = max_t() - 1;
603 } else {
604 result = first_free_;
605 first_free_ = tet_next(first_free_);
606 remove_tet_from_list(result);
607 }
608
609 cell_to_cell_store_[4 * result] = NO_INDEX;
610 cell_to_cell_store_[4 * result + 1] = NO_INDEX;
611 cell_to_cell_store_[4 * result + 2] = NO_INDEX;
612 cell_to_cell_store_[4 * result + 3] = NO_INDEX;
613
614 return result;
615 }
616
630 index_t v1, index_t v2,
631 index_t v3, index_t v4
632 ) {
633 index_t result = new_tetrahedron();
634 cell_to_v_store_[4 * result] = v1;
635 cell_to_v_store_[4 * result + 1] = v2;
636 cell_to_v_store_[4 * result + 2] = v3;
637 cell_to_v_store_[4 * result + 3] = v4;
638 return result;
639 }
640
649 cur_stamp_ = (stamp | NOT_IN_LIST_BIT);
650 }
651
667 bool tet_is_marked(index_t t) const {
668 return cell_next_[t] == cur_stamp_;
669 }
670
685 cell_next_[t] = cur_stamp_;
686 }
687
688 /********* Combinatorics ******************************************/
689
708 geo_debug_assert(f < 4);
709 geo_debug_assert(v < 3);
710 return index_t(tet_facet_vertex_[f][v]);
711 }
712
721 geo_debug_assert(t < max_t());
722 geo_debug_assert(lv < 4);
723 return cell_to_v_store_[4 * t + lv];
724 }
725
734 geo_debug_assert(t < max_t());
735 // Find local index of v in tetrahedron t vertices.
736 const index_t* T = &(cell_to_v_store_[4 * t]);
737 return find_4(T,v);
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] != NO_INDEX);
751 return 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 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] = t2;
791 }
792
802 geo_debug_assert(t1 < max_t());
803 geo_debug_assert(t2 < max_t());
804 geo_debug_assert(t1 != t2);
805
806 // Find local index of t2 in tetrahedron t1 adajcent tets.
807 const index_t* T = &(cell_to_cell_store_[4 * t1]);
808 index_t result = find_4(T,t2);
809
810 // Sanity check: make sure that t1 is adjacent to t2
811 // only once!
812 geo_debug_assert(tet_adjacent(t1,(result+1)%4) != t2);
813 geo_debug_assert(tet_adjacent(t1,(result+2)%4) != t2);
814 geo_debug_assert(tet_adjacent(t1,(result+3)%4) != t2);
815 return result;
816 }
817
832 index_t t,
833 index_t v0, index_t v1, index_t v2, index_t v3,
834 index_t a0, index_t a1, index_t a2, index_t a3
835 ) {
836 geo_debug_assert(t < max_t());
837 cell_to_v_store_[4 * t] = v0;
838 cell_to_v_store_[4 * t + 1] = v1;
839 cell_to_v_store_[4 * t + 2] = v2;
840 cell_to_v_store_[4 * t + 3] = v3;
841 cell_to_cell_store_[4 * t] = a0;
842 cell_to_cell_store_[4 * t + 1] = a1;
843 cell_to_cell_store_[4 * t + 2] = a2;
844 cell_to_cell_store_[4 * t + 3] = a3;
845 }
846
847 /****** Combinatorics - traversals ************************/
848
859 geo_debug_assert(t < max_t());
860 geo_debug_assert(v1 != v2);
861 // Find local index of v1 and v2 in tetrahedron t
862 const index_t* T = &(cell_to_v_store_[4 * t]);
863 index_t lv1 = find_4(T,v1);
864 index_t lv2 = find_4(T,v2);
865 geo_debug_assert(lv1 != lv2);
866 return index_t(halfedge_facet_[lv1][lv2]);
867 }
868
881 index_t t, index_t v1, index_t v2,
882 index_t& f12, index_t& f21
883 ) const {
884 geo_debug_assert(t < max_t());
885 geo_debug_assert(v1 != v2);
886
887 // Find local index of v1 and v2 in tetrahedron t
888 // The following expression is 10% faster than using
889 // if() statements (multiply by boolean result of test).
890 // Thank to Laurent Alonso for this idea.
891 const index_t* T = &(cell_to_v_store_[4 * t]);
892
893 index_t lv1 = index_t(
894 (T[1] == v1) | ((T[2] == v1) * 2) | ((T[3] == v1) * 3)
895 );
896
897 index_t lv2 = index_t(
898 (T[1] == v2) | ((T[2] == v2) * 2) | ((T[3] == v2) * 3)
899 );
900
901 geo_debug_assert(lv1 != 0 || T[0] == v1);
902 geo_debug_assert(lv2 != 0 || T[0] == v2);
903 geo_debug_assert(lv1 != lv2);
904
905 f12 = index_t(halfedge_facet_[lv1][lv2]);
906 f21 = index_t(halfedge_facet_[lv2][lv1]);
907 }
908
909
920 return (index_t)tet_adjacent(
921 t, get_facet_by_halfedge(t, v1, v2)
922 );
923 }
924
925 /****** Predicates **********************************************/
926
940 bool tet_is_conflict(index_t t, const double* p) const {
941
942 // Lookup tetrahedron vertices
943 const double* pv[4];
944 for(index_t i=0; i<4; ++i) {
945 index_t v = tet_vertex(t,i);
946 pv[i] = (v == NO_INDEX) ? nullptr : vertex_ptr(v);
947 }
948
949 // Check for virtual tetrahedra (then in_sphere()
950 // is replaced with orient3d())
951 for(index_t lf = 0; lf < 4; ++lf) {
952
953 if(pv[lf] == nullptr) {
954
955 // Facet of a virtual tetrahedron opposite to
956 // infinite vertex corresponds to
957 // the triangle on the convex hull of the points.
958 // Orientation is obtained by replacing vertex lf
959 // with p.
960 pv[lf] = p;
961 Sign sign = PCK::orient_3d(pv[0],pv[1],pv[2],pv[3]);
962
963 if(sign > 0) {
964 return true;
965 }
966
967 if(sign < 0) {
968 return false;
969 }
970
971 // If sign is zero, we check the real tetrahedron
972 // adjacent to the facet on the convex hull.
973 geo_debug_assert(tet_adjacent(t, lf) != NO_INDEX);
974 index_t t2 = tet_adjacent(t, lf);
975 geo_debug_assert(!tet_is_virtual(t2));
976
977 // If t2 is already chained in the conflict list,
978 // then it is conflict
979 if(tet_is_in_list(t2)) {
980 return true;
981 }
982
983 // If t2 is marked, then it is not in conflict.
984 if(tet_is_marked(t2)) {
985 return false;
986 }
987
988 return tet_is_conflict(t2, p);
989 }
990 }
991
992 // If the tetrahedron is a finite one, it is in conflict
993 // if its circumscribed sphere contains the point (this is
994 // the standard case).
995
996 if(weighted_) {
997 double h0 = heights_[finite_tet_vertex(t, 0)];
998 double h1 = heights_[finite_tet_vertex(t, 1)];
999 double h2 = heights_[finite_tet_vertex(t, 2)];
1000 double h3 = heights_[finite_tet_vertex(t, 3)];
1001 index_t pindex = index_t(
1002 (p - vertex_ptr(0)) / int(vertex_stride_)
1003 );
1004 double h = heights_[pindex];
1005 return (PCK::orient_3dlifted_SOS(
1006 pv[0],pv[1],pv[2],pv[3],p,h0,h1,h2,h3,h
1007 ) > 0) ;
1008 }
1009
1010 return (PCK::in_sphere_3d_SOS(pv[0], pv[1], pv[2], pv[3], p) > 0);
1011 }
1012
1013 protected:
1014
1023 static index_t find_4(const index_t* T, index_t v) {
1024 // The following expression is 10% faster than using
1025 // if() statements. This uses the C++ norm, that
1026 // ensures that the 'true' boolean value converted to
1027 // an int is always 1. With most compilers, this avoids
1028 // generating branching instructions.
1029 // Thank to Laurent Alonso for this idea.
1030 // Note: Laurent also has this version:
1031 // (T[0] != v)+(T[2]==v)+2*(T[3]==v)
1032 // that avoids a *3 multiply, but it is not faster in
1033 // practice.
1034 index_t result = index_t(
1035 (T[1] == v) | ((T[2] == v) * 2) | ((T[3] == v) * 3)
1036 );
1037 // Sanity check, important if it was T[0], not explicitly
1038 // tested (detects input that does not meet the precondition).
1039 geo_debug_assert(T[result] == v);
1040 return result;
1041 }
1042
1046 ~Delaunay3d() override;
1047
1052 void show_tet(index_t t) const;
1053
1061
1067 void show_list(index_t first, const std::string& list_name) const;
1068
1072 void check_combinatorics(bool verbose = false) const;
1073
1077 void check_geometry(bool verbose = false) const;
1078
1079 private:
1080 vector<index_t> cell_to_v_store_;
1081 vector<index_t> cell_to_cell_store_;
1082 vector<index_t> cell_next_;
1083 vector<index_t> reorder_;
1084 index_t cur_stamp_; // used for marking
1085 index_t first_free_;
1086 bool weighted_;
1087 vector<double> heights_; // only used in weighted mode
1088
1092 bool debug_mode_;
1093
1097 bool verbose_debug_mode_;
1098
1102 bool benchmark_mode_;
1103
1112 static char tet_facet_vertex_[4][3];
1113
1118 static char halfedge_facet_[4][4];
1119
1124 std::stack<index_t> S_;
1125
1131 class StellateConflictStack {
1132 public:
1133
1143 void push(index_t t1, index_t t1fbord, index_t t1fprev) {
1144 store_.resize(store_.size()+1);
1145 top().t1 = t1;
1146 top().t1fbord = Numeric::uint8(t1fbord);
1147 top().t1fprev = Numeric::uint8(t1fprev);
1148 }
1149
1156 void save_locals(index_t new_t, index_t t1ft2, index_t t2ft1) {
1157 geo_debug_assert(!empty());
1158 top().new_t = new_t;
1159 top().t1ft2 = Numeric::uint8(t1ft2);
1160 top().t2ft1 = Numeric::uint8(t2ft1);
1161 }
1162
1172 void get_parameters(
1173 index_t& t1, index_t& t1fbord, index_t& t1fprev
1174 ) const {
1175 geo_debug_assert(!empty());
1176 t1 = top().t1;
1177 t1fbord = index_t(top().t1fbord);
1178 t1fprev = index_t(top().t1fprev);
1179 }
1180
1181
1188 void get_locals(
1189 index_t& new_t, index_t& t1ft2, index_t& t2ft1
1190 ) const {
1191 geo_debug_assert(!empty());
1192 new_t = top().new_t;
1193 t1ft2 = index_t(top().t1ft2);
1194 t2ft1 = index_t(top().t2ft1);
1195 }
1196
1200 void pop() {
1201 geo_debug_assert(!empty());
1202 store_.pop_back();
1203 }
1204
1210 bool empty() const {
1211 return store_.empty();
1212 }
1213
1214 private:
1215
1220 struct Frame {
1221 // Parameters
1222 index_t t1;
1223 index_t new_t;
1224 Numeric::uint8 t1fbord ;
1225
1226 // Local variables
1227 Numeric::uint8 t1fprev ;
1228 Numeric::uint8 t1ft2 ;
1229 Numeric::uint8 t2ft1 ;
1230 };
1231
1238 Frame& top() {
1239 geo_debug_assert(!empty());
1240 return *store_.rbegin();
1241 }
1242
1249 const Frame& top() const {
1250 geo_debug_assert(!empty());
1251 return *store_.rbegin();
1252 }
1253
1254 std::vector<Frame> store_;
1255 };
1256
1261 StellateConflictStack S2_;
1262
1263 Cavity cavity_;
1264 };
1265
1266 /************************************************************************/
1267
1279 class GEOGRAM_API RegularWeightedDelaunay3d : public Delaunay3d {
1280 public:
1291
1292 protected:
1297 };
1298}
1299
1300#endif
#define geo_debug_assert(x)
Verifies that a condition is met.
Definition assert.h:196
Implementation of Delaunay in 3d.
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.
void set_tet_adjacent(index_t t1, index_t lf1, index_t t2)
Sets a tetrahedron-to-tetrahedron adjacency.
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.
void mark_tet(index_t t)
Marks a tetrahedron.
index_t stellate_conflict_zone_iterative(index_t v, index_t t_bndry, index_t f_bndry, index_t prev_f=NO_INDEX)
Creates a star of tetrahedra filling the conflict zone.
void get_facets_by_halfedge(index_t t, index_t v1, index_t v2, index_t &f12, index_t &f21) const
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.
index_t tet_next(index_t t) const
Gets the index of a successor of a tetrahedron.
bool tet_is_finite(index_t t) const
Tests whether a given tetrahedron is a finite one.
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.
bool tet_is_free(index_t t) const
Tests whether a tetrahedron is in the free list.
bool tet_is_virtual(index_t t) const
Tests whether a tetrahedron is a virtual one.
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.
index_t tet_vertex(index_t t, index_t lv) const
Gets the index of a vertex of a tetrahedron.
index_t get_facet_by_halfedge(index_t t, index_t v1, index_t v2) const
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.
index_t tet_adjacent(index_t t, index_t lf) const
Gets the index of a tetrahedron adjacent to another one.
void remove_tet_from_list(index_t t)
Removes a tetrahedron from the linked list it belongs to.
void add_tet_to_list(index_t t, index_t &first, index_t &last)
Adds a tetrahedron to a linked list.
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.
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.
index_t new_tetrahedron(index_t v1, index_t v2, index_t v3, index_t v4)
Creates a new tetrahedron.
index_t find_tet_adjacent(index_t t1, index_t t2) const
Finds the index of the facet accros which t1 is adjacent to t2_in.
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 next_around_halfedge(index_t &t, index_t v1, index_t v2) const
Gets the next tetrahedron around an oriented edge of a tetrahedron.
static index_t find_4(const index_t *T, index_t v)
Finds the index of an integer in an array of four integers.
index_t new_tetrahedron()
Creates a new tetrahedron.
void show_list(index_t first, const std::string &list_name) const
For debugging purposes, displays a tetrahedron.
bool tet_is_conflict(index_t t, const double *p) const
Tests whether a given tetrahedron is in conflict with a given 3d point.
index_t finite_tet_vertex(index_t t, index_t lv) const
Gets the index of a vertex of a tetrahedron.
void set_tet(index_t t, index_t v0, index_t v1, index_t v2, index_t v3, index_t a0, index_t a1, index_t a2, index_t a3)
Sets the vertices and adjacent tetrahedra of a tetrahedron.
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.
index_t max_t() const
Maximum valid index for a tetrahedron.
void set_tet_vertex(index_t t, index_t lv, index_t v)
Sets a tetrahedron-to-vertex adjacency.
index_t find_tet_vertex(index_t t, index_t v) const
Finds the index of the vertex in a tetrahedron.
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.
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.
RegularWeightedDelaunay3d(coord_index_t dimension=4)
Constructs a new Regular Delaunay3d triangulation.
~RegularWeightedDelaunay3d() override
RegularWeightedDelaunay3d destructor.
Vector with aligned memory allocation.
Definition memory.h:660
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
Global Vorpaline namespace.
Definition basic.h:55
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.