Geogram
Version 1.9.1-rc
A programming library of geometric algorithms
|
Abstract interface for Delaunay triangulation in Nd. More...
#include <geogram/delaunay/delaunay.h>
Classes | |
struct | InvalidDimension |
Invalid dimension exception. More... | |
struct | InvalidInput |
Invalid input exception. More... | |
Public Member Functions | |
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... | |
virtual void | set_vertices (index_t nb_vertices, const double *vertices) |
Sets the vertices of this Delaunay, and recomputes the cells. 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... | |
virtual index_t | nearest_vertex (const double *p) const |
Computes the nearest vertex from a query point. 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... | |
Static Public Member Functions | |
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 Member Functions | |
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... | |
Protected Attributes | |
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 | |
(Note that these are not member functions.) | |
typedef SmartPointer< Delaunay > | Delaunay_var |
Smart pointer that refers to a Delaunay object. | |
typedef Factory1< Delaunay, coord_index_t > | DelaunayFactory |
Delaunay Factory. More... | |
Abstract interface for Delaunay triangulation in Nd.
Delaunay objects are created using method create() which uses the Factory service. New Delaunay triangulations can be implemented and registered to the factory using geo_register_Delaunay_creator().
Definition at line 71 of file delaunay.h.
|
protected |
Creates a new Delaunay triangulation.
This creates a new Delaunay triangulation for the specified dimension
. Specific implementations of the Delaunay triangulation may not support the specified dimension
and will throw a InvalidDimension exception.
[in] | dimension | dimension of the triangulation |
InvalidDimension | This exception is thrown if the specified dimension is not supported by the Delaunay implementation. |
Retrieves a local facet index from two adacent cell global indices.
[in] | c1 | global index of first cell |
[in] | c2 | global index of second cell |
c2
is adjacent to c1
c1
and cell c2
are adjacent Definition at line 431 of file delaunay.h.
|
inline |
Gets an adjacent cell index by cell index and local facet index.
[in] | c | cell index |
[in] | lf | local facet index |
c
accros facet lf
if it exists, or -1 if on border Definition at line 379 of file delaunay.h.
|
inline |
Tests whether a cell is finite.
true | if cell c is finite |
false | otherwise |
Definition at line 399 of file delaunay.h.
bool GEO::Delaunay::cell_is_infinite | ( | index_t | c | ) | const |
Tests whether a cell is infinite.
true | if cell c is infinite |
false | otherwise |
|
inline |
Gets the number of vertices in each cell.
Cell_size = dimension + 1
Definition at line 180 of file delaunay.h.
|
inline |
Gets a pointer to the cell-to-cell adjacency array.
Definition at line 348 of file delaunay.h.
|
inline |
Gets a pointer to the cell-to-vertex incidence array.
Definition at line 340 of file delaunay.h.
|
inline |
Gets a vertex index by cell index and local vertex index.
[in] | c | cell index |
[in] | lv | local vertex index in cell c |
Definition at line 365 of file delaunay.h.
|
inline |
Gets the constraints.
Definition at line 311 of file delaunay.h.
|
static |
Creates a Delaunay triangulation of the specified dimension.
[in] | dim | dimension of the triangulation |
[in] | name | name of the implementation to use:
|
nullptr | if format is not a valid Delaunay algorithm name. |
otherwise,a | pointer to a Delaunay algorithm object. The returned pointer must be stored in an Delaunay_var that does automatic destruction: static Delaunay * create(coord_index_t dim, const std::string &name="default") Creates a Delaunay triangulation of the specified dimension. SmartPointer< Delaunay > Delaunay_var Smart pointer that refers to a Delaunay object. Definition: delaunay.h:792 |
|
inline |
Gets the default number of stored neighbors.
Storage of neighbors is optimized for a default neighborhood size.
Definition at line 602 of file delaunay.h.
|
inline |
Gets the dimension of this Delaunay.
Definition at line 171 of file delaunay.h.
Gets the one-ring neighbors of vertex v.
Depending on store_neighbors_ internal flag, the neighbors are computed or copied from the previously computed list.
[in] | v | vertex index |
[out] | neighbors | indices of the one-ring neighbors of vertex v |
Definition at line 481 of file delaunay.h.
|
protectedvirtual |
Internal implementation for get_neighbors (with vector).
[in] | v | index of the Delaunay vertex |
[in,out] | neighbors | the computed neighbors of vertex v . Its size is used to determine the number of queried neighbors. |
Reimplemented in GEO::Delaunay_NearestNeighbors.
|
inline |
Tests whether mesh refinement is selected.
true | if mesh refinement is selected |
false | otherwise |
Definition at line 287 of file delaunay.h.
|
inline |
Retrieves a local vertex index from cell index and global vertex index.
[in] | c | cell index |
[in] | v | global vertex index |
v
in cell c
c
is incident to vertex v
Definition at line 411 of file delaunay.h.
|
static |
This function needs to be called once before using the Delaunay class.
registers the factories.
|
inline |
Tests whether infinite elements are kept.
true | if infinite elements are kept |
false | otherwise |
Definition at line 551 of file delaunay.h.
|
inline |
Gets the number of cells.
Definition at line 319 of file delaunay.h.
|
inline |
Gets the number of finite cells.
Finite cells have indices 0..nb_finite_cells()-1 and infinite cells have indices nb_finite_cells()..nb_cells()-1
Definition at line 331 of file delaunay.h.
|
inline |
Gets the number of vertices.
Definition at line 242 of file delaunay.h.
|
virtual |
Computes the nearest vertex from a query point.
[in] | p | query point |
Reimplemented in GEO::PeriodicDelaunay3d, GEO::ParallelDelaunay3d, GEO::Delaunay_NearestNeighbors, GEO::Delaunay3d, and GEO::Delaunay2d.
|
inline |
Traverses the list of cells incident to a vertex.
Can only be used if set_stores_cicl(true) was called.
[in] | c | cell index |
[in] | lv | local vertex index |
c
or -1 if c
was the last one in the list Definition at line 465 of file delaunay.h.
Gets the region id associated with a tetrahedron.
Only callable if set_keep_region(true) was called before set_vertices() in constrained mode.
[in] | t | a tetrahedron index. |
t
. void GEO::Delaunay::save_histogram | ( | std::ostream & | out | ) | const |
Saves the histogram of vertex degree (can be visualized with gnuplot).
[out] | out | an ASCII stream where to output the histogram. |
|
protectedvirtual |
Sets the arrays that represent the combinatorics of this Delaunay.
[in] | nb_cells | number of cells |
[in] | cell_to_v | the cell-to-vertex incidence array |
[in] | cell_to_cell | the cell-to-cell adjacency array |
Specifies the bounds of each level to be used when hierarchic ordering is specified from outside.
This function is used by some implementation when set_reorder(false) was called.
[in] | levels | specifies the bounds of each level used by the hierarchical index. First level has indices between levels[0] ... levels[1]. |
Reimplemented in GEO::PeriodicDelaunay3d, and GEO::ParallelDelaunay3d.
|
inlinevirtual |
Defines the constraints.
The triangulation will be constrained to pass through the vertices and triangles of the mesh. This function should be called before set_vertices().
[in] | mesh | the definition of the constraints |
Definition at line 263 of file delaunay.h.
|
inline |
Sets the default number of stored neighbors.
Storage of neighbors is optimized for a default neighborhood size.
[in] | x | default number of stored neighbors |
Definition at line 591 of file delaunay.h.
|
inlineprotected |
Sets the dimension of this Delaunay.
Updates all the parameters related with the dimension. This includes vertex_stride (number of doubles between two consecutive vertices), cell size (number of vertices in a cell), cell_v_stride (number of integers between two consecutive cell vertex arrays), cell_neigh_stride (number of integers between two consecutive cell adjacency arrays).
[in] | dim | the dimension |
Definition at line 728 of file delaunay.h.
|
inline |
Specifies whether all internal regions should be kept.
Only relevant in constrained mode.
[in] | x | if true, all internal regions are kept, else only the outer most region is kept (default). |
Definition at line 619 of file delaunay.h.
|
inline |
Sets whether infinite elements should be kept.
Internally, Delaunay implementation uses an infinite vertex and infinite simplices indicent to it. By default they are discarded at the end of set_vertices().
[in] | x | true if infinite elements should be kept, false otherwise |
Definition at line 563 of file delaunay.h.
Sets the circular incident edge list.
[in] | c1 | index of a cell |
[in] | lv | local index of a vertex of c1 |
[in] | c2 | index of the next cell around c1's vertex lv |
Definition at line 696 of file delaunay.h.
|
inline |
Specifies the desired quality for mesh elements when refinement is enabled (.
Only taken into account after set_refine(true) is called. It is not taken into account by all implementations. This function should be called before set_vertices().
[in] | qual | typically in [1.0, 2.0], specifies the desired quality of mesh elements (1.0 means maximum quality, and generates a higher number of elements). |
Definition at line 302 of file delaunay.h.
|
inline |
Specifies whether the mesh should be refined.
If set, then the mesh elements are improved by inserting additional vertices in the mesh. It is not taken into account by all implementations. This function should be called before set_vertices().
[in] | x | true if the mesh should be refined, false otherwise. |
Definition at line 277 of file delaunay.h.
|
inline |
Specifies whether vertices should be reordered.
Reordering is activated by default. Some special usages of Delaunay3d may require to deactivate it (for instance if vertices are already known to be ordered).
[in] | x | if true, then vertices are reordered using BRIO-Hilbert ordering. This improves speed significantly (enabled by default). |
Definition at line 204 of file delaunay.h.
|
inline |
Specifies whether incident tetrahedra lists should be stored.
[in] | x | if true, incident trahedra lists are stored, else they are not. |
Definition at line 541 of file delaunay.h.
|
inline |
Specifies whether neighbors should be stored.
Vertices neighbors (i.e. Delaunay 1-skeleton) can be stored for faster access (used for instance by RestrictedVoronoiDiagram).
[in] | x | if true neighbors will be stored, else they will not |
Definition at line 518 of file delaunay.h.
|
inline |
Specifies whether thread-safe mode should be used.
[in] | x | if true then thread-safe mode will be used, else it will not. |
Definition at line 580 of file delaunay.h.
|
virtual |
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 in GEO::PeriodicDelaunay3d, GEO::ParallelDelaunay3d, GEO::Delaunay_NearestNeighbors, GEO::Delaunay3d, and GEO::Delaunay2d.
|
virtual |
Stores the neighbors of a vertex.
Used internally for parallel computation of the neighborhoods.
[in] | i | index of the vertex for which the neighbors should be stored. |
Reimplemented in GEO::Delaunay_NearestNeighbors.
|
inline |
Tests whether incident tetrahedra lists are stored.
true | if incident tetrahedra lists are stored. |
false | otherwise. |
Definition at line 531 of file delaunay.h.
|
inline |
Tests whether neighbors are stored.
Vertices neighbors (i.e. Delaunay 1-skeleton) can be stored for faster access (used for instance by RestrictedVoronoiDiagram).
true | if neighbors are stored. |
false | otherwise. |
Definition at line 507 of file delaunay.h.
|
virtual |
Tests whether constraints are supported by this Delaunay.
true | if constraints are supported |
false | otherwise |
|
inline |
Tests whether thread-safe mode is active.
Definition at line 571 of file delaunay.h.
|
protectedvirtual |
Updates the circular incident cell lists.
Used by next_around_vertex().
Reimplemented in GEO::PeriodicDelaunay3d.
|
inline |
Gets an incident cell index by a vertex index.
Can only be used if set_stores_cicl(true) was called.
[in] | v | a vertex index |
v
Definition at line 449 of file delaunay.h.
|
inline |
Gets a pointer to a vertex by its global index.
[in] | i | global index of the vertex |
i
Definition at line 233 of file delaunay.h.
|
inline |
Gets a pointer to the array of vertices.
Definition at line 224 of file delaunay.h.
|
related |
This Factory is used to create Delaunay objects. It can also be used to register new Delaunay implementations.
Definition at line 804 of file delaunay.h.