Geogram  Version 1.8.9-rc
A programming library of geometric algorithms
GEOGen::RestrictedVoronoiDiagram< DIM > Class Template Reference

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::Meshmesh () const
 Gets the input mesh.
 
GEO::Meshmesh ()
 Gets the input mesh.
 
const Delaunaydelaunay () const
 Gets the Delaunay triangulation.
 
Delaunaydelaunay ()
 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 Polygoncurrent_polygon () const
 Gets the current polygon. More...
 
index_t current_tet () const
 Gets the undex of the mesh tetrahedron currently processed. More...
 
const Polyhedroncurrent_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...
 
PointAllocatorpoint_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...
 
Polygonintersect_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::Meshmesh_
 
Delaunaydelaunay_
 
GEO::Delaunay_NearestNeighborsdelaunay_nn_
 
PointAllocator intersections_
 
Polygoncurrent_polygon_
 
Polygon P1
 
Polygon P2
 
GEO::vector< index_t > neighbors_
 
index_t current_facet_
 
index_t current_seed_
 
Polyhedroncurrent_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_
 
FacetSeedMarkingfacet_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...
 

Detailed Description

template<index_t DIM>
class GEOGen::RestrictedVoronoiDiagram< DIM >

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

Note
This is an internal implementation class, not meant to be used directly, use GEO::RestrictedVoronoiDiagram instead.

Definition at line 78 of file generic_RVD.h.

Constructor & Destructor Documentation

◆ RestrictedVoronoiDiagram()

template<index_t DIM>
GEOGen::RestrictedVoronoiDiagram< DIM >::RestrictedVoronoiDiagram ( Delaunay delaunay,
GEO::Mesh mesh 
)
inline

Constructs a new RestrictedVoronoiDiagram.

Parameters
[in]delaunaythe Delaunay triangulation
[in]meshthe input mesh

Definition at line 118 of file generic_RVD.h.

Member Function Documentation

◆ check_SR()

template<index_t DIM>
bool GEOGen::RestrictedVoronoiDiagram< DIM >::check_SR ( ) const
inline

Tests whether radius of security is enforced.

Return values
trueif radius of security test is used.
falseotherwise.

Definition at line 365 of file generic_RVD.h.

◆ clip_by_cell() [1/2]

template<index_t DIM>
void GEOGen::RestrictedVoronoiDiagram< DIM >::clip_by_cell ( index_t  i,
Polygon *&  ping,
Polygon *&  pong 
)
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).

Parameters
[in]iindex of the vertex that defines the Voronoi cell
[in,out]pingthe input polygon. On exit, contains the result.
[out]ponga buffer used to implement reentrant clipping. Its content is modified by the function.

Definition at line 2240 of file generic_RVD.h.

◆ clip_by_cell() [2/2]

template<index_t DIM>
void GEOGen::RestrictedVoronoiDiagram< DIM >::clip_by_cell ( index_t  seed,
Polyhedron C 
)
inlineprotected

Computes the intersection between a Voronoi cell and a cell in plain mode.

Parameters
[in]seedthe index of the seed that defines the Voronoi cell
[in,out]Cthe cell to be clipped

Definition at line 2378 of file generic_RVD.h.

◆ clip_by_cell_SR() [1/2]

template<index_t DIM>
void GEOGen::RestrictedVoronoiDiagram< DIM >::clip_by_cell_SR ( index_t  i,
Polygon *&  ping,
Polygon *&  pong 
)
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.

Parameters
[in]iindex of the vertex that defines the Voronoi cell
[in,out]pingthe input polygon. On exit, contains the result.
[out]ponga buffer used to implement reentrant clipping. Its content is modified by the function.

Definition at line 2157 of file generic_RVD.h.

◆ clip_by_cell_SR() [2/2]

template<index_t DIM>
void GEOGen::RestrictedVoronoiDiagram< DIM >::clip_by_cell_SR ( index_t  seed,
Polyhedron C 
)
inlineprotected

Computes the intersection between a Voronoi cell and a cell in radius-of-security mode.

Parameters
[in]seedthe index of the seed that defines the Voronoi cell
[in,out]Cthe cell to be clipped

Definition at line 2301 of file generic_RVD.h.

◆ clip_by_plane() [1/2]

template<index_t DIM>
void GEOGen::RestrictedVoronoiDiagram< DIM >::clip_by_plane ( Polygon ping,
Polygon pong,
index_t  i,
index_t  j 
)
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.

Parameters
[in]iindex of the first extremity of the bisector
[in]jindex of the second extremity of the bisector
[in]pingthe input polygon
[out]pongping clipped by the bisector

Definition at line 2263 of file generic_RVD.h.

◆ clip_by_plane() [2/2]

template<index_t DIM>
void GEOGen::RestrictedVoronoiDiagram< DIM >::clip_by_plane ( Polyhedron C,
index_t  i,
index_t  j 
)
inlineprotected

Computes the intersection between a Voronoi cell and a half-space determined by a bisector.

Parameters
[in,out]Ccell to be clipped
[in]iindex of the first extremity of the bisector
[in]jindex of the second extremity of the bisector

Definition at line 2398 of file generic_RVD.h.

◆ compute_surfacic()

template<index_t DIM>
template<class ACTION >
void GEOGen::RestrictedVoronoiDiagram< DIM >::compute_surfacic ( const ACTION &  action)
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.

Template Parameters
ACTIONneeds 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.

◆ compute_surfacic_with_cnx_priority()

template<index_t DIM>
template<class ACTION >
void GEOGen::RestrictedVoronoiDiagram< DIM >::compute_surfacic_with_cnx_priority ( const ACTION &  action)
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).

Note
This function is less efficient than compute_surfacic_with_seeds_priority() but is required by some traversals that need to be done in that order.
Template Parameters
ACTIONneeds 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.

◆ compute_surfacic_with_seeds_priority()

template<index_t DIM>
template<class ACTION >
void GEOGen::RestrictedVoronoiDiagram< DIM >::compute_surfacic_with_seeds_priority ( const ACTION &  action)
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.

Template Parameters
ACTIONneeds 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.

◆ compute_volumetric()

template<index_t DIM>
template<class ACTION >
void GEOGen::RestrictedVoronoiDiagram< DIM >::compute_volumetric ( const ACTION &  action)
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.

Template Parameters
ACTIONneeds 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.

◆ compute_volumetric_with_cnx_priority()

template<index_t DIM>
template<class ACTION >
void GEOGen::RestrictedVoronoiDiagram< DIM >::compute_volumetric_with_cnx_priority ( const ACTION &  action)
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).

Note
This function is less efficient than compute_volumetric_with_seeds_priority() but is required by some traversals that need to be done in that order.
Template Parameters
ACTIONneeds 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.

◆ compute_volumetric_with_seeds_priority()

template<index_t DIM>
template<class ACTION >
void GEOGen::RestrictedVoronoiDiagram< DIM >::compute_volumetric_with_seeds_priority ( const ACTION &  action)
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.

Template Parameters
ACTIONneeds 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.

◆ connected_components_priority()

template<index_t DIM>
bool GEOGen::RestrictedVoronoiDiagram< DIM >::connected_components_priority ( ) const
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.

Return values
trueif connected components priority is used.
falseotherwise.

Definition at line 171 of file generic_RVD.h.

◆ current_facet()

template<index_t DIM>
index_t GEOGen::RestrictedVoronoiDiagram< DIM >::current_facet ( ) const
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.

◆ current_polygon()

template<index_t DIM>
const Polygon& GEOGen::RestrictedVoronoiDiagram< DIM >::current_polygon ( ) const
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.

◆ current_polyhedron()

template<index_t DIM>
const Polyhedron& GEOGen::RestrictedVoronoiDiagram< DIM >::current_polyhedron ( ) const
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.

◆ current_seed()

template<index_t DIM>
index_t GEOGen::RestrictedVoronoiDiagram< DIM >::current_seed ( ) const
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.

◆ current_tet()

template<index_t DIM>
index_t GEOGen::RestrictedVoronoiDiagram< DIM >::current_tet ( ) const
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.

◆ facet_seed_is_visited()

template<index_t DIM>
bool GEOGen::RestrictedVoronoiDiagram< DIM >::facet_seed_is_visited ( index_t  f,
index_t  s 
) const
inline

Tests whether a (facet,seed) couple was visited.

Parameters
[in]findex of the facet
[in]sindex of the seed

Definition at line 1821 of file generic_RVD.h.

◆ find_seed_near_facet()

template<index_t DIM>
index_t GEOGen::RestrictedVoronoiDiagram< DIM >::find_seed_near_facet ( index_t  f)
inlineprotected

Finds a seed near a given facet.

Parameters
[in]findex of the facet in the mesh
Returns
the index of a Voronoi seed such that there is a non-empty intersection between the Voronoi cell of the seed and facet f.

Definition at line 2033 of file generic_RVD.h.

◆ find_seed_near_point()

template<index_t DIM>
index_t GEOGen::RestrictedVoronoiDiagram< DIM >::find_seed_near_point ( const double *  p)
inlineprotected

Finds a seed near a given point.

Parameters
[in]ppointer to the coordinates of the point
Returns
the index of a Voronoi seed such that its Voronoi cell contains the point p.

Definition at line 2059 of file generic_RVD.h.

◆ find_seed_near_tet()

template<index_t DIM>
index_t GEOGen::RestrictedVoronoiDiagram< DIM >::find_seed_near_tet ( index_t  t)
inlineprotected

Finds a seed near a given tetrahedron.

Parameters
[in]tindex of the tetrahedron in the mesh
Returns
the index of a Voronoi seed such that there is a non-empty intersection between the Voronoi cell of the seed and tetrahedron t.

Definition at line 2047 of file generic_RVD.h.

◆ for_each_border_halfedge()

template<index_t DIM>
template<class BOACTION >
void GEOGen::RestrictedVoronoiDiagram< DIM >::for_each_border_halfedge ( const BOACTION &  action)
inline

Iterates on the halfedges on the borders of the restricted Voronoi cells that are on the boundary of the input mesh.

Parameters
[in]actionthe user action object
Template Parameters
BOACTIONneeds 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.

◆ for_each_halfedge()

template<index_t DIM>
template<class HEACTION >
void GEOGen::RestrictedVoronoiDiagram< DIM >::for_each_halfedge ( const HEACTION &  action)
inline

Iterates on the halfedges on the borders of the restricted Voronoi cells.

Parameters
[in]actionthe user action object
Template Parameters
HEACTIONneeds 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.

◆ for_each_polygon()

template<index_t DIM>
template<class ACTION >
void GEOGen::RestrictedVoronoiDiagram< DIM >::for_each_polygon ( const ACTION &  action)
inline

Iterates on the facets of this RVD.

Parameters
[in]actionthe user action object
Template Parameters
ACTIONneeds 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.

◆ for_each_polyhedron()

template<index_t DIM>
template<class ACTION >
void GEOGen::RestrictedVoronoiDiagram< DIM >::for_each_polyhedron ( const ACTION &  action)
inline

Iterates on the polyhedra of this RVD.

Parameters
[in]actionthe user action object
Template Parameters
ACTIONneeds 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.

◆ for_each_primal_tetrahedron()

template<index_t DIM>
template<class ACTION >
void GEOGen::RestrictedVoronoiDiagram< DIM >::for_each_primal_tetrahedron ( const ACTION &  action)
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.

Parameters
[in]actionthe user action object
Template Parameters
ACTIONneeds 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.

◆ for_each_primal_triangle()

template<index_t DIM>
template<class PRIMTRIACTION >
void GEOGen::RestrictedVoronoiDiagram< DIM >::for_each_primal_triangle ( const PRIMTRIACTION &  action)
inline

Iterates on the triangles of the Restricted Delaunay Triangulation.

Parameters
[in]actionthe user action object
Template Parameters
PRIMTRIACTIONneeds 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.

◆ for_each_tetrahedron()

template<index_t DIM>
template<class ACTION >
void GEOGen::RestrictedVoronoiDiagram< DIM >::for_each_tetrahedron ( const ACTION &  action)
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.

Parameters
[in]actionthe user action object
Template Parameters
ACTIONneeds 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:
  • v is the index of the current Voronoi cell (or Delaunay vertex)
  • v_adj is the index of the Voronoi cell adjacent to t accros facet (v1, v2, v3) or -1 if it does not exists adjacent to v or -1 if current face is a tetrahedron facet
  • t is the index of the current tetrahedron
  • t_adj is the index of the tetrahedron adjacent to t accros facet (v1, v2, v3) or -1 if it does not exists
  • v0,v1,v2 and v3 are the four vertices of tetrahedron.

Definition at line 1259 of file generic_RVD.h.

◆ for_each_triangle()

template<index_t DIM>
template<class TRIACTION >
void GEOGen::RestrictedVoronoiDiagram< DIM >::for_each_triangle ( const TRIACTION &  action)
inline

Iterates on the facets of this RVD, triangulated on the fly.

Parameters
[in]actionthe user action object
Template Parameters
TRIACTIONneeds 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.

◆ for_each_volumetric_integration_simplex()

template<index_t DIM>
template<class ACTION >
void GEOGen::RestrictedVoronoiDiagram< DIM >::for_each_volumetric_integration_simplex ( const ACTION &  action,
bool  visit_inner_tets = false,
bool  coherent_triangles = false 
)
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.

Parameters
[in]actionthe user action object
[in]visit_inner_tetsif 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_trianglesif 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.
Template Parameters
ACTIONneeds 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:
  • v is the index of the current Voronoi cell (or Delaunay vertex)
  • v_adj is the index of the Voronoi cell adjacent to t accros facet (v1, v2, v3) or -1 if it does not exists adjacent to v or -1 if current face is a tetrahedron facet
  • t is the index of the current tetrahedron
  • t_adj is the index of the tetrahedron adjacent to t accros facet (v1, v2, v3) or -1 if it does not exists
  • v1,v2 and v3 are the three vertices of the facet on the border of the restricted Voronoi cell.

Definition at line 1222 of file generic_RVD.h.

◆ get_facet_seed_connected_component()

template<index_t DIM>
signed_index_t GEOGen::RestrictedVoronoiDiagram< DIM >::get_facet_seed_connected_component ( index_t  f,
index_t  s 
) const
inline

Gets the index of the connected component associated with a (facet,seed).

Parameters
[in]findex of the facet
[in]sindex of the seed
Returns
the index of the connected component or -1 if the (f, s) couple was not visited already

Definition at line 1834 of file generic_RVD.h.

◆ get_neighbors()

template<index_t DIM>
void GEOGen::RestrictedVoronoiDiagram< DIM >::get_neighbors ( index_t  v)
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.

◆ init_get_neighbors()

template<index_t DIM>
void GEOGen::RestrictedVoronoiDiagram< DIM >::init_get_neighbors ( )
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.

◆ intersect_cell_cell()

template<index_t DIM>
void GEOGen::RestrictedVoronoiDiagram< DIM >::intersect_cell_cell ( index_t  seed,
Polyhedron C 
)
inline

Computes the intersection between a Voronoi cell and a cell with radius of security or plain mode.

Parameters
[in]seedthe index of the seed that defines the Voronoi cell
[in,out]Cthe cell to be clipped

Definition at line 2285 of file generic_RVD.h.

◆ intersect_cell_facet()

template<index_t DIM>
Polygon* GEOGen::RestrictedVoronoiDiagram< DIM >::intersect_cell_facet ( index_t  seed,
Polygon F 
)
inlineprotected

Computes the intersection between the Voronoi cell of a seed and a facet.

Parameters
[in]seedthe index of the seed
[in]Fthe facet represented as a Polygon

The result is provided in current_polygon_

Definition at line 2123 of file generic_RVD.h.

◆ nb_facets_in_range()

template<index_t DIM>
index_t GEOGen::RestrictedVoronoiDiagram< DIM >::nb_facets_in_range ( ) const
inline

Gets the number of facets in the current range.

See also
set_facets_range()

Definition at line 250 of file generic_RVD.h.

◆ nb_tetrahedra_in_range()

template<index_t DIM>
index_t GEOGen::RestrictedVoronoiDiagram< DIM >::nb_tetrahedra_in_range ( ) const
inline

Gets the number of tetrahedra in the current range.

See also
set_tetrahedra_range()

Definition at line 258 of file generic_RVD.h.

◆ point_allocator()

template<index_t DIM>
PointAllocator* GEOGen::RestrictedVoronoiDiagram< DIM >::point_allocator ( )
inline

Gets the PointAllocator.

Returns
a pointer to the PointAllocator, used to create the new vertices generated by intersections.

Definition at line 375 of file generic_RVD.h.

◆ set_connected_components_priority()

template<index_t DIM>
void GEOGen::RestrictedVoronoiDiagram< DIM >::set_connected_components_priority ( bool  x)
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.

◆ set_exact_predicates()

template<index_t DIM>
void GEOGen::RestrictedVoronoiDiagram< DIM >::set_exact_predicates ( bool  x)
inline

Specifies whether exact predicates should be used.

If exact predicates are used, symbolic mode is ensured.

Parameters
[in]xif set, exact predicates are used.

Definition at line 338 of file generic_RVD.h.

◆ set_facets_range()

template<index_t DIM>
void GEOGen::RestrictedVoronoiDiagram< DIM >::set_facets_range ( index_t  facets_begin,
index_t  facets_end 
)
inline

Sets the facets range.

Computations can be restricted to a contiguous facet range.

Parameters
[in]facets_beginfirst facet in the range.
[in]facets_endone position past the last facet in the range.

Definition at line 226 of file generic_RVD.h.

◆ set_symbolic()

template<index_t DIM>
void GEOGen::RestrictedVoronoiDiagram< DIM >::set_symbolic ( bool  x)
inline

Sets symbolic mode.

If exact mode is active, symbolic mode is enforced.

Parameters
[in]xif set, the symbolic representation of the intersections are computed.

Definition at line 318 of file generic_RVD.h.

◆ set_tetrahedra_range()

template<index_t DIM>
void GEOGen::RestrictedVoronoiDiagram< DIM >::set_tetrahedra_range ( index_t  tets_begin,
index_t  tets_end 
)
inline

Sets the tetrahedra range.

Computations can be restricted to a contiguous tetrahedra range.

Parameters
[in]tets_beginfirst tetrahedron in the range.
[in]tets_endone position past the last tetrahedron in the range.

Definition at line 240 of file generic_RVD.h.

◆ swap_polygons()

template<index_t DIM>
void GEOGen::RestrictedVoronoiDiagram< DIM >::swap_polygons ( Polygon *&  ping,
Polygon *&  pong 
)
inlineprotected

Swaps two pointers between two polygons.

Used by re-entrant Sutherlang-Hogdman clipping.

Definition at line 2105 of file generic_RVD.h.


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