40 #ifndef GEOGRAM_VORONOI_GENERIC_RVD
41 #define GEOGRAM_VORONOI_GENERIC_RVD
77 template <index_t DIM>
133 facets_begin_ = UNSPECIFIED_RANGE;
134 facets_end_ = UNSPECIFIED_RANGE;
135 tets_begin_ = UNSPECIFIED_RANGE;
136 tets_end_ = UNSPECIFIED_RANGE;
137 connected_components_priority_ =
false;
138 facet_seed_marking_ =
nullptr;
139 connected_component_changed_ =
false;
140 current_connected_component_ = 0;
144 current_polygon_ =
nullptr;
146 current_polyhedron_ =
nullptr;
157 connected_components_priority_ = x;
172 return connected_components_priority_;
228 facets_begin_ = facets_begin;
229 facets_end_ = facets_end;
242 tets_begin_ = tets_begin;
243 tets_end_ = tets_end;
251 return facets_end_ - facets_begin_;
259 return tets_end_ - tets_begin_;
268 return current_facet_;
277 return current_seed_;
288 return *current_polygon_;
309 return *current_polyhedron_;
376 return &intersections_;
391 template <
class ACTION>
416 const_cast<ACTION&
> ( do_it_)(v, P);
420 const ACTION& do_it_;
430 template <
class ACTION>
459 const_cast<ACTION&
> (do_it_)(
466 const ACTION& do_it_;
476 template <
class ACTION>
507 const_cast<ACTION&
> (do_it_)(
515 const ACTION& do_it_;
525 template <
class ACTION>
557 const_cast<ACTION&
> (do_it_)(
566 const ACTION& do_it_;
575 template <
class ACTION>
610 if(iv1 < iv2 && iv1 < iv3) {
611 const_cast<ACTION&
> (do_it_)(iv1, iv2, iv3);
618 const ACTION& do_it_;
634 template <
class ACTION>
659 const_cast<ACTION&
> ( do_it_)(v, t, C);
663 const ACTION& do_it_;
688 template <
class ACTION>
709 bool visit_inner_tets =
false,
710 bool coherent_triangles =
false
713 visit_inner_tets_(visit_inner_tets),
714 coherent_triangles_(coherent_triangles)
747 if(!visit_inner_tets_) {
750 t_adj = -adjacent - 1;
751 }
else if(adjacent > 0) {
754 v_adj = adjacent - 1;
767 if(coherent_triangles_) {
783 const_cast<ACTION&
> (do_it_)(
784 v, v_adj, t, t_adj, v1, v2, v3
823 }
while(cur != first);
855 const ACTION& do_it_;
856 bool visit_inner_tets_;
857 bool coherent_triangles_;
882 template <
class ACTION>
914 const Vertex* v0 =
nullptr;
916 for(t0 = 0; t0 < C.
max_t(); ++t0) {
942 t_adj = -adjacent - 1;
943 }
else if(adjacent > 0) {
946 v_adj = adjacent - 1;
974 const_cast<ACTION&
> (do_it_)(
975 v, v_adj, t, t_adj, * v0, v1, v2, v3
1005 }
while(cur != first);
1010 const ACTION& do_it_;
1020 template <
class ACTION>
1054 if(v < v1 && v < v2 && v < v3) {
1055 const_cast<ACTION&
> (do_it_)(v, v1, v2, v3);
1063 const ACTION& do_it_;
1083 template <
class ACTION>
1085 const ACTION& action
1087 this->
template compute_surfacic<PolygonAction<ACTION> >(
1100 template <
class TRIACTION>
1102 const TRIACTION& action
1104 this->
template compute_surfacic<TriangleAction<TRIACTION> >(
1118 template <
class HEACTION>
1120 const HEACTION& action
1122 this->
template compute_surfacic<HalfedgeAction<HEACTION> >(
1136 template <
class BOACTION>
1138 const BOACTION& action
1140 this->
template compute_surfacic<BorderHalfedgeAction<BOACTION> >(
1154 template <
class PRIMTRIACTION>
1156 const PRIMTRIACTION& action
1178 template <
class ACTION>
1180 const ACTION& action
1182 this->
template compute_volumetric<PolyhedronAction<ACTION> >(
1221 template <
class ACTION>
1223 const ACTION& action,
1224 bool visit_inner_tets =
false,
bool coherent_triangles =
false
1230 action, visit_inner_tets, coherent_triangles
1258 template <
class ACTION>
1260 const ACTION& action
1262 this->
template compute_volumetric<TetrahedronAction<ACTION> >(
1281 template <
class ACTION>
1283 const ACTION& action
1287 this->
template compute_volumetric<PrimalTetrahedronAction<ACTION> >(
1313 template <
class ACTION>
1315 if(connected_components_priority_) {
1316 this->
template compute_surfacic_with_cnx_priority<ACTION>(
1320 this->
template compute_surfacic_with_seeds_priority<ACTION>(
1338 template <
class ACTION>
1340 const ACTION& action
1343 facets_begin_ == UNSPECIFIED_RANGE &&
1344 facets_end_ == UNSPECIFIED_RANGE
1347 facets_end_ = mesh_->facets.
nb();
1349 current_polygon_ =
nullptr;
1368 for(
index_t f = facets_begin_; f < facets_end_; f++) {
1369 if(!facet_is_marked[f-facets_begin_]) {
1371 facet_is_marked[f-facets_begin_] =
true;
1372 adjacent_facets.push(
1375 while(!adjacent_facets.empty()) {
1376 current_facet_ = adjacent_facets.top().f;
1377 current_seed_ = adjacent_facets.top().seed;
1378 adjacent_facets.pop();
1384 mesh_, current_facet_, symbolic_, vertex_weight
1391 seed_stamp[current_seed_] = current_facet_;
1392 adjacent_seeds.push(current_seed_);
1394 while(!adjacent_seeds.empty()) {
1395 current_seed_ = adjacent_seeds.top();
1396 adjacent_seeds.pop();
1417 if(!facet_is_marked[
1418 index_t(neigh_f)-facets_begin_
1421 index_t(neigh_f)-facets_begin_
1423 adjacent_facets.push(
1434 seed_stamp[neigh_s] != current_facet_
1436 seed_stamp[neigh_s] = current_facet_;
1437 adjacent_seeds.push(
index_t(neigh_s));
1445 current_polygon_ =
nullptr;
1461 template <
class ACTION>
1463 const ACTION& action
1465 if(connected_components_priority_) {
1466 this->
template compute_volumetric_with_cnx_priority<ACTION>(
1470 this->
template compute_volumetric_with_seeds_priority<ACTION>(
1488 template <
class ACTION>
1490 const ACTION& action
1493 tets_begin_ == UNSPECIFIED_RANGE &&
1494 tets_end_ == UNSPECIFIED_RANGE
1497 tets_end_ = mesh_->cells.
nb();
1500 geo_assert(tets_begin_ != UNSPECIFIED_RANGE);
1517 current_polyhedron_ = &C;
1522 for(
index_t t = tets_begin_; t < tets_end_; ++t) {
1523 if(!tet_is_marked[t-tets_begin_]) {
1525 tet_is_marked[t-tets_begin_] =
true;
1529 while(!adjacent_tets.empty()) {
1530 current_tet_ = adjacent_tets.top().f;
1531 current_seed_ = adjacent_tets.top().seed;
1532 adjacent_tets.pop();
1544 seed_stamp[current_seed_] = current_tet_;
1545 adjacent_seeds.push(current_seed_);
1547 while(!adjacent_seeds.empty()) {
1548 current_seed_ = adjacent_seeds.top();
1549 adjacent_seeds.pop();
1552 mesh_, current_tet_, symbolic_, vertex_weight
1560 current_seed_, current_tet_,
1585 if(seed_stamp[neigh_s] != current_tet_) {
1586 seed_stamp[neigh_s] = current_tet_;
1587 adjacent_seeds.push(neigh_s);
1622 current_polyhedron_ =
nullptr;
1645 template <
class ACTION>
1647 const ACTION& action
1651 tets_begin_ == UNSPECIFIED_RANGE &&
1652 tets_end_ == UNSPECIFIED_RANGE
1655 tets_end_ = mesh_->cells.
nb();
1658 geo_assert(tets_begin_ != UNSPECIFIED_RANGE);
1661 current_polyhedron_ =
nullptr;
1664 std::deque<TetSeed> adjacent_seeds;
1665 std::stack<index_t> adjacent_tets;
1669 tets_end_ - tets_begin_, NO_STAMP
1679 facet_seed_marking_ = &visited;
1681 current_polyhedron_ = &C;
1683 current_connected_component_ = 0;
1695 for(
index_t t = tets_begin_; t < tets_end_; ++t) {
1696 if(tet_stamp[t - tets_begin_] == NO_STAMP) {
1700 adjacent_seeds.push_back(
1701 TetSeed(current_tet_, current_seed_)
1705 while(!adjacent_seeds.empty()) {
1707 current_tet_ = adjacent_seeds.front().f;
1708 current_seed_ = adjacent_seeds.front().seed;
1709 adjacent_seeds.pop_front();
1711 tet_stamp[current_tet_ - tets_begin_] ==
1717 if(visited.
is_marked(current_tet_, current_seed_)) {
1721 connected_component_changed_ =
true;
1722 adjacent_tets.push(current_tet_);
1723 tet_stamp[current_tet_ - tets_begin_] =
1725 while(!adjacent_tets.empty()) {
1726 current_tet_ = adjacent_tets.top();
1727 adjacent_tets.pop();
1734 mesh_, current_tet_, symbolic_, vertex_weight
1748 connected_component_changed_ =
false;
1750 bool touches_RVC_border =
false;
1782 tet_stamp[neigh_t - tets_begin_] !=
1785 tet_stamp[neigh_t - tets_begin_] =
1787 adjacent_tets.push(neigh_t);
1792 touches_RVC_border =
true;
1793 TetSeed ts(current_tet_, neigh_s);
1795 adjacent_seeds.push_back(ts);
1800 if(touches_RVC_border) {
1802 TetSeed(current_tet_, current_seed_),
1803 current_connected_component_
1807 ++current_connected_component_;
1811 facet_seed_marking_ =
nullptr;
1847 return connected_component_changed_;
1854 return current_connected_component_;
1878 template <
class ACTION>
1880 const ACTION& action
1884 facets_begin_ == UNSPECIFIED_RANGE &&
1885 facets_end_ == UNSPECIFIED_RANGE
1888 facets_end_ = mesh_->facets.
nb();
1891 current_polygon_ =
nullptr;
1894 std::deque<FacetSeed> adjacent_seeds;
1895 std::stack<index_t> adjacent_facets;
1899 facets_end_ - facets_begin_, NO_STAMP
1903 facets_end_ - facets_begin_, delaunay_->
nb_vertices()
1906 facet_seed_marking_ = &visited;
1908 current_connected_component_ = 0;
1909 index_t F_index = facets_end_ + 1;
1920 for(
index_t f = facets_begin_; f < facets_end_; ++f) {
1921 if(facet_stamp[f - facets_begin_] == NO_STAMP) {
1925 adjacent_seeds.push_back(
1926 FacetSeed(current_facet_, current_seed_)
1931 !adjacent_seeds.empty()
1933 current_facet_ = adjacent_seeds.front().f;
1934 current_seed_ = adjacent_seeds.front().seed;
1935 adjacent_seeds.pop_front();
1937 facet_stamp[current_facet_ - facets_begin_] ==
1943 if(visited.
is_marked(current_facet_, current_seed_)) {
1947 connected_component_changed_ =
true;
1948 adjacent_facets.push(current_facet_);
1949 facet_stamp[current_facet_ - facets_begin_] =
1951 while(!adjacent_facets.empty()) {
1952 current_facet_ = adjacent_facets.top();
1953 adjacent_facets.pop();
1958 if(F_index != current_facet_) {
1960 mesh_, current_facet_, symbolic_,
1963 F_index = current_facet_;
1972 connected_component_changed_ =
false;
1974 bool touches_RVC_border =
false;
1993 facet_stamp[neigh_f - facets_begin_] !=
1996 facet_stamp[neigh_f - facets_begin_] =
1998 adjacent_facets.push(neigh_f);
2003 touches_RVC_border =
true;
2005 current_facet_,
index_t(neigh_s)
2008 adjacent_seeds.push_back(fs);
2012 if(touches_RVC_border) {
2014 FacetSeed(current_facet_, current_seed_),
2015 current_connected_component_
2019 ++current_connected_component_;
2023 facet_seed_marking_ =
nullptr;
2034 const double* p = mesh_->vertices.
point_ptr(
2035 mesh_->facets.
vertex(f,0)
2048 index_t v = mesh_->cells.tet_vertex(t, 0);
2049 const double* p = mesh_->vertices.
point_ptr(v);
2068 if(exact_ && delaunay_nn_ !=
nullptr) {
2071 double neighbors_sq_dist[10];
2077 nb, p, neighbors, neighbors_sq_dist
2079 index_t nearest = neighbors[0];
2080 double min_d = neighbors_sq_dist[0];
2081 for(
index_t i = 1; i < nb; ++i) {
2082 if(neighbors_sq_dist[i] != min_d) {
2085 if(neighbors[i] < nearest) {
2086 nearest = neighbors[i];
2106 if(ping != &P1 && ping != &P2) {
2112 std::swap(ping, pong);
2124 intersections_.
clear();
2132 if(delaunay_nn_ !=
nullptr) {
2162 const double* geo_restrict pi = delaunay_->
vertex_ptr(i);
2166 index_t prev_nb_neighbors = 0;
2167 neighbors_.resize(0);
2171 if(neighbors_.
size() == 0) {
2174 if(prev_nb_neighbors == neighbors_.
size()) {
2178 for(; jj < neighbors_.
size(); jj++) {
2182 geo_decl_aligned(
double dik);
2183 const double* geo_restrict pk = ping->
vertex(k).
point();
2186 R2 = std::max(R2, dik);
2188 geo_decl_aligned(
double dij);
2189 const double* geo_restrict pj = delaunay_->
vertex_ptr(j);
2195 if(dij > 4.1 * R2) {
2207 prev_nb_neighbors = nb_neighbors;
2209 if(nb_neighbors > 8) {
2210 nb_neighbors += nb_neighbors / 8;
2215 nb_neighbors = std::min(
2242 for(
index_t jj = 0; jj < neighbors_.
size(); jj++) {
2268 pong, intersections_, mesh_, delaunay_, i, j, exact_, symbolic_
2287 if(delaunay_nn_ !=
nullptr) {
2307 const double* geo_restrict pi = delaunay_->
vertex_ptr(seed);
2311 index_t prev_nb_neighbors = 0;
2312 neighbors_.resize(0);
2317 if(neighbors_.
size() == 0) {
2320 if(prev_nb_neighbors == neighbors_.
size()) {
2324 for(; jj < neighbors_.
size(); jj++) {
2331 geo_decl_aligned(
double dik);
2332 const double* geo_restrict pk =
2336 R2 = std::max(R2, dik);
2338 geo_decl_aligned(
double dij);
2339 const double* geo_restrict pj = delaunay_->
vertex_ptr(j);
2345 if(dij > 4.1 * R2) {
2356 prev_nb_neighbors = nb_neighbors;
2358 if(nb_neighbors > 8) {
2359 nb_neighbors += nb_neighbors / 8;
2364 nb_neighbors = std::min(
2382 if(neighbors_.
size() == 0) {
2385 for(
index_t jj = 0; jj < neighbors_.
size(); jj++) {
2400 mesh_, delaunay_, i, j, exact_, symbolic_
2436 if(stamp_.
size() == 0) {
2444 neighbors_.resize(0);
2451 if(stamp_[w] != cur_stamp_) {
2452 stamp_[w] = cur_stamp_;
2453 neighbors_.push_back(w);
2497 bool connected_components_priority_;
2499 bool connected_component_changed_;
2500 index_t current_connected_component_;
2511 thisclass& operator= (
const thisclass& rhs);
Declaration of base types for implementing user-defined code that queries the cells of a restricted V...
A function to suppress unused parameters compilation warnings.
#define geo_assert(x)
Verifies that a condition is met.
#define geo_debug_assert(x)
Verifies that a condition is met.
Generic mechanism for attributes.
A Corner corresponds to a vertex seen from a triangle.
Computes the intersection between a set of halfspaces.
index_t max_v() const
Gets the maximum valid vertex index plus one.
void initialize_from_mesh_tetrahedron(const Mesh *mesh, index_t t, bool symbolic, const GEO::Attribute< double > &vertex_weight)
Assigns a mesh tetrahedron to this ConvexCell.
signed_index_t clip_by_plane(const Mesh *mesh, const Delaunay *delaunay, index_t i, index_t j, bool exact, bool symbolic)
Clips this ConvexCell with a plane.
signed_index_t vertex_id(index_t v) const
Gets the id of a vertex.
void clear()
Clears this ConvexCell.
const GEOGen::Vertex & triangle_dual(index_t t) const
Gets the dual vertex that corresponds to a triangle.
index_t find_triangle_vertex(index_t t, index_t v) const
Finds the local index of a triangle vertex.
bool triangle_is_used(index_t t) const
Tests whether a given triangle is used.
void move_to_next_around_vertex(Corner &c) const
Replaces a corner by the next corner obtained by turing around the vertex.
signed_index_t vertex_triangle(index_t v) const
Gets one of the triangles incident to a vertex.
index_t max_t() const
Gets the maximum valid triangle index plus one.
Stores associations between (facet,seed) pairs and the index of a connected component.
signed_index_t get_connected_component(const FacetSeed &fs) const
Gets the index of the connected component associated with a given FacetSeed.
void mark(const FacetSeed &fs, index_t conn_comp)
Marks a FacetSeed and sets the associated connected component index.
bool is_marked(index_t facet, index_t seed) const
Tests whether a given facet,seed couple is marked.
An allocator for points that are created from intersections in GenericVoronoiDiagram.
void clear()
Clears this PointAllocator.
Internal representation of polygons for GenericVoronoiDiagram.
void initialize_from_mesh_facet(const Mesh *mesh, index_t f, bool symbolic, const GEO::Attribute< double > &vertex_weight)
Assigns a mesh facet to this Polygon.
index_t nb_vertices() const
Gets the number of vertices.
index_t next_vertex(index_t i) const
Gets the index of the successor of a Vertex.
const Vertex & vertex(index_t i) const
Gets a vertex by index.
void clip_by_plane(Polygon &target, PointAllocator &target_intersections, const Mesh *mesh, const Delaunay *delaunay, index_t i, index_t j, bool exact, bool symbolic)
Clips a polygon with a plane.
Adapter class used internally to implement for_each_border_halfedge().
BorderHalfedgeAction(const ACTION &do_it)
Creates a new BorderHalfedgeAction that wraps a user ACTION instance.
void operator()(index_t v, index_t f, const Polygon &P) const
Callback called for each integration simplex.
Adapter class used internally to implement for_each_halfedge().
void operator()(index_t v, index_t f, const Polygon &P) const
Callback called for each integration simplex.
HalfedgeAction(const ACTION &do_it)
Creates a new HalfedgeAction that wraps a user ACTION instance.
Adapter class used internally to implement for_each_polygon()
void operator()(index_t v, index_t f, const Polygon &P) const
Callback called for each polygon.
PolygonAction(const ACTION &do_it)
Creates a new PolygonAction around a user ACTION instance.
Adapter class used internally to implement for_each_polyhedron()
void operator()(index_t v, index_t t, const Polyhedron &C) const
Callback called for each polyhedron.
PolyhedronAction(const ACTION &do_it)
Creates a new PolyhedronAction that wraps a user ACTION instance.
Adapter class used internally to implement for_each_primal_tetrahedron()
PrimalTetrahedronAction(const ACTION &do_it)
Constructs a new PrimalTetrahedronAction.
void operator()(index_t v, index_t t, const Polyhedron &C) const
Callback called for each polyhedron.
Adapter class used internally to implement for_each_primal_triangle()
void operator()(index_t iv1, index_t f, const Polygon &P) const
Callback called for each primal triangle.
PrimalTriangleAction(const ACTION &do_it)
Creates a new PrimalTriangleAction that wraps a user ACTION instance.
Adapter class used internally to implement for_each_tetrahedron()
TetrahedronAction(const ACTION &do_it)
Creates a new TetrahedronAction that wraps a user ACTION instance.
bool facet_is_incident_to_vertex(const Polyhedron &C, Polyhedron::Corner &c, index_t t) const
Tests whether a Polyhedron facet is incident to a vertex.
void operator()(index_t v, index_t t, const Polyhedron &C) const
Callback called for each polyhedron.
Adapter class used internally to implement for_each_triangle().
void operator()(index_t v, index_t f, const Polygon &P) const
Callback called for each integration simplex.
TriangleAction(const ACTION &do_it)
Creates a new TriangleAction that wraps a user ACTION instance.
Adapter class used internally to implement for_each_volumetric_integration_simplex()
VolumetricIntegrationSimplexAction(const ACTION &do_it, bool visit_inner_tets=false, bool coherent_triangles=false)
Creates a new VolumetricIntegrationSimplexAction that wraps a user ACTION instance.
void operator()(index_t v, index_t t, const Polyhedron &C) const
Callback called for each polyhedron.
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.
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.
Computes the intersection between a surface (Mesh) and a Voronoi diagram (dual of a Delaunay).
GEOGen::Vertex Vertex
Internal representation of vertices.
index_t nb_facets_in_range() const
Gets the number of facets in the current range.
void clip_by_cell(index_t seed, Polyhedron &C)
Computes the intersection between a Voronoi cell and a cell in plain mode.
void set_exact_predicates(bool x)
Specifies whether exact predicates should be used.
GEO::Mesh * mesh()
Gets the input mesh.
void for_each_primal_triangle(const PRIMTRIACTION &action)
Iterates on the triangles of the Restricted Delaunay Triangulation.
static coord_index_t dimension()
Gets the dimension.
index_t find_seed_near_point(const double *p)
Finds a seed near a given point.
GEOGen::ConvexCell Polyhedron
Internal representation of volumetric cells.
GEOGen::Polygon Polygon
Internal representation of polygons.
void set_mesh(GEO::Mesh *mesh)
Sets the input mesh.
void for_each_polyhedron(const ACTION &action)
Iterates on the polyhedra of this RVD.
Polygon * intersect_cell_facet(index_t seed, Polygon &F)
Computes the intersection between the Voronoi cell of a seed and a facet.
void compute_volumetric_with_cnx_priority(const ACTION &action)
Low-level API of Restricted Voronoi Diagram traversal with connected components priority.
void for_each_triangle(const TRIACTION &action)
Iterates on the facets of this RVD, triangulated on the fly.
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.
void for_each_tetrahedron(const ACTION &action)
Iterates on the polyhedra of this RVD decomposed on the fly into tetrahedra.
index_t current_seed() const
Gets the index of the Delaunay vertex currently processed.
void clip_by_cell(index_t i, Polygon *&ping, Polygon *&pong)
Computes the intersection between a Voronoi cell and a polygon.
bool symbolic() const
Tests whether symbolic mode is active.
const GEO::Mesh * mesh() const
Gets the input mesh.
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 ...
void for_each_halfedge(const HEACTION &action)
Iterates on the halfedges on the borders of the restricted Voronoi cells.
const Polyhedron & current_polyhedron() const
Gets the current cell.
void compute_volumetric(const ACTION &action)
Low-level API of Restricted Voronoi Diagram traversal .
index_t current_connected_component() const
Gets the index of the current connected component.
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.
void for_each_polygon(const ACTION &action)
Iterates on the facets of this RVD.
index_t current_facet() const
Gets the index of the mesh facet currently processed.
const Polygon & current_polygon() const
Gets the current polygon.
index_t current_tet() const
Gets the undex of the mesh tetrahedron currently processed.
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.
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'.
void get_neighbors(index_t v)
Caches the neighbors of a Delaunay vertex.
GEOGen::PointAllocator PointAllocator
Used to allocate the generated points.
void compute_volumetric_with_seeds_priority(const ACTION &action)
Low-level API of Restricted Voronoi Diagram traversal with seeds priority in volumetric mode.
void set_check_SR(bool x)
Specifies whether radius of security should be enforced.
PointAllocator * point_allocator()
Gets the PointAllocator.
RestrictedVoronoiDiagram(Delaunay *delaunay, GEO::Mesh *mesh)
Constructs a new RestrictedVoronoiDiagram.
bool connected_component_changed() const
Tests whether the current connected component changed.
void set_delaunay(Delaunay *delaunay)
Sets the Delaunay triangulation.
bool connected_components_priority() const
Tests whether connected components priority is set.
bool exact_predicates() const
Tests whether exact predicates are used.
Delaunay * delaunay()
Gets the Delaunay triangulation.
index_t find_seed_near_facet(index_t f)
Finds a seed near a given facet.
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).
void swap_polygons(Polygon *&ping, Polygon *&pong)
Swaps two pointers between two polygons.
void set_connected_components_priority(bool x)
Sets traveral priority.
void for_each_primal_tetrahedron(const ACTION &action)
Iterates on the primal tetrahedra of this RVD.
void set_symbolic(bool x)
Sets symbolic mode.
void init_get_neighbors()
Creates the data structure for optimized get_neighbors() function.
bool facet_seed_is_visited(index_t f, index_t s) const
Tests whether a (facet,seed) couple was visited.
void compute_surfacic(const ACTION &action)
Low-level API of Restricted Voronoi Diagram traversal.
void set_facets_range(index_t facets_begin, index_t facets_end)
Sets the facets range.
void compute_surfacic_with_seeds_priority(const ACTION &action)
Low-level API of Restricted Voronoi Diagram traversal with seeds priority in surfacic mode.
index_t nb_tetrahedra_in_range() const
Gets the number of tetrahedra in the current range.
bool check_SR() const
Tests whether radius of security is enforced.
void compute_surfacic_with_cnx_priority(const ACTION &action)
Low-level API of Restricted Voronoi Diagram traversal with connected components priority.
index_t find_seed_near_tet(index_t t)
Finds a seed near a given tetrahedron.
void clip_by_plane(Polygon &ping, Polygon &pong, index_t i, index_t j)
Computes the intersection between a polygon and a half-space.
const Delaunay * delaunay() const
Gets the Delaunay triangulation.
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.
void set_tetrahedra_range(index_t tets_begin, index_t tets_end)
Sets the tetrahedra range.
A set of three integers that encodes the equation of a vertex in GenericVoronoiDiagram.
index_t nb_bisectors() const
Gets the number of bisectors in the symbolic representation.
index_t bisector(signed_index_t i) const
Gets a bisector.
Internal representation of vertices in GenericVoronoiDiagram.
const SymbolicVertex & sym() const
Gets the symbolic representation.
signed_index_t adjacent_seed() const
Gets the adjacent seed.
const double * point() const
Gets the geometric location at this Vertex.
signed_index_t adjacent_facet() const
Gets the adjacent facet.
bool check_flag(EdgeFlag f) const
Tests an EdgeFlag in this Vertex.
bool bind_if_is_defined(AttributesManager &manager, const std::string &name)
Binds this Attribute to an AttributesManager if it already exists in the AttributesManager.
Delaunay interface for NearestNeighbors search.
NearestNeighborSearch * nn_search()
Gets the NearestNeighborSearch used internally.
virtual void enlarge_neighborhood(index_t i, index_t nb)
Stores nb neighbors with vertex i.
Abstract interface for Delaunay triangulation in Nd.
signed_index_t cell_vertex(index_t c, index_t lv) const
Gets a vertex index by cell index and local vertex index.
virtual index_t nearest_vertex(const double *p) const
Computes the nearest vertex from a query point.
const double * vertex_ptr(index_t i) const
Gets a pointer to a vertex by its global index.
index_t nb_cells() const
Gets the number of cells.
void get_neighbors(index_t v, vector< index_t > &neighbors) const
Gets the one-ring neighbors of vertex v.
coord_index_t dimension() const
Gets the dimension of this Delaunay.
index_t cell_size() const
Gets the number of vertices in each cell.
index_t nb_vertices() const
Gets the number of vertices.
signed_index_t vertex_cell(index_t v) const
Gets an incident cell index by a vertex index.
index_t index(index_t c, signed_index_t v) const
Retrieves a local vertex index from cell index and global vertex index.
signed_index_t next_around_vertex(index_t c, index_t lv) const
Traverses the list of cells incident to a vertex.
index_t vertex(index_t f, index_t lv) const
Gets a vertex by facet and local vertex index.
index_t nb() const
Gets the number of (sub-)elements.
AttributesManager & attributes() const
Gets the attributes manager.
const double * point_ptr(index_t v) const
Gets a point.
virtual void get_nearest_neighbors(index_t nb_neighbors, const double *query_point, index_t *neighbors, double *neighbors_sq_dist) const =0
Finds the nearest neighbors of a point given by coordinates.
Specialization of vector for elements of type bool.
index_t size() const
Gets the number of elements.
Some utilities for implementing surfacic and volumetric restricted Voronoi diagrams.
std::stack< TetSeed > TetSeedStack
A stack of TetSeed.
std::stack< FacetSeed > FacetSeedStack
A stack of FacetSeed.
std::stack< index_t > SeedStack
A stack of seed indices (index_t).
Common include file, providing basic definitions. Should be included before anything else by all head...
Geometric functions in arbitrary dimension.
Classes for managing tuples of indices.
basic_quadindex< signed_index_t > signed_quadindex
A basic_quadindex made of 4 signed integers.
#define geo_dim_alignment(dim)
Gets a point alignment.
#define geo_assume_aligned(var, alignment)
Informs the compiler that a given pointer is memory-aligned.
double distance2(const COORD_T *p1, const COORD_T *p2, coord_index_t dim)
Computes the squared distance between two nd points.
Global Vorpaline namespace.
index_t max_index_t()
Gets the maximum positive value of type index_t.
void geo_argused(const T &)
Suppresses compiler warnings about unused parameters.
geo_signed_index_t signed_index_t
The type for storing and manipulating indices differences.
GEOGen::SymbolicVertex SymbolicVertex
Symbolic representation of a RestrictedVoronoiDiagram vertex.
geo_index_t index_t
The type for storing and manipulating indices.
geo_coord_index_t coord_index_t
The type for storing coordinate indices, and iterating on the coordinates of a point.
Types and functions for numbers manipulation.
Filtered exact predicates for restricted Voronoi diagrams.
Function and classes for process manipulation.