40#ifndef GEOGRAM_DELAUNAY_DELAUNAY_3D
41#define GEOGRAM_DELAUNAY_DELAUNAY_3D
45#include <geogram/delaunay/cavity.h>
142 static constexpr index_t NO_TETRAHEDRON = NO_INDEX;
180 const double* p,
index_t hint = NO_TETRAHEDRON,
181 bool thread_safe =
false,
182 Sign* orient =
nullptr
344 tet_vertex(t1,
index_t(halfedge_facet_[t1ft2][t1fborder]));
346 tet_vertex(t1,
index_t(halfedge_facet_[t1fborder][t1ft2]));
354 index_t next_t = tet_adjacent(cur_t,cur_f);
355 while(tet_is_in_list(next_t)) {
358 cur_f = get_facet_by_halfedge(cur_t,ev1,ev2);
359 next_t = tet_adjacent(cur_t, cur_f);
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);
387 return cell_to_v_store_.size() / 4;
430 static constexpr index_t END_OF_LIST = ~NOT_IN_LIST_BIT;
452 return (cell_next_[t] & NOT_IN_LIST_BIT) == 0;
468 return cell_next_[t];
484 if(last == END_OF_LIST) {
487 cell_next_[t] = END_OF_LIST;
489 cell_next_[t] = first;
506 cell_next_[t] = NOT_IN_LIST;
515 static constexpr index_t VERTEX_AT_INFINITY = NO_INDEX;
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;
547 return !tet_is_free(t) && tet_is_finite(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) ;
579 return tet_is_in_list(t);
591 if(first_free_ == END_OF_LIST) {
592 cell_to_v_store_.resize(
593 cell_to_v_store_.size() + 4, NO_INDEX
595 cell_to_cell_store_.resize(
596 cell_to_cell_store_.size() + 4, NO_INDEX
601 cell_next_.push_back(
index_t(NOT_IN_LIST));
602 result = max_t() - 1;
604 result = first_free_;
605 first_free_ = tet_next(first_free_);
606 remove_tet_from_list(result);
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;
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;
649 cur_stamp_ = (stamp | NOT_IN_LIST_BIT);
668 return cell_next_[t] == cur_stamp_;
685 cell_next_[t] = cur_stamp_;
710 return index_t(tet_facet_vertex_[f][v]);
723 return cell_to_v_store_[4 * t + lv];
736 const index_t* T = &(cell_to_v_store_[4 * t]);
751 return cell_to_v_store_[4 * t + lv];
763 cell_to_v_store_[4 * t + lv] = v;
775 index_t result = cell_to_cell_store_[4 * t + lf];
790 cell_to_cell_store_[4 * t1 + lf1] = t2;
807 const index_t* T = &(cell_to_cell_store_[4 * t1]);
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;
862 const index_t* T = &(cell_to_v_store_[4 * t]);
866 return index_t(halfedge_facet_[lv1][lv2]);
891 const index_t* T = &(cell_to_v_store_[4 * t]);
894 (T[1] == v1) | ((T[2] == v1) * 2) | ((T[3] == v1) * 3)
898 (T[1] == v2) | ((T[2] == v2) * 2) | ((T[3] == v2) * 3)
905 f12 =
index_t(halfedge_facet_[lv1][lv2]);
906 f21 =
index_t(halfedge_facet_[lv2][lv1]);
921 t, get_facet_by_halfedge(t, v1, v2)
946 pv[i] = (v == NO_INDEX) ?
nullptr : vertex_ptr(v);
951 for(
index_t lf = 0; lf < 4; ++lf) {
953 if(pv[lf] ==
nullptr) {
961 Sign sign = PCK::orient_3d(pv[0],pv[1],pv[2],pv[3]);
974 index_t t2 = tet_adjacent(t, lf);
979 if(tet_is_in_list(t2)) {
984 if(tet_is_marked(t2)) {
988 return tet_is_conflict(t2, p);
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)];
1002 (p - vertex_ptr(0)) /
int(vertex_stride_)
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
1010 return (PCK::in_sphere_3d_SOS(pv[0], pv[1], pv[2], pv[3], p) > 0);
1035 (T[1] == v) | ((T[2] == v) * 2) | ((T[3] == v) * 3)
1097 bool verbose_debug_mode_;
1102 bool benchmark_mode_;
1112 static char tet_facet_vertex_[4][3];
1118 static char halfedge_facet_[4][4];
1124 std::stack<index_t> S_;
1131 class StellateConflictStack {
1144 store_.resize(store_.size()+1);
1158 top().new_t = new_t;
1172 void get_parameters(
1173 index_t& t1, index_t& t1fbord, index_t& t1fprev
1177 t1fbord = index_t(top().t1fbord);
1178 t1fprev = index_t(top().t1fprev);
1189 index_t& new_t, index_t& t1ft2, index_t& t2ft1
1192 new_t = top().new_t;
1210 bool empty()
const {
1211 return store_.empty();
1224 Numeric::uint8 t1fbord ;
1227 Numeric::uint8 t1fprev ;
1228 Numeric::uint8 t1ft2 ;
1229 Numeric::uint8 t2ft1 ;
1240 return *store_.rbegin();
1249 const Frame& top()
const {
1251 return *store_.rbegin();
1254 std::vector<Frame> store_;
1261 StellateConflictStack S2_;
#define geo_debug_assert(x)
Verifies that a condition is met.
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.
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.
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.
Global Vorpaline namespace.
geo_index_t index_t
The type for storing and manipulating indices.
Sign
Integer constants that represent the sign of a value.
geo_coord_index_t coord_index_t
The type for storing coordinate indices, and iterating on the coordinates of a point.
Filtered exact predicates for restricted Voronoi diagrams.