Geogram  Version 1.8.9-rc
A programming library of geometric algorithms
GEO::Delaunay3d Class Reference

Implementation of Delaunay in 3d. More...

#include <geogram/delaunay/delaunay_3d.h>

Inheritance diagram for GEO::Delaunay3d:
GEO::Delaunay GEO::Counted GEO::RegularWeightedDelaunay3d

Public Member Functions

 Delaunay3d (coord_index_t dimension=3)
 Constructs a new Delaunay3d. More...
 
void set_vertices (index_t nb_vertices, const double *vertices) override
 Sets the vertices of this Delaunay, and recomputes the cells. More...
 
index_t nearest_vertex (const double *p) const override
 Computes the nearest vertex from a query point. More...
 
- Public Member Functions inherited from GEO::Delaunay
coord_index_t dimension () const
 Gets the dimension of this Delaunay. More...
 
index_t cell_size () const
 Gets the number of vertices in each cell. More...
 
void set_reorder (bool x)
 Specifies whether vertices should be reordered. More...
 
virtual void set_BRIO_levels (const vector< index_t > &levels)
 Specifies the bounds of each level to be used when hierarchic ordering is specified from outside. More...
 
const double * vertices_ptr () const
 Gets a pointer to the array of vertices. More...
 
const double * vertex_ptr (index_t i) const
 Gets a pointer to a vertex by its global index. More...
 
index_t nb_vertices () const
 Gets the number of vertices. More...
 
virtual bool supports_constraints () const
 Tests whether constraints are supported by this Delaunay. More...
 
virtual void set_constraints (const Mesh *mesh)
 Defines the constraints. More...
 
void set_refine (bool x)
 Specifies whether the mesh should be refined. More...
 
bool get_refine () const
 Tests whether mesh refinement is selected. More...
 
void set_quality (double qual)
 Specifies the desired quality for mesh elements when refinement is enabled (. More...
 
const Meshconstraints () const
 Gets the constraints. More...
 
index_t nb_cells () const
 Gets the number of cells. More...
 
index_t nb_finite_cells () const
 Gets the number of finite cells. More...
 
const signed_index_tcell_to_v () const
 Gets a pointer to the cell-to-vertex incidence array. More...
 
const signed_index_tcell_to_cell () const
 Gets a pointer to the cell-to-cell adjacency array. More...
 
signed_index_t cell_vertex (index_t c, index_t lv) const
 Gets a vertex index by cell index and local vertex index. More...
 
signed_index_t cell_adjacent (index_t c, index_t lf) const
 Gets an adjacent cell index by cell index and local facet index. More...
 
bool cell_is_infinite (index_t c) const
 Tests whether a cell is infinite. More...
 
bool cell_is_finite (index_t c) const
 Tests whether a cell is finite. More...
 
index_t index (index_t c, signed_index_t v) const
 Retrieves a local vertex index from cell index and global vertex index. More...
 
index_t adjacent_index (index_t c1, index_t c2) const
 Retrieves a local facet index from two adacent cell global indices. More...
 
signed_index_t vertex_cell (index_t v) const
 Gets an incident cell index by a vertex index. More...
 
signed_index_t next_around_vertex (index_t c, index_t lv) const
 Traverses the list of cells incident to a vertex. More...
 
void get_neighbors (index_t v, vector< index_t > &neighbors) const
 Gets the one-ring neighbors of vertex v. More...
 
void save_histogram (std::ostream &out) const
 Saves the histogram of vertex degree (can be visualized with gnuplot). More...
 
bool stores_neighbors () const
 Tests whether neighbors are stored. More...
 
void set_stores_neighbors (bool x)
 Specifies whether neighbors should be stored. More...
 
bool stores_cicl () const
 Tests whether incident tetrahedra lists are stored. More...
 
void set_stores_cicl (bool x)
 Specifies whether incident tetrahedra lists should be stored. More...
 
bool keeps_infinite () const
 Tests whether infinite elements are kept. More...
 
void set_keeps_infinite (bool x)
 Sets whether infinite elements should be kept. More...
 
bool thread_safe () const
 Tests whether thread-safe mode is active. More...
 
void set_thread_safe (bool x)
 Specifies whether thread-safe mode should be used. More...
 
void set_default_nb_neighbors (index_t x)
 Sets the default number of stored neighbors. More...
 
index_t default_nb_neighbors () const
 Gets the default number of stored neighbors. More...
 
void clear_neighbors ()
 Frees all memory used for neighbors storage.
 
void set_keep_regions (bool x)
 Specifies whether all internal regions should be kept. More...
 
virtual index_t region (index_t t) const
 Gets the region id associated with a tetrahedron. More...
 
virtual void store_neighbors_CB (index_t i)
 Stores the neighbors of a vertex. More...
 
- Public Member Functions inherited from GEO::Counted
void ref () const
 Increments the reference count. More...
 
void unref () const
 Decrements the reference count. More...
 
bool is_shared () const
 Check if the object is shared. More...
 
int nb_refs () const
 Gets the number of references that point to this object. More...
 

Protected Member Functions

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. More...
 
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. More...
 
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. More...
 
index_t insert (index_t v, index_t hint=NO_TETRAHEDRON)
 Inserts a point in the triangulation. More...
 
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. More...
 
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. More...
 
index_t stellate_cavity (index_t v)
 Creates a star of tetrahedra filling the conflict zone. More...
 
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. More...
 
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. More...
 
index_t max_t () const
 Maximum valid index for a tetrahedron. More...
 
bool tet_is_in_list (index_t t) const
 Tests whether a tetrahedron belongs to a linked list. More...
 
index_t tet_next (index_t t) const
 Gets the index of a successor of a tetrahedron. More...
 
void add_tet_to_list (index_t t, index_t &first, index_t &last)
 Adds a tetrahedron to a linked list. More...
 
void remove_tet_from_list (index_t t)
 Removes a tetrahedron from the linked list it belongs to. More...
 
bool tet_is_finite (index_t t) const
 Tests whether a given tetrahedron is a finite one. More...
 
bool tet_is_real (index_t t) const
 Tests whether a tetrahedron is a real one. More...
 
bool tet_is_virtual (index_t t) const
 Tests whether a tetrahedron is a virtual one. More...
 
bool tet_is_free (index_t t) const
 Tests whether a tetrahedron is in the free list. More...
 
index_t new_tetrahedron ()
 Creates a new tetrahedron. More...
 
index_t new_tetrahedron (signed_index_t v1, signed_index_t v2, signed_index_t v3, signed_index_t v4)
 Creates a new tetrahedron. More...
 
void set_tet_mark_stamp (index_t stamp)
 Generates a unique stamp for marking tets. More...
 
bool tet_is_marked (index_t t) const
 Tests whether a tetrahedron is marked. More...
 
void mark_tet (index_t t)
 Marks a tetrahedron. More...
 
signed_index_t tet_vertex (index_t t, index_t lv) const
 Gets the index of a vertex of a tetrahedron. More...
 
index_t find_tet_vertex (index_t t, signed_index_t v) const
 Finds the index of the vertex in a tetrahedron. More...
 
index_t finite_tet_vertex (index_t t, index_t lv) const
 Gets the index of a vertex of a tetrahedron. More...
 
void set_tet_vertex (index_t t, index_t lv, signed_index_t v)
 Sets a tetrahedron-to-vertex adjacency. More...
 
signed_index_t tet_adjacent (index_t t, index_t lf) const
 Gets the index of a tetrahedron adjacent to another one. More...
 
void set_tet_adjacent (index_t t1, index_t lf1, index_t t2)
 Sets a tetrahedron-to-tetrahedron adjacency. More...
 
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. More...
 
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. More...
 
index_t get_facet_by_halfedge (index_t t, signed_index_t v1, signed_index_t v2) const
 
void get_facets_by_halfedge (index_t t, signed_index_t v1, signed_index_t v2, index_t &f12, index_t &f21) const
 
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. More...
 
bool tet_is_conflict (index_t t, const double *p) const
 Tests whether a given tetrahedron is in conflict with a given 3d point. More...
 
 ~Delaunay3d () override
 Delaunay3d destructor.
 
void show_tet (index_t t) const
 For debugging purposes, displays a tetrahedron. More...
 
void show_tet_adjacent (index_t t, index_t lf) const
 For debugging purposes, displays a tetrahedron adjacency. More...
 
void show_list (index_t first, const std::string &list_name) const
 For debugging purposes, displays a tetrahedron. More...
 
void check_combinatorics (bool verbose=false) const
 For debugging purposes, tests some combinatorial properties.
 
void check_geometry (bool verbose=false) const
 For debugging purposes, test some geometrical properties.
 
- Protected Member Functions inherited from GEO::Delaunay
 Delaunay (coord_index_t dimension)
 Creates a new Delaunay triangulation. More...
 
 ~Delaunay () override
 Delaunay destructor.
 
virtual void get_neighbors_internal (index_t v, vector< index_t > &neighbors) const
 Internal implementation for get_neighbors (with vector). More...
 
virtual void set_arrays (index_t nb_cells, const signed_index_t *cell_to_v, const signed_index_t *cell_to_cell)
 Sets the arrays that represent the combinatorics of this Delaunay. More...
 
virtual void update_v_to_cell ()
 Stores for each vertex v a cell incident to v.
 
virtual void update_cicl ()
 Updates the circular incident cell lists. More...
 
virtual void update_neighbors ()
 Computes the stored neighbor lists.
 
void set_next_around_vertex (index_t c1, index_t lv, index_t c2)
 Sets the circular incident edge list. More...
 
void set_dimension (coord_index_t dim)
 Sets the dimension of this Delaunay. More...
 
- Protected Member Functions inherited from GEO::Counted
 Counted ()
 Creates a reference counted object. More...
 
virtual ~Counted ()
 Destroys a reference counted object. More...
 

Static Protected Member Functions

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. More...
 
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. More...
 

Static Protected Attributes

static const index_t NO_TETRAHEDRON = index_t(-1)
 Symbolic constant for uninitialized hint. More...
 
static const index_t NOT_IN_LIST = index_t(~0)
 Default symbolic value of the cell_next_ field that indicates that a tetrahedron is not in a linked list. More...
 
static const index_t NOT_IN_LIST_BIT = index_t(1u << 31)
 If cell_next_[t] & NOT_IN_LIST_BIT != 0, then t is not in a linked list. More...
 
static const index_t END_OF_LIST = ~(NOT_IN_LIST_BIT)
 Symbolic value of the cell_next_ field that indicates the end of list in a linked list of tetrahedra.
 
static const signed_index_t VERTEX_AT_INFINITY = -1
 Symbolic value for a vertex of a tetrahedron that indicates a virtual tetrahedron. More...
 

Additional Inherited Members

- Static Public Member Functions inherited from GEO::Delaunay
static Delaunaycreate (coord_index_t dim, const std::string &name="default")
 Creates a Delaunay triangulation of the specified dimension. More...
 
static void initialize ()
 This function needs to be called once before using the Delaunay class. More...
 
- Static Public Member Functions inherited from GEO::Counted
static void ref (const Counted *counted)
 Increments the reference count. More...
 
static void unref (const Counted *counted)
 Decrements the reference count. More...
 
- Protected Attributes inherited from GEO::Delaunay
coord_index_t dimension_
 
index_t vertex_stride_
 
index_t cell_size_
 
index_t cell_v_stride_
 
index_t cell_neigh_stride_
 
const double * vertices_
 
index_t nb_vertices_
 
index_t nb_cells_
 
const signed_index_tcell_to_v_
 
const signed_index_tcell_to_cell_
 
vector< signed_index_tv_to_cell_
 
vector< signed_index_tcicl_
 
bool is_locked_
 
PackedArrays neighbors_
 
bool store_neighbors_
 
index_t default_nb_neighbors_
 
bool do_reorder_
 If true, uses BRIO reordering (in some implementations)

 
const Meshconstraints_
 
bool refine_
 
double quality_
 
bool store_cicl_
 It true, circular incident tet lists are stored.
 
bool keep_infinite_
 If true, infinite vertex and infinite simplices are kept.
 
index_t nb_finite_cells_
 If keep_infinite_ is true, then finite cells are 0..nb_finite_cells_-1 and infinite cells are nb_finite_cells_ ... nb_cells_.
 
bool keep_regions_
 

Detailed Description

Implementation of Delaunay in 3d.

This package uses concepts inspired by two triangulation softwares, CGAL and tetgen, described in the following references. This package follows the idea used in CGAL of traversing the cavity from inside, since it traverses less tetrahedra than when traversing from outside.

  • Jean-Daniel Boissonnat, Olivier Devillers, Monique Teillaud, and Mariette Yvinec. Triangulations in CGAL. In Proc. 16th Annu. ACM Sympos. Comput. Geom., pages 11-18, 2000.
  • Hang Si, Constrained Delaunay tetrahedral mesh generation and refinement. Finite elements in Analysis and Design, 46 (1-2):33–46, 2010.

Note that the algorithm here does not support vertex deletion nor degenerate input with all coplanar or all colinear points (use CGAL instead if you have these requirements).

The core algorithm used in this code, CGAL and tetgen was independently and simultaneously discovered by Bowyer and Watson:

  • Adrian Bowyer, "Computing Dirichlet tessellations", Comput. J., vol. 24, no 2, 1981, p. 162-166
  • David F. Watson, "Computing the n-dimensional Delaunay tessellation with application to Voronoi polytopes", Comput. J., vol. 24, no 2, 1981, p. 167-172

The spatial reordering method, that dramatically increases the performances, also used in this code, CGAL and tetgen was introduced in the following references. The second one is a smart implementation based on the std::nth_element() function of the STL, that inspired the compute_BRIO_ordering() function of this package.

  • Nina Amenta, Sunghee Choi and Gunter Rote, "Incremental constructions con brio", ACM Symposium on Computational Geometry 2003.
  • Christophe Delage and Olivier Devillers. Spatial Sorting. In CGAL User and Reference Manual. CGAL Editorial Board, 3.9 edition, 2011

The locate() function is based on the following two references. The first one randomizes the choice of the next tetrahedron. The second one uses an inexact locate() function to initialize the exact one (it is called "structural filtering"). The first idea is used in both CGAL and tetgen, and the second one is used in CGAL.

  • Walking in a triangulation, O Devillers, S Pion, M Teillaud 17th Annual Symposium on Computational geometry, 106-114
  • Stefan Funke , Kurt Mehlhorn and Stefan Naher, "Structural filtering, a paradigm for efficient and exact geometric programs", Comput. Geom., 1999

Definition at line 107 of file delaunay_3d.h.

Constructor & Destructor Documentation

◆ Delaunay3d()

GEO::Delaunay3d::Delaunay3d ( coord_index_t  dimension = 3)

Constructs a new Delaunay3d.

Parameters
[in]dimensiondimension of the triangulation (3 or 4). If dimension = 4, this creates a regular triangulation (dual of a power diagram). In this case:
  • the input points are 4d points, were the fourth coordinate of point \( i \) is \( \sqrt{W - w_i} \) where \( W \) is the maximum of the weights of all the points and \d$ w_i $ is the weight associated with vertex \( i \).
  • the constructed combinatorics is a tetrahedralized volume (3d and not 4d although dimension() returns 4). This tetrahedralized volume corresponds to the regular triangulation of the weighted points.

Member Function Documentation

◆ add_tet_to_list()

void GEO::Delaunay3d::add_tet_to_list ( index_t  t,
index_t first,
index_t last 
)
inlineprotected

Adds a tetrahedron to a linked list.

Tetrahedra can be linked, it is used to manage the free list that recycles deleted tetrahedra.

Parameters
[in]tthe index of the tetrahedron
[in,out]firstfirst item of the list or END_OF_LIST if the list is empty
[in,out]lastlast item of the list or END_OF_LIST if the list is empty

Definition at line 484 of file delaunay_3d.h.

◆ create_first_tetrahedron()

bool GEO::Delaunay3d::create_first_tetrahedron ( index_t iv0,
index_t iv1,
index_t iv2,
index_t iv3 
)
protected

Finds in the pointset a set of four non-coplanar points.

This function is used to initiate the incremental Delaunay construction.

Parameters
[out]iv0index of the first vertex
[out]iv1index of the second vertex
[out]iv2index of the third vertex
[out]iv3index of the fourth vertex
Return values
trueif a set of four non-coplanar points was found
falseif all the points are coplanar

◆ find_4()

static index_t GEO::Delaunay3d::find_4 ( const signed_index_t T,
signed_index_t  v 
)
inlinestaticprotected

Finds the index of an integer in an array of four integers.

Parameters
[in]Ta const pointer to an array of four integers
[in]vthe integer to retrieve in T
Returns
the index (0,1,2 or 3) of v in T
Precondition
The four entries of T are different and one of them is equal to v.

Definition at line 1035 of file delaunay_3d.h.

◆ find_conflict_zone()

void GEO::Delaunay3d::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 
)
protected

Determines the list of tetrahedra in conflict with a given point.

Parameters
[in]vthe index of the point to be inserted
[in]tthe index of a tetrahedron that contains p, as returned by locate()
[in]orientan array of four signs indicating the orientation of p with respect to the four faces of t, as returned by locate()
[out]t_bndrya tetrahedron adjacent to the boundary of the conflict zone
[out]f_bndrythe facet along which t_bndry is adjacent to the boundary of the conflict zone
[out]firstthe index of the first tetrahedron in conflict
[out]lastthe index of the last tetrahedron in conflict The other tetrahedra are linked, and can be traversed from first by using tet_next() until last or END_OF_LIST is reached. The conflict zone can be empty under two circumstances:
  • the vertex v already exists in the triangulation
  • the triangulation is weighted and v is not visible in either cases, both first and last contain END_OF_LIST

◆ find_conflict_zone_iterative()

void GEO::Delaunay3d::find_conflict_zone_iterative ( const double *  p,
index_t  t,
index_t t_bndry,
index_t f_bndry,
index_t first,
index_t last 
)
protected

This function is used to implement find_conflict_zone.

This function detects the neighbors of t that are in the conflict zone and calls itself recursively on them.

Parameters
[in]pthe point to be inserted
[in]tindex of a tetrahedron in the fonflict zone
[out]t_bndrya tetrahedron adjacent to the boundary of the conflict zone
[out]f_bndrythe facet along which t_bndry is adjacent to the boundary of the conflict zone
[out]firstthe index of the first tetrahedron in conflict
[out]lastthe index of the last tetrahedron in conflict
Precondition
The tetrahedron t was alredy marked as conflict (tet_is_in_list(t))

◆ find_tet_adjacent()

index_t GEO::Delaunay3d::find_tet_adjacent ( index_t  t1,
index_t  t2_in 
) const
inlineprotected

Finds the index of the facet accros which t1 is adjacent to t2_in.

Parameters
[in]t1first tetrahedron
[in]t2_insecond tetrahedron
Returns
f such that tet_adjacent(t1,f)==t2_in
Precondition
t1 and t2_in are adjacent

Definition at line 801 of file delaunay_3d.h.

◆ find_tet_vertex()

index_t GEO::Delaunay3d::find_tet_vertex ( index_t  t,
signed_index_t  v 
) const
inlineprotected

Finds the index of the vertex in a tetrahedron.

Parameters
[in]tthe tetrahedron
[in]vthe vertex
Returns
iv such that tet_vertex(t,v)==iv
Precondition
t is incident to v

Definition at line 732 of file delaunay_3d.h.

◆ finite_tet_vertex()

index_t GEO::Delaunay3d::finite_tet_vertex ( index_t  t,
index_t  lv 
) const
inlineprotected

Gets the index of a vertex of a tetrahedron.

Parameters
[in]tindex of the tetrahedron
[in]lvlocal vertex (0,1,2 or 3) index in t
Returns
the global index of the lvth vertex of tetrahedron t
Precondition
Vertex lv of tetrahedron t is not at infinity

Definition at line 747 of file delaunay_3d.h.

◆ get_facet_by_halfedge()

index_t GEO::Delaunay3d::get_facet_by_halfedge ( index_t  t,
signed_index_t  v1,
signed_index_t  v2 
) const
inlineprotected

Gets the local facet index incident to an oriented halfedge.

Parameters
[in]tindex of the tetrahedron
[in]v1global index of the first extremity
[in]v2global index of the second extremity
Returns
the local index of the facet incident to the oriented edge v1, v2.

Definition at line 865 of file delaunay_3d.h.

◆ get_facets_by_halfedge()

void GEO::Delaunay3d::get_facets_by_halfedge ( index_t  t,
signed_index_t  v1,
signed_index_t  v2,
index_t f12,
index_t f21 
) const
inlineprotected

Gets the local facet indices incident to an oriented halfedge.

Parameters
[in]tindex of the tetrahedron
[in]v1global index of the first extremity
[in]v2global index of the second extremity
[out]f12the local index of the facet indicent to the halfedge [v1,v2]
[out]f21the local index of the facet indicent to the halfedge [v2,v1]

Definition at line 890 of file delaunay_3d.h.

◆ get_neighbor_along_conflict_zone_border()

bool GEO::Delaunay3d::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
inlineprotected

Finds the neighbor of a tetrahedron on the border of the conflict zone.

This function is used by stellate_conflict_zone_iterative()

Parameters
[in]t1a tetrahedron on the border of the conflict zone
[in]t1fborderthe local facet index of t1 along which it is on the border of the conflict zone
[in]t1ft2the local facet index of t1 that will be traversed
[out]t2a tetrahedron on the border of the conflict zone, with an edge common to facets t1fborder and t1ft2 of tetrahedron t1
[out]t2fborderthe local facet index of t2 along which it is on the border of the conflict zone
[out]t2ft1the local index of the facet of t2 that has a common edge with facets t1fborder and t1ft2 of tetrahedron t1
Return values
trueif t2 is a newly created tetrahedron
falseif t2 is an old tetrahedron in conflict

Definition at line 325 of file delaunay_3d.h.

◆ insert()

index_t GEO::Delaunay3d::insert ( index_t  v,
index_t  hint = NO_TETRAHEDRON 
)
protected

Inserts a point in the triangulation.

Parameters
[in]vthe index of the point to be inserted
[in]hintthe index of a tetrahedron as near as possible to v, or -1 if unspecified
Returns
the index of one of the tetrahedra incident to point v

◆ locate()

index_t GEO::Delaunay3d::locate ( const double *  p,
index_t  hint = NO_TETRAHEDRON,
bool  thread_safe = false,
Sign orient = nullptr 
) const
protected

Finds the tetrahedron that contains a point.

If the point is on a face, edge or vertex, the function returns one of the tetrahedra incident to that face, edge or vertex.

Parameters
[in]pa pointer to the coordinates of the point
[in]thread_safeif true, a global spinlock is used to protect the calls to random(), this is necessary if multiple threads use locate() simultaneously
[out]orienta pointer to an array of four Signs or nullptr. If non-nullptr, returns the orientation with respect to the four facets of the tetrahedron that contains p.
Returns
the index of a tetrahedron that contains p. If the point is outside the convex hull of the inserted so-far points, then the returned tetrahedron is a virtual one (first vertex is the "vertex at infinity" of index -1) or NO_TETRAHEDRON if the virtual tetrahedra were previously removed.

◆ locate_inexact()

index_t GEO::Delaunay3d::locate_inexact ( const double *  p,
index_t  hint,
index_t  max_iter 
) const
protected

Finds the tetrahedron that (approximately) contains a point using inexact predicates.

The result of this function can be used as a hint for locate(). It accelerates locate as compared to calling it directly. This technique is referred to as "structural filtering".

Parameters
[in]pa pointer to the coordinates of the point
[in]max_itermaximum number of traversed tets
Returns
the index of a tetrahedron that (approximately) contains p. If the point is outside the convex hull of the inserted so-far points, then the returned tetrahedron is a virtual one (first vertex is the "vertex at infinity" of index -1) or NO_TETRAHEDRON if the virtual tetrahedra were previously removed.

◆ mark_tet()

void GEO::Delaunay3d::mark_tet ( index_t  t)
inlineprotected

Marks a tetrahedron.

A tetrahedron is marked whenever it is detected as non-conflict. The same space is used for marking and for chaining the conflict list. The index of the point being inserted is used as a time-stamp for marking tetrahedra. A tetrahedron can be in the following states:

  • in list
  • not in list and marked
  • not in list and not marked
    Parameters
    [in]tindex of the tetrahedron to be marked

Definition at line 683 of file delaunay_3d.h.

◆ max_t()

index_t GEO::Delaunay3d::max_t ( ) const
inlineprotected

Maximum valid index for a tetrahedron.

This includes not only real tetrahedra, but also the virtual ones on the border, the conflict list and the free list.

Returns
the maximum valid index for a tetrahedron

Definition at line 389 of file delaunay_3d.h.

◆ nearest_vertex()

index_t GEO::Delaunay3d::nearest_vertex ( const double *  p) const
overridevirtual

Computes the nearest vertex from a query point.

Parameters
[in]pquery point
Returns
the index of the nearest vertex

Reimplemented from GEO::Delaunay.

◆ new_tetrahedron() [1/2]

index_t GEO::Delaunay3d::new_tetrahedron ( )
inlineprotected

Creates a new tetrahedron.

Uses either a tetrahedron recycled from the free list, or creates a new one by expanding the two indices arrays.

Returns
the index of the newly created tetrahedron

Definition at line 592 of file delaunay_3d.h.

◆ new_tetrahedron() [2/2]

index_t GEO::Delaunay3d::new_tetrahedron ( signed_index_t  v1,
signed_index_t  v2,
signed_index_t  v3,
signed_index_t  v4 
)
inlineprotected

Creates a new tetrahedron.

Sets the vertices. Adjacent tetrahedra index are left uninitialized. Uses either a tetrahedron recycled from the free list, or creates a new one by expanding the two indices arrays.

Parameters
[in]v1index of the first vertex
[in]v2index of the second vertex
[in]v3index of the third vertex
[in]v4index of the fourth vertex
Returns
the index of the newly created tetrahedron

Definition at line 628 of file delaunay_3d.h.

◆ next_around_halfedge()

index_t GEO::Delaunay3d::next_around_halfedge ( index_t t,
signed_index_t  v1,
signed_index_t  v2 
) const
inlineprotected

Gets the next tetrahedron around an oriented edge of a tetrahedron.

Parameters
[in,out]tthe tetrahedron
[in]v1global index of the first extremity of the edge
[in]v2global index of the second extremity of the edge
Returns
the next tetrahedron from t around the oriented edge (v1 v2).

Definition at line 929 of file delaunay_3d.h.

◆ remove_tet_from_list()

void GEO::Delaunay3d::remove_tet_from_list ( index_t  t)
inlineprotected

Removes a tetrahedron from the linked list it belongs to.

Tetrahedra can be linked, it is used to manage both the free list that recycles deleted tetrahedra and the list of tetrahedra in conflict with the inserted point.

Parameters
[in]tthe index of the tetrahedron

Definition at line 506 of file delaunay_3d.h.

◆ set_tet()

void GEO::Delaunay3d::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 
)
inlineprotected

Sets the vertices and adjacent tetrahedra of a tetrahedron.

Parameters
[in]tindex of the tetrahedron
[in]v0index of the first vertex
[in]v1index of the second vertex
[in]v2index of the third vertex
[in]v3index of the fourth vertex
[in]a0index of the adjacent tetrahedron opposite to v0
[in]a1index of the adjacent tetrahedron opposite to v1
[in]a2index of the adjacent tetrahedron opposite to v2
[in]a3index of the adjacent tetrahedron opposite to v3

Definition at line 837 of file delaunay_3d.h.

◆ set_tet_adjacent()

void GEO::Delaunay3d::set_tet_adjacent ( index_t  t1,
index_t  lf1,
index_t  t2 
)
inlineprotected

Sets a tetrahedron-to-tetrahedron adjacency.

Parameters
[in]t1index of the first tetrahedron
[in]lf1local facet index (0,1,2 or 3) in t1
[in]t2index of the tetrahedron adjacent to t1 accros lf1

Definition at line 786 of file delaunay_3d.h.

◆ set_tet_mark_stamp()

void GEO::Delaunay3d::set_tet_mark_stamp ( index_t  stamp)
inlineprotected

Generates a unique stamp for marking tets.

Storage is shared for list-chaining and stamp-marking (both are mutually exclusive), therefore the stamp has the NOT_IN_LIST_BIT set.

Parameters
[in]stampthe unique stamp for marking tets

Definition at line 647 of file delaunay_3d.h.

◆ set_tet_vertex()

void GEO::Delaunay3d::set_tet_vertex ( index_t  t,
index_t  lv,
signed_index_t  v 
)
inlineprotected

Sets a tetrahedron-to-vertex adjacency.

Parameters
[in]tindex of the tetrahedron
[in]lvlocal vertex index (0,1,2 or 3) in t
[in]vglobal index of the vertex

Definition at line 760 of file delaunay_3d.h.

◆ set_vertices()

void GEO::Delaunay3d::set_vertices ( index_t  nb_vertices,
const double *  vertices 
)
overridevirtual

Sets the vertices of this Delaunay, and recomputes the cells.

Parameters
[in]nb_verticesnumber of vertices
[in]verticesa pointer to the coordinates of the vertices, as a contiguous array of doubles

Reimplemented from GEO::Delaunay.

◆ show_list()

void GEO::Delaunay3d::show_list ( index_t  first,
const std::string &  list_name 
) const
protected

For debugging purposes, displays a tetrahedron.

Parameters
[in]firstindex of the first tetrahedron in the list
[in]list_namename of the list, will be displayed as well

◆ show_tet()

void GEO::Delaunay3d::show_tet ( index_t  t) const
protected

For debugging purposes, displays a tetrahedron.

Parameters
[in]tindex of the tetrahedron to display.

◆ show_tet_adjacent()

void GEO::Delaunay3d::show_tet_adjacent ( index_t  t,
index_t  lf 
) const
protected

For debugging purposes, displays a tetrahedron adjacency.

Parameters
[in]tindex of the tetrahedron to display.
[in]lflocal index (0,1,2 or 3) of the tetrahedron facet adjacenty to display.

◆ stellate_cavity()

index_t GEO::Delaunay3d::stellate_cavity ( index_t  v)
protected

Creates a star of tetrahedra filling the conflict zone.

Parameters
[in]vthe index of the point to be inserted

This function is used when the Cavity computed when traversing the conflict zone is OK, that is to say when its array sizes were not exceeded.

Returns
the index of one the newly created tetrahedron

◆ stellate_conflict_zone_iterative()

index_t GEO::Delaunay3d::stellate_conflict_zone_iterative ( index_t  v,
index_t  t_bndry,
index_t  f_bndry,
index_t  prev_f = index_t(-1) 
)
protected

Creates a star of tetrahedra filling the conflict zone.

For each tetrahedron facet on the border of the conflict zone, a new tetrahedron is created, resting on the facet and incident to vertex v. The function is called recursively until the entire conflict zone is filled.

Parameters
[in]vthe index of the point to be inserted
[in]t_bndryindex of a tetrahedron on the border of the conflict zone.
[in]f_bndryindex of the facet along which t_bndry is incident to the border of the conflict zone
[in]prev_fthe facet of t_bndry connected to the tetrahedron that t_bndry was reached from, or index_t(-1) if it is the first tetrahedron.
Returns
the index of one the newly created tetrahedron

◆ tet_adjacent()

signed_index_t GEO::Delaunay3d::tet_adjacent ( index_t  t,
index_t  lf 
) const
inlineprotected

Gets the index of a tetrahedron adjacent to another one.

Parameters
[in]tindex of the tetrahedron
[in]lflocal facet (0,1,2 or 3) index in t
Returns
the tetrahedron adjacent to t accorss facet lf

Definition at line 772 of file delaunay_3d.h.

◆ tet_facet_vertex()

static index_t GEO::Delaunay3d::tet_facet_vertex ( index_t  f,
index_t  v 
)
inlinestaticprotected

Returns the local index of a vertex by facet and by local vertex index in the facet.

tet facet vertex is such that the tetrahedron formed with:

  • vertex lv
  • tet_facet_vertex(lv,0)
  • tet_facet_vertex(lv,1)
  • tet_facet_vertex(lv,2) has the same orientation as the original tetrahedron for any vertex lv.
    Parameters
    [in]flocal facet index, in (0,1,2,3)
    [in]vlocal vertex index, in (0,1,2)
    Returns
    the local tetrahedron vertex index of vertex v in facet f

Definition at line 706 of file delaunay_3d.h.

◆ tet_is_conflict()

bool GEO::Delaunay3d::tet_is_conflict ( index_t  t,
const double *  p 
) const
inlineprotected

Tests whether a given tetrahedron is in conflict with a given 3d point.

A real tetrahedron is in conflict with a point whenever the point is contained by its circumscribed sphere, and a virtual tetrahedron is in conflict with a point whenever the tetrahedron formed by its real face and with the point has positive orientation.

Parameters
[in]tthe index of the tetrahedron
[in]pa pointer to the coordinates of the point
Return values
trueif point p is in conflict with tetrahedron t
falseotherwise

Definition at line 952 of file delaunay_3d.h.

◆ tet_is_finite()

bool GEO::Delaunay3d::tet_is_finite ( index_t  t) const
inlineprotected

Tests whether a given tetrahedron is a finite one.

Infinite tetrahedra are the ones that are incident to the infinite vertex (index -1)

Parameters
[in]tthe index of the tetrahedron
Return values
trueif t is finite
falseotherwise

Definition at line 530 of file delaunay_3d.h.

◆ tet_is_free()

bool GEO::Delaunay3d::tet_is_free ( index_t  t) const
inlineprotected

Tests whether a tetrahedron is in the free list.

Deleted tetrahedra are recycled in a free list.

Parameters
[in]tindex of the tetrahedron
Return values
trueif tetrahedron t is in the free list
falseotherwise

Definition at line 581 of file delaunay_3d.h.

◆ tet_is_in_list()

bool GEO::Delaunay3d::tet_is_in_list ( index_t  t) const
inlineprotected

Tests whether a tetrahedron belongs to a linked list.

Tetrahedra can be linked, it is used to manage both the free list that recycles deleted tetrahedra, the conflict region and the list of newly created tetrahedra. In addition, a tetrahedron that is not in a list can be marked. The same space is used for marking and chaining tetrahedra in lists. A tetrahedron can be in the following states:

  • in list
  • not in list and marked
  • not in list and not marked
    Parameters
    [in]tthe index of the tetrahedron
    Return values
    trueif tetrahedron t belongs to a linked list
    falseotherwise

Definition at line 453 of file delaunay_3d.h.

◆ tet_is_marked()

bool GEO::Delaunay3d::tet_is_marked ( index_t  t) const
inlineprotected

Tests whether a tetrahedron is marked.

A tetrahedron is marked whenever it is detected as non-conflict. The index of the point being inserted is used as a time-stamp for marking tetrahedra. The same space is used for marking and for chaining the conflict list. A tetrahedron can be in the following states:

  • in list
  • not in list and marked
  • not in list and not marked
    Parameters
    [in]tindex of the tetrahedron
    Return values
    trueif tetrahedron t is marked
    falseotherwise

Definition at line 666 of file delaunay_3d.h.

◆ tet_is_real()

bool GEO::Delaunay3d::tet_is_real ( index_t  t) const
inlineprotected

Tests whether a tetrahedron is a real one.

Real tetrahedra are incident to four user-specified vertices (there are also virtual tetrahedra that are incident to the vertex at infinity, with index -1)

Parameters
[in]tindex of the tetrahedron
Return values
trueif tetrahedron t is a real one
falseotherwise

Definition at line 549 of file delaunay_3d.h.

◆ tet_is_virtual()

bool GEO::Delaunay3d::tet_is_virtual ( index_t  t) const
inlineprotected

Tests whether a tetrahedron is a virtual one.

Virtual tetrahedra are tetrahedra incident to the vertex at infinity.

Parameters
[in]tindex of the tetrahedron
Return values
trueif tetrahedron t is virtual
falseotherwise

Definition at line 562 of file delaunay_3d.h.

◆ tet_next()

index_t GEO::Delaunay3d::tet_next ( index_t  t) const
inlineprotected

Gets the index of a successor of a tetrahedron.

Tetrahedra can be linked, it is used to manage both the free list that recycles deleted tetrahedra.

Parameters
[in]tthe index of the tetrahedron
Return values
END_OF_LISTif the end of the list is reached
theindex of the successor of tetrahedron \t otherwise
Precondition
tet_is_in_list(t)

Definition at line 468 of file delaunay_3d.h.

◆ tet_vertex()

signed_index_t GEO::Delaunay3d::tet_vertex ( index_t  t,
index_t  lv 
) const
inlineprotected

Gets the index of a vertex of a tetrahedron.

Parameters
[in]tindex of the tetrahedron
[in]lvlocal vertex (0,1,2 or 3) index in t
Returns
the global index of the lvth vertex of tetrahedron t or -1 if the vertex is at infinity

Definition at line 719 of file delaunay_3d.h.

Member Data Documentation

◆ NO_TETRAHEDRON

const index_t GEO::Delaunay3d::NO_TETRAHEDRON = index_t(-1)
staticprotected

Symbolic constant for uninitialized hint.

Locate functions can be accelerated by specifying a hint. This constant indicates that no hint is given.

Definition at line 145 of file delaunay_3d.h.

◆ NOT_IN_LIST

const index_t GEO::Delaunay3d::NOT_IN_LIST = index_t(~0)
staticprotected

Default symbolic value of the cell_next_ field that indicates that a tetrahedron is not in a linked list.

This is the default value. Note that it suffices that NOT_IN_LIST_BIT is set for a tetrahedron to be not in any list. A tetrahedron can be:

  • in a list (cell_next_[t] & NOT_IN_LIST_BIT == 0)
  • not in a list and not marked (cell_next_[t] & NOT_IN_LIST_BIT != 0) && (cell_next_[t] != cur_stamp_)
  • not in a list and marked (cell_next_[t] == cur_stamp_)

Definition at line 409 of file delaunay_3d.h.

◆ NOT_IN_LIST_BIT

const index_t GEO::Delaunay3d::NOT_IN_LIST_BIT = index_t(1u << 31)
staticprotected

If cell_next_[t] & NOT_IN_LIST_BIT != 0, then t is not in a linked list.

The other bits of cell_next_[t] are used to store the stamp (i.e. index of the current point being inserted). The stamp is used for marking tetrahedra that were detected as non-conflict when inserting a point. A tetrahedron can be:

  • in a list (cell_next_[t] & NOT_IN_LIST_BIT == 0)
  • not in a list and not marked (cell_next_[t] & NOT_IN_LIST_BIT != 0) && (cell_next_[t] != cur_stamp_)
  • not in a list and marked (cell_next_[t] == cur_stamp_)

Definition at line 426 of file delaunay_3d.h.

◆ VERTEX_AT_INFINITY

const signed_index_t GEO::Delaunay3d::VERTEX_AT_INFINITY = -1
staticprotected

Symbolic value for a vertex of a tetrahedron that indicates a virtual tetrahedron.

The three other vertices then correspond to a facet on the convex hull of the points.

Definition at line 518 of file delaunay_3d.h.


The documentation for this class was generated from the following file: