\( \newcommand{\E}{\mathrm{E}} \) \( \newcommand{\A}{\mathrm{A}} \) \( \newcommand{\R}{\mathrm{R}} \) \( \newcommand{\N}{\mathrm{N}} \) \( \newcommand{\Q}{\mathrm{Q}} \) \( \newcommand{\Z}{\mathrm{Z}} \) \( \def\ccSum #1#2#3{ \sum_{#1}^{#2}{#3} } \def\ccProd #1#2#3{ \sum_{#1}^{#2}{#3} }\)
CGAL 4.7 - dD Convex Hulls and Delaunay Triangulations
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Groups Pages
CGAL::Delaunay_d< R, Lifted_R > Class Template Reference

#include <CGAL/Delaunay_d.h>

Inherits from

CGAL::Convex_hull_d< Lifted_R >.

Definition

Deprecated:
This package is deprecated since the version 4.6 of CGAL. The package dD Triangulations should be used instead.

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.

Template Parameters
Rmust be a model of the concept DelaunayTraits_d.
Lifted_Rmust be a model of the concept DelaunayLiftedTraits_d.

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 Delaunay_d::Point_d Point;
int main()
{
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:

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
  • Construct_sphere_d
  • Squared_distance_d
  • Point_of_sphere_d
  • Affinely_independent_d
  • Contained_in_simplex_d

Related Functions

(Note that these are not member functions.)

template<typename R , typename Lifted_R >
void d2_map (const 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. More...
 

Types

enum  Delaunay_voronoi_kind { NEAREST, FURTHEST }
 interface flags More...
 
typedef unspecified_type Simplex_handle
 handles to the simplices of the complex.
 
typedef unspecified_type Vertex_handle
 handles to vertices of the complex.
 
typedef unspecified_type Point_d
 the point type More...
 
typedef unspecified_type Sphere_d
 the sphere type More...
 
typedef unspecified_type Point_const_iterator
 the iterator for points.
 
typedef unspecified_type Vertex_iterator
 the iterator for vertices.
 
typedef unspecified_type Simplex_iterator
 the iterator for simplices.
 

Creation

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

 Delaunay_d (int d, R k1=R(), Lifted_R k2=Lifted_R())
 creates an instance DT of type Delaunay_d. More...
 

Operations

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

int dimension ()
 returns the dimension of ambient space.
 
int 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, etc.
 
bool is_simplex_of_nearest (Simplex_handle s)
 returns true if s is a simplex of the nearest site triangulation.
 
bool is_simplex_of_furthest (Simplex_handle s)
 returns true if s is a simplex of the furthest site triangulation.
 
Vertex_handle vertex_of_simplex (Simplex_handle s, int i)
 returns the vertex associated with the i-th node of s. More...
 
Point_d associated_point (Vertex_handle v)
 returns the point associated with vertex v.
 
Point_d point_of_simplex (Simplex_handle s, int i)
 returns the point associated with the i-th vertex of s. More...
 
Simplex_handle 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). More...
 
int 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. More...
 
Simplex_handle simplex (Vertex_handle v)
 returns a simplex of the nearest site triangulation incident to v.
 
int index (Vertex_handle v)
 returns the index of v in DT.simplex(v).
 
bool contains (Simplex_handle s, const Point_d &x)
 returns true if x is contained in the closure of simplex s.
 
bool empty ()
 decides whether DT is empty.
 
void clear ()
 re-initializes DT to the empty Delaunay triangulation.
 
Vertex_handle insert (const Point_d &x)
 inserts point x into DT and returns the corresponding Vertex_handle. More...
 
Simplex_handle locate (const 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 lookup (const 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 nearest_neighbor (const 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 \(\in S\) }.
 
std::list< Vertex_handlerange_search (const Sphere_d &C)
 returns the list of all vertices contained in the closure of sphere \( C\).
 
std::list< Vertex_handlerange_search (const std::vector< Point_d > &A)
 returns the list of all vertices contained in the closure of the simplex whose corners are given by A. More...
 
std::list< Simplex_handleall_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_handleall_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_dall_points ()
 returns \( S\).
 
Point_const_iterator points_begin ()
 returns the start iterator for points in DT.
 
Point_const_iterator points_end ()
 returns the past the end iterator for points in DT.
 
Simplex_iterator simplices_begin (Delaunay_voronoi_kind k=NEAREST)
 returns the start iterator for simplices of DT.
 
Simplex_iterator simplices_end ()
 returns the past the end iterator for simplices of DT.
 

Additional Inherited Members

- Private Types inherited from CGAL::Convex_hull_d< Lifted_R >
typedef unspecified_type R
 the representation class.
 
typedef unspecified_type Point_d
 the point type.
 
typedef unspecified_type Hyperplane_d
 the hyperplane type.
 
typedef unspecified_type Simplex_handle
 handle for simplices.
 
typedef unspecified_type Facet_handle
 handle for facets.
 
typedef unspecified_type Vertex_handle
 handle for vertices.
 
typedef unspecified_type Simplex_iterator
 iterator for simplices.
 
typedef unspecified_type Facet_iterator
 iterator for facets.
 
typedef unspecified_type Vertex_iterator
 iterator for vertices.
 
typedef unspecified_type Hull_vertex_iterator
 iterator for vertices that are part of the convex hull.
 
typedef unspecified_type Point_const_iterator
 const iterator for all inserted points.
 
typedef unspecified_type Hull_point_const_iterator
 const iterator for all points that are part of the convex hull.
 
- Private Member Functions inherited from CGAL::Convex_hull_d< Lifted_R >
 Convex_hull_d (int d, Lifted_RKernel=Lifted_R())
 creates an instance C of type Convex_hull_d. More...
 
int dimension ()
 returns the dimension of ambient space.
 
int current_dimension ()
 returns the affine dimension dcur of \( S\).
 
Point_d associated_point (Vertex_handle v)
 returns the point associated with vertex \( v\).
 
Vertex_handle vertex_of_simplex (Simplex_handle s, int i)
 returns the vertex corresponding to the \( i\)-th vertex of \( s\). More...
 
Point_d point_of_simplex (Simplex_handle s, int i)
 same as C.associated_point(C.vertex_of_simplex(s,i)).
 
Simplex_handle 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). More...
 
int 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\). More...
 
Simplex_handle simplex (Vertex_handle v)
 returns a simplex of which \( v\) is a node. More...
 
int index (Vertex_handle v)
 returns the index of \( v\) in simplex(v).
 
Vertex_handle vertex_of_facet (Facet_handle f, int i)
 returns the vertex corresponding to the \( i\)-th vertex of \( f\). More...
 
Point_d point_of_facet (Facet_handle f, int i)
 same as C.associated_point(C.vertex_of_facet(f,i)).
 
Facet_handle opposite_facet (Facet_handle f, int i)
 returns the facet opposite to the \( i\)-th vertex of \( f\) (Facet_handle() if there is no such facet). More...
 
int index_of_vertex_in_opposite_facet (Facet_handle f, int i)
 returns the index of the vertex opposite to the \( i\)-th vertex of \( f\). More...
 
Hyperplane_d hyperplane_supporting (Facet_handle f)
 returns a hyperplane supporting facet f. More...
 
Vertex_handle insert (const Point_d &x)
 adds point x to the underlying set of points. More...
 
void insert (Forward_iterator first, Forward_iterator last)
 adds S = set [first,last) to the underlying set of points. More...
 
bool is_dimension_jump (const Point_d &x)
 returns true if \( x\) is not contained in the affine hull of S.
 
std::list< Facet_handlefacets_visible_from (const Point_d &x)
 returns the list of all facets that are visible from x. More...
 
Bounded_side bounded_side (const Point_d &x)
 returns ON_BOUNDED_SIDE (ON_BOUNDARY,ON_UNBOUNDED_SIDE) if x is contained in the interior (lies on the boundary, is contained in the exterior) of C. More...
 
void clear (int d)
 re-initializes C to an empty hull in \( d\)-dimensional space.
 
int number_of_vertices ()
 returns the number of vertices of C.
 
int number_of_facets ()
 returns the number of facets of C.
 
int number_of_simplices ()
 returns the number of bounded simplices of C.
 
void print_statistics ()
 gives information about the size of the current hull and the number of visibility tests performed.
 
bool is_valid (bool throw_exceptions=false)
 checks the validity of the data structure. More...
 
Vertex_iterator vertices_begin ()
 an iterator of C to start the iteration over all vertices of the complex.
 
Vertex_iterator vertices_end ()
 the past the end iterator for vertices.
 
Simplex_iterator simplices_begin ()
 an iterator of C to start the iteration over all simplices of the complex.
 
Simplex_iterator simplices_end ()
 the past the end iterator for simplices.
 
Facet_iterator facets_begin ()
 an iterator of C to start the iteration over all facets of the complex.
 
Facet_iterator facets_end ()
 the past the end iterator for facets.
 
Hull_vertex_iterator hull_vertices_begin ()
 an iterator to start the iteration over all vertices of C that are part of the convex hull.
 
Hull_vertex_iterator hull_vertices_end ()
 the past the end iterator for hull vertices.
 
Point_const_iterator points_begin ()
 returns the start iterator for all points that have been inserted to construct C.
 
Point_const_iterator points_end ()
 returns the past the end iterator for points.
 
Hull_point_const_iterator hull_points_begin ()
 returns an iterator to start the iteration over all points in the convex hull C. More...
 
Hull_point_const_iterator hull_points_end ()
 returns the past the end iterator for points in the convex hull.
 
void visit_all_facets (const Visitor &V)
 each facet of C is visited by the visitor object V. More...
 
const std::list< Point_d > & all_points ()
 returns a list of all points that have been inserted to construct C.
 
std::list< Vertex_handleall_vertices ()
 returns a list of all vertices of C (also interior ones).
 
std::list< Simplex_handleall_simplices ()
 returns a list of all simplices in C.
 
std::list< Facet_handleall_facets ()
 returns a list of all facets of C.
 

Member Typedef Documentation

template<typename R , typename Lifted_R >
typedef unspecified_type CGAL::Delaunay_d< R, Lifted_R >::Point_d

the point type

template<typename R , typename Lifted_R >
typedef unspecified_type CGAL::Delaunay_d< R, Lifted_R >::Sphere_d

the sphere type

Member Enumeration Documentation

template<typename R , typename Lifted_R >
enum CGAL::Delaunay_d::Delaunay_voronoi_kind

interface flags

Enumerator
NEAREST 
FURTHEST 

Constructor & Destructor Documentation

template<typename R , typename Lifted_R >
CGAL::Delaunay_d< R, Lifted_R >::Delaunay_d ( 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 >.

Member Function Documentation

template<typename R , typename Lifted_R >
int CGAL::Delaunay_d< R, Lifted_R >::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.
template<typename R , typename Lifted_R >
Vertex_handle CGAL::Delaunay_d< R, Lifted_R >::insert ( const 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.

template<typename R , typename Lifted_R >
Simplex_handle CGAL::Delaunay_d< R, Lifted_R >::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.
template<typename R , typename Lifted_R >
Point_d CGAL::Delaunay_d< R, Lifted_R >::point_of_simplex ( Simplex_handle  s,
int  i 
)

returns the point associated with the i-th vertex of s.

Precondition
0 <= i <= dcur.
template<typename R , typename Lifted_R >
std::list<Vertex_handle> CGAL::Delaunay_d< R, Lifted_R >::range_search ( const 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.
template<typename R , typename Lifted_R >
Vertex_handle CGAL::Delaunay_d< R, Lifted_R >::vertex_of_simplex ( Simplex_handle  s,
int  i 
)

returns the vertex associated with the i-th node of s.

Precondition
0 <= i <= dcur.

Friends And Related Function Documentation

template<typename R , typename Lifted_R >
template<typename R , typename Lifted_R >
void d2_map ( const 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_dR, Lifted_R >::NEAREST 
)
related

constructs a LEDA graph representation of the nearest (kind = NEAREST or the furthest (kind = FURTHEST) site Delaunay triangulation.

Precondition
dim() == 2.