CGAL::Delaunay_triangulation_3<DelaunayTriangulationTraits_3,TriangulationDataStructure_3>

Definition

The class Delaunay_triangulation_3 represents a three-dimensional Delaunay triangulation.

The user is advised to use the class Triangulation_hierarchy_3 rather than this basic Delaunay triangulation class: it offers the same functionalities but is much more efficient for large data sets.

#include <CGAL/Delaunay_triangulation_3.h>

Parameters

The first template argument must be a model of the DelaunayTriangulationTraits_3 concept.

The second template argument must be a model of the TriangulationDataStructure_3 concept. It has the default value Triangulation_data_structure_3<Triangulation_vertex_base_3<DelaunayTriangulationTraits_3>, Triangulation_cell_base_3<DelaunayTriangulationTraits_3> >.

Inherits From

Triangulation_3<DelaunayTriangulationTraits_3,TriangulationDataStructure_3>

Types

In addition to those inherited, the following types are defined, for use by the construction of the Voronoi diagram:

typedef DelaunayTriangulationTraits_3::Line_3
Line;
typedef DelaunayTriangulationTraits_3::Ray_3
Ray;
typedef DelaunayTriangulationTraits_3::Plane_3
Plane;
typedef DelaunayTriangulationTraits_3::Object_3
Object;

Creation

Delaunay_triangulation_3<DelaunayTriangulationTraits_3,TriangulationDataStructure_3> dt ( DelaunayTriangulationTraits_3 traits = DelaunayTriangulationTraits_3());
Creates an empty Delaunay triangulation, possibly specifying a traits class traits.


Delaunay_triangulation_3<DelaunayTriangulationTraits_3,TriangulationDataStructure_3> dt ( Delaunay_triangulation_3 dt1);
Copy constructor.


template < class InputIterator >
Delaunay_triangulation_3<DelaunayTriangulationTraits_3,TriangulationDataStructure_3> dt ( InputIterator first,
InputIterator last,
DelaunayTriangulationTraits_3 traits = DelaunayTriangulationTraits_3());
Creates a Delaunay triangulation of the points specified by the iterator range [first,last) of value type Point, possibly specifying a traits class traits.

Operations

Insertion

The following methods overload the corresponding methods of triangulations to ensure the empty sphere property of Delaunay triangulations.

In the degenerate case when there are cospherical points, the Delaunay triangulation is known not to be uniquely defined. In this case, CGAL chooses a particular Delaunay triangulation using a symbolic perturbation scheme [DT03].

Vertex_handle
dt.insert ( Point p,
Cell_handle start = Cell_handle())
Inserts point p in the triangulation and returns the corresponding vertex. Similar to the insertion in a triangulation, but ensures in addition the empty sphere property of all the created faces. The optional argument start is used as a starting place for the search.

Vertex_handle
dt.insert ( Point p,
Locate_type lt,
Cell_handle loc,
int li,
int lj)
Inserts point p in the triangulation and returns the corresponding vertex. Similar to the above insert() function, but takes as additional parameter the return values of a previous location query. See description of Triangulation_3::locate().

The following method allows one to insert several points. It returns the number of inserted points.

template < class InputIterator >
int dt.insert ( InputIterator first, InputIterator last)
Inserts the points in the iterator range [.first, last.), of value type Point.

Point moving

Vertex_handle dt.move_point ( Vertex_handle v, Point p)
Moves the point stored in v to p, while preserving the Delaunay property. This performs an action semantically equivalent to remove(v) followed by insert(p), but is supposedly faster when the point has not moved much. Returns the handle to the new vertex.
Precondition: v is a finite vertex of the triangulation.

Removal

When a vertex v is removed from a triangulation, all the cells incident to v must be removed, and the polyhedral region consisting of all the tetrahedra that are incident to v must be retriangulated. So, the problem reduces to triangulating a polyhedral region, while preserving its boundary, or to compute a constrained triangulation. This is known to be sometimes impossible: the Schönhardt polyhedron cannot be triangulated [She98].

However, when dealing with Delaunay triangulations, the case of such polyhedra that cannot be retriangulated cannot happen, so CGAL proposes a vertex removal.

void dt.remove ( Vertex_handle v)
Removes the vertex v from the triangulation.
Note that in CGAL 3.0 we have implemented a new algorithm for retriangulating the hole after the removal. In case that you experience problems when removing a vertex, the old code can be enabled with a #define CGAL_DELAUNAY_3_OLD_REMOVE 1.
Precondition: v is a finite vertex of the triangulation.

template < typename InputIterator >
int dt.remove ( InputIterator first, InputIterator beyond)
Removes the vertices specified by the iterator range [first, beyond) of value type Vertex_handle. remove() is called over each element of the range. The number of vertices removed is returned.
Precondition: All vertices of the range are finite vertices of the triangulation.

Queries

Bounded_side dt.side_of_sphere ( Cell_handle c, Point p)
Returns a value indicating on which side of the circumscribed sphere of c the point p lies. More precisely, it returns:
- ON_BOUNDED_SIDE if p is inside the sphere. For an infinite cell this means that p lies strictly either in the half space limited by its finite facet and not containing any other point of the triangulation, or in the interior of the disk circumscribing the finite facet.
- ON_BOUNDARY if p on the boundary of the sphere. For an infinite cell this means that p lies on the circle circumscribing the finite facet.
- ON_UNBOUNDED_SIDE if p lies outside the sphere. For an infinite cell this means that p does not satisfy either of the two previous conditions.
Precondition: dt.dimension() =3.

Bounded_side dt.side_of_circle ( Facet f, Point p)
Returns a value indicating on which side of the circumscribed circle of f the point p lies. More precisely, it returns:
- in dimension 3:
- For a finite facet, ON_BOUNDARY if p lies on the circle, ON_UNBOUNDED_SIDE when it lies in the exterior of the disk, ON_BOUNDED_SIDE when it lies in its interior.
- For an infinite facet, it considers the plane defined by the finite facet of the same cell, and does the same as in dimension 2 in this plane.
- in dimension 2:
- For a finite facet, ON_BOUNDARY if p lies on the circle, ON_UNBOUNDED_SIDE when it lies in the exterior of the disk, ON_BOUNDED_SIDE when it lies in its interior.
- For an infinite facet, ON_BOUNDARY if the point lies on the finite edge of f (endpoints included), ON_BOUNDED_SIDE for a point in the open half plane defined by f and not containing any other point of the triangulation, ON_UNBOUNDED_SIDE elsewhere.
Precondition: dt.dimension() 2 and in dimension 3, p is coplanar with f.

Bounded_side dt.side_of_circle ( Cell_handle c, int i, Point p)
Same as the previous method for facet i of cell c.

Vertex_handle
dt.nearest_vertex ( Point p,
Cell_handle c = Cell_handle())
Returns any nearest vertex to the point p, or the default constructed handle if the triangulation is empty. The optional argument c is a hint specifying where to start the search.
Precondition: c is a cell of dt.

Vertex_handle dt.nearest_vertex_in_cell ( Point p, Cell_handle c)
Returns the vertex of the cell c that is nearest to p.

A point p is said to be in conflict with a cell c in dimension 3 (resp. a facet f in dimension 2) iff dt.side_of_sphere(c, p) (resp. dt.side_of_circle(f, p)) returns ON_BOUNDED_SIDE. The set of cells (resp. facets in dimension 2) which are in conflict with p is connected, and it forms a hole.

template <class OutputIteratorBoundaryFacets, class OutputIteratorCells>
std::pair<OutputIteratorBoundaryFacets, OutputIteratorCells>
dt.find_conflicts ( Point p,
Cell_handle c,
OutputIteratorBoundaryFacets bfit,
OutputIteratorCells cit)
Computes the conflict hole induced by p. The starting cell (resp. facet) c must be in conflict. Then this function returns respectively in the output iterators:
- cit: the cells (resp. facets) in conflict.
- bfit: the facets (resp. edges) on the boundary, that is, the facets (resp. edges) (t, i) where the cell (resp. facet) t is in conflict, but t->neighbor(i) is not.
This function can be used in conjunction with insert_in_hole() in order to decide the insertion of a point after seeing which elements of the triangulation are affected. Returns the pair composed of the resulting output iterators.
Precondition: dt.dimension() 2, and c is in conflict with p.

template <class OutputIteratorBoundaryFacets, class OutputIteratorCells, class OutputIteratorInternalFacets>
Triple<OutputIteratorBoundaryFacets, OutputIteratorCells, OutputIteratorInternalFacets>
dt.find_conflicts ( Point p,
Cell_handle c,
OutputIteratorBoundaryFacets bfit,
OutputIteratorCells cit,
OutputIteratorInternalFacets ifit)
Same as the other find_conflicts() function, except that it also computes the internal facets, i.e. the facets common to two cells which are in conflict with p. Then this function returns respectively in the output iterators:
- cit: the cells (resp. facets) in conflict.
- bfit: the facets (resp. edges) on the boundary, that is, the facets (resp. edges) (t, i) where the cell (resp. facet) t is in conflict, but t->neighbor(i) is not.
- ifit: the facets (resp. edges) inside the hole, that is, delimiting two cells (resp facets) in conflict.
Returns the Triple composed of the resulting output iterators.
Precondition: dt.dimension() 2, and c is in conflict with p.

template <class OutputIterator>
OutputIterator
dt.vertices_in_conflict ( Point p,
Cell_handle c,
OutputIterator res)
Similar to find_conflicts(), but reports the vertices which are on the boundary of the conflict hole of p, in the output iterator res. Returns the resulting output iterator.
Precondition: dt.dimension() 2, and c is in conflict with p.

Voronoi diagram

CGAL offers several functionalities to display the Voronoi diagram of a set of points in 3D.

Note that the user should use a kernel with exact constructions in order to guarantee the computation of the Voronoi diagram (as opposed to computing the triangulation only, which requires only exact predicates).

Point dt.dual ( Cell_handle c)
Returns the circumcenter of the four vertices of c.
Precondition: dt.dimension()=3 and c is not infinite.

Object dt.dual ( Facet f) Returns the dual of facet f, which is
in dimension 3: either a segment, if the two cells incident to f are finite, or a ray, if one of them is infinite;
in dimension 2: a point.
Precondition: dt.dimension() 2 and f is not infinite.

Object dt.dual ( Cell_handle c, int i)
same as the previous method for facet (c,i).

template <class Stream>
Stream & dt.draw_dual ( Stream & os)
Sends the set of duals to all the facets of dt into os.


begin of advanced section  advanced  begin of advanced section

Checking

bool dt.is_valid ( bool verbose = false)
Checks the combinatorial validity of the triangulation and the validity of its geometric embedding (see Section reference). Also checks that all the circumscribing spheres (resp. circles in dimension 2) of cells (resp. facets in dimension 2) are empty.
When verbose is set to true, messages describing the first invalidity encountered are printed.

bool dt.is_valid ( Cell_handle c, bool verbose = false)
Checks the combinatorial and geometric validity of the cell (see Section reference). Also checks that the circumscribing sphere (resp. circle in dimension 2) of cells (resp. facet in dimension 2) is empty.
When verbose is set to true, messages are printed to give a precise indication of the kind of invalidity encountered.

These methods are mainly a debugging help for the users of advanced features.
end of advanced section  advanced  end of advanced section

See Also

CGAL::Triangulation_hierarchy_3.