Geogram
Version 1.9.1-rc
A programming library of geometric algorithms
|
Computes the intersection between a surface (Mesh) and a Voronoi diagram (dual of a Delaunay). More...
#include <geogram/voronoi/generic_RVD.h>
Classes | |
class | BorderHalfedgeAction |
Adapter class used internally to implement for_each_border_halfedge(). More... | |
class | HalfedgeAction |
Adapter class used internally to implement for_each_halfedge(). More... | |
class | PolygonAction |
Adapter class used internally to implement for_each_polygon() More... | |
class | PolyhedronAction |
Adapter class used internally to implement for_each_polyhedron() More... | |
class | PrimalTetrahedronAction |
Adapter class used internally to implement for_each_primal_tetrahedron() More... | |
class | PrimalTriangleAction |
Adapter class used internally to implement for_each_primal_triangle() More... | |
class | TetrahedronAction |
Adapter class used internally to implement for_each_tetrahedron() More... | |
class | TriangleAction |
Adapter class used internally to implement for_each_triangle(). More... | |
class | VolumetricIntegrationSimplexAction |
Adapter class used internally to implement for_each_volumetric_integration_simplex() More... | |
Public Types | |
typedef GEOGen::PointAllocator | PointAllocator |
Used to allocate the generated points. | |
typedef GEOGen::Vertex | Vertex |
Internal representation of vertices. | |
typedef GEOGen::Polygon | Polygon |
Internal representation of polygons. | |
typedef GEOGen::ConvexCell | Polyhedron |
Internal representation of volumetric cells. | |
Public Member Functions | |
RestrictedVoronoiDiagram (Delaunay *delaunay, GEO::Mesh *mesh) | |
Constructs a new RestrictedVoronoiDiagram. More... | |
void | set_connected_components_priority (bool x) |
Sets traveral priority. More... | |
bool | connected_components_priority () const |
Tests whether connected components priority is set. More... | |
const GEO::Mesh * | mesh () const |
Gets the input mesh. | |
GEO::Mesh * | mesh () |
Gets the input mesh. | |
const Delaunay * | delaunay () const |
Gets the Delaunay triangulation. | |
Delaunay * | delaunay () |
Gets the Delaunay triangulation. | |
void | set_delaunay (Delaunay *delaunay) |
Sets the Delaunay triangulation. | |
void | set_mesh (GEO::Mesh *mesh) |
Sets the input mesh. | |
void | set_facets_range (index_t facets_begin, index_t facets_end) |
Sets the facets range. More... | |
void | set_tetrahedra_range (index_t tets_begin, index_t tets_end) |
Sets the tetrahedra range. More... | |
index_t | nb_facets_in_range () const |
Gets the number of facets in the current range. More... | |
index_t | nb_tetrahedra_in_range () const |
Gets the number of tetrahedra in the current range. More... | |
index_t | current_facet () const |
Gets the index of the mesh facet currently processed. More... | |
index_t | current_seed () const |
Gets the index of the Delaunay vertex currently processed. More... | |
const Polygon & | current_polygon () const |
Gets the current polygon. More... | |
index_t | current_tet () const |
Gets the undex of the mesh tetrahedron currently processed. More... | |
const Polyhedron & | current_polyhedron () const |
Gets the current cell. More... | |
void | set_symbolic (bool x) |
Sets symbolic mode. More... | |
bool | symbolic () const |
Tests whether symbolic mode is active. | |
void | set_exact_predicates (bool x) |
Specifies whether exact predicates should be used. More... | |
bool | exact_predicates () const |
Tests whether exact predicates are used. | |
void | set_check_SR (bool x) |
Specifies whether radius of security should be enforced. | |
bool | check_SR () const |
Tests whether radius of security is enforced. More... | |
PointAllocator * | point_allocator () |
Gets the PointAllocator. More... | |
Public interface for computation/iteration | |
template<class ACTION > | |
void | for_each_polygon (const ACTION &action) |
Iterates on the facets of this RVD. More... | |
template<class TRIACTION > | |
void | for_each_triangle (const TRIACTION &action) |
Iterates on the facets of this RVD, triangulated on the fly. More... | |
template<class HEACTION > | |
void | for_each_halfedge (const HEACTION &action) |
Iterates on the halfedges on the borders of the restricted Voronoi cells. More... | |
template<class BOACTION > | |
void | for_each_border_halfedge (const BOACTION &action) |
Iterates on the halfedges on the borders of the restricted Voronoi cells that are on the boundary of the input mesh. More... | |
template<class PRIMTRIACTION > | |
void | for_each_primal_triangle (const PRIMTRIACTION &action) |
Iterates on the triangles of the Restricted Delaunay Triangulation. More... | |
template<class ACTION > | |
void | for_each_polyhedron (const ACTION &action) |
Iterates on the polyhedra of this RVD. More... | |
template<class ACTION > | |
void | for_each_volumetric_integration_simplex (const ACTION &action, bool visit_inner_tets=false, bool coherent_triangles=false) |
Iterates on the polyhedra of this RVD decomposed on the fly into tetrahedra. More... | |
template<class ACTION > | |
void | for_each_tetrahedron (const ACTION &action) |
Iterates on the polyhedra of this RVD decomposed on the fly into tetrahedra. More... | |
template<class ACTION > | |
void | for_each_primal_tetrahedron (const ACTION &action) |
Iterates on the primal tetrahedra of this RVD. More... | |
Static Public Member Functions | |
static coord_index_t | dimension () |
Gets the dimension. | |
Protected Member Functions | |
Clipping for surfacic mode | |
void | swap_polygons (Polygon *&ping, Polygon *&pong) |
Swaps two pointers between two polygons. More... | |
Polygon * | intersect_cell_facet (index_t seed, Polygon &F) |
Computes the intersection between the Voronoi cell of a seed and a facet. More... | |
void | clip_by_cell_SR (index_t i, Polygon *&ping, Polygon *&pong) |
Computes the intersection between the Voronoi cell of a vertex and the Mesh 'ping'. More... | |
void | clip_by_cell (index_t i, Polygon *&ping, Polygon *&pong) |
Computes the intersection between a Voronoi cell and a polygon. More... | |
void | clip_by_plane (Polygon &ping, Polygon &pong, index_t i, index_t j) |
Computes the intersection between a polygon and a half-space. More... | |
Optimized get neighbors | |
void | init_get_neighbors () |
Creates the data structure for optimized get_neighbors() function. More... | |
void | get_neighbors (index_t v) |
Caches the neighbors of a Delaunay vertex. More... | |
Protected Attributes | |
GEO::Mesh * | mesh_ |
Delaunay * | delaunay_ |
GEO::Delaunay_NearestNeighbors * | delaunay_nn_ |
PointAllocator | intersections_ |
Polygon * | current_polygon_ |
Polygon | P1 |
Polygon | P2 |
GEO::vector< index_t > | neighbors_ |
index_t | current_facet_ |
index_t | current_seed_ |
Polyhedron * | current_polyhedron_ |
index_t | current_tet_ |
signed_index_t | cur_stamp_ |
GEO::vector< signed_index_t > | stamp_ |
bool | symbolic_ |
bool | check_SR_ |
bool | exact_ |
coord_index_t | dimension_ |
index_t | facets_begin_ |
index_t | facets_end_ |
index_t | tets_begin_ |
index_t | tets_end_ |
bool | connected_components_priority_ |
FacetSeedMarking * | facet_seed_marking_ |
bool | connected_component_changed_ |
index_t | current_connected_component_ |
Static Protected Attributes | |
static const index_t | UNSPECIFIED_RANGE = index_t(-1) |
Computation | |
bool | facet_seed_is_visited (index_t f, index_t s) const |
Tests whether a (facet,seed) couple was visited. More... | |
signed_index_t | get_facet_seed_connected_component (index_t f, index_t s) const |
Gets the index of the connected component associated with a (facet,seed). More... | |
bool | connected_component_changed () const |
Tests whether the current connected component changed. | |
index_t | current_connected_component () const |
Gets the index of the current connected component. | |
template<class ACTION > | |
void | compute_surfacic (const ACTION &action) |
Low-level API of Restricted Voronoi Diagram traversal. More... | |
template<class ACTION > | |
void | compute_surfacic_with_seeds_priority (const ACTION &action) |
Low-level API of Restricted Voronoi Diagram traversal with seeds priority in surfacic mode. More... | |
template<class ACTION > | |
void | compute_volumetric (const ACTION &action) |
Low-level API of Restricted Voronoi Diagram traversal . More... | |
template<class ACTION > | |
void | compute_volumetric_with_seeds_priority (const ACTION &action) |
Low-level API of Restricted Voronoi Diagram traversal with seeds priority in volumetric mode. More... | |
template<class ACTION > | |
void | compute_volumetric_with_cnx_priority (const ACTION &action) |
Low-level API of Restricted Voronoi Diagram traversal with connected components priority. More... | |
template<class ACTION > | |
void | compute_surfacic_with_cnx_priority (const ACTION &action) |
Low-level API of Restricted Voronoi Diagram traversal with connected components priority. More... | |
index_t | find_seed_near_facet (index_t f) |
Finds a seed near a given facet. More... | |
index_t | find_seed_near_tet (index_t t) |
Finds a seed near a given tetrahedron. More... | |
index_t | find_seed_near_point (const double *p) |
Finds a seed near a given point. More... | |
Clipping for volumetric mode | |
void | intersect_cell_cell (index_t seed, Polyhedron &C) |
Computes the intersection between a Voronoi cell and a cell with radius of security or plain mode. More... | |
void | clip_by_cell_SR (index_t seed, Polyhedron &C) |
Computes the intersection between a Voronoi cell and a cell in radius-of-security mode. More... | |
void | clip_by_cell (index_t seed, Polyhedron &C) |
Computes the intersection between a Voronoi cell and a cell in plain mode. More... | |
void | clip_by_plane (Polyhedron &C, index_t i, index_t j) |
Computes the intersection between a Voronoi cell and a half-space determined by a bisector. More... | |
Computes the intersection between a surface (Mesh) and a Voronoi diagram (dual of a Delaunay).
The surface may be embedded in nD (the Voronoi diagram is then of dimension n).
Definition at line 78 of file generic_RVD.h.
|
inline |
Constructs a new RestrictedVoronoiDiagram.
[in] | delaunay | the Delaunay triangulation |
[in] | mesh | the input mesh |
Definition at line 118 of file generic_RVD.h.
|
inline |
Tests whether radius of security is enforced.
true | if radius of security test is used. |
false | otherwise. |
Definition at line 365 of file generic_RVD.h.
|
inlineprotected |
Computes the intersection between a Voronoi cell and a polygon.
The Voronoi cell is determined by vertex i
and the input polygon is in ping
. The result is returned in ping
(Note that ping
and pong
are references, and that they are swapped after each bisector clipping, this is why the final result is in ping
(and not in pong
).
[in] | i | index of the vertex that defines the Voronoi cell |
[in,out] | ping | the input polygon. On exit, contains the result. |
[out] | pong | a buffer used to implement reentrant clipping. Its content is modified by the function. |
Definition at line 2240 of file generic_RVD.h.
|
inlineprotected |
Computes the intersection between a Voronoi cell and a cell in plain mode.
[in] | seed | the index of the seed that defines the Voronoi cell |
[in,out] | C | the cell to be clipped |
Definition at line 2378 of file generic_RVD.h.
|
inlineprotected |
Computes the intersection between the Voronoi cell of a vertex and the Mesh 'ping'.
The result is returned in ping
(Note that ping
and pong
are references, and that they are swapped after each bisector clipping, this is why the final result is in ping
(and not in pong
). This version uses the Security Radius algorithm.
[in] | i | index of the vertex that defines the Voronoi cell |
[in,out] | ping | the input polygon. On exit, contains the result. |
[out] | pong | a buffer used to implement reentrant clipping. Its content is modified by the function. |
Definition at line 2157 of file generic_RVD.h.
|
inlineprotected |
Computes the intersection between a Voronoi cell and a cell in radius-of-security mode.
[in] | seed | the index of the seed that defines the Voronoi cell |
[in,out] | C | the cell to be clipped |
Definition at line 2301 of file generic_RVD.h.
|
inlineprotected |
Computes the intersection between a polygon and a half-space.
The input polygon is in ping
and the half-space is determined by the positive side of the bisector of segment [i
,j
] (the side of i
). The result is stored into the Polygon pong
.
[in] | i | index of the first extremity of the bisector |
[in] | j | index of the second extremity of the bisector |
[in] | ping | the input polygon |
[out] | pong | ping clipped by the bisector |
Definition at line 2263 of file generic_RVD.h.
|
inlineprotected |
Computes the intersection between a Voronoi cell and a half-space determined by a bisector.
[in,out] | C | cell to be clipped |
[in] | i | index of the first extremity of the bisector |
[in] | j | index of the second extremity of the bisector |
Definition at line 2398 of file generic_RVD.h.
|
inlineprotected |
Low-level API of Restricted Voronoi Diagram traversal.
Client code may use for_each_facet(),for_each_triangle() or for_each_primal_triangle() instead.
ACTION | needs to implement: operator()(index_t v, index_t f, const Polygon& P) const where v denotes the index of the current Voronoi cell (or Delaunay vertex), f the index of the current facet and P the computed intersection between the Voronoi cell of v and facet f. |
Definition at line 1314 of file generic_RVD.h.
|
inlineprotected |
Low-level API of Restricted Voronoi Diagram traversal with connected components priority.
Client code may use for_each_facet(),for_each_triangle() or for_each_primal_triangle() instead. This version of the algorithm traverses the RVD and ensures that the group of subfacets that belong to the same restricted Voronoi cell will be traversed consecutively. It is used by the algorithm that computes the final surface in CVT (i.e., the dual of the connected components).
ACTION | needs to implement: operator()(index_t v, index_t f, const Polygon& P) const where v denotes the index of the current Voronoi cell (or Delaunay vertex), f the index of the current facet and P the computed intersection between the Voronoi cell of v and facet f. |
Definition at line 1879 of file generic_RVD.h.
|
inlineprotected |
Low-level API of Restricted Voronoi Diagram traversal with seeds priority in surfacic mode.
Client code may use for_each_facet(),for_each_triangle() or for_each_primal_triangle() instead.
ACTION | needs to implement: operator()(index_t v, index_t f, const Polygon& P) const where v denotes the index of the current Voronoi cell (or Delaunay vertex), f the index of the current facet and P the computed intersection between the Voronoi cell of v and facet f. |
Definition at line 1339 of file generic_RVD.h.
|
inlineprotected |
Low-level API of Restricted Voronoi Diagram traversal .
Selects seed-priority or tetrahedron-priority modes according to connected_components_priority mode. Client code may use for_each_polyhedron() or for_each_volumetric_integration_simplex() instead of this function.
ACTION | needs to implement: operator()(index_t v, index_t t, const Polyhedron& C) const where v denotes the index of the current Voronoi cell (or Delaunay vertex), t the index of the current tetrahedron and C the computed intersection between the Voronoi cell of v and tetrahedron t |
Definition at line 1462 of file generic_RVD.h.
|
inlineprotected |
Low-level API of Restricted Voronoi Diagram traversal with connected components priority.
Client code may use for_each_cell() instead. This version of the algorithm traverses the RVD and ensures that the group of subfacets that belong to the same restricted Voronoi cell will be traversed consecutively. It is used by the algorithm that computes the final surface in CVT (i.e., the dual of the connected components).
ACTION | needs to implement: operator()(index_t v, index_t t, const Polyhedron& C) const where v denotes the index of the current Voronoi cell (or Delaunay vertex), c the index of the current tetrahedron and C the computed intersection between the Voronoi cell of v and tetrahedron t. |
Definition at line 1646 of file generic_RVD.h.
|
inlineprotected |
Low-level API of Restricted Voronoi Diagram traversal with seeds priority in volumetric mode.
Client code may use for_each_polyhedron() or for_each_volumetric_integration_simplex() instead.
ACTION | needs to implement: operator()(index_t v, index_t t, const Polyhedron& C) const where v denotes the index of the current Voronoi cell (or Delaunay vertex), t the index of the current tetrahedron and C the computed intersection between the Voronoi cell of v and tetrahedron t |
Definition at line 1489 of file generic_RVD.h.
|
inline |
Tests whether connected components priority is set.
If connected_components_priority is set, then the connected components of the restricted Voronoi cells will be traversed one by one.
true | if connected components priority is used. |
false | otherwise. |
Definition at line 171 of file generic_RVD.h.
|
inline |
Gets the index of the mesh facet currently processed.
Can be used in surfacic traversals (and not volumetric traversals).
Definition at line 267 of file generic_RVD.h.
|
inline |
Gets the current polygon.
The current polygon corresponds to the intersection between the current facet and the Voronoi cell of the current seed. Can be used in surfacic traversals (and not volumetric traversals).
Definition at line 287 of file generic_RVD.h.
|
inline |
Gets the current cell.
The current cell corresponds to the intersection between the current tetrahedron and the Voronoi cell of the current seed. Can be used in volumetric traversals (and not in surfacic traversals).
Definition at line 308 of file generic_RVD.h.
|
inline |
Gets the index of the Delaunay vertex currently processed.
Can be used in both surfacic traversals and volumetric traversals.
Definition at line 276 of file generic_RVD.h.
|
inline |
Gets the undex of the mesh tetrahedron currently processed.
Can be used in volumetric traversals (and not in surfacic traversals).
Definition at line 296 of file generic_RVD.h.
|
inline |
Tests whether a (facet,seed) couple was visited.
[in] | f | index of the facet |
[in] | s | index of the seed |
Definition at line 1821 of file generic_RVD.h.
|
inlineprotected |
Finds a seed near a given facet.
[in] | f | index of the facet in the mesh |
f
. Definition at line 2033 of file generic_RVD.h.
|
inlineprotected |
Finds a seed near a given point.
[in] | p | pointer to the coordinates of the point |
p
. Definition at line 2059 of file generic_RVD.h.
|
inlineprotected |
Finds a seed near a given tetrahedron.
[in] | t | index of the tetrahedron in the mesh |
t
. Definition at line 2047 of file generic_RVD.h.
|
inline |
Iterates on the halfedges on the borders of the restricted Voronoi cells that are on the boundary of the input mesh.
[in] | action | the user action object |
BOACTION | needs to implement: operator()(index_t c, const TopoPolyVertexEdge& v1, v2) const where c denotes the index of the current Voronoi cell (or Delaunay vertex). |
Definition at line 1137 of file generic_RVD.h.
|
inline |
Iterates on the halfedges on the borders of the restricted Voronoi cells.
[in] | action | the user action object |
HEACTION | needs to implement: operator()(index_t c, const TopoPolyVertex& v1, v2) const where c denotes the index of the current Voronoi cell (or Delaunay vertex). |
Definition at line 1119 of file generic_RVD.h.
|
inline |
Iterates on the facets of this RVD.
[in] | action | the user action object |
ACTION | needs to implement: operator()(index_t v, index_t f, const Polygon& P) const where v denotes the index of the current Voronoi cell (or Delaunay vertex), f the index of the current facet and P the computed intersection between facet f and the Voronoi cell of v. |
Definition at line 1084 of file generic_RVD.h.
|
inline |
Iterates on the polyhedra of this RVD.
[in] | action | the user action object |
ACTION | needs to implement: operator()(index_t v, index_t t, const Polyhedron& C) const where v denotes the index of the current Voronoi cell (or Delaunay vertex), t the index of the current tetrahedron and C the computed intersection between tetrahedron t and the Voronoi cell of v. |
Definition at line 1179 of file generic_RVD.h.
|
inline |
Iterates on the primal tetrahedra of this RVD.
The tetrahedra are not coherently oriented, and need a subsequent traversal operation to reorient them. They can be also reoriented geometrically using the orient3d() predicate.
[in] | action | the user action object |
ACTION | needs to implement: operator()(index_t v0, index_t v1, index_t v2, index_t v3) where v0,v1,v2 and v3 are the indices of the four vertices of tetrahedron. |
Definition at line 1282 of file generic_RVD.h.
|
inline |
Iterates on the triangles of the Restricted Delaunay Triangulation.
[in] | action | the user action object |
PRIMTRIACTION | needs to implement: operator()(index_t i, unsigned j, index_t k) const where i,j,k denote the three indices of the Delaunay vertices that define the primal triangle. |
Definition at line 1155 of file generic_RVD.h.
|
inline |
Iterates on the polyhedra of this RVD decomposed on the fly into tetrahedra.
The tetrahedra are generated by connecting one of the vertices of the cell to the other ones.
[in] | action | the user action object |
ACTION | needs to implement: operator()(index_t v, signed_index_t v_adj, index_t t, index_t t_adj, const Vertex& v0, const Vertex& v1, const Vertex& v2, const Vertex& v3 ) where the parameters are as follows:
|
Definition at line 1259 of file generic_RVD.h.
|
inline |
Iterates on the facets of this RVD, triangulated on the fly.
[in] | action | the user action object |
TRIACTION | needs to implement: operator()(index_t c, const TopoPolyVertex& v1, v2, v3) const where c denotes the index of the current Voronoi cell (or Delaunay vertex). |
Definition at line 1101 of file generic_RVD.h.
|
inline |
Iterates on the polyhedra of this RVD decomposed on the fly into tetrahedra.
The generated tetrahedra may be geometrically incorrect, but they are algebraically correct. In other word, their signed volumes sum as the volume of the restricted Voronoi cell.
[in] | action | the user action object |
[in] | visit_inner_tets | if set, all the tetrahedron-cell intersections are visited, else only tetrahedra on the border of the restricted Voronoi cell are visited. Since all the visited triangles are connected to the current Voronoi seed by a tetrahedron, the computed volume is the same in both cases. |
[in] | coherent_triangles | if set, this ensures that the polygonal facets of the cells are always triangulated in a coherent manner when seen from two different cells. For instance, it is required if a tetrahedral mesh is reconstructed. |
ACTION | needs to implement: operator()(index_t v, signed_index_t v_adj, index_t t, index_t t_adj, const Vertex& v1, const Vertex& v2, const Vertex& v3 ) where the parameters are as follows:
|
Definition at line 1222 of file generic_RVD.h.
|
inline |
Gets the index of the connected component associated with a (facet,seed).
[in] | f | index of the facet |
[in] | s | index of the seed |
f
, s
) couple was not visited already Definition at line 1834 of file generic_RVD.h.
|
inlineprotected |
Caches the neighbors of a Delaunay vertex.
This function is only used when the stored delaunay triangulation is a traditional one. When the stored delaunay triangulation is represented by a KdTree, function is not used.
Definition at line 2435 of file generic_RVD.h.
|
inlineprotected |
Creates the data structure for optimized get_neighbors() function.
This function is only used when the stored delaunay triangulation is a traditional one. When the stored delaunay triangulation is represented by a KdTree, function is not used.
Definition at line 2418 of file generic_RVD.h.
|
inline |
Computes the intersection between a Voronoi cell and a cell with radius of security or plain mode.
[in] | seed | the index of the seed that defines the Voronoi cell |
[in,out] | C | the cell to be clipped |
Definition at line 2285 of file generic_RVD.h.
|
inlineprotected |
Computes the intersection between the Voronoi cell of a seed and a facet.
[in] | seed | the index of the seed |
[in] | F | the facet represented as a Polygon |
The result is provided in current_polygon_
Definition at line 2123 of file generic_RVD.h.
|
inline |
Gets the number of facets in the current range.
Definition at line 250 of file generic_RVD.h.
|
inline |
Gets the number of tetrahedra in the current range.
Definition at line 258 of file generic_RVD.h.
|
inline |
Gets the PointAllocator.
Definition at line 375 of file generic_RVD.h.
|
inline |
Sets traveral priority.
If connected_components_priority is set, then the connected components of the restricted Voronoi cells will be traversed one by one.
Definition at line 156 of file generic_RVD.h.
|
inline |
Specifies whether exact predicates should be used.
If exact predicates are used, symbolic mode is ensured.
[in] | x | if set, exact predicates are used. |
Definition at line 338 of file generic_RVD.h.
|
inline |
Sets the facets range.
Computations can be restricted to a contiguous facet range.
[in] | facets_begin | first facet in the range. |
[in] | facets_end | one position past the last facet in the range. |
Definition at line 226 of file generic_RVD.h.
|
inline |
Sets symbolic mode.
If exact mode is active, symbolic mode is enforced.
[in] | x | if set, the symbolic representation of the intersections are computed. |
Definition at line 318 of file generic_RVD.h.
|
inline |
Sets the tetrahedra range.
Computations can be restricted to a contiguous tetrahedra range.
[in] | tets_begin | first tetrahedron in the range. |
[in] | tets_end | one position past the last tetrahedron in the range. |
Definition at line 240 of file generic_RVD.h.
|
inlineprotected |
Swaps two pointers between two polygons.
Used by re-entrant Sutherlang-Hogdman clipping.
Definition at line 2105 of file generic_RVD.h.