CGAL 5.5.1 - 2D Hyperbolic Delaunay Triangulations
User Manual

Author
Mikhail Bogdanov, Iordan Iordanov, and Monique Teillaud

This package enables the computation of Delaunay triangulations of point sets in the Poincaré disk model of the hyperbolic plane.

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.

Figure 43.1 The Poincaré disk model for the hyperbolic plane. The figure shows two hyperbolic lines and three concentric hyperbolic circles with different radii.

Euclidean and Hyperbolic Delaunay Triangulations

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 43.2 - Left). More precisely, the hyperbolic Delaunay triangulation of \(\mathcal P\) only contains the simplices of the Euclidean Delaunay triangulation that are hyperbolic:

  • A Euclidean Delaunay face is hyperbolic if its circumscribing circle is contained in \(\mathbb H^2\).
  • A Euclidean Delaunay edge is hyperbolic if one of the empty disks (i.e., not containing any point of \(\mathcal P\)) passing through its endpoints is contained in \(\mathbb H^2\).

In the Euclidean Delaunay triangulation, there is a bijection between non-hyperbolic faces and non-hyperbolic edges [1], illustrated by Figure 43.2 - Right.

Figure 43.2 Left: The Euclidean (red) and hyperbolic (black) Delaunay triangulations of a given set of points in the unit disk. Only the colored faces are faces of the hyperbolic Delaunay triangulation. The hyperbolic and Euclidean geometric embeddings of a Delaunay face that exists in both triangulations are different. Right: The shaded face is non-hyperbolic. Its dashed edge is non-hyperbolic, as no empty circle through its endpoints is contained in \(\mathbb H^2\). Its other two edges are hyperbolic.

The hyperbolic Delaunay triangulation is a simplicial complex, i.e., a set of simplices such that

  • any face of a simplex is a simplex,
  • two simplices either are disjoint or share a common face.

Moreover, it is connected [1].

Software Design

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:

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.

Examples

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

#include <CGAL/Hyperbolic_Delaunay_triangulation_2.h>
#include <CGAL/Hyperbolic_Delaunay_triangulation_traits_2.h>
#include <CGAL/point_generators_2.h>
#include <CGAL/Timer.h>
#include <iostream>
#include <vector>
typedef Gt::Point_2 Point_2;
typedef Gt::Circle_2 Circle_2;
typedef Gt::FT FT;
int main(int argc, char** argv)
{
int N;
if(argc < 2) {
std::cout << "usage: " << argv[0] << " [number_of_points]" << std::endl;
std::cout << "Defaulting to 100k points..." << std::endl;
N = 100000;
} else {
N = atoi(argv[1]);
}
CGAL::Timer timer;
CGAL::Random_points_in_disc_2<Point_2, Creator> in_disc;
std::vector<Point_2> pts;
std::vector<Point_2>::iterator ip;
for(int i=0; i<N; ++i)
pts.push_back(*(in_disc++));
Dt dt_during;
std::cout << "Insertion of points one by one (hyperbolic filtering at each step)" << std::endl;
std::cout << "===================================================================" << std::endl;
timer.start();
for(ip = pts.begin(); ip != pts.end(); ++ip)
dt_during.insert(*ip);
timer.stop();
std::cout << "Number of vertices: " << dt_during.number_of_vertices() << std::endl;
std::cout << "Number of hyperbolic faces: " << dt_during.number_of_hyperbolic_faces() << std::endl;
std::cout << "Number of hyperbolic edges: " << dt_during.number_of_hyperbolic_edges() << std::endl;
std::cout << "Time: " << timer.time() << std::endl << std::endl;
Dt dt_end;
std::cout << "Insertion of point set (hyperbolic filtering only once at the end)" << std::endl;
std::cout << "===================================================================" << std::endl;
timer.reset();
timer.start();
dt_end.insert(pts.begin(),pts.end());
timer.stop();
std::cout << "Number of vertices: " << dt_end.number_of_vertices() << std::endl;
std::cout << "Number of hyperbolic faces: " << dt_end.number_of_hyperbolic_faces() << std::endl;
std::cout << "Number of hyperbolic edges: " << dt_end.number_of_hyperbolic_edges() << std::endl;
std::cout << "Time: " << timer.time() << std::endl;
return EXIT_SUCCESS;
}

The example below shows how user-defined info can be added to the hyperbolic faces.
File Hyperbolic_triangulation_2/ht2_example_color.cpp

#include <CGAL/Hyperbolic_triangulation_face_base_2.h>
#include <CGAL/Hyperbolic_Delaunay_triangulation_2.h>
#include <CGAL/Hyperbolic_Delaunay_triangulation_CK_traits_2.h>
#include <CGAL/Triangulation_face_base_with_info_2.h>
#include <CGAL/Triangulation_vertex_base_2.h>
#include <CGAL/point_generators_2.h>
#include <CGAL/IO/Color.h>
#include <iostream>
#include <vector>
Hyperbolic_face_base> Face_base_with_info;
typedef CGAL::Triangulation_vertex_base_2<Gt> Vertex_base;
Face_base_with_info> TDS;
typedef Dt::Point Point_2;
int main(int argc, char** argv)
{
int N;
if(argc < 2) {
std::cout << "usage: " << argv[0] << " [number_of_points]" << std::endl;
std::cout << "Defaulting to 100k points..." << std::endl;
N = 100000;
} else {
N = atoi(argv[1]);
}
std::cout << "Number of points: " << N << std::endl;
CGAL::Random_points_in_disc_2<Point_2, Creator> in_disc;
std::vector<Point_2> pts;
for(int i=0; i<N; ++i)
pts.push_back(*(in_disc++));
Dt dt;
dt.insert(pts.begin(), pts.end());
Dt::Vertex_handle vo = dt.insert(Point_2(0,0));
int origin_faces = 0;
Dt::Hyperbolic_faces_iterator fit;
for(fit = dt.hyperbolic_faces_begin(); fit != dt.hyperbolic_faces_end(); ++fit)
{
if(fit->has_vertex(vo))
{
fit->info() = CGAL::IO::red();
origin_faces++;
}
}
int red_faces = 0;
for(fit = dt.hyperbolic_faces_begin(); fit != dt.hyperbolic_faces_end(); ++fit)
{
if(fit->info() == CGAL::IO::red())
{
red_faces++;
}
}
assert(red_faces == origin_faces);
std::cout << "Number of points " << N << std::endl;
std::cout << "Number of vertices: " << dt.number_of_vertices() << std::endl;
std::cout << "Number of hyperbolic faces: " << dt.number_of_hyperbolic_faces() << std::endl;
std::cout << "Number of faces incident to the origin: " << origin_faces << std::endl;
return EXIT_SUCCESS;
}

Performance

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:

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:

  • Machine 1: MacBook Pro (2015), CPU: Intel Core i5 @ 2.9 GHz, RAM: 16 GB @ 1867 MHz, OS: Mac OS X (10.10.5), Compiler: gcc version 7.3.0;
  • Machine 2: Dell Vostro 5471 (2018), CPU: Intel Core i5 @ 1.6 GHZ (up to 3.4 GHz in TurboMode), RAM: 8GB @ 2400 MHz, OS: Ubuntu 18.04 (kernel 4.15.0), Compiler: gcc version 7.3.0.
Table 1: Comparison of insertion times of 1 million random points
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.

Design and Implementation History

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.