CGAL 5.1 - 2D Hyperbolic Delaunay Triangulations
|
This package enables the computation of Delaunay triangulations of point sets in the Poincaré disk model of the hyperbolic plane.
The Poincaré disk model represents the hyperbolic plane \(\mathbb H^2\) in the unit disk centered at the origin in the Euclidean plane: points of \(\mathbb H^2\) lie in the interior of the disk, while its boundary, the unit circle, represents the set \(\mathcal H_\infty\) of points at infinity.
In this model, a hyperbolic line is either an arc of circle perpendicular to the unit circle or, if it passes through the origin, a diameter of the unit disk. A hyperbolic circle is a Euclidean circle contained in the unit disk; however, its hyperbolic center and radius are not the same as its Euclidean center and radius.
As hyperbolic circles coincide with Euclidean circles contained in the unit disk, the combinatorial structure of the hyperbolic Delaunay triangulation of a set \(\mathcal P\) of points in \(\mathbb H^2\) is a subset of the Euclidean Delaunay triangulation of \(\mathcal P\) (see Figure 42.2 - Left). More precisely, the hyperbolic Delaunay triangulation of \(\mathcal P\) only contains the simplices of the Euclidean Delaunay triangulation that are hyperbolic:
In the Euclidean Delaunay triangulation, there is a bijection between non-hyperbolic faces and non-hyperbolic edges [1], illustrated by Figure 42.2 - Right.
The hyperbolic Delaunay triangulation is a simplicial complex, i.e., a set of simplices such that
Moreover, it is connected [1].
From what was said above, it is natural that the class Hyperbolic_Delaunay_triangulation_2
privately inherits from the class Delaunay_triangulation_2
. Consequently, users are encouraged to look at Chapter 2D Triangulation of the CGAL manual to know more in particular about the representation of triangulations in CGAL and the flexibility of the design.
The class Hyperbolic_Delaunay_triangulation_2
has two template parameters:
Gt
, which provides geometric primitives. The requirements on this first template parameter are described by the concept HyperbolicDelaunayTriangulationTraits_2
, which refines DelaunayTriangulationTraits_2
. TriangulationDataStructure_2
. The default for this second template parameter is Triangulation_data_structure_2< Triangulation_vertex_base_2<Gt>, Hyperbolic_triangulation_face_base_2<Gt> >
. Two models of the concept HyperbolicDelaunayTriangulationTraits_2
are proposed for the geometric traits class. The first one, CGAL::Hyperbolic_Delaunay_triangulation_CK_traits_2
, is based upon CGAL::Circular_kernel_2
and guarantees exact constructions of Delaunay triangulations and dual objects when the input points have rational coordinates. The second one, CGAL::Hyperbolic_Delaunay_triangulation_traits_2
, is more general, as it guarantees exact constructions even for input points with algebraic coordinates; however the first model is more efficient for rational points.
The example below shows insertion of random points in a hyperbolic Delaunay triangulation. The same set of points is inserted twice. The first time points are inserted one by one, which causes Euclidean faces to be filtered at each insertion. The second time, all points are inserted and the filtering is done once at the end.
File Hyperbolic_triangulation_2/ht2_example.cpp
The example below shows how user-defined info can be added to the hyperbolic faces.
File Hyperbolic_triangulation_2/ht2_example_color.cpp
We have measured the insertion execution time of our implementation with both traits classes CGAL::Hyperbolic_Delaunay_triangulation_CK_traits_2
and CGAL::Hyperbolic_Delaunay_triangulation_traits_2
with their default template parameters against the insertion time in a Euclidean CGAL triangulation. We generate 1 million random points, uniformly distributed in the unit disk with respect to the Euclidean metric. We insert the same set of points in three triangulations:
CGAL::Hyperbolic_Delaunay_triangulation_traits_2
(CORE traits) as traits class; CGAL::Hyperbolic_Delaunay_triangulation_CK_traits_2
(CK traits) as traits class; CGAL::Exact_predicates_inexact_constructions_kernel
(EPICK) as traits class. We create two instances of each type of triangulation. In one instance we insert the points one by one, which causes non-hyperbolic faces to be filtered out at each insertion. In the other instance we insert the points via iterator input, which causes non-hyperbolic faces to be filtered only once, after all points have been inserted. We report results averaged over 10 executions of this experiment in Table 1 below. The experiments have been executed on two machines:
Triangulation type | Machine 1 | Machine 2 | ||
---|---|---|---|---|
Sequential insertion | Iterator insertion | Sequential insertion | Iterator insertion | |
Hyperbolic (CORE traits) | 955 sec. | 23 sec. | 884 sec. | 20 sec. |
Hyperbolic (CK traits) | 330 sec. | 1 sec. | 289 sec. | 1 sec. |
Euclidean (EPICK) | 131 sec. | 0.71 sec. | 114 sec. | 0.68 sec. |
This package implements the algorithms for computing Delaunay complexes in the hyperbolic plane, described by Mikhail Bogdanov, Olivier Devillers, and Monique Teillaud [1].
Mikhail Bogdanov wrote most of the code. Iordan Iordanov added the traits class CGAL::Hyperbolic_Delaunay_triangulation_traits_2
and worked on the documentations. Both were PhD candidates advised by Monique Teillaud.
Authors acknowledge partial support from ANR SoS.