Geogram  Version 1.9.1-rc
A programming library of geometric algorithms
generic_RVD.h
Go to the documentation of this file.
1 /*
2  * Copyright (c) 2000-2022 Inria
3  * All rights reserved.
4  *
5  * Redistribution and use in source and binary forms, with or without
6  * modification, are permitted provided that the following conditions are met:
7  *
8  * * Redistributions of source code must retain the above copyright notice,
9  * this list of conditions and the following disclaimer.
10  * * Redistributions in binary form must reproduce the above copyright notice,
11  * this list of conditions and the following disclaimer in the documentation
12  * and/or other materials provided with the distribution.
13  * * Neither the name of the ALICE Project-Team nor the names of its
14  * contributors may be used to endorse or promote products derived from this
15  * software without specific prior written permission.
16  *
17  * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
18  * AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
19  * IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
20  * ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
21  * LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
22  * CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
23  * SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
24  * INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
25  * CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
26  * ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
27  * POSSIBILITY OF SUCH DAMAGE.
28  *
29  * Contact: Bruno Levy
30  *
31  * https://www.inria.fr/fr/bruno-levy
32  *
33  * Inria,
34  * Domaine de Voluceau,
35  * 78150 Le Chesnay - Rocquencourt
36  * FRANCE
37  *
38  */
39 
40 #ifndef GEOGRAM_VORONOI_GENERIC_RVD
41 #define GEOGRAM_VORONOI_GENERIC_RVD
42 
43 #include <geogram/basic/common.h>
44 #include <geogram/basic/numeric.h>
48 #include <geogram/mesh/index.h>
50 #include <geogram/basic/process.h>
52 #include <geogram/basic/argused.h>
53 
54 #include <deque>
55 #include <algorithm>
56 #include <iostream>
57 
67 namespace GEOGen {
68 
77  template <index_t DIM>
79 
82 
83  public:
88  return DIM;
89  }
90 
95 
100 
105 
110 
111  /********************************************************************/
112 
120  GEO::Mesh* mesh
121  ) :
122  mesh_(mesh),
123  delaunay_(delaunay),
124  intersections_(DIM),
125  symbolic_(false),
126  check_SR_(true),
127  exact_(false)
128  {
129  delaunay_nn_ = dynamic_cast<GEO::Delaunay_NearestNeighbors*>(
130  delaunay_
131  );
132  dimension_ = 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;
141  cur_stamp_ = -1;
142  current_facet_ = GEO::max_index_t();
143  current_seed_ = GEO::max_index_t();
144  current_polygon_ = nullptr;
145  current_tet_ = GEO::max_index_t();
146  current_polyhedron_ = nullptr;
147  }
148 
157  connected_components_priority_ = x;
158  }
159 
160 
172  return connected_components_priority_;
173  }
174 
178  const GEO::Mesh* mesh() const {
179  return mesh_;
180  }
181 
186  return mesh_;
187  }
188 
192  const Delaunay* delaunay() const {
193  return delaunay_;
194  }
195 
200  return delaunay_;
201  }
202 
207  delaunay_ = delaunay;
208  delaunay_nn_ = dynamic_cast<GEO::Delaunay_NearestNeighbors*>(
209  delaunay_
210  );
211  }
212 
217  mesh_ = mesh;
218  }
219 
226  void set_facets_range(index_t facets_begin, index_t facets_end) {
227  geo_debug_assert(facets_end >= facets_begin);
228  facets_begin_ = facets_begin;
229  facets_end_ = facets_end;
230  }
231 
240  void set_tetrahedra_range(index_t tets_begin, index_t tets_end) {
241  geo_debug_assert(tets_end >= tets_begin);
242  tets_begin_ = tets_begin;
243  tets_end_ = tets_end;
244  }
245 
251  return facets_end_ - facets_begin_;
252  }
253 
259  return tets_end_ - tets_begin_;
260  }
261 
268  return current_facet_;
269  }
270 
277  return current_seed_;
278  }
279 
287  const Polygon& current_polygon() const {
288  return *current_polygon_;
289  }
290 
297  return current_tet_;
298  }
299 
309  return *current_polyhedron_;
310  }
311 
318  void set_symbolic(bool x) {
319  symbolic_ = x;
320  // exact mode requires symbolic mode.
321  if(exact_) {
322  symbolic_ = true;
323  }
324  }
325 
329  bool symbolic() const {
330  return symbolic_;
331  }
332 
338  void set_exact_predicates(bool x) {
339  exact_ = x;
340  // exact mode requires symbolic mode.
341  if(exact_) {
342  symbolic_ = true;
343  }
344  }
345 
349  bool exact_predicates() const {
350  return exact_;
351  }
352 
356  void set_check_SR(bool x) {
357  check_SR_ = x;
358  }
359 
365  bool check_SR() const {
366  return check_SR_;
367  }
368 
376  return &intersections_;
377  }
378 
379  protected:
391  template <class ACTION>
393  public:
398  PolygonAction(const ACTION& do_it) :
399  do_it_(do_it) {
400  }
401 
411  index_t v,
412  index_t f,
413  const Polygon& P
414  ) const {
415  GEO::geo_argused(f);
416  const_cast<ACTION&> ( do_it_)(v, P);
417  }
418 
419  protected:
420  const ACTION& do_it_;
421  };
422 
430  template <class ACTION>
432  public:
438  TriangleAction(const ACTION& do_it) :
439  do_it_(do_it) {
440  }
441 
453  index_t v,
454  index_t f,
455  const Polygon& P
456  ) const {
457  GEO::geo_argused(f);
458  for(index_t i = 1; i + 1 < P.nb_vertices(); i++) {
459  const_cast<ACTION&> (do_it_)(
460  v, P.vertex(0), P.vertex(i), P.vertex(i + 1)
461  );
462  }
463  }
464 
465  protected:
466  const ACTION& do_it_;
467  };
468 
476  template <class ACTION>
478  public:
484  HalfedgeAction(const ACTION& do_it) :
485  do_it_(do_it) {
486  }
487 
498  index_t v,
499  index_t f,
500  const Polygon& P
501  ) const {
502  GEO::geo_argused(v);
503  GEO::geo_argused(f);
504  for(index_t i = 0; i < P.nb_vertices(); i++) {
505  if(P.vertex(i).check_flag(INTERSECT)) {
506  index_t j = P.next_vertex(i);
507  const_cast<ACTION&> (do_it_)(
508  P.vertex(i), P.vertex(j)
509  );
510  }
511  }
512  }
513 
514  protected:
515  const ACTION& do_it_;
516  };
517 
525  template <class ACTION>
527  public:
533  BorderHalfedgeAction(const ACTION& do_it) :
534  do_it_(do_it) {
535  }
536 
548  index_t v,
549  index_t f,
550  const Polygon& P
551  ) const {
552  GEO::geo_argused(f);
553  for(index_t i = 0; i < P.nb_vertices(); i++) {
554  if(P.vertex(i).check_flag(ORIGINAL)) {
555  if(P.vertex(i).adjacent_facet() == -1) {
556  index_t j = P.next_vertex(i);
557  const_cast<ACTION&> (do_it_)(
558  v, P.vertex(i), P.vertex(j)
559  );
560  }
561  }
562  }
563  }
564 
565  private:
566  const ACTION& do_it_;
567  };
568 
575  template <class ACTION>
577  public:
583  PrimalTriangleAction(const ACTION& do_it) :
584  do_it_(do_it) {
585  }
586 
595  index_t iv1,
596  index_t f,
597  const Polygon& P
598  ) const {
599  GEO::geo_argused(f);
600  for(index_t i = 0; i < P.nb_vertices(); i++) {
601  const Vertex& ve = P.vertex(i);
602  // Primal triangles correspond to vertices of
603  // the RVD that are on two bisectors.
604  if(ve.sym().nb_bisectors() == 2) {
605  index_t iv2 = ve.sym().bisector(0);
606  index_t iv3 = ve.sym().bisector(1);
607  // This test generates triangle (iv1,iv2,iv3)
608  // only once (i.e. if iv1 is the vertex with
609  // the smallest index).
610  if(iv1 < iv2 && iv1 < iv3) {
611  const_cast<ACTION&> (do_it_)(iv1, iv2, iv3);
612  }
613  }
614  }
615  }
616 
617  protected:
618  const ACTION& do_it_;
619  };
620 
634  template <class ACTION>
636  public:
642  PolyhedronAction(const ACTION& do_it) :
643  do_it_(do_it) {
644  }
645 
655  index_t v,
656  index_t t,
657  const Polyhedron& C
658  ) const {
659  const_cast<ACTION&> ( do_it_)(v, t, C);
660  }
661 
662  protected:
663  const ACTION& do_it_;
664  };
665 
688  template <class ACTION>
690  public:
708  const ACTION& do_it,
709  bool visit_inner_tets = false,
710  bool coherent_triangles = false
711  ) :
712  do_it_(do_it),
713  visit_inner_tets_(visit_inner_tets),
714  coherent_triangles_(coherent_triangles)
715  {
716  }
717 
727  index_t v,
728  index_t t,
729  const Polyhedron& C
730  ) const {
731  for(index_t cv = 0; cv < C.max_v(); ++cv) {
732  signed_index_t ct = C.vertex_triangle(cv);
733  if(ct == -1) {
734  continue;
735  }
737 
738  signed_index_t adjacent = C.vertex_id(cv);
739  signed_index_t v_adj = -1;
740  signed_index_t t_adj = -1;
741 
742  if(adjacent < 0) {
743  // Negative adjacent indices correspond to
744  // tet-tet links (ignored when we want to triangulate
745  // the border of the restricted Voronoi cell while
746  // ignoring internal structures).
747  if(!visit_inner_tets_) {
748  continue;
749  }
750  t_adj = -adjacent - 1;
751  } else if(adjacent > 0) {
752  // Positive adjacent indices correspond to
753  // Voronoi seed - Voronoi seed link
754  v_adj = adjacent - 1;
755  }
756  // and adjacent indicex equal to zero corresponds
757  // to tet on border.
758 
760  index_t(ct),
762  );
763 
764  // If required, ensure that two polygonal facets
765  // seen from two different volumetric cells will
766  // be triangulated coherently.
767  if(coherent_triangles_) {
769  }
770 
771  const Vertex& v1 = C.triangle_dual(c1.t);
772 
773  Polyhedron::Corner c2 = c1;
775  geo_debug_assert(c2 != c1);
776 
777  Polyhedron::Corner c3 = c2;
779  geo_debug_assert(c3 != c1);
780  do {
781  const Vertex& v2 = C.triangle_dual(c2.t);
782  const Vertex& v3 = C.triangle_dual(c3.t);
783  const_cast<ACTION&> (do_it_)(
784  v, v_adj, t, t_adj, v1, v2, v3
785  );
786  c2 = c3;
788  } while(c3 != c1);
789  }
790  }
791 
809  const Polyhedron& C, Polyhedron::Corner& c,
810  index_t center_vertex_id
811  ) const {
812  Polyhedron::Corner first = c;
813  Polyhedron::Corner cur = c;
814  do {
815  if(symbolic_compare(
816  C.triangle_dual(cur.t),
817  C.triangle_dual(c.t),
818  center_vertex_id
819  )) {
820  c = cur;
821  }
823  } while(cur != first);
824  }
825 
840  static bool symbolic_compare(
841  const Vertex& p1, const Vertex& p2, index_t center_vertex_id
842  ) {
844  signed_index_t(center_vertex_id),
845  p1.sym()[0], p1.sym()[1], p1.sym()[2]
846  );
848  signed_index_t(center_vertex_id),
849  p2.sym()[0], p2.sym()[1], p2.sym()[2]
850  );
851  return K1 < K2;
852  }
853 
854  protected:
855  const ACTION& do_it_;
856  bool visit_inner_tets_;
857  bool coherent_triangles_;
858  };
859 
882  template <class ACTION>
884  public:
891  const ACTION& do_it
892  ) :
893  do_it_(do_it)
894  {
895  }
896 
906  index_t v,
907  index_t t,
908  const Polyhedron& C
909  ) const {
910 
911  // Find a vertex of the current cell,
912  // that will be used as the 'origin'
913  // vertex
914  const Vertex* v0 = nullptr;
915  index_t t0;
916  for(t0 = 0; t0 < C.max_t(); ++t0) {
917  if(C.triangle_is_used(t0)) {
918  v0 = &C.triangle_dual(t0);
919  break;
920  }
921  }
922 
923  // If current cell is empty, return
924  if(v0 == nullptr) {
925  return;
926  }
927 
928  for(index_t cv = 0; cv < C.max_v(); ++cv) {
929  signed_index_t ct = C.vertex_triangle(cv);
930  if(ct == -1) {
931  continue;
932  }
934 
935  signed_index_t adjacent = C.vertex_id(cv);
936  signed_index_t v_adj = -1;
937  signed_index_t t_adj = -1;
938 
939  if(adjacent < 0) {
940  // Negative adjacent indices correspond to
941  // tet-tet links
942  t_adj = -adjacent - 1;
943  } else if(adjacent > 0) {
944  // Positive adjacent indices correspond to
945  // Voronoi seed - Voroni seed link
946  v_adj = adjacent - 1;
947  }
948  // and adjacent indicex equal to zero corresponds
949  // to tet on border.
950 
952  index_t(ct), C.find_triangle_vertex(index_t(ct), cv)
953  );
954 
955  // If the current facet is incident to
956  // the origin vertex, then skip it (else
957  // it would generate flat tetrahedra)
958  if(facet_is_incident_to_vertex(C, c1, t0)) {
959  continue;
960  }
961 
962  const Vertex& v1 = C.triangle_dual(c1.t);
963 
964  Polyhedron::Corner c2 = c1;
966  geo_debug_assert(c2 != c1);
967 
968  Polyhedron::Corner c3 = c2;
970  geo_debug_assert(c3 != c1);
971  do {
972  const Vertex& v2 = C.triangle_dual(c2.t);
973  const Vertex& v3 = C.triangle_dual(c3.t);
974  const_cast<ACTION&> (do_it_)(
975  v, v_adj, t, t_adj, * v0, v1, v2, v3
976  );
977  c2 = c3;
979  } while(c3 != c1);
980  }
981  }
982 
983  protected:
995  const Polyhedron& C, Polyhedron::Corner& c,
996  index_t t
997  ) const {
998  Polyhedron::Corner first = c;
999  Polyhedron::Corner cur = c;
1000  do {
1001  if(cur.t == t) {
1002  return true;
1003  }
1005  } while(cur != first);
1006  return false;
1007  }
1008 
1009  protected:
1010  const ACTION& do_it_;
1011  };
1012 
1020  template <class ACTION>
1022  public:
1027  PrimalTetrahedronAction(const ACTION& do_it) :
1028  do_it_(do_it) {
1029  }
1030 
1040  index_t v,
1041  index_t t,
1042  const Polyhedron& C
1043  ) const {
1044  GEO::geo_argused(t);
1045  for(index_t it = 0; it < C.max_t(); ++it) {
1046  if(C.triangle_is_used(it)) {
1047  const SymbolicVertex& sym = C.triangle_dual(it).sym();
1048  if(sym.nb_bisectors() == 3) {
1049  index_t v1 = sym.bisector(0);
1050  index_t v2 = sym.bisector(1);
1051  index_t v3 = sym.bisector(2);
1052  // This test ensures that the tet (v,v1,v2,v3)
1053  // is generated only once.
1054  if(v < v1 && v < v2 && v < v3) {
1055  const_cast<ACTION&> (do_it_)(v, v1, v2, v3);
1056  }
1057  }
1058  }
1059  }
1060  }
1061 
1062  protected:
1063  const ACTION& do_it_;
1064  };
1065 
1066  public:
1083  template <class ACTION>
1084  inline void for_each_polygon(
1085  const ACTION& action
1086  ) {
1087  this->template compute_surfacic<PolygonAction<ACTION> >(
1088  PolygonAction<ACTION>(action)
1089  );
1090  }
1091 
1100  template <class TRIACTION>
1101  inline void for_each_triangle(
1102  const TRIACTION& action
1103  ) {
1104  this->template compute_surfacic<TriangleAction<TRIACTION> >(
1106  );
1107  }
1108 
1118  template <class HEACTION>
1119  inline void for_each_halfedge(
1120  const HEACTION& action
1121  ) {
1122  this->template compute_surfacic<HalfedgeAction<HEACTION> >(
1123  HalfedgeAction<HEACTION>(action)
1124  );
1125  }
1126 
1136  template <class BOACTION>
1138  const BOACTION& action
1139  ) {
1140  this->template compute_surfacic<BorderHalfedgeAction<BOACTION> >(
1142  );
1143  }
1144 
1154  template <class PRIMTRIACTION>
1156  const PRIMTRIACTION& action
1157  ) {
1158  bool sym_backup = symbolic();
1159  set_symbolic(true);
1160  this->template compute_surfacic<
1161  PrimalTriangleAction<PRIMTRIACTION
1162  > >(
1164  );
1165  set_symbolic(sym_backup);
1166  }
1167 
1178  template <class ACTION>
1179  inline void for_each_polyhedron(
1180  const ACTION& action
1181  ) {
1182  this->template compute_volumetric<PolyhedronAction<ACTION> >(
1183  PolyhedronAction<ACTION>(action)
1184  );
1185  }
1186 
1221  template <class ACTION>
1223  const ACTION& action,
1224  bool visit_inner_tets = false, bool coherent_triangles = false
1225  ) {
1226  this->template compute_volumetric<
1228  >(
1230  action, visit_inner_tets, coherent_triangles
1231  )
1232  );
1233  }
1234 
1258  template <class ACTION>
1260  const ACTION& action
1261  ) {
1262  this->template compute_volumetric<TetrahedronAction<ACTION> >(
1264  action
1265  )
1266  );
1267  }
1268 
1281  template <class ACTION>
1283  const ACTION& action
1284  ) {
1285  bool sym_backup = symbolic();
1286  set_symbolic(true);
1287  this->template compute_volumetric<PrimalTetrahedronAction<ACTION> >(
1289  action
1290  )
1291  );
1292  set_symbolic(sym_backup);
1293  }
1294 
1295  protected:
1313  template <class ACTION>
1314  inline void compute_surfacic(const ACTION& action) {
1315  if(connected_components_priority_) {
1316  this->template compute_surfacic_with_cnx_priority<ACTION>(
1317  action
1318  );
1319  } else {
1320  this->template compute_surfacic_with_seeds_priority<ACTION>(
1321  action
1322  );
1323  }
1324  }
1325 
1338  template <class ACTION>
1340  const ACTION& action
1341  ) {
1342  if(
1343  facets_begin_ == UNSPECIFIED_RANGE &&
1344  facets_end_ == UNSPECIFIED_RANGE
1345  ) {
1346  facets_begin_ = 0;
1347  facets_end_ = mesh_->facets.nb();
1348  }
1349  current_polygon_ = nullptr;
1350  GEO::vector<index_t> seed_stamp(
1351  delaunay_->nb_vertices(), index_t(-1)
1352  );
1353  GEO::vector<bool> facet_is_marked(facets_end_-facets_begin_, false);
1355 
1356  FacetSeedStack adjacent_facets;
1357  SeedStack adjacent_seeds;
1358  Polygon F;
1359  GEO::Attribute<double> vertex_weight;
1360  vertex_weight.bind_if_is_defined(
1361  mesh_->vertices.attributes(), "weight"
1362  );
1363 
1364  // The algorithm propagates along both the facet-graph of
1365  // the surface and the 1-skeleton of the Delaunay triangulation,
1366  // and computes all the relevant intersections between
1367  // each Voronoi cell and facet.
1368  for(index_t f = facets_begin_; f < facets_end_; f++) {
1369  if(!facet_is_marked[f-facets_begin_]) {
1370  // Propagate along the facet-graph.
1371  facet_is_marked[f-facets_begin_] = true;
1372  adjacent_facets.push(
1374  );
1375  while(!adjacent_facets.empty()) {
1376  current_facet_ = adjacent_facets.top().f;
1377  current_seed_ = adjacent_facets.top().seed;
1378  adjacent_facets.pop();
1379 
1380  // Copy the current facet from the Mesh into
1381  // RestrictedVoronoiDiagram's Polygon data structure
1382  // (gathers all the necessary information)
1384  mesh_, current_facet_, symbolic_, vertex_weight
1385  );
1386 
1387  // Propagate along the Delaunay 1-skeleton
1388  // This will traverse all the seeds such that their
1389  // Voronoi cell has a non-empty intersection with
1390  // the current facet.
1391  seed_stamp[current_seed_] = current_facet_;
1392  adjacent_seeds.push(current_seed_);
1393 
1394  while(!adjacent_seeds.empty()) {
1395  current_seed_ = adjacent_seeds.top();
1396  adjacent_seeds.pop();
1397 
1398  current_polygon_ = intersect_cell_facet(
1399  current_seed_, F
1400  );
1401 
1402  action(
1403  current_seed_, current_facet_, current_polygon()
1404  );
1405 
1406  // Propagate to adjacent facets and adjacent seeds
1407  for(index_t v = 0;
1408  v < current_polygon().nb_vertices(); v++
1409  ) {
1410  const Vertex& ve = current_polygon().vertex(v);
1411  signed_index_t neigh_f = ve.adjacent_facet();
1412  if(
1413  neigh_f >= signed_index_t(facets_begin_) &&
1414  neigh_f < signed_index_t(facets_end_) &&
1415  neigh_f != signed_index_t(current_facet_)
1416  ) {
1417  if(!facet_is_marked[
1418  index_t(neigh_f)-facets_begin_
1419  ]) {
1420  facet_is_marked[
1421  index_t(neigh_f)-facets_begin_
1422  ] = true;
1423  adjacent_facets.push(
1424  FacetSeed(
1425  index_t(neigh_f),
1426  current_seed_
1427  )
1428  );
1429  }
1430  }
1431  signed_index_t neigh_s = ve.adjacent_seed();
1432  if(neigh_s != -1) {
1433  if(
1434  seed_stamp[neigh_s] != current_facet_
1435  ) {
1436  seed_stamp[neigh_s] = current_facet_;
1437  adjacent_seeds.push(index_t(neigh_s));
1438  }
1439  }
1440  }
1441  }
1442  }
1443  }
1444  }
1445  current_polygon_ = nullptr;
1446  }
1447 
1461  template <class ACTION>
1462  inline void compute_volumetric(
1463  const ACTION& action
1464  ) {
1465  if(connected_components_priority_) {
1466  this->template compute_volumetric_with_cnx_priority<ACTION>(
1467  action
1468  );
1469  } else {
1470  this->template compute_volumetric_with_seeds_priority<ACTION>(
1471  action
1472  );
1473  }
1474  }
1475 
1488  template <class ACTION>
1490  const ACTION& action
1491  ) {
1492  if(
1493  tets_begin_ == UNSPECIFIED_RANGE &&
1494  tets_end_ == UNSPECIFIED_RANGE
1495  ) {
1496  tets_begin_ = 0;
1497  tets_end_ = mesh_->cells.nb();
1498  }
1499 
1500  geo_assert(tets_begin_ != UNSPECIFIED_RANGE);
1501  geo_assert(tets_end_ != UNSPECIFIED_RANGE);
1502 
1503  GEO::vector<index_t> seed_stamp(
1504  delaunay_->nb_vertices(), index_t(-1)
1505  );
1506  GEO::vector<bool> tet_is_marked(tets_end_-tets_begin_, false);
1508 
1509  TetSeedStack adjacent_tets;
1510  SeedStack adjacent_seeds;
1511  Polyhedron C(dimension());
1512  GEO::Attribute<double> vertex_weight;
1513  vertex_weight.bind_if_is_defined(
1514  mesh_->vertices.attributes(), "weight"
1515  );
1516 
1517  current_polyhedron_ = &C;
1518  // The algorithm propagates along both the facet-graph of
1519  // the surface and the 1-skeleton of the Delaunay triangulation,
1520  // and computes all the relevant intersections between
1521  // each Voronoi cell and facet.
1522  for(index_t t = tets_begin_; t < tets_end_; ++t) {
1523  if(!tet_is_marked[t-tets_begin_]) {
1524  // Propagate along the tet-graph.
1525  tet_is_marked[t-tets_begin_] = true;
1526  adjacent_tets.push(
1528  );
1529  while(!adjacent_tets.empty()) {
1530  current_tet_ = adjacent_tets.top().f;
1531  current_seed_ = adjacent_tets.top().seed;
1532  adjacent_tets.pop();
1533 
1534  // Note: current cell could be looked up here,
1535  // (from current_tet_) if we chose to keep it
1536  // and copy it right before clipping (I am
1537  // not sure that it is worth it, lookup time
1538  // will be probably fast enough)
1539 
1540  // Propagate along the Delaunay 1-skeleton
1541  // This will traverse all the seeds such that their
1542  // Voronoi cell has a non-empty intersection with
1543  // the current facet.
1544  seed_stamp[current_seed_] = current_tet_;
1545  adjacent_seeds.push(current_seed_);
1546 
1547  while(!adjacent_seeds.empty()) {
1548  current_seed_ = adjacent_seeds.top();
1549  adjacent_seeds.pop();
1550 
1552  mesh_, current_tet_, symbolic_, vertex_weight
1553  );
1554 
1556  current_seed_, C
1557  );
1558 
1559  action(
1560  current_seed_, current_tet_,
1562  );
1563 
1564  // Propagate to adjacent tets and adjacent seeds
1565  // Iterate on the vertices of the cell (remember:
1566  // the cell is represented in dual form)
1567  for(index_t v = 0;
1568  v < current_polyhedron().max_v(); ++v
1569  ) {
1570 
1571  // Skip clipping planes that are no longer
1572  // connected to a cell facet.
1573  if(
1575  == -1
1576  ) {
1577  continue;
1578  }
1579 
1580  signed_index_t id =
1582  if(id > 0) {
1583  // Propagate to adjacent seed
1584  index_t neigh_s = index_t(id - 1);
1585  if(seed_stamp[neigh_s] != current_tet_) {
1586  seed_stamp[neigh_s] = current_tet_;
1587  adjacent_seeds.push(neigh_s);
1588  }
1589  } else if(id < 0) {
1590  // id==0 corresponds to facet on boundary
1591  // (skipped)
1592  // id<0 corresponds to adjacent tet index
1593 
1594  // Propagate to adjacent tet
1595  signed_index_t neigh_t = -id - 1;
1596  if(
1597  neigh_t >=
1598  signed_index_t(tets_begin_) &&
1599  neigh_t < signed_index_t(tets_end_) &&
1600  neigh_t != signed_index_t(current_tet_)
1601  ) {
1602  if(!tet_is_marked[
1603  index_t(neigh_t)-tets_begin_
1604  ]) {
1605  tet_is_marked[
1606  index_t(neigh_t)-tets_begin_
1607  ] = true;
1608  adjacent_tets.push(
1609  TetSeed(
1610  index_t(neigh_t),
1611  current_seed_
1612  )
1613  );
1614  }
1615  }
1616  }
1617  }
1618  }
1619  }
1620  }
1621  }
1622  current_polyhedron_ = nullptr;
1623  }
1624 
1625 
1645  template <class ACTION>
1647  const ACTION& action
1648  ) {
1649 
1650  if(
1651  tets_begin_ == UNSPECIFIED_RANGE &&
1652  tets_end_ == UNSPECIFIED_RANGE
1653  ) {
1654  tets_begin_ = 0;
1655  tets_end_ = mesh_->cells.nb();
1656  }
1657 
1658  geo_assert(tets_begin_ != UNSPECIFIED_RANGE);
1659  geo_assert(tets_end_ != UNSPECIFIED_RANGE);
1660 
1661  current_polyhedron_ = nullptr;
1663 
1664  std::deque<TetSeed> adjacent_seeds;
1665  std::stack<index_t> adjacent_tets;
1666 
1667  static const index_t NO_STAMP = index_t(-1);
1668  GEO::vector<index_t> tet_stamp(
1669  tets_end_ - tets_begin_, NO_STAMP
1670  );
1671 
1672  TetSeedMarking visited(
1673  tets_end_ - tets_begin_, delaunay_->nb_vertices()
1674  );
1675 
1676  // Yes, facet_seed_marking_ points to the TetSeedMarking,
1677  // (TetSeedMarking is typedef-ed as FacetSeedMarking),
1678  // ugly I know... to be revised.
1679  facet_seed_marking_ = &visited;
1680  Polyhedron C(dimension());
1681  current_polyhedron_ = &C;
1682 
1683  current_connected_component_ = 0;
1684  // index_t C_index = tets_end_ + 1; // Unused (see comment later)
1685 
1686  GEO::Attribute<double> vertex_weight;
1687  vertex_weight.bind_if_is_defined(
1688  mesh_->vertices.attributes(),"weight"
1689  );
1690 
1691  // The algorithm propagates along both the facet-graph of
1692  // the surface and the 1-skeleton of the Delaunay triangulation,
1693  // and computes all the relevant intersections between
1694  // each Voronoi cell and facet.
1695  for(index_t t = tets_begin_; t < tets_end_; ++t) {
1696  if(tet_stamp[t - tets_begin_] == NO_STAMP) {
1697  current_tet_ = t;
1698  current_seed_ = find_seed_near_tet(t);
1699 
1700  adjacent_seeds.push_back(
1701  TetSeed(current_tet_, current_seed_)
1702  );
1703 
1704  // Propagate along the Delaunay-graph.
1705  while(!adjacent_seeds.empty()) {
1706  // Yes, f, because TetSeed is typedef-ed as FacetSeed
1707  current_tet_ = adjacent_seeds.front().f;
1708  current_seed_ = adjacent_seeds.front().seed;
1709  adjacent_seeds.pop_front();
1710  if(
1711  tet_stamp[current_tet_ - tets_begin_] ==
1712  current_seed_
1713  ) {
1714  continue;
1715  }
1716 
1717  if(visited.is_marked(current_tet_, current_seed_)) {
1718  continue;
1719  }
1720 
1721  connected_component_changed_ = true;
1722  adjacent_tets.push(current_tet_);
1723  tet_stamp[current_tet_ - tets_begin_] =
1724  current_seed_;
1725  while(!adjacent_tets.empty()) {
1726  current_tet_ = adjacent_tets.top();
1727  adjacent_tets.pop();
1728 
1729  // Copy the current tet from the Mesh into
1730  // RestrictedVoronoiDiagram's Polyhedron data structure
1731  // (gathers all the necessary information)
1732 
1734  mesh_, current_tet_, symbolic_, vertex_weight
1735  );
1736 
1737  // Note: difference with compute_surfacic_with_cnx_priority():
1738  // Since intersect_cell_cell() overwrites C, we
1739  // need to initialize C from the mesh for each visited
1740  // (tet,seed) pair (and the test for current_tet_ change
1741  // with C_index is not used here).
1742  // C_index = current_tet_;
1743 
1744  intersect_cell_cell(current_seed_, C);
1745  action(
1746  current_seed_, current_tet_, current_polyhedron()
1747  );
1748  connected_component_changed_ = false;
1749 
1750  bool touches_RVC_border = false;
1751 
1752  // Propagate to adjacent tets and adjacent seeds
1753  for(index_t v = 0;
1754  v < current_polyhedron().max_v(); ++v
1755  ) {
1756 
1757  // Skip clipping planes that are no longer
1758  // connected to a cell facet.
1759  if(
1761  == -1
1762  ) {
1763  continue;
1764  }
1765 
1767  if(id < 0) {
1768  // id == 0 corresponds to facet on boundary
1769  // (skipped)
1770  // id < 0 corresponds to adjacent tet index
1771  signed_index_t s_neigh_t = -id-1;
1772  if(
1773  s_neigh_t >= signed_index_t(tets_begin_) &&
1774  s_neigh_t < signed_index_t(tets_end_)
1775  ) {
1777  s_neigh_t !=
1778  signed_index_t(current_tet_)
1779  );
1780  index_t neigh_t = index_t(s_neigh_t);
1781  if(
1782  tet_stamp[neigh_t - tets_begin_] !=
1783  current_seed_
1784  ) {
1785  tet_stamp[neigh_t - tets_begin_] =
1786  current_seed_;
1787  adjacent_tets.push(neigh_t);
1788  }
1789  }
1790  } else if(id > 0) {
1791  index_t neigh_s = index_t(id-1);
1792  touches_RVC_border = true;
1793  TetSeed ts(current_tet_, neigh_s);
1794  if(!visited.is_marked(ts)) {
1795  adjacent_seeds.push_back(ts);
1796  }
1797  }
1798 
1799  }
1800  if(touches_RVC_border) {
1801  visited.mark(
1802  TetSeed(current_tet_, current_seed_),
1803  current_connected_component_
1804  );
1805  }
1806  }
1807  ++current_connected_component_;
1808  }
1809  }
1810  }
1811  facet_seed_marking_ = nullptr;
1812  }
1813 
1814 
1815  public:
1822  geo_debug_assert(facet_seed_marking_ != nullptr);
1823  return facet_seed_marking_->is_marked(FacetSeed(f, s));
1824  }
1825 
1835  index_t f, index_t s
1836  ) const {
1837  geo_debug_assert(facet_seed_marking_ != nullptr);
1838  return facet_seed_marking_->get_connected_component(
1839  FacetSeed(f, s)
1840  );
1841  }
1842 
1847  return connected_component_changed_;
1848  }
1849 
1854  return current_connected_component_;
1855  }
1856 
1857  protected:
1878  template <class ACTION>
1880  const ACTION& action
1881  ) {
1882 
1883  if(
1884  facets_begin_ == UNSPECIFIED_RANGE &&
1885  facets_end_ == UNSPECIFIED_RANGE
1886  ) {
1887  facets_begin_ = 0;
1888  facets_end_ = mesh_->facets.nb();
1889  }
1890 
1891  current_polygon_ = nullptr;
1893 
1894  std::deque<FacetSeed> adjacent_seeds;
1895  std::stack<index_t> adjacent_facets;
1896 
1897  static const index_t NO_STAMP = index_t(-1);
1898  GEO::vector<index_t> facet_stamp(
1899  facets_end_ - facets_begin_, NO_STAMP
1900  );
1901 
1902  FacetSeedMarking visited(
1903  facets_end_ - facets_begin_, delaunay_->nb_vertices()
1904  );
1905 
1906  facet_seed_marking_ = &visited;
1907  Polygon F;
1908  current_connected_component_ = 0;
1909  index_t F_index = facets_end_ + 1;
1910 
1911  GEO::Attribute<double> vertex_weight;
1912  vertex_weight.bind_if_is_defined(
1913  mesh_->vertices.attributes(),"weight"
1914  );
1915 
1916  // The algorithm propagates along both the facet-graph of
1917  // the surface and the 1-skeleton of the Delaunay triangulation,
1918  // and computes all the relevant intersections between
1919  // each Voronoi cell and facet.
1920  for(index_t f = facets_begin_; f < facets_end_; ++f) {
1921  if(facet_stamp[f - facets_begin_] == NO_STAMP) {
1922  current_facet_ = f;
1923  current_seed_ = find_seed_near_facet(f);
1924 
1925  adjacent_seeds.push_back(
1926  FacetSeed(current_facet_, current_seed_)
1927  );
1928 
1929  // Propagate along the Delaunay-graph.
1930  while(
1931  !adjacent_seeds.empty()
1932  ) {
1933  current_facet_ = adjacent_seeds.front().f;
1934  current_seed_ = adjacent_seeds.front().seed;
1935  adjacent_seeds.pop_front();
1936  if(
1937  facet_stamp[current_facet_ - facets_begin_] ==
1938  current_seed_
1939  ) {
1940  continue;
1941  }
1942 
1943  if(visited.is_marked(current_facet_, current_seed_)) {
1944  continue;
1945  }
1946 
1947  connected_component_changed_ = true;
1948  adjacent_facets.push(current_facet_);
1949  facet_stamp[current_facet_ - facets_begin_] =
1950  current_seed_;
1951  while(!adjacent_facets.empty()) {
1952  current_facet_ = adjacent_facets.top();
1953  adjacent_facets.pop();
1954 
1955  // Copy the current facet from the Mesh into
1956  // RestrictedVoronoiDiagram's Polygon data structure
1957  // (gathers all the necessary information)
1958  if(F_index != current_facet_) {
1960  mesh_, current_facet_, symbolic_,
1961  vertex_weight
1962  );
1963  F_index = current_facet_;
1964  }
1965 
1966  current_polygon_ = intersect_cell_facet(
1967  current_seed_, F
1968  );
1969  action(
1970  current_seed_, current_facet_, current_polygon()
1971  );
1972  connected_component_changed_ = false;
1973 
1974  bool touches_RVC_border = false;
1975 
1976  // Propagate to adjacent facets and adjacent seeds
1977  for(index_t v = 0;
1978  v < current_polygon().nb_vertices(); v++
1979  ) {
1980  const Vertex& ve = current_polygon().vertex(v);
1981  signed_index_t s_neigh_f = ve.adjacent_facet();
1982  if(
1983  s_neigh_f >= signed_index_t(facets_begin_)
1984  &&
1985  s_neigh_f < signed_index_t(facets_end_)
1986  ) {
1988  s_neigh_f !=
1989  signed_index_t(current_facet_)
1990  );
1991  index_t neigh_f = index_t(s_neigh_f);
1992  if(
1993  facet_stamp[neigh_f - facets_begin_] !=
1994  current_seed_
1995  ) {
1996  facet_stamp[neigh_f - facets_begin_] =
1997  current_seed_;
1998  adjacent_facets.push(neigh_f);
1999  }
2000  }
2001  signed_index_t neigh_s = ve.adjacent_seed();
2002  if(neigh_s != -1) {
2003  touches_RVC_border = true;
2004  FacetSeed fs(
2005  current_facet_, index_t(neigh_s)
2006  );
2007  if(!visited.is_marked(fs)) {
2008  adjacent_seeds.push_back(fs);
2009  }
2010  }
2011  }
2012  if(touches_RVC_border) {
2013  visited.mark(
2014  FacetSeed(current_facet_, current_seed_),
2015  current_connected_component_
2016  );
2017  }
2018  }
2019  ++current_connected_component_;
2020  }
2021  }
2022  }
2023  facet_seed_marking_ = nullptr;
2024  }
2025 
2034  const double* p = mesh_->vertices.point_ptr(
2035  mesh_->facets.vertex(f,0)
2036  );
2037  return find_seed_near_point(p);
2038  }
2039 
2048  index_t v = mesh_->cells.tet_vertex(t, 0);
2049  const double* p = mesh_->vertices.point_ptr(v);
2050  return find_seed_near_point(p);
2051  }
2052 
2059  index_t find_seed_near_point(const double* p) {
2060  // In order to be compatible with the symbolic
2061  // perturbation, if the nearest neighbor is
2062  // non-unique, we need to return the one of
2063  // lowest index (because in case of several seeds
2064  // at equal distance, the one of lowest index
2065  // is guaranteed to have the facet in its Voronoi
2066  // cell from the point of view of symbolic
2067  // perturbation).
2068  if(exact_ && delaunay_nn_ != nullptr) {
2069  // TODO: may need more than 10
2070  index_t neighbors[10];
2071  double neighbors_sq_dist[10];
2072  index_t nb = 10;
2073  if(delaunay_nn_->nb_vertices() < nb) {
2074  nb = delaunay_nn_->nb_vertices();
2075  }
2076  delaunay_nn_->nn_search()->get_nearest_neighbors(
2077  nb, p, neighbors, neighbors_sq_dist
2078  );
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) {
2083  break;
2084  }
2085  if(neighbors[i] < nearest) {
2086  nearest = neighbors[i];
2087  }
2088  }
2089  return nearest;
2090  }
2091 
2092  return delaunay_->nearest_vertex(p);
2093  }
2094 
2105  void swap_polygons(Polygon*& ping, Polygon*& pong) {
2106  if(ping != &P1 && ping != &P2) {
2107  // First clipping operation, ping points to F
2108  // (current facet copied)
2109  ping = &P2;
2110  pong = &P1;
2111  } else {
2112  std::swap(ping, pong);
2113  }
2114  }
2115 
2124  intersections_.clear();
2125 
2126  // Initialize ping-pong pointers for Sutherland-Hodgman
2127  // re-entrant clipping and copy current facet into 'ping' buffer.
2128  Polygon* ping = &F;
2129  Polygon* pong = &P2;
2130 
2131  // Clip current facet by current Voronoi cell (associated with seed)
2132  if(delaunay_nn_ != nullptr) {
2133  clip_by_cell_SR(seed, ping, pong); // "Security Radius" mode.
2134  } else {
2135  clip_by_cell(seed, ping, pong); // Standard mode.
2136  }
2137 
2138  return ping; // Yes, 'ping', and not 'pong'
2139  // see comments in clip_by_cell()
2140  }
2141 
2157  void clip_by_cell_SR(index_t i, Polygon*& ping, Polygon*& pong) {
2158  // 'Security radius' mode.
2159  // Note: the vertices of the neighborhood are returned in
2160  // increasing distance to pi. We stop the clippings as soon as the
2161  // 'security radius' is reached.
2162  const double* geo_restrict pi = delaunay_->vertex_ptr(i);
2164 
2165  index_t jj = 0;
2166  index_t prev_nb_neighbors = 0;
2167  neighbors_.resize(0);
2168  while(neighbors_.size() < delaunay_nn_->nb_vertices() - 1) {
2169 
2170  delaunay_nn_->get_neighbors(i, neighbors_);
2171  if(neighbors_.size() == 0) {
2172  return;
2173  }
2174  if(prev_nb_neighbors == neighbors_.size()) {
2175  return;
2176  }
2177 
2178  for(; jj < neighbors_.size(); jj++) {
2179  index_t j = neighbors_[jj];
2180  double R2 = 0.0;
2181  for(index_t k = 0; k < ping->nb_vertices(); k++) {
2182  geo_decl_aligned(double dik);
2183  const double* geo_restrict pk = ping->vertex(k).point();
2185  dik = GEO::Geom::distance2(pi, pk, dimension());
2186  R2 = std::max(R2, dik);
2187  }
2188  geo_decl_aligned(double dij);
2189  const double* geo_restrict pj = delaunay_->vertex_ptr(j);
2191  dij = GEO::Geom::distance2(pi, pj, dimension());
2192  // A little bit more than 4, because when
2193  // exact predicates are used, we need to
2194  // include tangent bisectors in the computation.
2195  if(dij > 4.1 * R2) {
2196  return;
2197  }
2198  clip_by_plane(*ping, *pong, i, j);
2199  swap_polygons(ping, pong);
2200  }
2201 
2202  if(!check_SR_) {
2203  return;
2204  }
2205 
2206  index_t nb_neighbors = neighbors_.size();
2207  prev_nb_neighbors = nb_neighbors;
2208 
2209  if(nb_neighbors > 8) {
2210  nb_neighbors += nb_neighbors / 8;
2211  } else {
2212  nb_neighbors++;
2213  }
2214 
2215  nb_neighbors = std::min(
2216  nb_neighbors,
2217  delaunay_nn_->nb_vertices() - 1
2218  );
2219 
2220  delaunay_nn_->enlarge_neighborhood(i, nb_neighbors);
2221  }
2222  }
2223 
2240  void clip_by_cell(index_t i, Polygon*& ping, Polygon*& pong) {
2241  get_neighbors(i);
2242  for(index_t jj = 0; jj < neighbors_.size(); jj++) {
2243  index_t j = neighbors_[jj];
2244  clip_by_plane(*ping, *pong, i, j);
2245  swap_polygons(ping, pong);
2246  }
2247  }
2248 
2264  Polygon& ping, Polygon& pong,
2265  index_t i, index_t j
2266  ) {
2267  ping.clip_by_plane<DIM>(
2268  pong, intersections_, mesh_, delaunay_, i, j, exact_, symbolic_
2269  );
2270  }
2271 
2278  public:
2286  // Clip current facet by current Voronoi cell (associated with seed)
2287  if(delaunay_nn_ != nullptr) {
2288  clip_by_cell_SR(seed, C); // "Security Radius" mode.
2289  } else {
2290  clip_by_cell(seed, C); // Standard mode.
2291  }
2292  }
2293 
2294  protected:
2302 
2303  // 'Security radius' mode.
2304  // Note: the vertices of the neighborhood are returned in
2305  // increasing distance to pi. We stop the clippings as soon as the
2306  // 'security radius' is reached.
2307  const double* geo_restrict pi = delaunay_->vertex_ptr(seed);
2309 
2310  index_t jj = 0;
2311  index_t prev_nb_neighbors = 0;
2312  neighbors_.resize(0);
2313 
2314  while(neighbors_.size() < delaunay_nn_->nb_vertices() - 1) {
2315 
2316  delaunay_nn_->get_neighbors(seed, neighbors_);
2317  if(neighbors_.size() == 0) {
2318  return;
2319  }
2320  if(prev_nb_neighbors == neighbors_.size()) {
2321  return;
2322  }
2323 
2324  for(; jj < neighbors_.size(); jj++) {
2325  index_t j = neighbors_[jj];
2326  double R2 = 0.0;
2327  for(index_t k = 0; k < C.max_t(); ++k) {
2328  if(!C.triangle_is_used(k)) {
2329  continue;
2330  }
2331  geo_decl_aligned(double dik);
2332  const double* geo_restrict pk =
2333  C.triangle_dual(k).point();
2335  dik = GEO::Geom::distance2(pi, pk, dimension());
2336  R2 = std::max(R2, dik);
2337  }
2338  geo_decl_aligned(double dij);
2339  const double* geo_restrict pj = delaunay_->vertex_ptr(j);
2341  dij = GEO::Geom::distance2(pi, pj, dimension());
2342  // A little bit more than 4, because when
2343  // exact predicates are used, we need to
2344  // include tangent bisectors in the computation.
2345  if(dij > 4.1 * R2) {
2346  return;
2347  }
2348  clip_by_plane(C, seed, j);
2349  }
2350 
2351  if(!check_SR_) {
2352  return;
2353  }
2354 
2355  index_t nb_neighbors = neighbors_.size();
2356  prev_nb_neighbors = nb_neighbors;
2357 
2358  if(nb_neighbors > 8) {
2359  nb_neighbors += nb_neighbors / 8;
2360  } else {
2361  nb_neighbors++;
2362  }
2363 
2364  nb_neighbors = std::min(
2365  nb_neighbors,
2366  delaunay_nn_->nb_vertices() - 1
2367  );
2368  delaunay_nn_->enlarge_neighborhood(seed, nb_neighbors);
2369  }
2370  }
2371 
2379  get_neighbors(seed);
2380  // Check whether cell is empty (may happen with
2381  // power diagrams)
2382  if(neighbors_.size() == 0) {
2383  C.clear();
2384  }
2385  for(index_t jj = 0; jj < neighbors_.size(); jj++) {
2386  index_t j = neighbors_[jj];
2387  clip_by_plane(C, seed, j);
2388  }
2389  }
2390 
2399  C.clip_by_plane<DIM>(
2400  mesh_, delaunay_, i, j, exact_, symbolic_
2401  );
2402  }
2403 
2419  // In dimension 3 (and if we do not used the ANN-based algorithm),
2420  // we can use the faster 'stamp-based' algorithm for finding the
2421  // neighbors.
2422  if(delaunay_->dimension() == 3 && delaunay_->nb_cells() != 0) {
2423  cur_stamp_ = 0;
2424  stamp_.assign(delaunay_->nb_vertices(), -1);
2425  }
2426  }
2427 
2436  if(stamp_.size() == 0) {
2437  // Used in ANN mode and with higher dimensions.
2438  delaunay_->get_neighbors(v, neighbors_);
2439  } else {
2440  // Used in 3D mode with standard Delaunay.
2441  // The following loop replaces
2442  // delaunay_->get_neighbors(v,neighbors_) ;
2443  // (and makes the overall algorithm 10 to 30% more efficient)
2444  neighbors_.resize(0);
2445  index_t t = index_t(delaunay_->vertex_cell(v));
2446  do {
2447  index_t lv = delaunay_->index(t, signed_index_t(v));
2448  for(index_t lw = 0; lw < delaunay_->cell_size(); lw++) {
2449  if(lw != lv) {
2450  index_t w = index_t(delaunay_->cell_vertex(t, lw));
2451  if(stamp_[w] != cur_stamp_) {
2452  stamp_[w] = cur_stamp_;
2453  neighbors_.push_back(w);
2454  }
2455  }
2456  }
2457  t = index_t(delaunay_->next_around_vertex(t, lv));
2458  } while(t != index_t(delaunay_->vertex_cell(v)));
2459  cur_stamp_++;
2460  }
2461  }
2462 
2465  protected:
2466  GEO::Mesh* mesh_;
2467  Delaunay* delaunay_;
2468  GEO::Delaunay_NearestNeighbors* delaunay_nn_;
2469 
2470  PointAllocator intersections_;
2471  Polygon* current_polygon_;
2472  Polygon P1, P2;
2473  GEO::vector<index_t> neighbors_;
2474  index_t current_facet_;
2475  index_t current_seed_;
2476  Polyhedron* current_polyhedron_;
2477  index_t current_tet_;
2478 
2479  // For optimized get_neighbors().
2480  signed_index_t cur_stamp_;
2482 
2483  bool symbolic_;
2484  bool check_SR_;
2485  bool exact_;
2486 
2487  coord_index_t dimension_;
2488 
2489  static const index_t UNSPECIFIED_RANGE = index_t(-1);
2490 
2491  index_t facets_begin_;
2492  index_t facets_end_;
2493 
2494  index_t tets_begin_;
2495  index_t tets_end_;
2496 
2497  bool connected_components_priority_;
2498  FacetSeedMarking* facet_seed_marking_;
2499  bool connected_component_changed_;
2500  index_t current_connected_component_;
2501 
2502  private:
2506  RestrictedVoronoiDiagram(const thisclass& rhs);
2507 
2511  thisclass& operator= (const thisclass& rhs);
2512  };
2513 }
2514 
2515 namespace GEO {
2516 
2521 }
2522 
2523 #endif
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.
Definition: assert.h:149
#define geo_debug_assert(x)
Verifies that a condition is met.
Definition: assert.h:196
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().
Definition: generic_RVD.h:526
BorderHalfedgeAction(const ACTION &do_it)
Creates a new BorderHalfedgeAction that wraps a user ACTION instance.
Definition: generic_RVD.h:533
void operator()(index_t v, index_t f, const Polygon &P) const
Callback called for each integration simplex.
Definition: generic_RVD.h:547
Adapter class used internally to implement for_each_halfedge().
Definition: generic_RVD.h:477
void operator()(index_t v, index_t f, const Polygon &P) const
Callback called for each integration simplex.
Definition: generic_RVD.h:497
HalfedgeAction(const ACTION &do_it)
Creates a new HalfedgeAction that wraps a user ACTION instance.
Definition: generic_RVD.h:484
Adapter class used internally to implement for_each_polygon()
Definition: generic_RVD.h:392
void operator()(index_t v, index_t f, const Polygon &P) const
Callback called for each polygon.
Definition: generic_RVD.h:410
PolygonAction(const ACTION &do_it)
Creates a new PolygonAction around a user ACTION instance.
Definition: generic_RVD.h:398
Adapter class used internally to implement for_each_polyhedron()
Definition: generic_RVD.h:635
void operator()(index_t v, index_t t, const Polyhedron &C) const
Callback called for each polyhedron.
Definition: generic_RVD.h:654
PolyhedronAction(const ACTION &do_it)
Creates a new PolyhedronAction that wraps a user ACTION instance.
Definition: generic_RVD.h:642
Adapter class used internally to implement for_each_primal_tetrahedron()
Definition: generic_RVD.h:1021
PrimalTetrahedronAction(const ACTION &do_it)
Constructs a new PrimalTetrahedronAction.
Definition: generic_RVD.h:1027
void operator()(index_t v, index_t t, const Polyhedron &C) const
Callback called for each polyhedron.
Definition: generic_RVD.h:1039
Adapter class used internally to implement for_each_primal_triangle()
Definition: generic_RVD.h:576
void operator()(index_t iv1, index_t f, const Polygon &P) const
Callback called for each primal triangle.
Definition: generic_RVD.h:594
PrimalTriangleAction(const ACTION &do_it)
Creates a new PrimalTriangleAction that wraps a user ACTION instance.
Definition: generic_RVD.h:583
Adapter class used internally to implement for_each_tetrahedron()
Definition: generic_RVD.h:883
TetrahedronAction(const ACTION &do_it)
Creates a new TetrahedronAction that wraps a user ACTION instance.
Definition: generic_RVD.h:890
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.
Definition: generic_RVD.h:994
void operator()(index_t v, index_t t, const Polyhedron &C) const
Callback called for each polyhedron.
Definition: generic_RVD.h:905
Adapter class used internally to implement for_each_triangle().
Definition: generic_RVD.h:431
void operator()(index_t v, index_t f, const Polygon &P) const
Callback called for each integration simplex.
Definition: generic_RVD.h:452
TriangleAction(const ACTION &do_it)
Creates a new TriangleAction that wraps a user ACTION instance.
Definition: generic_RVD.h:438
Adapter class used internally to implement for_each_volumetric_integration_simplex()
Definition: generic_RVD.h:689
VolumetricIntegrationSimplexAction(const ACTION &do_it, bool visit_inner_tets=false, bool coherent_triangles=false)
Creates a new VolumetricIntegrationSimplexAction that wraps a user ACTION instance.
Definition: generic_RVD.h:707
void operator()(index_t v, index_t t, const Polyhedron &C) const
Callback called for each polyhedron.
Definition: generic_RVD.h:726
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.
Definition: generic_RVD.h:808
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.
Definition: generic_RVD.h:840
Computes the intersection between a surface (Mesh) and a Voronoi diagram (dual of a Delaunay).
Definition: generic_RVD.h:78
GEOGen::Vertex Vertex
Internal representation of vertices.
Definition: generic_RVD.h:99
index_t nb_facets_in_range() const
Gets the number of facets in the current range.
Definition: generic_RVD.h:250
void clip_by_cell(index_t seed, Polyhedron &C)
Computes the intersection between a Voronoi cell and a cell in plain mode.
Definition: generic_RVD.h:2378
void set_exact_predicates(bool x)
Specifies whether exact predicates should be used.
Definition: generic_RVD.h:338
GEO::Mesh * mesh()
Gets the input mesh.
Definition: generic_RVD.h:185
void for_each_primal_triangle(const PRIMTRIACTION &action)
Iterates on the triangles of the Restricted Delaunay Triangulation.
Definition: generic_RVD.h:1155
static coord_index_t dimension()
Gets the dimension.
Definition: generic_RVD.h:87
index_t find_seed_near_point(const double *p)
Finds a seed near a given point.
Definition: generic_RVD.h:2059
GEOGen::ConvexCell Polyhedron
Internal representation of volumetric cells.
Definition: generic_RVD.h:109
GEOGen::Polygon Polygon
Internal representation of polygons.
Definition: generic_RVD.h:104
void set_mesh(GEO::Mesh *mesh)
Sets the input mesh.
Definition: generic_RVD.h:216
void for_each_polyhedron(const ACTION &action)
Iterates on the polyhedra of this RVD.
Definition: generic_RVD.h:1179
Polygon * intersect_cell_facet(index_t seed, Polygon &F)
Computes the intersection between the Voronoi cell of a seed and a facet.
Definition: generic_RVD.h:2123
void compute_volumetric_with_cnx_priority(const ACTION &action)
Low-level API of Restricted Voronoi Diagram traversal with connected components priority.
Definition: generic_RVD.h:1646
void for_each_triangle(const TRIACTION &action)
Iterates on the facets of this RVD, triangulated on the fly.
Definition: generic_RVD.h:1101
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.
Definition: generic_RVD.h:2398
void for_each_tetrahedron(const ACTION &action)
Iterates on the polyhedra of this RVD decomposed on the fly into tetrahedra.
Definition: generic_RVD.h:1259
index_t current_seed() const
Gets the index of the Delaunay vertex currently processed.
Definition: generic_RVD.h:276
void clip_by_cell(index_t i, Polygon *&ping, Polygon *&pong)
Computes the intersection between a Voronoi cell and a polygon.
Definition: generic_RVD.h:2240
bool symbolic() const
Tests whether symbolic mode is active.
Definition: generic_RVD.h:329
const GEO::Mesh * mesh() const
Gets the input mesh.
Definition: generic_RVD.h:178
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 ...
Definition: generic_RVD.h:1137
void for_each_halfedge(const HEACTION &action)
Iterates on the halfedges on the borders of the restricted Voronoi cells.
Definition: generic_RVD.h:1119
const Polyhedron & current_polyhedron() const
Gets the current cell.
Definition: generic_RVD.h:308
void compute_volumetric(const ACTION &action)
Low-level API of Restricted Voronoi Diagram traversal .
Definition: generic_RVD.h:1462
index_t current_connected_component() const
Gets the index of the current connected component.
Definition: generic_RVD.h:1853
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.
Definition: generic_RVD.h:2301
void for_each_polygon(const ACTION &action)
Iterates on the facets of this RVD.
Definition: generic_RVD.h:1084
index_t current_facet() const
Gets the index of the mesh facet currently processed.
Definition: generic_RVD.h:267
const Polygon & current_polygon() const
Gets the current polygon.
Definition: generic_RVD.h:287
index_t current_tet() const
Gets the undex of the mesh tetrahedron currently processed.
Definition: generic_RVD.h:296
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.
Definition: generic_RVD.h:2285
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'.
Definition: generic_RVD.h:2157
void get_neighbors(index_t v)
Caches the neighbors of a Delaunay vertex.
Definition: generic_RVD.h:2435
GEOGen::PointAllocator PointAllocator
Used to allocate the generated points.
Definition: generic_RVD.h:94
void compute_volumetric_with_seeds_priority(const ACTION &action)
Low-level API of Restricted Voronoi Diagram traversal with seeds priority in volumetric mode.
Definition: generic_RVD.h:1489
void set_check_SR(bool x)
Specifies whether radius of security should be enforced.
Definition: generic_RVD.h:356
PointAllocator * point_allocator()
Gets the PointAllocator.
Definition: generic_RVD.h:375
RestrictedVoronoiDiagram(Delaunay *delaunay, GEO::Mesh *mesh)
Constructs a new RestrictedVoronoiDiagram.
Definition: generic_RVD.h:118
bool connected_component_changed() const
Tests whether the current connected component changed.
Definition: generic_RVD.h:1846
void set_delaunay(Delaunay *delaunay)
Sets the Delaunay triangulation.
Definition: generic_RVD.h:206
bool connected_components_priority() const
Tests whether connected components priority is set.
Definition: generic_RVD.h:171
bool exact_predicates() const
Tests whether exact predicates are used.
Definition: generic_RVD.h:349
Delaunay * delaunay()
Gets the Delaunay triangulation.
Definition: generic_RVD.h:199
index_t find_seed_near_facet(index_t f)
Finds a seed near a given facet.
Definition: generic_RVD.h:2033
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).
Definition: generic_RVD.h:1834
void swap_polygons(Polygon *&ping, Polygon *&pong)
Swaps two pointers between two polygons.
Definition: generic_RVD.h:2105
void set_connected_components_priority(bool x)
Sets traveral priority.
Definition: generic_RVD.h:156
void for_each_primal_tetrahedron(const ACTION &action)
Iterates on the primal tetrahedra of this RVD.
Definition: generic_RVD.h:1282
void set_symbolic(bool x)
Sets symbolic mode.
Definition: generic_RVD.h:318
void init_get_neighbors()
Creates the data structure for optimized get_neighbors() function.
Definition: generic_RVD.h:2418
bool facet_seed_is_visited(index_t f, index_t s) const
Tests whether a (facet,seed) couple was visited.
Definition: generic_RVD.h:1821
void compute_surfacic(const ACTION &action)
Low-level API of Restricted Voronoi Diagram traversal.
Definition: generic_RVD.h:1314
void set_facets_range(index_t facets_begin, index_t facets_end)
Sets the facets range.
Definition: generic_RVD.h:226
void compute_surfacic_with_seeds_priority(const ACTION &action)
Low-level API of Restricted Voronoi Diagram traversal with seeds priority in surfacic mode.
Definition: generic_RVD.h:1339
index_t nb_tetrahedra_in_range() const
Gets the number of tetrahedra in the current range.
Definition: generic_RVD.h:258
bool check_SR() const
Tests whether radius of security is enforced.
Definition: generic_RVD.h:365
void compute_surfacic_with_cnx_priority(const ACTION &action)
Low-level API of Restricted Voronoi Diagram traversal with connected components priority.
Definition: generic_RVD.h:1879
index_t find_seed_near_tet(index_t t)
Finds a seed near a given tetrahedron.
Definition: generic_RVD.h:2047
void clip_by_plane(Polygon &ping, Polygon &pong, index_t i, index_t j)
Computes the intersection between a polygon and a half-space.
Definition: generic_RVD.h:2263
const Delaunay * delaunay() const
Gets the Delaunay triangulation.
Definition: generic_RVD.h:192
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.
Definition: generic_RVD.h:1222
void set_tetrahedra_range(index_t tets_begin, index_t tets_end)
Sets the tetrahedra range.
Definition: generic_RVD.h:240
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.
Definition: attributes.h:1190
Delaunay interface for NearestNeighbors search.
Definition: delaunay_nn.h:63
NearestNeighborSearch * nn_search()
Gets the NearestNeighborSearch used internally.
Definition: delaunay_nn.h:92
virtual void enlarge_neighborhood(index_t i, index_t nb)
Stores nb neighbors with vertex i.
Abstract interface for Delaunay triangulation in Nd.
Definition: delaunay.h:71
signed_index_t cell_vertex(index_t c, index_t lv) const
Gets a vertex index by cell index and local vertex index.
Definition: delaunay.h:365
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.
Definition: delaunay.h:233
index_t nb_cells() const
Gets the number of cells.
Definition: delaunay.h:319
void get_neighbors(index_t v, vector< index_t > &neighbors) const
Gets the one-ring neighbors of vertex v.
Definition: delaunay.h:481
coord_index_t dimension() const
Gets the dimension of this Delaunay.
Definition: delaunay.h:171
index_t cell_size() const
Gets the number of vertices in each cell.
Definition: delaunay.h:180
index_t nb_vertices() const
Gets the number of vertices.
Definition: delaunay.h:242
signed_index_t vertex_cell(index_t v) const
Gets an incident cell index by a vertex index.
Definition: delaunay.h:449
index_t index(index_t c, signed_index_t v) const
Retrieves a local vertex index from cell index and global vertex index.
Definition: delaunay.h:411
signed_index_t next_around_vertex(index_t c, index_t lv) const
Traverses the list of cells incident to a vertex.
Definition: delaunay.h:465
index_t vertex(index_t f, index_t lv) const
Gets a vertex by facet and local vertex index.
Definition: mesh.h:1054
index_t nb() const
Gets the number of (sub-)elements.
Definition: mesh.h:89
AttributesManager & attributes() const
Gets the attributes manager.
Definition: mesh.h:100
const double * point_ptr(index_t v) const
Gets a point.
Definition: mesh.h:455
Represents a mesh.
Definition: mesh.h:2693
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.
Definition: memory.h:795
index_t size() const
Gets the number of elements.
Definition: memory.h:674
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.
Definition: index.h:605
#define geo_dim_alignment(dim)
Gets a point alignment.
Definition: memory.h:270
#define geo_assume_aligned(var, alignment)
Informs the compiler that a given pointer is memory-aligned.
Definition: memory.h:374
double distance2(const COORD_T *p1, const COORD_T *p2, coord_index_t dim)
Computes the squared distance between two nd points.
Definition: geometry_nd.h:65
Global Vorpaline namespace.
Definition: basic.h:55
index_t max_index_t()
Gets the maximum positive value of type index_t.
Definition: numeric.h:334
void geo_argused(const T &)
Suppresses compiler warnings about unused parameters.
Definition: argused.h:60
geo_signed_index_t signed_index_t
The type for storing and manipulating indices differences.
Definition: numeric.h:343
GEOGen::SymbolicVertex SymbolicVertex
Symbolic representation of a RestrictedVoronoiDiagram vertex.
Definition: generic_RVD.h:2520
geo_index_t index_t
The type for storing and manipulating indices.
Definition: numeric.h:329
geo_coord_index_t coord_index_t
The type for storing coordinate indices, and iterating on the coordinates of a point.
Definition: numeric.h:363
Types and functions for numbers manipulation.
Filtered exact predicates for restricted Voronoi diagrams.
Function and classes for process manipulation.
A (facet,seed) pair.