CGAL::Delaunay_d< R, Lifted_R >

Definition

An instance DT of type Delaunay_d< R, Lifted_R > is the nearest and furthest site Delaunay triangulation of a set S of points in some d-dimensional space. We call S the underlying point set and d or dim the dimension of the underlying space. We use dcur to denote the affine dimension of S. The data type supports incremental construction of Delaunay triangulations and various kind of query operations (in particular, nearest and furthest neighbor queries and range queries with spheres and simplices).

A Delaunay triangulation is a simplicial complex. All simplices in the Delaunay triangulation have dimension dcur. In the nearest site Delaunay triangulation the circumsphere of any simplex in the triangulation contains no point of S in its interior. In the furthest site Delaunay triangulation the circumsphere of any simplex contains no point of S in its exterior. If the points in S are co-circular then any triangulation of S is a nearest as well as a furthest site Delaunay triangulation of S. If the points in S are not co-circular then no simplex can be a simplex of both triangulations. Accordingly, we view DT as either one or two collection(s) of simplices. If the points in S are co-circular there is just one collection: the set of simplices of some triangulation. If the points in S are not co-circular there are two collections. One collection consists of the simplices of a nearest site Delaunay triangulation and the other collection consists of the simplices of a furthest site Delaunay triangulation.

For each simplex of maximal dimension there is a handle of type Simplex_handle and for each vertex of the triangulation there is a handle of type Vertex_handle. Each simplex has 1 + dcur vertices indexed from 0 to dcur. For any simplex s and any index i, DT.vertex_of(s,i) returns the i-th vertex of s. There may or may not be a simplex t opposite to the vertex of s with index i. The function DT.opposite_simplex(s,i) returns t if it exists and returns Simplex_handle() otherwise. If t exists then s and t share dcur vertices, namely all but the vertex with index i of s and the vertex with index DT.index_of_vertex_in_opposite_simplex(s,i) of t. Assume that t = DT.opposite_simplex(s,i) exists and let j = DT.index_of_vertex_in_opposite_simplex(s,i). Then s = DT.opposite_simplex(t,j) and i = DT.index_of_vertex_in_opposite_simplex(t,j). In general, a vertex belongs to many simplices.

Any simplex of DT belongs either to the nearest or to the furthest site Delaunay triangulation or both. The test DT.simplex_of_nearest(dt_simplex s) returns true if s belongs to the nearest site triangulation and the test DT.simplex_of_furthest(dt_simplex s) returns true if s belongs to the furthest site triangulation.

Inherits From

Convex_hull_d<Lifted_R>

Types

Delaunay_d< R, Lifted_R >::Simplex_handle
handles to the simplices of the complex.


Delaunay_d< R, Lifted_R >::Vertex_handle
handles to vertices of the complex.


Delaunay_d< R, Lifted_R >::Point_d
the point type


Delaunay_d< R, Lifted_R >::Sphere_d
the sphere type


enum Delaunay_voronoi_kind { NEAREST, FURTHEST};
interface flags

To use these types you can typedef them into the global scope after instantiation of the class. We use Vertex_handle instead of Delaunay_d< R, Lifted_R >::Vertex_handle from now on. Similarly we use Simplex_handle.

Delaunay_d< R, Lifted_R >::Point_const_iterator
the iterator for points.


Delaunay_d< R, Lifted_R >::Vertex_iterator
the iterator for vertices.


Delaunay_d< R, Lifted_R >::Simplex_iterator
the iterator for simplices.

Creation

Delaunay_d< R, Lifted_R > DT ( int d,
R k1 = R(),
Lifted_R k2 = Lifted_R());
creates an instance DT of type Delaunay_d. The dimension of the underlying space is d and S is initialized to the empty point set. The traits class R specifies the models of all types and the implementations of all geometric primitives used by the Delaunay class. The traits class Lifted_R specifies the models of all types and the implementations of all geometric primitives used by the base class of Delaunay_d< R, Lifted_R >. The second template parameter defaults to the first: Delaunay_d<R> = Delaunay_d<R, Lifted_R = R >.

The data type Delaunay_d offers neither copy constructor nor assignment operator.

Requirements

R is a model of the concept DelaunayTraits_d . Lifted_R is a model of the concept DelaunayLiftedTraits_d .

Operations

All operations below that take a point x as an argument have the common precondition that x.dimension() = DT.dimension().

int DT.dimension () returns the dimension of ambient space

int DT.current_dimension ()
returns the affine dimension of the current point set, i.e., -1 is S is empty, 0 if S consists of a single point, 1 if all points of S lie on a common line, etcetera.

bool DT.is_simplex_of_nearest ( Simplex_handle s)
returns true if s is a simplex of the nearest site triangulation.

bool DT.is_simplex_of_furthest ( Simplex_handle s)
returns true if s is a simplex of the furthest site triangulation.

Vertex_handle DT.vertex_of_simplex ( Simplex_handle s, int i)
returns the vertex associated with the i-th node of s.
Precondition: 0 i dcur.

Point_d DT.associated_point ( Vertex_handle v)
returns the point associated with vertex v.

Point_d DT.point_of_simplex ( Simplex_handle s, int i)
returns the point associated with the i-th vertex of s.
Precondition: 0 i dcur.

Simplex_handle DT.opposite_simplex ( Simplex_handle s, int i)
returns the simplex opposite to the i-th vertex of s (Simplex_handle() if there is no such simplex).
Precondition: 0 i dcur.

int
DT.index_of_vertex_in_opposite_simplex ( Simplex_handle s,
int i)
returns the index of the vertex opposite to the i-th vertex of s.
Precondition: 0 i dcur.

Simplex_handle DT.simplex ( Vertex_handle v)
returns a simplex of the nearest site triangulation incident to v.

int DT.index ( Vertex_handle v)
returns the index of v in DT.simplex(v).

bool DT.contains ( Simplex_handle s, Point_d x)
returns true if x is contained in the closure of simplex s.

bool DT.empty () decides whether DT is empty.

void DT.clear () reinitializes DT to the empty Delaunay triangulation.

Vertex_handle DT.insert ( Point_d x)
inserts point x into DT and returns the corresponding Vertex_handle. More precisely, if there is already a vertex v in DT positioned at x (i.e., associated_point(v) is equal to x) then associated_point(v) is changed to x (i.e., associated_point(v) is made identical to x) and if there is no such vertex then a new vertex v with associated_point(v) = x is added to DT. In either case, v is returned.

Simplex_handle DT.locate ( Point_d x)
returns a simplex of the nearest site triangulation containing x in its closure (returns Simplex_handle() if x lies outside the convex hull of S).

Vertex_handle DT.lookup ( Point_d x)
if DT contains a vertex v with associated_point(v) = x the result is v otherwise the result is Vertex_handle().

Vertex_handle DT.nearest_neighbor ( Point_d x)
computes a vertex v of DT that is closest to x, i.e.,
dist(x,associated_point(v)) = min{ dist(x, associated_point(u)) u S }.

std::list<Vertex_handle>
DT.range_search ( Sphere_d C)
returns the list of all vertices contained in the closure of sphere C.

std::list<Vertex_handle>
DT.range_search ( std::vector<Point_d> A)
returns the list of all vertices contained in the closure of the simplex whose corners are given by A.
Precondition: A must consist of d+1 affinely independent points in base space.

std::list<Simplex_handle>
DT.all_simplices ( Delaunay_voronoi_kind k = NEAREST)
returns a list of all simplices of either the nearest or the furthest site Delaunay triangulation of S.

std::list<Vertex_handle>
DT.all_vertices ( Delaunay_voronoi_kind k = NEAREST)
returns a list of all vertices of either the nearest or the furthest site Delaunay triangulation of S.

std::list<Point_d> DT.all_points () returns S.

Point_const_iterator
DT.points_begin () returns the start iterator for points in DT.

Point_const_iterator
DT.points_end () returns the past the end iterator for points in DT.

Simplex_iterator DT.simplices_begin ( Delaunay_voronoi_kind k = NEAREST)
returns the start iterator for simplices of DT.

Simplex_iterator DT.simplices_end ()
returns the past the end iterator for simplices of DT.

Implementation

The data type is derived from Convex_hull_d via the lifting map. For a point x in d-dimensional space let lift(x) be its lifting to the unit paraboloid of revolution. There is an intimate relationship between the Delaunay triangulation of a point set S and the convex hull of lift(S): The nearest site Delaunay triangulation is the projection of the lower hull and the furthest site Delaunay triangulation is the upper hull. For implementation details we refer the reader to the implementation report available from the CGAL server.

The space requirement is the same as for convex hulls. The time requirement for an insert is the time to insert the lifted point into the convex hull of the lifted points.

Example

The abstract data type Delaunay_d has a default instantiation by means of the d-dimensional geometric kernel.

#include <CGAL/Homogeneous_d.h>
#include <CGAL/leda_integer.h>
#include <CGAL/Delaunay_d.h>

typedef leda_integer RT;
typedef CGAL::Homogeneous_d<RT> Kernel;
typedef CGAL::Delaunay_d<Kernel> Delaunay_d;
typedef Delaunay_d::Point_d Point;
typedef Delaunay_d::Simplex_handle Simplex_handle;
typedef Delaunay_d::Vertex_handle Vertex_handle;

int main()
{
  Delaunay_d T(2);
  Vertex_handle v1 = T.insert(Point_d(2,11));
  ...
}

Traits requirements

Delaunay_d< R, Lifted_R > requires the following types from the kernel traits Lifted_R:

  RT Point_d Vector_d Ray_d Hyperplane_d 
and uses the following function objects from the kernel traits:
  Construct_hyperplane_d
  Construct_vector_d
  Vector_to_point_d / Point_to_vector_d
  Orientation_d
  Orthogonal_vector_d
  Oriented_side_d / Has_on_positive_side_d
  Affinely_independent_d
  Contained_in_simplex_d
  Contained_in_affine_hull_d
  Intersect_d
  Lift_to_paraboloid_d / Project_along_d_axis_d
  Component_accessor_d
Delaunay_d< R, Lifted_R > requires the following types from the kernel traits R:
  FT Point_d Sphere_d 
and uses the following function objects from the kernel traits R:
  Construct_sphere_d
  Squared_distance_d
  Point_of_sphere_d
  Affinely_independent_d
  Contained_in_simplex_d

Low Dimensional Output Routines

include <CGAL/IO/Delaunay_d_window_stream.h>

template <typename R, typename Lifted_R>
template <typename R, typename Lifted_R> void
d2_show ( Delaunay_d<R,Lifted_R> D,
CGAL::Window_stream& W,
typename Delaunay_d<R,Lifted_R>::Delaunay_voronoi_kind k = Delaunay_d<R,Lifted_R>::NEAREST)
draws the underlying simplicial complex D into window W.

Precondition: dim == 2.

template <typename R, typename Lifted_R>
template <typename R, typename Lifted_R> void
d2_map ( Delaunay_d<R,Lifted_R> D,
GRAPH< typename Delaunay_d<R,Lifted_R>::Point_d, int >& DTG,
typename Delaunay_d<R,Lifted_R>::Delaunay_voronoi_kind k = Delaunay_d<R,Lifted_R>::NEAREST)
constructs a LEDA graph representation of the nearest (kind = NEAREST or the furthest (kind = FURTHEST) site Delaunay triangulation.

Precondition: dim() == 2.