Graphite  Version 3
An experimental 3D geometry processing program
GEOGen::RestrictedVoronoiDiagram< DIM >::VolumetricIntegrationSimplexAction< ACTION > Class Template Reference

Adapter class used internally to implement for_each_volumetric_integration_simplex() More...

#include <geogram/voronoi/generic_RVD.h>

Public Member Functions

 VolumetricIntegrationSimplexAction (const ACTION &do_it, bool visit_inner_tets=false, bool coherent_triangles=false)
 Creates a new VolumetricIntegrationSimplexAction that wraps a user ACTION instance. More...
 
void operator() (index_t v, index_t t, const Polyhedron &C) const
 Callback called for each polyhedron. More...
 
void move_to_first_corner_of_facet (const Polyhedron &C, Polyhedron::Corner &c, index_t center_vertex_id) const
 Finds the first corner of a facet in a Polyhedron. More...
 

Static Public Member Functions

static bool symbolic_compare (const Vertex &p1, const Vertex &p2, index_t center_vertex_id)
 Compares the symbolic information of two vertices in such a way that a global order is defined. More...
 

Protected Attributes

const ACTION & do_it_
 
bool visit_inner_tets_
 
bool coherent_triangles_
 

Detailed Description

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

Adapter class used internally to implement for_each_volumetric_integration_simplex()

Overrides constness checks, to allow using temporaries as argument of for_each_xxx().

Template Parameters
ACTIONthe user action class. It needs to implement: operator()(index_t v, signed_index_t v_adj, index_t t, signed_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 689 of file generic_RVD.h.

Constructor & Destructor Documentation

◆ VolumetricIntegrationSimplexAction()

template<index_t DIM>
template<class ACTION >
GEOGen::RestrictedVoronoiDiagram< DIM >::VolumetricIntegrationSimplexAction< ACTION >::VolumetricIntegrationSimplexAction ( const ACTION &  do_it,
bool  visit_inner_tets = false,
bool  coherent_triangles = false 
)
inline

Creates a new VolumetricIntegrationSimplexAction that wraps a user ACTION instance.

Parameters
[in]do_itthe user ACTION instance
[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.

Definition at line 707 of file generic_RVD.h.

Member Function Documentation

◆ move_to_first_corner_of_facet()

template<index_t DIM>
template<class ACTION >
void GEOGen::RestrictedVoronoiDiagram< DIM >::VolumetricIntegrationSimplexAction< ACTION >::move_to_first_corner_of_facet ( const Polyhedron C,
Polyhedron::Corner c,
index_t  center_vertex_id 
) const
inline

Finds the first corner of a facet in a Polyhedron.

This function is used to ensure that a facet is triangulated coherently when seen from two different volumetric cells, by generating a fan of triangles that radiates from the first corner. The global order used to find the first corner is defined by the function symbolic_compare().

Parameters
[in]Cthe Polyhedron
[in,out]ca corner of the facet, replaced by the first corner of the facet on exit.
[in]center_vertex_idindex of the current Voronoi seed (needed to determine the full symbolic information in the vertices).

Definition at line 808 of file generic_RVD.h.

◆ operator()()

template<index_t DIM>
template<class ACTION >
void GEOGen::RestrictedVoronoiDiagram< DIM >::VolumetricIntegrationSimplexAction< ACTION >::operator() ( index_t  v,
index_t  t,
const Polyhedron C 
) const
inline

Callback called for each polyhedron.

Routes the callback to the wrapped user action class.

Parameters
[in]vindex of current Delaunay seed
[in]tindex of current mesh tetrahedron
[in]Cintersection between current mesh tetrahedron and the Voronoi cell of v

Definition at line 726 of file generic_RVD.h.

◆ symbolic_compare()

template<index_t DIM>
template<class ACTION >
static bool GEOGen::RestrictedVoronoiDiagram< DIM >::VolumetricIntegrationSimplexAction< ACTION >::symbolic_compare ( const Vertex p1,
const Vertex p2,
index_t  center_vertex_id 
)
inlinestatic

Compares the symbolic information of two vertices in such a way that a global order is defined.

This function is used to ensure that a facet is triangulated coherently when seen from two different volumetric cells (it uniquely determines the "first" vertex).

Parameters
[in]p1first vertex to compare
[in]p2second vertex to compare
[in]center_vertex_idindex of the current Voronoi seed (needed to determine the full symbolic information in p1 and p2).
Returns
true if p1 is before p2 in the global order, false otherwise.

Definition at line 840 of file generic_RVD.h.


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