Graphite Version 3
An experimental 3D geometry processing program
|
Multithreaded implementation of Delaunay in 3d with optional periodic boundary conditions. More...
#include <geogram/delaunay/periodic_delaunay_3d.h>
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. | |
PeriodicDelaunay3d (const vec3 &period) | |
Constructs a new PeriodicDelaunay3d. | |
void | set_vertices (index_t nb_vertices, const double *vertices) override |
Sets the vertices of this Delaunay, and recomputes the cells. | |
void | set_weights (const double *weights) |
Sets the weights. | |
void | compute () |
Computes the Delaunay triangulation. | |
void | use_exact_predicates_for_convex_cell (bool x) |
Use exact predicates in convex cell computations. | |
vec3 | vertex (index_t v) const |
Gets a vertex by index. | |
double | weight (index_t v) const |
Gets a weight by index. | |
index_t | nearest_vertex (const double *p) const override |
Computes the nearest vertex from a query point. | |
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. | |
void | get_incident_tets (index_t v, IncidentTetrahedra &W) const |
computes the set of tetrahedra that are incident to a vertex. | |
void | copy_Laguerre_cell_from_Delaunay (GEO::index_t i, ConvexCell &C, IncidentTetrahedra &W) const |
Copies a Laguerre cell from the triangulation. | |
void | copy_Laguerre_cell_from_Delaunay (GEO::index_t i, ConvexCell &C) const |
Copies a Laguerre cell from the triangulation. | |
bool | has_empty_cells () const |
Tests whether the Laguerre diagram has empty cells. | |
void | save_cells (const std::string &basename, bool clipped) |
Saves the cells in an Alias-Wavefront file. | |
![]() | |
coord_index_t | dimension () const |
Gets the dimension of this Delaunay. | |
index_t | cell_size () const |
Gets the number of vertices in each cell. | |
void | set_reorder (bool x) |
Specifies whether vertices should be reordered. | |
const double * | vertices_ptr () const |
Gets a pointer to the array of vertices. | |
const double * | vertex_ptr (index_t i) const |
Gets a pointer to a vertex by its global index. | |
index_t | nb_vertices () const |
Gets the number of vertices. | |
virtual bool | supports_constraints () const |
Tests whether constraints are supported by this Delaunay. | |
virtual void | set_constraints (const Mesh *mesh) |
Defines the constraints. | |
void | set_refine (bool x) |
Specifies whether the mesh should be refined. | |
bool | get_refine () const |
Tests whether mesh refinement is selected. | |
void | set_quality (double qual) |
Specifies the desired quality for mesh elements when refinement is enabled (. | |
const Mesh * | constraints () const |
Gets the constraints. | |
index_t | nb_cells () const |
Gets the number of cells. | |
index_t | nb_finite_cells () const |
Gets the number of finite cells. | |
const signed_index_t * | cell_to_v () const |
Gets a pointer to the cell-to-vertex incidence array. | |
const signed_index_t * | cell_to_cell () const |
Gets a pointer to the cell-to-cell adjacency array. | |
signed_index_t | cell_vertex (index_t c, index_t lv) const |
Gets a vertex index by cell index and local vertex index. | |
signed_index_t | cell_adjacent (index_t c, index_t lf) const |
Gets an adjacent cell index by cell index and local facet index. | |
bool | cell_is_infinite (index_t c) const |
Tests whether a cell is infinite. | |
bool | cell_is_finite (index_t c) const |
Tests whether a cell is finite. | |
index_t | index (index_t c, signed_index_t v) const |
Retrieves a local vertex index from cell index and global vertex index. | |
index_t | adjacent_index (index_t c1, index_t c2) const |
Retrieves a local facet index from two adacent cell global indices. | |
signed_index_t | vertex_cell (index_t v) const |
Gets an incident cell index by a vertex index. | |
signed_index_t | next_around_vertex (index_t c, index_t lv) const |
Traverses the list of cells incident to a vertex. | |
void | get_neighbors (index_t v, vector< index_t > &neighbors) const |
Gets the one-ring neighbors of vertex v. | |
void | save_histogram (std::ostream &out) const |
Saves the histogram of vertex degree (can be visualized with gnuplot). | |
bool | stores_neighbors () const |
Tests whether neighbors are stored. | |
void | set_stores_neighbors (bool x) |
Specifies whether neighbors should be stored. | |
bool | stores_cicl () const |
Tests whether incident tetrahedra lists are stored. | |
void | set_stores_cicl (bool x) |
Specifies whether incident tetrahedra lists should be stored. | |
bool | keeps_infinite () const |
Tests whether infinite elements are kept. | |
void | set_keeps_infinite (bool x) |
Sets whether infinite elements should be kept. | |
bool | thread_safe () const |
Tests whether thread-safe mode is active. | |
void | set_thread_safe (bool x) |
Specifies whether thread-safe mode should be used. | |
void | set_default_nb_neighbors (index_t x) |
Sets the default number of stored neighbors. | |
index_t | default_nb_neighbors () const |
Gets the default number of stored neighbors. | |
void | clear_neighbors () |
Frees all memory used for neighbors storage. | |
void | set_keep_regions (bool x) |
Specifies whether all internal regions should be kept. | |
virtual index_t | region (index_t t) const |
Gets the region id associated with a tetrahedron. | |
virtual void | store_neighbors_CB (index_t i) |
Stores the neighbors of a vertex. | |
![]() | |
void | ref () const |
Increments the reference count. | |
void | unref () const |
Decrements the reference count. | |
bool | is_shared () const |
Check if the object is shared. | |
int | nb_refs () const |
Gets the number of references that point to this object. | |
![]() | |
index_t | nb_vertices_non_periodic () const |
Gets the number of non-periodic vertices. | |
index_t | periodic_vertex_instance (index_t pv) const |
Gets the instance from a periodic vertex. | |
index_t | periodic_vertex_real (index_t pv) const |
Gets the real vertex from a periodic vertex. | |
signed_index_t | periodic_vertex_real (signed_index_t pv) const |
Gets the real vertex from a periodic vertex. | |
index_t | make_periodic_vertex (index_t real, index_t instance) const |
Makes a periodic vertex from a real vertex and instance. | |
void | periodic_vertex_get_T (index_t pv, int &Tx, int &Ty, int &Tz) const |
Gets the translation from a periodic vertex. | |
void | periodic_vertex_set_T (index_t &pv, int Tx, int Ty, int Tz) const |
Sets the translation in a periodic vertex. | |
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. | |
index_t | compress (bool shrink=true) |
Removes unused tetrahedra. | |
void | update_v_to_cell () override |
Stores for each vertex v a cell incident to v. | |
void | update_cicl () override |
Updates the circular incident cell lists. | |
void | handle_periodic_boundaries () |
Duplicates the points with Voronoi cells that cross the boundary. | |
void | handle_periodic_boundaries_phase_I () |
Phase I of periodic boundaries handling. | |
bool | Laguerre_vertex_is_in_conflict_with_plane (index_t t, vec4 P) const |
Tests the position of a Laguerre vertex w.r.t. a plane. | |
void | handle_periodic_boundaries_phase_II () |
Phase II of periodic boundaries handling. | |
void | insert_vertices (const char *phase, index_t b, index_t e) |
Inserts vertices from reorder_[b] to reorder_[e-1] using multithreaded Delaunay. Called by insert_vertices() if there are many vertices to insert. | |
void | insert_vertices_with_BRIO (const char *phase, const vector< index_t > &levels) |
Inserts vertices as indicated by a reordering vector and a vector of BRIO levels, as obtained using compute_BRIO_order_periodic() (internal function, in the .cpp) | |
void | check_volume () |
Checks the volume of Laguerre cells. | |
PeriodicDelaunay3dThread * | thread (index_t t) |
Gets a thread by index. | |
index_t | nb_threads () const |
Gets the number of threads. | |
void | check_max_t () |
![]() | |
Delaunay (coord_index_t dimension) | |
Creates a new Delaunay triangulation. | |
~Delaunay () override | |
Delaunay destructor. | |
virtual void | get_neighbors_internal (index_t v, vector< index_t > &neighbors) const |
Internal implementation for get_neighbors (with vector). | |
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. | |
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. | |
void | set_dimension (coord_index_t dim) |
Sets the dimension of this Delaunay. | |
![]() | |
Counted () | |
Creates a reference counted object. | |
virtual | ~Counted () |
Destroys a reference counted object. | |
Friends | |
class | PeriodicDelaunay3dThread |
class | LaguerreDiagramOmegaSimple3d |
Additional Inherited Members | |
![]() | |
static Delaunay * | create (coord_index_t dim, const std::string &name="default") |
Creates a Delaunay triangulation of the specified dimension. | |
static void | initialize () |
This function needs to be called once before using the Delaunay class. | |
![]() | |
static void | ref (const Counted *counted) |
Increments the reference count. | |
static void | unref (const Counted *counted) |
Decrements the reference count. | |
![]() | |
static index_t | T_to_instance (int Tx, int Ty, int Tz) |
Gets the instance from a translation. | |
![]() | |
index_t | nb_vertices_non_periodic_ |
Number of real vertices. | |
![]() | |
static int | translation [27][3] |
Gives for each instance the integer translation coordinates in {-1,0,1}. | |
static int | reorder_instances [27] |
Used to back-map an integer translation to an instance. | |
static bool | instance_is_positive [27] |
Tests whether all the coordinates of the translation vector associated with an instance are 0 or 1. | |
![]() | |
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_ |
![]() | |
typedef SmartPointer< Delaunay > | Delaunay_var |
Smart pointer that refers to a Delaunay object. | |
typedef Factory1< Delaunay, coord_index_t > | DelaunayFactory |
Delaunay Factory. | |
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).
Definition at line 67 of file periodic_delaunay_3d.h.
GEO::PeriodicDelaunay3d::PeriodicDelaunay3d | ( | bool | periodic, |
double | period = 1.0 |
||
) |
Constructs a new PeriodicDelaunay3d.
[in] | periodic | if true, constructs a periodic triangulation. |
[in] | period | the edge length of the periodic domain. |
GEO::PeriodicDelaunay3d::PeriodicDelaunay3d | ( | const vec3 & | period | ) |
Constructs a new PeriodicDelaunay3d.
[in] | period | the edge lengths along x,y,z in the periodic domain |
|
protected |
Removes unused tetrahedra.
[in] | shrink | if true, then array space is shrunk to fit the new number of tetrahedra. |
|
protected |
Copies a Laguerre cell facet from the triangulation.
[in] | i | the index of the vertex of which the Laguerre cell should be computed. |
[in] | Pi | the coordinates of vertex i |
[in] | wi | the weight associated to vertex i |
[in] | Pi_len2 | the squared length of vertex i (considered as a vector). |
[in] | t | a tetrahedron of the Delaunay triangulation, incident to vertex i |
[out] | C | the Laguerre cell. |
[in,out] | W | a reference to a PeriodicDelaunay3d::IncidentTetrahedra |
i
within tetrahedron t
|
inline |
Copies a Laguerre cell from the triangulation.
Delaunay neigbhors are stored in ConvexCell vertex global indices.
[in] | i | the index of the vertex of which the Laguerre cell should be computed. |
[out] | C | the Laguerre cell. |
Definition at line 248 of file periodic_delaunay_3d.h.
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.
[in] | i | the index of the vertex of which the Laguerre cell should be computed. |
[out] | C | the Laguerre cell. |
[in,out] | W | a reference to a PeriodicDelaunay3d::IncidentTetrahedra |
void GEO::PeriodicDelaunay3d::get_incident_tets | ( | index_t | v, |
IncidentTetrahedra & | W | ||
) | const |
computes the set of tetrahedra that are incident to a vertex.
[in] | v | the index of the vertex. |
[in,out] | W | a reference to a PeriodicDelaunay3d::IncidentTetrahedra. On exit it contains the list of incident tets. |
|
protected |
Phase I of periodic boundaries handling.
For each cell that traverses the boundary, generates the periodic vertices instances corresponding to the crossed faces of the domain, and inserts them in the list of vertices to be inserted (the reorder_ vector).
|
protected |
Phase II of periodic boundaries handling.
Adds the newly discovered neighbors of the vertices inserted during phase I, back-translated to their original domain instances in the list of vertices to be inserted (in the reorder_ member).
|
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).
true | if the Laguerre diagram has empty cells. |
false | otherwise. |
Definition at line 264 of file periodic_delaunay_3d.h.
|
protected |
Inserts vertices from reorder_[b] to reorder_[e-1] using multithreaded Delaunay. Called by insert_vertices() if there are many vertices to insert.
If an empty cells is detected, has_empty_cells_ is set and the function exits.
|
protected |
Inserts vertices as indicated by a reordering vector and a vector of BRIO levels, as obtained using compute_BRIO_order_periodic() (internal function, in the .cpp)
used by insert_vertices() (and maybe also compute(), we shall see if we can). Returns immediatly if an empty cell is encountered, then it sets has_empty_cells_.
[in] | levels | a const reference to the vector of BRIO levels, as offsets in the member reordering vector reorder_. |
|
protected |
Tests the position of a Laguerre vertex w.r.t. a plane.
The positive side of the plane equation corresponds to what is kept. In other words, the normal vector P.x, P.y, P.z points towards the interior of this ConvexCell.
[in] | t | a tetrahedron index. The considered Laguerre vertex is the dual of this tetrahedron. |
[in] | P | the plane equation. |
true | if the Laguerre vertex would be clipped by plane P |
false | otherwise |
|
inlineprotected |
Gets the number of threads.
Definition at line 409 of file periodic_delaunay_3d.h.
|
overridevirtual |
Computes the nearest vertex from a query point.
[in] | p | query point |
Reimplemented from GEO::Delaunay.
void GEO::PeriodicDelaunay3d::save_cells | ( | const std::string & | basename, |
bool | clipped | ||
) |
Saves the cells in an Alias-Wavefront file.
[in] | basename | the basename of the file. Filename is basename_callid.obj, where callid is the number of times the function was invoked. |
[in] | clipped | if true, clip the cells by the domain without saving, else keep the cells as is. |
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 from GEO::Delaunay.
|
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.
void GEO::PeriodicDelaunay3d::set_weights | ( | const double * | weights | ) |
Sets the weights.
[in] | weights | pointer to the array of weights. Size is the number of real vertices, i.e., the parameter nb_vertices passed to set_vertices(). |
|
inlineprotected |
Gets a thread by index.
[in] | t | the index of the thread, in 0..nb_threads()-1 |
Definition at line 398 of file periodic_delaunay_3d.h.
|
overrideprotectedvirtual |
Updates the circular incident cell lists.
Used by next_around_vertex().
Reimplemented from GEO::Delaunay.
|
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.
|
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.
[in] | x | true if exact predicates should be used (default), false otherwise. |
Definition at line 165 of file periodic_delaunay_3d.h.
Gets a vertex by index.
[in] | v | a vertex index. Can be a virtual vertex index when in periodic mode. |
v
is a virtual vertex, then the translation is applied to the real vertex. Definition at line 177 of file periodic_delaunay_3d.h.
|
inline |
Gets a weight by index.
[in] | v | a vertex index. Can be a virtual vertex index in periodic mode. |
Definition at line 197 of file periodic_delaunay_3d.h.
|
friend |
Definition at line 515 of file periodic_delaunay_3d.h.
|
friend |
Definition at line 416 of file periodic_delaunay_3d.h.