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

Multithreaded implementation of Delaunay in 3d with optional periodic boundary conditions. More...

#include <geogram/delaunay/periodic_delaunay_3d.h>

Inheritance diagram for GEO::PeriodicDelaunay3d:
GEO::Delaunay GEO::Periodic GEO::Counted

Classes

struct  IncidentTetrahedra
 Gathers some structures used by some algorithms, makes multithreading more efficient by avoiding dynamic reallocations. More...
 

Public Member Functions

 PeriodicDelaunay3d (bool periodic, double period=1.0)
 Constructs a new PeriodicDelaunay3d. More...
 
 PeriodicDelaunay3d (const vec3 &period)
 Constructs a new PeriodicDelaunay3d. More...
 
void set_vertices (index_t nb_vertices, const double *vertices) override
 Sets the vertices of this Delaunay, and recomputes the cells. More...
 
void set_weights (const double *weights)
 Sets the weights. More...
 
void compute ()
 Computes the Delaunay triangulation.
 
void use_exact_predicates_for_convex_cell (bool x)
 Use exact predicates in convex cell computations. More...
 
vec3 vertex (index_t v) const
 Gets a vertex by index. More...
 
double weight (index_t v) const
 Gets a weight by index. More...
 
index_t nearest_vertex (const double *p) const override
 Computes the nearest vertex from a query point. More...
 
void set_BRIO_levels (const vector< index_t > &levels) override
 Specifies the bounds of each level to be used when hierarchic ordering is specified from outside. More...
 
void get_incident_tets (index_t v, IncidentTetrahedra &W) const
 computes the set of tetrahedra that are incident to a vertex. More...
 
void copy_Laguerre_cell_from_Delaunay (GEO::index_t i, ConvexCell &C, IncidentTetrahedra &W) const
 Copies a Laguerre cell from the triangulation. More...
 
void copy_Laguerre_cell_from_Delaunay (GEO::index_t i, ConvexCell &C) const
 Copies a Laguerre cell from the triangulation. More...
 
bool has_empty_cells () const
 Tests whether the Laguerre diagram has empty cells. More...
 
void save_cells (const std::string &basename, bool clipped)
 Saves the cells in an Alias-Wavefront file. 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...
 
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...
 
- Public Member Functions inherited from GEO::Periodic
index_t periodic_vertex_instance (index_t pv) const
 Gets the instance from a periodic vertex. More...
 
index_t periodic_vertex_real (index_t pv) const
 Gets the real vertex from a periodic vertex. More...
 
signed_index_t periodic_vertex_real (signed_index_t pv) const
 Gets the real vertex from a periodic vertex. More...
 
index_t make_periodic_vertex (index_t real, index_t instance) const
 Makes a periodic vertex from a real vertex and instance. More...
 
void periodic_vertex_get_T (index_t pv, int &Tx, int &Ty, int &Tz) const
 Gets the translation from a periodic vertex. More...
 
void periodic_vertex_set_T (index_t &pv, int Tx, int Ty, int Tz) const
 Sets the translation in a periodic vertex. More...
 
std::string periodic_vertex_to_string (index_t v) const
 
std::string binary_to_string (Numeric::uint32 m) const
 

Protected Member Functions

GEO::index_t copy_Laguerre_cell_facet_from_Delaunay (GEO::index_t i, const GEO::vec3 &Pi, double wi, double Pi_len2, GEO::index_t t, ConvexCell &C, IncidentTetrahedra &W) const
 Copies a Laguerre cell facet from the triangulation. More...
 
index_t compress (bool shrink=true)
 Removes unused tetrahedra. More...
 
void update_v_to_cell () override
 Stores for each vertex v a cell incident to v. More...
 
void update_cicl () override
 Updates the circular incident cell lists. More...
 
void handle_periodic_boundaries ()
 Duplicates the points with Voronoi cells that cross the boundary.
 
index_t get_periodic_vertex_instances_to_create (index_t v, ConvexCell &C, bool use_instance[27], bool &cell_is_on_boundary, bool &cell_is_outside_cube, IncidentTetrahedra &W)
 Computes the periodic vertex instances that should be generated. More...
 
void insert_vertices (index_t b, index_t e)
 Insert vertices from reorder_[b] to reorder_[e-1]. More...
 
void check_volume ()
 Checks the volume of Laguerre cells.
 
PeriodicDelaunay3dThread * thread (index_t t)
 Gets a thread by index. More...
 
index_t nb_threads () const
 Gets the number of threads. More...
 
- 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_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...
 

Friends

class PeriodicDelaunay3dThread
 

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...
 
- Static Public Member Functions inherited from GEO::Periodic
static index_t T_to_instance (int Tx, int Ty, int Tz)
 Gets the instance from a translation. More...
 
- Public Attributes inherited from GEO::Periodic
index_t nb_vertices_non_periodic_
 Number of real vertices.
 
- Static Public Attributes inherited from GEO::Periodic
static int translation [27][3]
 Gives for each instance the integer translation coordinates in {-1,0,1}. More...
 
static int reorder_instances [27]
 Used to back-map an integer translation to an instance. More...
 
static bool instance_is_positive [27]
 Tests whether all the coordinates of the translation vector associated with an instance are 0 or 1.
 
- 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

Multithreaded implementation of Delaunay in 3d with optional periodic boundary conditions.

Periodicity is taken into account by detecting the vertices with Voronoi cells that straddle the boundary and duplicating them as need be (with an additional propagation to have the correct neighborhoods).

See also
Delaunay3d, ParallelDelaunay3d

Definition at line 68 of file periodic_delaunay_3d.h.

Constructor & Destructor Documentation

◆ PeriodicDelaunay3d() [1/2]

GEO::PeriodicDelaunay3d::PeriodicDelaunay3d ( bool  periodic,
double  period = 1.0 
)

Constructs a new PeriodicDelaunay3d.

Parameters
[in]periodicif true, constructs a periodic triangulation.
[in]periodthe edge length of the periodic domain.

◆ PeriodicDelaunay3d() [2/2]

GEO::PeriodicDelaunay3d::PeriodicDelaunay3d ( const vec3 period)

Constructs a new PeriodicDelaunay3d.

Parameters
[in]periodthe edge lengths along x,y,z in the periodic domain

Member Function Documentation

◆ compress()

index_t GEO::PeriodicDelaunay3d::compress ( bool  shrink = true)
protected

Removes unused tetrahedra.

Returns
the final number of tetrahedra.
Parameters
[in]shrinkif true, then array space is shrunk to fit the new number of tetrahedra.

◆ copy_Laguerre_cell_facet_from_Delaunay()

GEO::index_t GEO::PeriodicDelaunay3d::copy_Laguerre_cell_facet_from_Delaunay ( GEO::index_t  i,
const GEO::vec3 Pi,
double  wi,
double  Pi_len2,
GEO::index_t  t,
ConvexCell C,
IncidentTetrahedra W 
) const
protected

Copies a Laguerre cell facet from the triangulation.

Parameters
[in]ithe index of the vertex of which the Laguerre cell should be computed.
[in]Pithe coordinates of vertex i
[in]withe weight associated to vertex i
[in]Pi_len2the squared length of vertex i (considered as a vector).
[in]ta tetrahedron of the Delaunay triangulation, incident to vertex i
[out]Cthe Laguerre cell.
[in,out]Wa reference to a PeriodicDelaunay3d::IncidentTetrahedra
Returns
the local index of vertex i within tetrahedron t

◆ copy_Laguerre_cell_from_Delaunay() [1/2]

void GEO::PeriodicDelaunay3d::copy_Laguerre_cell_from_Delaunay ( GEO::index_t  i,
ConvexCell C 
) const
inline

Copies a Laguerre cell from the triangulation.

Delaunay neigbhors are stored in ConvexCell vertex global indices.

Parameters
[in]ithe index of the vertex of which the Laguerre cell should be computed.
[out]Cthe Laguerre cell.

Definition at line 249 of file periodic_delaunay_3d.h.

◆ copy_Laguerre_cell_from_Delaunay() [2/2]

void GEO::PeriodicDelaunay3d::copy_Laguerre_cell_from_Delaunay ( GEO::index_t  i,
ConvexCell C,
IncidentTetrahedra W 
) const

Copies a Laguerre cell from the triangulation.

Delaunay neigbhors are stored in ConvexCell vertex global indices.

Parameters
[in]ithe index of the vertex of which the Laguerre cell should be computed.
[out]Cthe Laguerre cell.
[in,out]Wa reference to a PeriodicDelaunay3d::IncidentTetrahedra

◆ get_incident_tets()

void GEO::PeriodicDelaunay3d::get_incident_tets ( index_t  v,
IncidentTetrahedra W 
) const

computes the set of tetrahedra that are incident to a vertex.

Parameters
[in]vthe index of the vertex.
[in,out]Wa reference to a PeriodicDelaunay3d::IncidentTetrahedra. On exit it contains the list of incident tets.

◆ get_periodic_vertex_instances_to_create()

index_t GEO::PeriodicDelaunay3d::get_periodic_vertex_instances_to_create ( index_t  v,
ConvexCell C,
bool  use_instance[27],
bool &  cell_is_on_boundary,
bool &  cell_is_outside_cube,
IncidentTetrahedra W 
)
protected

Computes the periodic vertex instances that should be generated.

Parameters
[in]vvertex index, in 0..nb_vertices_non_periodic_-1
[out]Cthe clipped Laguerre cell
[out]use_instancethe array of booleans that indicates which instance should be generated.
[out]cell_is_on_boundarytrue if the cell has an intersection with the cube, false otherwise.
[out]cell_is_outside_cubetrue if the cell is completely outside the cube, false otherwise.
[in,out]Wa reference to a PeriodicDelaunay3d::IncidentTetrahedra
Returns
the number of instances to generate.

◆ has_empty_cells()

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

Tests whether the Laguerre diagram has empty cells.

If the Laguerre diagram has empty cells, 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 265 of file periodic_delaunay_3d.h.

◆ insert_vertices()

void GEO::PeriodicDelaunay3d::insert_vertices ( index_t  b,
index_t  e 
)
protected

Insert vertices from reorder_[b] to reorder_[e-1].

If an empty cells is detected, has_empty_cells_ is set and the function exits.

◆ nb_threads()

index_t GEO::PeriodicDelaunay3d::nb_threads ( ) const
inlineprotected

Gets the number of threads.

Returns
the number of threads.

Definition at line 387 of file periodic_delaunay_3d.h.

◆ nearest_vertex()

index_t GEO::PeriodicDelaunay3d::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.

◆ save_cells()

void GEO::PeriodicDelaunay3d::save_cells ( const std::string &  basename,
bool  clipped 
)

Saves the cells in an Alias-Wavefront file.

Parameters
[in]basenamethe basename of the file. Filename is basename_callid.obj, where callid is the number of times the function was invoked.
[in]clippedif true, clip the cells by the domain without saving, else keep the cells as is.

◆ set_BRIO_levels()

void GEO::PeriodicDelaunay3d::set_BRIO_levels ( const vector< index_t > &  levels)
overridevirtual

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 from GEO::Delaunay.

◆ set_vertices()

void GEO::PeriodicDelaunay3d::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
Note
compute() needs to be called after.

Reimplemented from GEO::Delaunay.

◆ set_weights()

void GEO::PeriodicDelaunay3d::set_weights ( const double *  weights)

Sets the weights.

Parameters
[in]weightspointer to the array of weights. Size is the number of real vertices, i.e., the parameter nb_vertices passed to set_vertices().
Note
compute() needs to be called after.

◆ thread()

PeriodicDelaunay3dThread* GEO::PeriodicDelaunay3d::thread ( index_t  t)
inlineprotected

Gets a thread by index.

Parameters
[in]tthe index of the thread, in 0..nb_threads()-1
Returns
a pointer to the t-th thread.

Definition at line 376 of file periodic_delaunay_3d.h.

◆ update_cicl()

void GEO::PeriodicDelaunay3d::update_cicl ( )
overrideprotectedvirtual

Updates the circular incident cell lists.

Used by next_around_vertex().

Reimplemented from GEO::Delaunay.

◆ update_v_to_cell()

void GEO::PeriodicDelaunay3d::update_v_to_cell ( )
overrideprotectedvirtual

Stores for each vertex v a cell incident to v.

if update_periodic_v_to_cell_ is set to true, also updates the map periodic_v_to_cell_ that maps each virtual vertex to a tet incident to it.

Reimplemented from GEO::Delaunay.

◆ use_exact_predicates_for_convex_cell()

void GEO::PeriodicDelaunay3d::use_exact_predicates_for_convex_cell ( bool  x)
inline

Use exact predicates in convex cell computations.

Convex cell computations are used in periodic mode for determining the cells that straddle the domain boundary.

Parameters
[in]xtrue if exact predicates should be used (default), false otherwise.

Definition at line 166 of file periodic_delaunay_3d.h.

◆ vertex()

vec3 GEO::PeriodicDelaunay3d::vertex ( index_t  v) const
inline

Gets a vertex by index.

Parameters
[in]va vertex index. Can be a virtual vertex index when in periodic mode.
Returns
the 3d point associated with the vertex. In periodic mode, if v is a virtual vertex, then the translation is applied to the real vertex.

Definition at line 178 of file periodic_delaunay_3d.h.

◆ weight()

double GEO::PeriodicDelaunay3d::weight ( index_t  v) const
inline

Gets a weight by index.

Parameters
[in]va vertex index. Can be a virtual vertex index in periodic mode.
Returns
the weight associated with the vertex.

Definition at line 198 of file periodic_delaunay_3d.h.


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