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

Implementation of Delaunay in 2d. More...

#include <geogram/delaunay/delaunay_2d.h>

Inheritance diagram for GEO::Delaunay2d:
GEO::Delaunay GEO::Counted GEO::RegularWeightedDelaunay2d

Public Member Functions

 Delaunay2d (coord_index_t dimension=2)
 Constructs a new Delaunay2d. 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...
 
bool has_empty_cells () const
 Tests whether the Laguerre diagram has empty cells. More...
 
void abort_if_empty_cell (bool x)
 Specifies behavior if an empty cell is detected. 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_triangle (index_t &iv0, index_t &iv1, index_t &iv2)
 Finds in the pointset a set of three non-colinear points. More...
 
index_t locate (const double *p, index_t hint=NO_TRIANGLE, bool thread_safe=false, Sign *orient=nullptr) const
 Finds the triangle that contains a point. More...
 
index_t locate_inexact (const double *p, index_t hint, index_t max_iter) const
 Finds the triangle that (approximately) contains a point using inexact predicates. More...
 
index_t insert (index_t v, index_t hint=NO_TRIANGLE)
 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 &e_bndry, index_t &first, index_t &last)
 Determines the list of triangles in conflict with a given point. More...
 
void find_conflict_zone_iterative (const double *p, index_t t, index_t &t_bndry, index_t &e_bndry, index_t &first, index_t &last)
 This function is used to implement find_conflict_zone. More...
 
index_t stellate_conflict_zone (index_t v, index_t t_bndry, index_t e_bndry)
 Creates a star of triangles filling the conflict zone. More...
 
index_t max_t () const
 Maximum valid index for a triangle. More...
 
bool triangle_is_in_list (index_t t) const
 Tests whether a triangle belongs to a linked list. More...
 
index_t triangle_next (index_t t) const
 Gets the index of a successor of a triangle. More...
 
void add_triangle_to_list (index_t t, index_t &first, index_t &last)
 Adds a triangle to a linked list. More...
 
void remove_triangle_from_list (index_t t)
 Removes a triangle from the linked list it belongs to. More...
 
bool triangle_is_finite (index_t t) const
 Tests whether a given triangle is a finite one. More...
 
bool triangle_is_real (index_t t) const
 Tests whether a triangle is a real one. More...
 
bool triangle_is_virtual (index_t t) const
 Tests whether a triangle is a virtual one. More...
 
bool triangle_is_free (index_t t) const
 Tests whether a triangle is in the free list. More...
 
index_t new_triangle ()
 Creates a new triangle. More...
 
index_t new_triangle (signed_index_t v1, signed_index_t v2, signed_index_t v3)
 Creates a new triangle. More...
 
void set_triangle_mark_stamp (index_t stamp)
 Generates a unique stamp for marking tets. More...
 
bool triangle_is_marked (index_t t) const
 Tests whether a triangle is marked. More...
 
void mark_triangle (index_t t)
 Marks a triangle. More...
 
signed_index_t triangle_vertex (index_t t, index_t lv) const
 Gets the index of a vertex of a triangle. More...
 
index_t find_triangle_vertex (index_t t, signed_index_t v) const
 Finds the index of the vertex in a triangle. More...
 
index_t finite_triangle_vertex (index_t t, index_t lv) const
 Gets the index of a vertex of a triangle. More...
 
void set_triangle_vertex (index_t t, index_t lv, signed_index_t v)
 Sets a triangle-to-vertex adjacency. More...
 
signed_index_t triangle_adjacent (index_t t, index_t le) const
 Gets the index of a triangle adjacent to another one. More...
 
void set_triangle_adjacent (index_t t1, index_t le1, index_t t2)
 Sets a triangle-to-triangle adjacency. More...
 
index_t find_triangle_adjacent (index_t t1, index_t t2_in) const
 Finds the index of the edge 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, index_t a0, index_t a1, index_t a2)
 Sets the vertices and adjacent triangles of a triangle. More...
 
bool triangle_is_conflict (index_t t, const double *p) const
 Tests whether a given triangle is in conflict with a given 3d point. More...
 
 ~Delaunay2d () override
 Delaunay2d destructor.
 
void show_triangle (index_t t) const
 For debugging purposes, displays a triangle. More...
 
void show_triangle_adjacent (index_t t, index_t le) const
 For debugging purposes, displays a triangle adjacency. More...
 
void show_list (index_t first, const std::string &list_name) const
 For debugging purposes, displays a triangle. 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 triangle_edge_vertex (index_t e, index_t v)
 Returns the local index of a vertex by edge and by local vertex index in the edge. More...
 
static index_t find_3 (const signed_index_t *T, signed_index_t v)
 Finds the index of an integer in an array of three integers. More...
 

Static Protected Attributes

static const index_t NO_TRIANGLE = 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 triangle 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 triangles.
 
static const signed_index_t VERTEX_AT_INFINITY = -1
 Symbolic value for a vertex of a triangle that indicates a virtual triangle. 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 2d.

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 triangles 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 trianglesl 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 triangle. 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 106 of file delaunay_2d.h.

Constructor & Destructor Documentation

◆ Delaunay2d()

GEO::Delaunay2d::Delaunay2d ( coord_index_t  dimension = 2)

Constructs a new Delaunay2d.

Parameters
[in]dimensiondimension of the triangulation (2 or 3). If dimension = 3, this creates a regular triangulation (dual of a power diagram). In this case:
  • the input points are 3d points, were the third 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 triangulated surface (2d and not 3d although dimension() returns 3). This triangulated surface corresponds to the regular triangulation of the weighted points.

Member Function Documentation

◆ abort_if_empty_cell()

void GEO::Delaunay2d::abort_if_empty_cell ( bool  x)
inline

Specifies behavior if an empty cell is detected.

Parameters
[in]xif set, then computation is aborted as soon as an empty cell is detected.

only happens in RegularTriangulation/Laguerre diagram.

Definition at line 155 of file delaunay_2d.h.

◆ add_triangle_to_list()

void GEO::Delaunay2d::add_triangle_to_list ( index_t  t,
index_t first,
index_t last 
)
inlineprotected

Adds a triangle to a linked list.

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

Parameters
[in]tthe index of the triangle
[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 415 of file delaunay_2d.h.

◆ create_first_triangle()

bool GEO::Delaunay2d::create_first_triangle ( index_t iv0,
index_t iv1,
index_t iv2 
)
protected

Finds in the pointset a set of three non-colinear 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
Return values
trueif a set of three non-colinear points was found
falseif all the points are colinear

◆ find_3()

static index_t GEO::Delaunay2d::find_3 ( const signed_index_t T,
signed_index_t  v 
)
inlinestaticprotected

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

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

Definition at line 873 of file delaunay_2d.h.

◆ find_conflict_zone()

void GEO::Delaunay2d::find_conflict_zone ( index_t  v,
index_t  t,
const Sign orient,
index_t t_bndry,
index_t e_bndry,
index_t first,
index_t last 
)
protected

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

Parameters
[in]vthe index of the point to be inserted
[in]tthe index of a triangle that contains p, as returned by locate()
[in]orientan array of three signs indicating the orientation of p with respect to the three edges of t, as returned by locate()
[out]t_bndrya triangle adjacent to the boundary of the conflict zone
[out]e_bndrythe edge along which t_bndry is adjacent to the boundary of the conflict zone
[out]firstthe index of the first triangle in conflict
[out]lastthe index of the last triangle in conflict The other triangles are linked, and can be traversed from first by using triangle_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::Delaunay2d::find_conflict_zone_iterative ( const double *  p,
index_t  t,
index_t t_bndry,
index_t e_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 triangle in the fonflict zone
[out]t_bndrya triangle adjacent to the boundary of the conflict zone
[out]e_bndrythe edge along which t_bndry is adjacent to the boundary of the conflict zone
[out]firstthe index of the first triangle in conflict
[out]lastthe index of the last triangle in conflict
Precondition
The triangle t was already marked as conflict (triangle_is_in_list(t))

◆ find_triangle_adjacent()

index_t GEO::Delaunay2d::find_triangle_adjacent ( index_t  t1,
index_t  t2_in 
) const
inlineprotected

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

Parameters
[in]t1first triangle
[in]t2_insecond triangle
Returns
e such that triangle_adjacent(t1,e)==t2_in
Precondition
t1 and t2_in are adjacent

Definition at line 728 of file delaunay_2d.h.

◆ find_triangle_vertex()

index_t GEO::Delaunay2d::find_triangle_vertex ( index_t  t,
signed_index_t  v 
) const
inlineprotected

Finds the index of the vertex in a triangle.

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

Definition at line 658 of file delaunay_2d.h.

◆ finite_triangle_vertex()

index_t GEO::Delaunay2d::finite_triangle_vertex ( index_t  t,
index_t  lv 
) const
inlineprotected

Gets the index of a vertex of a triangle.

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

Definition at line 673 of file delaunay_2d.h.

◆ has_empty_cells()

bool GEO::Delaunay2d::has_empty_cells ( ) const
inline

Tests whether the Laguerre diagram has empty cells.

If the Laguerre diagram has empty cells and abort_if_empty_cell is set, then computation is stopped, and all the queries on the Laguerre diagram will not work (including the non-empty cells).

Return values
trueif the Laguerre diagram has empty cells.
falseotherwise.

Definition at line 144 of file delaunay_2d.h.

◆ insert()

index_t GEO::Delaunay2d::insert ( index_t  v,
index_t  hint = NO_TRIANGLE 
)
protected

Inserts a point in the triangulation.

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

◆ locate()

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

Finds the triangle that contains a point.

If the point is on an edge or vertex, the function returns one of the triangles incident to that 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 three Signs or nullptr. If non-nullptr, returns the orientation with respect to the three edges of the triangle that contains p.
Returns
the index of a triangle that contains p. If the point is outside the convex hull of the inserted so-far points, then the returned triangle is a virtual one (first vertex is the "vertex at infinity" of index -1) or NO_TRIANGLE if the virtual triangles were previously removed.

◆ locate_inexact()

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

Finds the triangle 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 triangle that (approximately) contains p. If the point is outside the convex hull of the inserted so-far points, then the returned triangle is a virtual one (first vertex is the "vertex at infinity" of index -1) or NO_TRIANGLE if the virtual triangles were previously removed.

◆ mark_triangle()

void GEO::Delaunay2d::mark_triangle ( index_t  t)
inlineprotected

Marks a triangle.

A triangle 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 triangles. A triangle can be in the following states:

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

Definition at line 610 of file delaunay_2d.h.

◆ max_t()

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

Maximum valid index for a triangle.

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

Returns
the maximum valid index for a triangle.

Definition at line 320 of file delaunay_2d.h.

◆ nearest_vertex()

index_t GEO::Delaunay2d::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_triangle() [1/2]

index_t GEO::Delaunay2d::new_triangle ( )
inlineprotected

Creates a new triangle.

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

Returns
the index of the newly created triangle

Definition at line 522 of file delaunay_2d.h.

◆ new_triangle() [2/2]

index_t GEO::Delaunay2d::new_triangle ( signed_index_t  v1,
signed_index_t  v2,
signed_index_t  v3 
)
inlineprotected

Creates a new triangle.

Sets the vertices. Adjacent triangles index are left uninitialized. Uses either a triangle 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
Returns
the index of the newly created triangle

Definition at line 556 of file delaunay_2d.h.

◆ remove_triangle_from_list()

void GEO::Delaunay2d::remove_triangle_from_list ( index_t  t)
inlineprotected

Removes a triangle from the linked list it belongs to.

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

Parameters
[in]tthe index of the triangle

Definition at line 437 of file delaunay_2d.h.

◆ set_tet()

void GEO::Delaunay2d::set_tet ( index_t  t,
signed_index_t  v0,
signed_index_t  v1,
signed_index_t  v2,
index_t  a0,
index_t  a1,
index_t  a2 
)
inlineprotected

Sets the vertices and adjacent triangles of a triangle.

Parameters
[in]tindex of the triangle
[in]v0index of the first vertex
[in]v1index of the second vertex
[in]v2index of the third vertex
[in]a0index of the adjacent triangle opposite to v0
[in]a1index of the adjacent triangle opposite to v1
[in]a2index of the adjacent triangle opposite to v2

Definition at line 761 of file delaunay_2d.h.

◆ set_triangle_adjacent()

void GEO::Delaunay2d::set_triangle_adjacent ( index_t  t1,
index_t  le1,
index_t  t2 
)
inlineprotected

Sets a triangle-to-triangle adjacency.

Parameters
[in]t1index of the first triangle
[in]le1local facet index (0,1,2) in t1
[in]t2index of the triangle adjacent to t1 accros lf1

Definition at line 712 of file delaunay_2d.h.

◆ set_triangle_mark_stamp()

void GEO::Delaunay2d::set_triangle_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 574 of file delaunay_2d.h.

◆ set_triangle_vertex()

void GEO::Delaunay2d::set_triangle_vertex ( index_t  t,
index_t  lv,
signed_index_t  v 
)
inlineprotected

Sets a triangle-to-vertex adjacency.

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

Definition at line 686 of file delaunay_2d.h.

◆ set_vertices()

void GEO::Delaunay2d::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::Delaunay2d::show_list ( index_t  first,
const std::string &  list_name 
) const
protected

For debugging purposes, displays a triangle.

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

◆ show_triangle()

void GEO::Delaunay2d::show_triangle ( index_t  t) const
protected

For debugging purposes, displays a triangle.

Parameters
[in]tindex of the triangle to display.

◆ show_triangle_adjacent()

void GEO::Delaunay2d::show_triangle_adjacent ( index_t  t,
index_t  le 
) const
protected

For debugging purposes, displays a triangle adjacency.

Parameters
[in]tindex of the triangle to display.
[in]lelocal index (0,1,2) of the triangle facet adjacenty to display.

◆ stellate_conflict_zone()

index_t GEO::Delaunay2d::stellate_conflict_zone ( index_t  v,
index_t  t_bndry,
index_t  e_bndry 
)
protected

Creates a star of triangles filling the conflict zone.

For each triangle edge on the border of the conflict zone, a new triangle is created, resting on the edge 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 triangle on the border of the conflict zone.
[in]e_bndryindex of the facet along which t_bndry is incident to the border of the conflict zone
Returns
the index of one the newly created triangles

◆ triangle_adjacent()

signed_index_t GEO::Delaunay2d::triangle_adjacent ( index_t  t,
index_t  le 
) const
inlineprotected

Gets the index of a triangle adjacent to another one.

Parameters
[in]tindex of the triangle
[in]lelocal edge (0,1,2) index in t
Returns
the triangle adjacent to t accros edge le

Definition at line 698 of file delaunay_2d.h.

◆ triangle_edge_vertex()

static index_t GEO::Delaunay2d::triangle_edge_vertex ( index_t  e,
index_t  v 
)
inlinestaticprotected

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

tri edge vertex is such that the triangle formed with:

  • vertex lv
  • triangle_edge_vertex(lv,0)
  • triangle_edge_vertex(lv,1) has the same orientation as the original triangle for any vertex lv.
    Parameters
    [in]elocal facet index, in (0,1,2)
    [in]vlocal vertex index, in (0,1)
    Returns
    the local triangle vertex index of vertex v in edge f

Definition at line 632 of file delaunay_2d.h.

◆ triangle_is_conflict()

bool GEO::Delaunay2d::triangle_is_conflict ( index_t  t,
const double *  p 
) const
inlineprotected

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

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

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

Definition at line 791 of file delaunay_2d.h.

◆ triangle_is_finite()

bool GEO::Delaunay2d::triangle_is_finite ( index_t  t) const
inlineprotected

Tests whether a given triangle is a finite one.

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

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

Definition at line 461 of file delaunay_2d.h.

◆ triangle_is_free()

bool GEO::Delaunay2d::triangle_is_free ( index_t  t) const
inlineprotected

Tests whether a triangle is in the free list.

Deleted triangles are recycled in a free list.

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

Definition at line 511 of file delaunay_2d.h.

◆ triangle_is_in_list()

bool GEO::Delaunay2d::triangle_is_in_list ( index_t  t) const
inlineprotected

Tests whether a triangle belongs to a linked list.

Triangles can be linked, it is used to manage both the free list that recycles deleted triangles, the conflict region and the list of newly created triangles. In addition, a triangle that is not in a list can be marked. The same space is used for marking and chaining triangles in lists. A triangle 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 triangle
    Return values
    trueif triangle t belongs to a linked list
    falseotherwise

Definition at line 384 of file delaunay_2d.h.

◆ triangle_is_marked()

bool GEO::Delaunay2d::triangle_is_marked ( index_t  t) const
inlineprotected

Tests whether a triangle is marked.

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

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

Definition at line 593 of file delaunay_2d.h.

◆ triangle_is_real()

bool GEO::Delaunay2d::triangle_is_real ( index_t  t) const
inlineprotected

Tests whether a triangle is a real one.

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

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

Definition at line 479 of file delaunay_2d.h.

◆ triangle_is_virtual()

bool GEO::Delaunay2d::triangle_is_virtual ( index_t  t) const
inlineprotected

Tests whether a triangle is a virtual one.

Virtual triangles are triangles incident to the vertex at infinity.

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

Definition at line 492 of file delaunay_2d.h.

◆ triangle_next()

index_t GEO::Delaunay2d::triangle_next ( index_t  t) const
inlineprotected

Gets the index of a successor of a triangle.

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

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

Definition at line 399 of file delaunay_2d.h.

◆ triangle_vertex()

signed_index_t GEO::Delaunay2d::triangle_vertex ( index_t  t,
index_t  lv 
) const
inlineprotected

Gets the index of a vertex of a triangle.

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

Definition at line 645 of file delaunay_2d.h.

Member Data Documentation

◆ NO_TRIANGLE

const index_t GEO::Delaunay2d::NO_TRIANGLE = 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 167 of file delaunay_2d.h.

◆ NOT_IN_LIST

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

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

This is the default value. Note that it suffices that NOT_IN_LIST_BIT is set for a triangle to be not in any list. A triangle 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 340 of file delaunay_2d.h.

◆ NOT_IN_LIST_BIT

const index_t GEO::Delaunay2d::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 triangles that were detected as non-conflict when inserting a point. A triangle 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 357 of file delaunay_2d.h.

◆ VERTEX_AT_INFINITY

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

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

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

Definition at line 449 of file delaunay_2d.h.


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