Geogram
Version 1.9.1-rc
A programming library of geometric algorithms
|
Implementation of Delaunay in 3d. More...
#include <geogram/delaunay/delaunay_3d.h>
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 Mesh * | constraints () 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_t * | cell_to_v () const |
Gets a pointer to the cell-to-vertex incidence array. More... | |
const signed_index_t * | cell_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 Delaunay * | create (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_t * | cell_to_v_ |
const signed_index_t * | cell_to_cell_ |
vector< signed_index_t > | v_to_cell_ |
vector< signed_index_t > | cicl_ |
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 Mesh * | constraints_ |
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_ |
Related Functions inherited from GEO::Delaunay | |
typedef SmartPointer< Delaunay > | Delaunay_var |
Smart pointer that refers to a Delaunay object. | |
typedef Factory1< Delaunay, coord_index_t > | DelaunayFactory |
Delaunay Factory. More... | |
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.
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:
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.
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.
Definition at line 107 of file delaunay_3d.h.
GEO::Delaunay3d::Delaunay3d | ( | coord_index_t | dimension = 3 | ) |
Constructs a new Delaunay3d.
[in] | dimension | dimension of the triangulation (3 or 4). If dimension = 4, this creates a regular triangulation (dual of a power diagram). In this case:
|
|
inlineprotected |
Adds a tetrahedron to a linked list.
Tetrahedra can be linked, it is used to manage the free list that recycles deleted tetrahedra.
[in] | t | the index of the tetrahedron |
[in,out] | first | first item of the list or END_OF_LIST if the list is empty |
[in,out] | last | last item of the list or END_OF_LIST if the list is empty |
Definition at line 484 of file delaunay_3d.h.
|
protected |
Finds in the pointset a set of four non-coplanar points.
This function is used to initiate the incremental Delaunay construction.
[out] | iv0 | index of the first vertex |
[out] | iv1 | index of the second vertex |
[out] | iv2 | index of the third vertex |
[out] | iv3 | index of the fourth vertex |
true | if a set of four non-coplanar points was found |
false | if all the points are coplanar |
|
inlinestaticprotected |
Finds the index of an integer in an array of four integers.
[in] | T | a const pointer to an array of four integers |
[in] | v | the integer to retrieve in T |
v
in T
T
are different and one of them is equal to v
. Definition at line 1035 of file delaunay_3d.h.
|
protected |
Determines the list of tetrahedra in conflict with a given point.
[in] | v | the index of the point to be inserted |
[in] | t | the index of a tetrahedron that contains p , as returned by locate() |
[in] | orient | an array of four signs indicating the orientation of p with respect to the four faces of t , as returned by locate() |
[out] | t_bndry | a tetrahedron adjacent to the boundary of the conflict zone |
[out] | f_bndry | the facet along which t_bndry is adjacent to the boundary of the conflict zone |
[out] | first | the index of the first tetrahedron in conflict |
[out] | last | the 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:
|
|
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.
[in] | p | the point to be inserted |
[in] | t | index of a tetrahedron in the fonflict zone |
[out] | t_bndry | a tetrahedron adjacent to the boundary of the conflict zone |
[out] | f_bndry | the facet along which t_bndry is adjacent to the boundary of the conflict zone |
[out] | first | the index of the first tetrahedron in conflict |
[out] | last | the index of the last tetrahedron in conflict |
t
was alredy marked as conflict (tet_is_in_list(t)) Finds the index of the facet accros which t1 is adjacent to t2_in.
[in] | t1 | first tetrahedron |
[in] | t2_in | second tetrahedron |
t1
and t2_in
are adjacent Definition at line 801 of file delaunay_3d.h.
|
inlineprotected |
Finds the index of the vertex in a tetrahedron.
[in] | t | the tetrahedron |
[in] | v | the vertex |
t
is incident to v
Definition at line 732 of file delaunay_3d.h.
Gets the index of a vertex of a tetrahedron.
[in] | t | index of the tetrahedron |
[in] | lv | local vertex (0,1,2 or 3) index in t |
lvth
vertex of tetrahedron t
lv
of tetrahedron t
is not at infinity Definition at line 747 of file delaunay_3d.h.
|
inlineprotected |
Gets the local facet index incident to an oriented halfedge.
[in] | t | index of the tetrahedron |
[in] | v1 | global index of the first extremity |
[in] | v2 | global index of the second extremity |
v1
, v2
. Definition at line 865 of file delaunay_3d.h.
|
inlineprotected |
Gets the local facet indices incident to an oriented halfedge.
[in] | t | index of the tetrahedron |
[in] | v1 | global index of the first extremity |
[in] | v2 | global index of the second extremity |
[out] | f12 | the local index of the facet indicent to the halfedge [v1,v2] |
[out] | f21 | the local index of the facet indicent to the halfedge [v2,v1] |
Definition at line 890 of file delaunay_3d.h.
|
inlineprotected |
Finds the neighbor of a tetrahedron on the border of the conflict zone.
This function is used by stellate_conflict_zone_iterative()
[in] | t1 | a tetrahedron on the border of the conflict zone |
[in] | t1fborder | the local facet index of t1 along which it is on the border of the conflict zone |
[in] | t1ft2 | the local facet index of t1 that will be traversed |
[out] | t2 | a tetrahedron on the border of the conflict zone, with an edge common to facets t1fborder and t1ft2 of tetrahedron t1 |
[out] | t2fborder | the local facet index of t2 along which it is on the border of the conflict zone |
[out] | t2ft1 | the local index of the facet of t2 that has a common edge with facets t1fborder and t1ft2 of tetrahedron t1 |
true | if t2 is a newly created tetrahedron |
false | if t2 is an old tetrahedron in conflict |
Definition at line 325 of file delaunay_3d.h.
|
protected |
Inserts a point in the triangulation.
[in] | v | the index of the point to be inserted |
[in] | hint | the index of a tetrahedron as near as possible to v , or -1 if unspecified |
v
|
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.
[in] | p | a pointer to the coordinates of the point |
[in] | thread_safe | if true, a global spinlock is used to protect the calls to random(), this is necessary if multiple threads use locate() simultaneously |
[out] | orient | a 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 . |
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.
|
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".
[in] | p | a pointer to the coordinates of the point |
[in] | max_iter | maximum number of traversed tets |
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.
|
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] | t | index of the tetrahedron to be marked |
Definition at line 683 of file delaunay_3d.h.
|
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.
Definition at line 389 of file delaunay_3d.h.
|
overridevirtual |
Computes the nearest vertex from a query point.
[in] | p | query point |
Reimplemented from GEO::Delaunay.
|
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.
Definition at line 592 of file delaunay_3d.h.
|
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.
[in] | v1 | index of the first vertex |
[in] | v2 | index of the second vertex |
[in] | v3 | index of the third vertex |
[in] | v4 | index of the fourth vertex |
Definition at line 628 of file delaunay_3d.h.
|
inlineprotected |
Gets the next tetrahedron around an oriented edge of a tetrahedron.
[in,out] | t | the tetrahedron |
[in] | v1 | global index of the first extremity of the edge |
[in] | v2 | global index of the second extremity of the edge |
t
around the oriented edge (v1
v2
). Definition at line 929 of file delaunay_3d.h.
|
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.
[in] | t | the index of the tetrahedron |
Definition at line 506 of file delaunay_3d.h.
|
inlineprotected |
Sets the vertices and adjacent tetrahedra of a tetrahedron.
[in] | t | index of the tetrahedron |
[in] | v0 | index of the first vertex |
[in] | v1 | index of the second vertex |
[in] | v2 | index of the third vertex |
[in] | v3 | index of the fourth vertex |
[in] | a0 | index of the adjacent tetrahedron opposite to v0 |
[in] | a1 | index of the adjacent tetrahedron opposite to v1 |
[in] | a2 | index of the adjacent tetrahedron opposite to v2 |
[in] | a3 | index of the adjacent tetrahedron opposite to v3 |
Definition at line 837 of file delaunay_3d.h.
Sets a tetrahedron-to-tetrahedron adjacency.
[in] | t1 | index of the first tetrahedron |
[in] | lf1 | local facet index (0,1,2 or 3) in t1 |
[in] | t2 | index of the tetrahedron adjacent to t1 accros lf1 |
Definition at line 786 of file delaunay_3d.h.
|
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.
[in] | stamp | the unique stamp for marking tets |
Definition at line 647 of file delaunay_3d.h.
|
inlineprotected |
Sets a tetrahedron-to-vertex adjacency.
[in] | t | index of the tetrahedron |
[in] | lv | local vertex index (0,1,2 or 3) in t |
[in] | v | global index of the vertex |
Definition at line 760 of file delaunay_3d.h.
|
overridevirtual |
Sets the vertices of this Delaunay, and recomputes the cells.
[in] | nb_vertices | number of vertices |
[in] | vertices | a pointer to the coordinates of the vertices, as a contiguous array of doubles |
Reimplemented from GEO::Delaunay.
|
protected |
For debugging purposes, displays a tetrahedron.
[in] | first | index of the first tetrahedron in the list |
[in] | list_name | name of the list, will be displayed as well |
|
protected |
For debugging purposes, displays a tetrahedron.
[in] | t | index of the tetrahedron to display. |
For debugging purposes, displays a tetrahedron adjacency.
[in] | t | index of the tetrahedron to display. |
[in] | lf | local index (0,1,2 or 3) of the tetrahedron facet adjacenty to display. |
Creates a star of tetrahedra filling the conflict zone.
[in] | v | the 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.
|
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.
[in] | v | the index of the point to be inserted |
[in] | t_bndry | index of a tetrahedron on the border of the conflict zone. |
[in] | f_bndry | index of the facet along which t_bndry is incident to the border of the conflict zone |
[in] | prev_f | the facet of t_bndry connected to the tetrahedron that t_bndry was reached from, or index_t(-1) if it is the first tetrahedron. |
|
inlineprotected |
Gets the index of a tetrahedron adjacent to another one.
[in] | t | index of the tetrahedron |
[in] | lf | local facet (0,1,2 or 3) index in t |
t
accorss facet lf
Definition at line 772 of file delaunay_3d.h.
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:
[in] | f | local facet index, in (0,1,2,3) |
[in] | v | local vertex index, in (0,1,2) |
v
in facet f
Definition at line 706 of file delaunay_3d.h.
|
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.
[in] | t | the index of the tetrahedron |
[in] | p | a pointer to the coordinates of the point |
true | if point p is in conflict with tetrahedron t |
false | otherwise |
Definition at line 952 of file delaunay_3d.h.
|
inlineprotected |
Tests whether a given tetrahedron is a finite one.
Infinite tetrahedra are the ones that are incident to the infinite vertex (index -1)
[in] | t | the index of the tetrahedron |
true | if t is finite |
false | otherwise |
Definition at line 530 of file delaunay_3d.h.
|
inlineprotected |
Tests whether a tetrahedron is in the free list.
Deleted tetrahedra are recycled in a free list.
[in] | t | index of the tetrahedron |
true | if tetrahedron t is in the free list |
false | otherwise |
Definition at line 581 of file delaunay_3d.h.
|
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] | t | the index of the tetrahedron |
true | if tetrahedron t belongs to a linked list |
false | otherwise |
Definition at line 453 of file delaunay_3d.h.
|
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] | t | index of the tetrahedron |
true | if tetrahedron t is marked |
false | otherwise |
Definition at line 666 of file delaunay_3d.h.
|
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)
[in] | t | index of the tetrahedron |
true | if tetrahedron t is a real one |
false | otherwise |
Definition at line 549 of file delaunay_3d.h.
|
inlineprotected |
Tests whether a tetrahedron is a virtual one.
Virtual tetrahedra are tetrahedra incident to the vertex at infinity.
[in] | t | index of the tetrahedron |
true | if tetrahedron t is virtual |
false | otherwise |
Definition at line 562 of file delaunay_3d.h.
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.
[in] | t | the index of the tetrahedron |
END_OF_LIST | if the end of the list is reached |
the | index of the successor of tetrahedron \t otherwise |
Definition at line 468 of file delaunay_3d.h.
|
inlineprotected |
Gets the index of a vertex of a tetrahedron.
[in] | t | index of the tetrahedron |
[in] | lv | local vertex (0,1,2 or 3) index in t |
lvth
vertex of tetrahedron t
or -1 if the vertex is at infinity Definition at line 719 of file delaunay_3d.h.
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.
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:
Definition at line 409 of file delaunay_3d.h.
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:
Definition at line 426 of file delaunay_3d.h.
|
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.