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

Abstract interface for Delaunay triangulation in Nd. More...

#include <geogram/delaunay/delaunay.h>

Inheritance diagram for GEO::Delaunay:
GEO::Counted GEO::Delaunay2d GEO::Delaunay3d GEO::Delaunay_NearestNeighbors GEO::ParallelDelaunay3d GEO::PeriodicDelaunay3d GEO::RegularWeightedDelaunay2d GEO::RegularWeightedDelaunay3d

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

Related Functions

(Note that these are not member functions.)

typedef SmartPointer< DelaunayDelaunay_var
 Smart pointer that refers to a Delaunay object.
 
typedef Factory1< Delaunay, coord_index_tDelaunayFactory
 Delaunay Factory. More...
 

Detailed Description

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().

See also
DelaunayFactory
geo_register_Delaunay_creator

Definition at line 71 of file delaunay.h.

Constructor & Destructor Documentation

◆ Delaunay()

GEO::Delaunay::Delaunay ( coord_index_t  dimension)
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.

Parameters
[in]dimensiondimension of the triangulation
Exceptions
InvalidDimensionThis exception is thrown if the specified dimension is not supported by the Delaunay implementation.
Note
This function is never called directly, use create()

Member Function Documentation

◆ adjacent_index()

index_t GEO::Delaunay::adjacent_index ( index_t  c1,
index_t  c2 
) const
inline

Retrieves a local facet index from two adacent cell global indices.

Parameters
[in]c1global index of first cell
[in]c2global index of second cell
Returns
the local index of the face accros which c2 is adjacent to c1
Precondition
cell c1 and cell c2 are adjacent

Definition at line 431 of file delaunay.h.

◆ cell_adjacent()

signed_index_t GEO::Delaunay::cell_adjacent ( index_t  c,
index_t  lf 
) const
inline

Gets an adjacent cell index by cell index and local facet index.

Parameters
[in]ccell index
[in]lflocal facet index
Returns
the index of the cell adjacent to c accros facet lf if it exists, or -1 if on border

Definition at line 379 of file delaunay.h.

◆ cell_is_finite()

bool GEO::Delaunay::cell_is_finite ( index_t  c) const
inline

Tests whether a cell is finite.

Return values
trueif cell c is finite
falseotherwise
See also
keeps_infinite(), set_keeps_infinite()

Definition at line 399 of file delaunay.h.

◆ cell_is_infinite()

bool GEO::Delaunay::cell_is_infinite ( index_t  c) const

Tests whether a cell is infinite.

Return values
trueif cell c is infinite
falseotherwise
See also
keeps_infinite(), set_keeps_infinite()

◆ cell_size()

index_t GEO::Delaunay::cell_size ( ) const
inline

Gets the number of vertices in each cell.

Cell_size = dimension + 1

Returns
the number of vertices in each cell

Definition at line 180 of file delaunay.h.

◆ cell_to_cell()

const signed_index_t* GEO::Delaunay::cell_to_cell ( ) const
inline

Gets a pointer to the cell-to-cell adjacency array.

Returns
a const pointer to the cell-to-cell adjacency array

Definition at line 348 of file delaunay.h.

◆ cell_to_v()

const signed_index_t* GEO::Delaunay::cell_to_v ( ) const
inline

Gets a pointer to the cell-to-vertex incidence array.

Returns
a const pointer to the cell-to-vertex incidence array

Definition at line 340 of file delaunay.h.

◆ cell_vertex()

signed_index_t GEO::Delaunay::cell_vertex ( index_t  c,
index_t  lv 
) const
inline

Gets a vertex index by cell index and local vertex index.

Parameters
[in]ccell index
[in]lvlocal vertex index in cell c
Returns
the index of the lv-th vertex of cell c.

Definition at line 365 of file delaunay.h.

◆ constraints()

const Mesh* GEO::Delaunay::constraints ( ) const
inline

Gets the constraints.

Returns
the constraints or nullptr if no constraints were definied.

Definition at line 311 of file delaunay.h.

◆ create()

static Delaunay* GEO::Delaunay::create ( coord_index_t  dim,
const std::string &  name = "default" 
)
static

Creates a Delaunay triangulation of the specified dimension.

Parameters
[in]dimdimension of the triangulation
[in]namename of the implementation to use:
  • "tetgen" - Delaunay with the Tetgen library (dimension 3 only)
  • "BDEL" - Delaunay in 3D (dimension 3 only)
  • "BPOW" - Weighted regular 3D triangulation (dimension 4 only)
  • "NN" - Delaunay with NearestNeighborSearch (any dimension)
  • "default" - uses the command line argument "algo:delaunay"
Return values
nullptrif format is not a valid Delaunay algorithm name.
otherwise,apointer to a Delaunay algorithm object. The returned pointer must be stored in an Delaunay_var that does automatic destruction:
Delaunay_var handler = Delaunay::create(3, "default");
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

◆ default_nb_neighbors()

index_t GEO::Delaunay::default_nb_neighbors ( ) const
inline

Gets the default number of stored neighbors.

Storage of neighbors is optimized for a default neighborhood size.

See also
store_neighbors()
Returns
The default number of stored neighbors.

Definition at line 602 of file delaunay.h.

◆ dimension()

coord_index_t GEO::Delaunay::dimension ( ) const
inline

Gets the dimension of this Delaunay.

Returns
the dimension of this Delauna

Definition at line 171 of file delaunay.h.

◆ get_neighbors()

void GEO::Delaunay::get_neighbors ( index_t  v,
vector< index_t > &  neighbors 
) const
inline

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.

Parameters
[in]vvertex index
[out]neighborsindices of the one-ring neighbors of vertex v
See also
stores_neighbors(), set_stores_neighbors()

Definition at line 481 of file delaunay.h.

◆ get_neighbors_internal()

virtual void GEO::Delaunay::get_neighbors_internal ( index_t  v,
vector< index_t > &  neighbors 
) const
protectedvirtual

Internal implementation for get_neighbors (with vector).

Parameters
[in]vindex of the Delaunay vertex
[in,out]neighborsthe computed neighbors of vertex v. Its size is used to determine the number of queried neighbors.

Reimplemented in GEO::Delaunay_NearestNeighbors.

◆ get_refine()

bool GEO::Delaunay::get_refine ( ) const
inline

Tests whether mesh refinement is selected.

Return values
trueif mesh refinement is selected
falseotherwise
See also
set_refine()

Definition at line 287 of file delaunay.h.

◆ index()

index_t GEO::Delaunay::index ( index_t  c,
signed_index_t  v 
) const
inline

Retrieves a local vertex index from cell index and global vertex index.

Parameters
[in]ccell index
[in]vglobal vertex index
Returns
the local index of vertex v in cell c
Precondition
cell c is incident to vertex v

Definition at line 411 of file delaunay.h.

◆ initialize()

static void GEO::Delaunay::initialize ( )
static

This function needs to be called once before using the Delaunay class.

registers the factories.

◆ keeps_infinite()

bool GEO::Delaunay::keeps_infinite ( ) const
inline

Tests whether infinite elements are kept.

Return values
trueif infinite elements are kept
falseotherwise

Definition at line 551 of file delaunay.h.

◆ nb_cells()

index_t GEO::Delaunay::nb_cells ( ) const
inline

Gets the number of cells.

Returns
the number of cells in this Delaunay

Definition at line 319 of file delaunay.h.

◆ nb_finite_cells()

index_t GEO::Delaunay::nb_finite_cells ( ) const
inline

Gets the number of finite cells.

Precondition
this function can only be called if keep_finite is set

Finite cells have indices 0..nb_finite_cells()-1 and infinite cells have indices nb_finite_cells()..nb_cells()-1

See also
set_keeps_infinite(), keeps_infinite()

Definition at line 331 of file delaunay.h.

◆ nb_vertices()

index_t GEO::Delaunay::nb_vertices ( ) const
inline

Gets the number of vertices.

Returns
the number of vertices in this Delaunay

Definition at line 242 of file delaunay.h.

◆ nearest_vertex()

virtual index_t GEO::Delaunay::nearest_vertex ( const double *  p) const
virtual

Computes the nearest vertex from a query point.

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

Reimplemented in GEO::PeriodicDelaunay3d, GEO::ParallelDelaunay3d, GEO::Delaunay_NearestNeighbors, GEO::Delaunay3d, and GEO::Delaunay2d.

◆ next_around_vertex()

signed_index_t GEO::Delaunay::next_around_vertex ( index_t  c,
index_t  lv 
) const
inline

Traverses the list of cells incident to a vertex.

Can only be used if set_stores_cicl(true) was called.

Parameters
[in]ccell index
[in]lvlocal vertex index
Returns
the index of the next cell around vertex c or -1 if c was the last one in the list
See also
stores_cicl(), set_store_cicl()

Definition at line 465 of file delaunay.h.

◆ region()

virtual index_t GEO::Delaunay::region ( index_t  t) const
virtual

Gets the region id associated with a tetrahedron.

Only callable if set_keep_region(true) was called before set_vertices() in constrained mode.

Parameters
[in]ta tetrahedron index.
Returns
the region associated with t.

◆ save_histogram()

void GEO::Delaunay::save_histogram ( std::ostream &  out) const

Saves the histogram of vertex degree (can be visualized with gnuplot).

Parameters
[out]outan ASCII stream where to output the histogram.

◆ set_arrays()

virtual void GEO::Delaunay::set_arrays ( index_t  nb_cells,
const signed_index_t cell_to_v,
const signed_index_t cell_to_cell 
)
protectedvirtual

Sets the arrays that represent the combinatorics of this Delaunay.

Parameters
[in]nb_cellsnumber of cells
[in]cell_to_vthe cell-to-vertex incidence array
[in]cell_to_cellthe cell-to-cell adjacency array

◆ set_BRIO_levels()

virtual void GEO::Delaunay::set_BRIO_levels ( const vector< index_t > &  levels)
virtual

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.

Parameters
[in]levelsspecifies 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.

◆ set_constraints()

virtual void GEO::Delaunay::set_constraints ( const Mesh mesh)
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().

Parameters
[in]meshthe definition of the constraints
Precondition
constraints_supported()

Definition at line 263 of file delaunay.h.

◆ set_default_nb_neighbors()

void GEO::Delaunay::set_default_nb_neighbors ( index_t  x)
inline

Sets the default number of stored neighbors.

Storage of neighbors is optimized for a default neighborhood size.

See also
store_neighbors()
Parameters
[in]xdefault number of stored neighbors

Definition at line 591 of file delaunay.h.

◆ set_dimension()

void GEO::Delaunay::set_dimension ( coord_index_t  dim)
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).

Parameters
[in]dimthe dimension

Definition at line 728 of file delaunay.h.

◆ set_keep_regions()

void GEO::Delaunay::set_keep_regions ( bool  x)
inline

Specifies whether all internal regions should be kept.

Only relevant in constrained mode.

Parameters
[in]xif true, all internal regions are kept, else only the outer most region is kept (default).

Definition at line 619 of file delaunay.h.

◆ set_keeps_infinite()

void GEO::Delaunay::set_keeps_infinite ( bool  x)
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().

Parameters
[in]xtrue if infinite elements should be kept, false otherwise

Definition at line 563 of file delaunay.h.

◆ set_next_around_vertex()

void GEO::Delaunay::set_next_around_vertex ( index_t  c1,
index_t  lv,
index_t  c2 
)
inlineprotected

Sets the circular incident edge list.

Parameters
[in]c1index of a cell
[in]lvlocal index of a vertex of c1
[in]c2index of the next cell around c1's vertex lv

Definition at line 696 of file delaunay.h.

◆ set_quality()

void GEO::Delaunay::set_quality ( double  qual)
inline

Specifies the desired quality for mesh elements when refinement is enabled (.

See also
set_refine).

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().

Parameters
[in]qualtypically 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.

◆ set_refine()

void GEO::Delaunay::set_refine ( bool  x)
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().

Parameters
[in]xtrue if the mesh should be refined, false otherwise.

Definition at line 277 of file delaunay.h.

◆ set_reorder()

void GEO::Delaunay::set_reorder ( bool  x)
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).

Parameters
[in]xif true, then vertices are reordered using BRIO-Hilbert ordering. This improves speed significantly (enabled by default).

Definition at line 204 of file delaunay.h.

◆ set_stores_cicl()

void GEO::Delaunay::set_stores_cicl ( bool  x)
inline

Specifies whether incident tetrahedra lists should be stored.

Parameters
[in]xif true, incident trahedra lists are stored, else they are not.

Definition at line 541 of file delaunay.h.

◆ set_stores_neighbors()

void GEO::Delaunay::set_stores_neighbors ( bool  x)
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).

Parameters
[in]xif true neighbors will be stored, else they will not

Definition at line 518 of file delaunay.h.

◆ set_thread_safe()

void GEO::Delaunay::set_thread_safe ( bool  x)
inline

Specifies whether thread-safe mode should be used.

Parameters
[in]xif true then thread-safe mode will be used, else it will not.

Definition at line 580 of file delaunay.h.

◆ set_vertices()

virtual void GEO::Delaunay::set_vertices ( index_t  nb_vertices,
const double *  vertices 
)
virtual

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 in GEO::PeriodicDelaunay3d, GEO::ParallelDelaunay3d, GEO::Delaunay_NearestNeighbors, GEO::Delaunay3d, and GEO::Delaunay2d.

◆ store_neighbors_CB()

virtual void GEO::Delaunay::store_neighbors_CB ( index_t  i)
virtual

Stores the neighbors of a vertex.

Used internally for parallel computation of the neighborhoods.

Parameters
[in]iindex of the vertex for which the neighbors should be stored.

Reimplemented in GEO::Delaunay_NearestNeighbors.

◆ stores_cicl()

bool GEO::Delaunay::stores_cicl ( ) const
inline

Tests whether incident tetrahedra lists are stored.

Return values
trueif incident tetrahedra lists are stored.
falseotherwise.

Definition at line 531 of file delaunay.h.

◆ stores_neighbors()

bool GEO::Delaunay::stores_neighbors ( ) const
inline

Tests whether neighbors are stored.

Vertices neighbors (i.e. Delaunay 1-skeleton) can be stored for faster access (used for instance by RestrictedVoronoiDiagram).

Return values
trueif neighbors are stored.
falseotherwise.

Definition at line 507 of file delaunay.h.

◆ supports_constraints()

virtual bool GEO::Delaunay::supports_constraints ( ) const
virtual

Tests whether constraints are supported by this Delaunay.

Return values
trueif constraints are supported
falseotherwise

◆ thread_safe()

bool GEO::Delaunay::thread_safe ( ) const
inline

Tests whether thread-safe mode is active.

Returns
true if thread-safe mode is active, false otherwise.

Definition at line 571 of file delaunay.h.

◆ update_cicl()

virtual void GEO::Delaunay::update_cicl ( )
protectedvirtual

Updates the circular incident cell lists.

Used by next_around_vertex().

Reimplemented in GEO::PeriodicDelaunay3d.

◆ vertex_cell()

signed_index_t GEO::Delaunay::vertex_cell ( index_t  v) const
inline

Gets an incident cell index by a vertex index.

Can only be used if set_stores_cicl(true) was called.

Parameters
[in]va vertex index
Returns
the index of a cell incident to vertex v
See also
stores_cicl(), set_store_cicl()

Definition at line 449 of file delaunay.h.

◆ vertex_ptr()

const double* GEO::Delaunay::vertex_ptr ( index_t  i) const
inline

Gets a pointer to a vertex by its global index.

Parameters
[in]iglobal index of the vertex
Returns
a pointer to vertex i

Definition at line 233 of file delaunay.h.

◆ vertices_ptr()

const double* GEO::Delaunay::vertices_ptr ( ) const
inline

Gets a pointer to the array of vertices.

Returns
A const pointer to the array of vertices.

Definition at line 224 of file delaunay.h.

Friends And Related Function Documentation

◆ DelaunayFactory

Delaunay Factory.

This Factory is used to create Delaunay objects. It can also be used to register new Delaunay implementations.

See also
geo_register_Delaunay_creator
Factory

Definition at line 804 of file delaunay.h.


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