\( \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.14 - 2D Triangulation
CGAL::Delaunay_triangulation_2< Traits, Tds > Class Template Reference

#include <CGAL/Delaunay_triangulation_2.h>

Inherits from

CGAL::Triangulation_2< Traits, Tds >.

Definition

The class Delaunay_triangulation_2 is designed to represent the Delaunay triangulation of a set of points in a plane.

A Delaunay triangulation of a set of points is a triangulation of the sets of points that fulfills the following empty circle property (also called Delaunay property): the circumscribing circle of any facet of the triangulation contains no point of the set in its interior. For a point set with no case of co-circularity of more than three points, the Delaunay triangulation is unique, it is the dual of the Voronoi diagram of the points.

Template Parameters
Tdsmust be a model of TriangulationDataStructure_2. CGAL provides a default instantiation for this parameter, which is the class CGAL::Triangulation_data_structure_2 < CGAL::Triangulation_vertex_base_2<Traits>, CGAL::Triangulation_face_base_2<Traits> >.
Traitsmust be a model of DelaunayTriangulationTraits_2. The concept DelaunayTriangulationTraits_2 refines the concept TriangulationTraits_2, providing a predicate type to check the empty circle property.

Changing this predicate type allows the user to build Delaunay triangulations for different metrics such that \( L_1\) or \( L_{\infty}\) or any metric defined by a convex object. However, the user of an exotic metric must be careful that the constructed triangulation has to be a triangulation of the convex hull which means that convex hull edges have to be Delaunay edges. This is granted for any smooth convex metric (like \( L_2\)) and can be ensured for other metrics (like \( L_{\infty}\)) by the addition to the point set of well chosen sentinel points. The concept of DelaunayTriangulationTraits_2 is described DelaunayTriangulationTraits_2.

When dealing with a large triangulations, the user is advised to encapsulate the Delaunay triangulation class into a triangulation hierarchy, which means to use the class Triangulation_hierarchy_2<Tr> with the template parameter instantiated with Delaunay_triangulation_2. The triangulation hierarchy will then offer the same functionalities but be much more for efficient for locations and insertions.

Types

All the types defined in Triangulation_2<Traits,Tds> are inherited.

Implementation

Insertion is implemented by inserting in the triangulation, then performing a sequence of Delaunay flips. The number of flips is \( O(d)\) if the new vertex is of degree \( d\) in the new triangulation. For points distributed uniformly at random, insertion takes time \( O(1)\) on average.

Removal calls the removal in the triangulation and then re-triangulates the hole in such a way that the Delaunay criterion is satisfied. Removal of a vertex of degree \( d\) takes time \( O(d^2)\). The degree \( d\) is \( O(1)\) for a random vertex in the triangulation.

After a point location step, the nearest neighbor is found in time \( O(n)\) in the worst case, but in time \( O(1)\) for vertices distributed uniformly at random and any query point.

See also
CGAL::Triangulation_2<Traits,Tds>
TriangulationDataStructure_2
DelaunayTriangulationTraits_2
Triangulation_hierarchy_2<Tr>
Examples:
Triangulation_2/copy_triangulation_2.cpp, Triangulation_2/hierarchy.cpp, Triangulation_2/info_insert_with_pair_iterator_2.cpp, Triangulation_2/info_insert_with_transform_iterator_2.cpp, Triangulation_2/info_insert_with_zip_iterator_2.cpp, Triangulation_2/print_cropped_voronoi.cpp, Triangulation_2/terrain.cpp, Triangulation_2/terrain_with_info.cpp, and Triangulation_2/voronoi.cpp.

Creation

 Delaunay_triangulation_2 (const Traits &gt=Traits())
 default constructor.
 
 Delaunay_triangulation_2 (const Delaunay_triangulation_2< Traits, Tds > &tr)
 copy constructor. More...
 
template<class InputIterator >
Delaunay_triangulation_2< Traits, Tds > dt (InputIterator first, InputIterator last, Traits gt=Traits())
 Equivalent to constructing an empty triangulation with the optional traits class argument and calling insert(first,last).
 

Insertion and Removal

The following insertion and removal functions overwrite the functions inherited from the class Triangulation_2<Traits,Tds> to maintain the Delaunay property.

In the degenerate case when there are co-circular 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 [2]. Note that the other modifier functions of Triangulation_2<Traits,Tds> are not overwritten. Thus a call to Triangulation_2::insert_in_face() insert_in_face() insert_in_edge(), insert_outside_convex_hull(), insert_outside_affine_hull() or flip() on a valid Delaunay triangulation might lead to a triangulation which is no longer a Delaunay triangulation.

Vertex_handle insert (const Point &p, Face_handle f=Face_handle())
 inserts point p. More...
 
Vertex_handle insert (const Point &p, Locate_type &lt, Face_handle loc, int li)
 inserts a point p at the location given by (lt,loc,li). More...
 
Vertex_handle push_back (const Point &p)
 equivalent to insert(p).
 
template<class PointInputIterator >
std::ptrdiff_t insert (PointInputIterator first, PointInputIterator last)
 inserts the points in the range [first,last). More...
 
template<class PointWithInfoInputIterator >
std::ptrdiff_t insert (PointWithInfoInputIterator first, PointWithInfoInputIterator last)
 inserts the points in the iterator range [first,last). More...
 
void remove (Vertex_handle v)
 removes the vertex from the triangulation.
 

Displacement

Vertex_handle move_if_no_collision (Vertex_handle v, const Point &p)
 if there is not already another vertex placed on p, the triangulation is modified such that the new position of vertex v is p, and v is returned. More...
 
Vertex_handle move (Vertex_handle v, const Point &p)
 same as move_if_no_collision(), if there is no collision. More...
 

Queries

Vertex_handle nearest_vertex (const Point &p, Face_handle f=Face_handle()) const
 returns any nearest vertex of p. More...
 
template<class OutputItFaces , class OutputItBoundaryEdges >
std::pair< OutputItFaces, OutputItBoundaryEdges > get_conflicts_and_boundary (const Point &p, OutputItFaces fit, OutputItBoundaryEdges eit, Face_handle start) const
 outputs the faces and boundary edges of the conflict zone of point p into output iterators. More...
 
template<class OutputItFaces >
OutputItFaces get_conflicts (const Point &p, OutputItFaces fit, Face_handle start) const
 outputs the faces of the conflict zone of point p into an output iterator. More...
 
template<class OutputItBoundaryEdges >
OutputItBoundaryEdges get_boundary_of_conflicts (const Point &p, OutputItBoundaryEdges eit, Face_handle start) const
 outputs the boundary edges of the conflict zone of point p into an output iterator. More...
 

Voronoi Diagram

The following member functions provide the elements of the dual Voronoi diagram.

Point dual (const Face_handle &f) const
 Returns the center of the circle circumscribed to face f. More...
 
Object dual (const Edge &e) const
 returns a segment, a ray or a line supported by the bisector of the endpoints of e. More...
 
Object dual (const Edge_circulator &ec) const
 Idem.
 
Object dual (const Edge_iterator &ei) const
 Idem.
 
template<class Stream >
Stream & draw_dual (Stream &ps)
 output the dual Voronoi diagram to stream ps.
 

Predicates

Oriented_side side_of_oriented_circle (Face_handle f, const Point &p) const
 Returns the side of p with respect to the circle circumscribing the triangle associated with f.
 

Miscellaneous

The checking function is_valid() is also overwritten to additionally test the empty circle property.

bool is_valid (bool verbose=false, int level=0) const
 tests the validity of the triangulation as a Triangulation_2 and additionally tests the Delaunay property. More...
 

Additional Inherited Members

- Public Types inherited from CGAL::Triangulation_2< Traits, Tds >
enum  Locate_type {
  VERTEX =0, EDGE, FACE, OUTSIDE_CONVEX_HULL,
  OUTSIDE_AFFINE_HULL
}
 specifies which case occurs when locating a point in the triangulation. More...
 
typedef Tds::Vertex_handle Vertex_handle
 handle to a vertex.
 
typedef Tds::Face_handle Face_handle
 handle to a face.
 
typedef Tds::Face_iterator All_faces_iterator
 iterator over all faces.
 
typedef Tds::Edge_iterator All_edges_iterator
 iterator over all edges.
 
typedef Tds::Vertex_iterator All_vertices_iterator
 iterator over all vertices.
 
typedef unspecified_type Finite_faces_iterator
 iterator over finite faces.
 
typedef unspecified_type Finite_edges_iterator
 iterator over finite edges.
 
typedef unspecified_type Finite_vertices_iterator
 iterator over finite vertices.
 
typedef unspecified_type Point_iterator
 iterator over the points corresponding the finite vertices of the triangulation.
 
typedef unspecified_type Line_face_circulator
 circulator over all faces intersected by a line.
 
typedef unspecified_type Face_circulator
 circulator over all faces incident to a given vertex.
 
typedef unspecified_type Edge_circulator
 circulator over all edges incident to a given vertex.
 
typedef unspecified_type Vertex_circulator
 circulator over all vertices incident to a given vertex.
 
typedef Traits Geom_traits
 the traits class.
 
typedef Tds Triangulation_data_structure
 the triangulation data structure type.
 
typedef Traits::Point_2 Point
 the point type.
 
typedef Traits::Segment_2 Segment
 the segment type.
 
typedef Traits::Triangle_2 Triangle
 the triangle type.
 
typedef Tds::Vertex Vertex
 the vertex type.
 
typedef Tds::Face Face
 the face type.
 
typedef Tds::Edge Edge
 the edge type.
 
typedef Tds::size_type size_type
 Size type (an unsigned integral type).
 
typedef Tds::difference_type difference_type
 Difference type (a signed integral type).
 
- Public Member Functions inherited from CGAL::Triangulation_2< Traits, Tds >
 Triangulation_2 (const Traits &gt=Traits())
 Introduces an empty triangulation.
 
 Triangulation_2 (const Triangulation_2 &tr)
 Copy constructor. More...
 
Triangulation_2 operator= (const Triangulation_2< Traits, Tds > &tr)
 Assignment. More...
 
void swap (Triangulation_2 &tr)
 The triangulations tr and *this are swapped. More...
 
void clear ()
 Deletes all faces and finite vertices resulting in an empty triangulation.
 
int dimension () const
 Returns the dimension of the convex hull.
 
size_type number_of_vertices () const
 Returns the number of finite vertices.
 
size_type number_of_faces () const
 Returns the number of finite faces.
 
Face_handle infinite_face () const
 a face incident to the infinite vertex.
 
Vertex_handle infinite_vertex () const
 the infinite vertex.
 
Vertex_handle finite_vertex () const
 a vertex distinct from the infinite vertex.
 
const Geom_traitsgeom_traits () const
 Returns a const reference to the triangulation traits object.
 
const TriangulationDataStructure_2tds () const
 Returns a const reference to the triangulation data structure.
 
TriangulationDataStructure_2tds ()
 Returns a reference to the triangulation data structure.
 
bool is_infinite (Vertex_handle v) const
 true iff v is the infinite vertex.
 
bool is_infinite (Face_handle f) const
 true iff face f is infinite.
 
bool is_infinite (Face_handle f, int i) const
 true iff edge (f,i) is infinite.
 
bool is_infinite (Edge e) const
 true iff edge e is infinite.
 
bool is_infinite (Edge_circulator ec) const
 true iff edge *ec is infinite.
 
bool is_infinite (All_edges_iterator ei) const
 true iff edge *ei is infinite.
 
bool is_edge (Vertex_handle va, Vertex_handle vb)
 true if there is an edge having va and vb as vertices.
 
bool is_edge (Vertex_handle va, Vertex_handle vb, Face_handle &fr, int &i)
 as above. More...
 
bool includes_edge (Vertex_handle va, Vertex_handle vb, Vertex_handle &vbr, Face_handle &fr, int &i)
 true if the line segment from va to vb includes an edge e incident to va. More...
 
bool is_face (Vertex_handle v1, Vertex_handle v2, Vertex_handle v3)
 true if there is a face having v1, v2 and v3 as vertices.
 
bool is_face (Vertex_handle v1, Vertex_handle v2, Vertex_handle v3, Face_handle &fr)
 as above. More...
 
Face_handle locate (const Point &query, Face_handle f=Face_handle()) const
 If the point query lies inside the convex hull of the points, a face that contains the query in its interior or on its boundary is returned. More...
 
Face_handle inexact_locate (const Point &query, Face_handle start=Face_handle()) const
 Same as locate() but uses inexact predicates. More...
 
Face_handle locate (const Point &query, Locate_type &lt, int &li, Face_handle h=Face_handle()) const
 Same as above. More...
 
Oriented_side oriented_side (Face_handle f, const Point &p) const
 Returns on which side of the oriented boundary of f lies the point p. More...
 
Oriented_side side_of_oriented_circle (Face_handle f, const Point &p)
 Returns on which side of the circumcircle of face f lies the point p. More...
 
void flip (Face_handle f, int i)
 Exchanges the edge incident to f and f->neighbor(i) with the other diagonal of the quadrilateral formed by f and f->neighbor(i). More...
 
Vertex_handle insert (const Point &p, Face_handle f=Face_handle())
 Inserts point p in the triangulation and returns the corresponding vertex. More...
 
Vertex_handle insert (const Point &p, Locate_type lt, Face_handle loc, int li)
 Same as above except that the location of the point p to be inserted is assumed to be given by (lt,loc,i) (see the description of the locate method above.)
 
Vertex_handle push_back (const Point &p)
 Equivalent to insert(p).
 
template<class InputIterator >
std::ptrdiff_t insert (InputIterator first, InputIterator last)
 Inserts the points in the range [first,last). More...
 
void remove (Vertex_handle v)
 Removes the vertex from the triangulation. More...
 
Vertex_handle move_if_no_collision (Vertex_handle v, const Point &p)
 If there is not already another vertex placed on p, the triangulation is modified such that the new position of vertex v is p, and v is returned. More...
 
Vertex_handle move (Vertex_handle v, const Point &p)
 If there is no collision during the move, this function is the same as move_if_no_collision . More...
 
Vertex_handle insert_first (const Point &p)
 Inserts the first finite vertex .
 
Vertex_handle insert_second (const Point &p)
 Inserts the second finite vertex .
 
Vertex_handle insert_in_face (const Point &p, Face_handle f)
 Inserts vertex v in face f. More...
 
Vertex_handle insert_in_edge (const Point &p, Face_handle f, int i)
 Inserts vertex v in edge i of f. More...
 
Vertex_handle insert_outside_convex_hull (const Point &p, Face_handle f)
 Inserts a point which is outside the convex hull but in the affine hull. More...
 
Vertex_handle insert_outside_affine_hull (const Point &p)
 Inserts a point which is outside the affine hull.
 
void remove_degree_3 (Vertex_handle v)
 Removes a vertex of degree three. More...
 
void remove_second (Vertex_handle v)
 Removes the before last finite vertex.
 
void remove_first (Vertex_handle v)
 Removes the last finite vertex.
 
template<class EdgeIt >
Vertex_handle star_hole (Point p, EdgeIt edge_begin, EdgeIt edge_end)
 creates a new vertex v and use it to star the hole whose boundary is described by the sequence of edges [edge_begin, edge_end). More...
 
template<class EdgeIt , class FaceIt >
Vertex_handle star_hole (Point p, EdgeIt edge_begin, EdgeIt edge_end, FaceIt face_begin, FaceIt face_end)
 same as above, except that the algorithm first recycles faces in the sequence [face_begin, face_end) and create new ones only when the sequence is exhausted. More...
 
Finite_vertices_iterator finite_vertices_begin () const
 Starts at an arbitrary finite vertex.
 
Finite_vertices_iterator finite_vertices_end () const
 Past-the-end iterator.
 
Finite_edges_iterator finite_edges_begin () const
 Starts at an arbitrary finite edge.
 
Finite_edges_iterator finite_edges_end () const
 Past-the-end iterator.
 
Finite_faces_iterator finite_faces_begin () const
 Starts at an arbitrary finite face.
 
Finite_faces_iterator finite_faces_end () const
 Past-the-end iterator.
 
Point_iterator points_begin () const
 
Point_iterator points_end () const
 Past-the-end iterator.
 
All_vertices_iterator all_vertices_begin () const
 Starts at an arbitrary vertex.
 
All_vertices_iterator all_vertices_end () const
 Past-the-end iterator.
 
All_edges_iterator all_edges_begin () const
 Starts at an arbitrary edge.
 
All_edges_iterator all_edges_end () const
 Past-the-end iterator.
 
All_faces_iterator all_faces_begin () const
 Starts at an arbitrary face.
 
All_faces_iterator all_faces_end () const
 Past-the-end iterator.
 
Line_face_circulator line_walk (const Point &p, const Point &q, Face_handle f=Face_handle()) const
 This function returns a circulator that allows to visit the faces intersected by the line pq. More...
 
Face_circulator incident_faces (Vertex_handle v) const
 Starts at an arbitrary face incident to v.
 
Face_circulator incident_faces (Vertex_handle v, Face_handle f) const
 Starts at face f. More...
 
Edge_circulator incident_edges (Vertex_handle v) const
 Starts at an arbitrary edge incident to v.
 
Edge_circulator incident_edges (Vertex_handle v, Face_handle f) const
 Starts at the first edge of f incident to v, in counterclockwise order around v. More...
 
Vertex_circulator incident_vertices (Vertex_handle v) const
 Starts at an arbitrary vertex incident to v.
 
Vertex_circulator incident_vertices (Vertex_handle v, Face_handle f)
 Starts at the first vertex of f adjacent to v in counterclockwise order around v. More...
 
Vertex_handle mirror_vertex (Face_handle f, int i) const
 returns the vertex of the \( i^{th}\) neighbor of f that is opposite to f. More...
 
int mirror_index (Face_handle f, int i) const
 returns the index of f in its \( i^{th}\) neighbor. More...
 
Edge mirror_edge (Edge e) const
 returns the same edge seen from the other adjacent face. More...
 
int ccw (int i) const
 Returns \( i+1\) modulo 3. More...
 
int cw (int i) const
 Returns \( i+2\) modulo 3. More...
 
Triangle triangle (Face_handle f) const
 Returns the triangle formed by the three vertices of f. More...
 
Segment segment (Face_handle f, int i) const
 Returns the line segment formed by the vertices ccw(i) and cw(i) of face f. More...
 
Segment segment (const Edge &e) const
 Returns the line segment corresponding to edge e. More...
 
Segment segment (const Edge_circulator &ec) const
 Returns the line segment corresponding to edge *ec. More...
 
Segment segment (const Edge_iterator &ei) const
 Returns the line segment corresponding to edge *ei. More...
 
Point circumcenter (Face_handle f) const
 Compute the circumcenter of the face pointed to by f. More...
 
void set_infinite_vertex (const Vertex_handle &v)
 This is an advanced function. More...
 
bool is_valid (bool verbose=false, int level=0) const
 Checks the combinatorial validity of the triangulation and also the validity of its geometric embedding. More...
 
- Public Member Functions inherited from CGAL::Triangulation_cw_ccw_2
 Triangulation_cw_ccw_2 ()
 default constructor.
 
int ccw (const int i) const
 returns the index of the neighbor or vertex that is next to the neighbor or vertex with index i in counterclockwise order around a face.
 
int cw (const int i) const
 returns the index of the neighbor or vertex that is next to the neighbor or vertex with index i in counterclockwise order around a face.
 

Constructor & Destructor Documentation

◆ Delaunay_triangulation_2()

template<typename Traits , typename Tds >
CGAL::Delaunay_triangulation_2< Traits, Tds >::Delaunay_triangulation_2 ( const Delaunay_triangulation_2< Traits, Tds > &  tr)

copy constructor.

All the vertices and faces are duplicated.

Member Function Documentation

◆ dual() [1/2]

template<typename Traits , typename Tds >
Point CGAL::Delaunay_triangulation_2< Traits, Tds >::dual ( const Face_handle f) const

Returns the center of the circle circumscribed to face f.

Precondition
f is not infinite.

◆ dual() [2/2]

template<typename Traits , typename Tds >
Object CGAL::Delaunay_triangulation_2< Traits, Tds >::dual ( const Edge e) const

returns a segment, a ray or a line supported by the bisector of the endpoints of e.

If faces incident to e are both finite, a segment whose endpoints are the duals of each incident face is returned. If only one incident face is finite, a ray whose endpoint is the dual of the finite incident face is returned. Otherwise both incident faces are infinite and the bisector line is returned.

◆ get_boundary_of_conflicts()

template<typename Traits , typename Tds >
template<class OutputItBoundaryEdges >
OutputItBoundaryEdges CGAL::Delaunay_triangulation_2< Traits, Tds >::get_boundary_of_conflicts ( const Point p,
OutputItBoundaryEdges  eit,
Face_handle  start 
) const

outputs the boundary edges of the conflict zone of point p into an output iterator.

This function outputs in the container pointed to by eit, the boundary of the zone in conflict with p. The boundary edges of the conflict zone are output in counterclockwise order and each edge is described through the incident face which is not in conflict with p. The function returns the resulting output iterator.

Template Parameters
OutputItBoundaryEdgesis an output iterator with Edge as value type.

◆ get_conflicts()

template<typename Traits , typename Tds >
template<class OutputItFaces >
OutputItFaces CGAL::Delaunay_triangulation_2< Traits, Tds >::get_conflicts ( const Point p,
OutputItFaces  fit,
Face_handle  start 
) const

outputs the faces of the conflict zone of point p into an output iterator.

same as get_conflicts_and_boundary() except that only the faces in conflict with p are output. The function returns the resulting output iterator.

Precondition
dimension()==2.

◆ get_conflicts_and_boundary()

template<typename Traits , typename Tds >
template<class OutputItFaces , class OutputItBoundaryEdges >
std::pair<OutputItFaces,OutputItBoundaryEdges> CGAL::Delaunay_triangulation_2< Traits, Tds >::get_conflicts_and_boundary ( const Point p,
OutputItFaces  fit,
OutputItBoundaryEdges  eit,
Face_handle  start 
) const

outputs the faces and boundary edges of the conflict zone of point p into output iterators.

This function outputs in the container pointed to by fit the faces which are in conflict with point p, i. e., the faces whose circumcircle contains p. It outputs in the container pointed to by eit the the boundary of the zone in conflict with p. The boundary edges of the conflict zone are output in counter-clockwise order and each edge is described through its incident face which is not in conflict with p. The function returns in a std::pair the resulting output iterators.

Template Parameters
OutItFacesis an output iterator with Face_handle as value type.
OutItBoundaryEdgesis an output iterator with Edge as value type.
Precondition
dimension()==2.

◆ insert() [1/4]

template<typename Traits , typename Tds >
Vertex_handle CGAL::Delaunay_triangulation_2< Traits, Tds >::insert ( const Point p,
Face_handle  f = Face_handle() 
)

inserts point p.

If point p coincides with an already existing vertex, this vertex is returned and the triangulation is not updated. Optional parameter f is used to initialize the location of p.

Examples:
Triangulation_2/terrain_with_info.cpp.

◆ insert() [2/4]

template<typename Traits , typename Tds >
Vertex_handle CGAL::Delaunay_triangulation_2< Traits, Tds >::insert ( const Point p,
Locate_type lt,
Face_handle  loc,
int  li 
)

inserts a point p at the location given by (lt,loc,li).

See also
Triangulation_2::locate()

◆ insert() [3/4]

template<typename Traits , typename Tds >
template<class PointInputIterator >
std::ptrdiff_t CGAL::Delaunay_triangulation_2< Traits, Tds >::insert ( PointInputIterator  first,
PointInputIterator  last 
)

inserts the points in the range [first,last).

Returns the number of inserted points. Note that this function is not guaranteed to insert the points following the order of PointInputIterator, as spatial_sort() is used to improve efficiency.

Template Parameters
PointInputIteratormust be an input iterator with the value type Point.

◆ insert() [4/4]

template<typename Traits , typename Tds >
template<class PointWithInfoInputIterator >
std::ptrdiff_t CGAL::Delaunay_triangulation_2< Traits, Tds >::insert ( PointWithInfoInputIterator  first,
PointWithInfoInputIterator  last 
)

inserts the points in the iterator range [first,last).

Returns the number of inserted points. Note that this function is not guaranteed to insert the points following the order of PointWithInfoInputIterator, as spatial_sort() is used to improve efficiency. Given a pair (p,i), the vertex v storing p also stores i, that is v.point() == p and v.info() == i. If several pairs have the same point, only one vertex is created, and one of the objects of type Vertex::Info will be stored in the vertex.

Precondition
Vertex must be model of the concept TriangulationVertexBaseWithInfo_2.
Template Parameters
PointWithInfoInputIteratormust be an input iterator with the value type std::pair<Point,Vertex::Info>.

◆ is_valid()

template<typename Traits , typename Tds >
bool CGAL::Delaunay_triangulation_2< Traits, Tds >::is_valid ( bool  verbose = false,
int  level = 0 
) const

tests the validity of the triangulation as a Triangulation_2 and additionally tests the Delaunay property.

This method is mainly useful for debugging Delaunay triangulation algorithms designed by the user.

◆ move()

template<typename Traits , typename Tds >
Vertex_handle CGAL::Delaunay_triangulation_2< Traits, Tds >::move ( Vertex_handle  v,
const Point p 
)

same as move_if_no_collision(), if there is no collision.

Otherwise, v is deleted and the vertex placed on p is returned.

Precondition
Vertex v must be finite.

◆ move_if_no_collision()

template<typename Traits , typename Tds >
Vertex_handle CGAL::Delaunay_triangulation_2< Traits, Tds >::move_if_no_collision ( Vertex_handle  v,
const Point p 
)

if there is not already another vertex placed on p, the triangulation is modified such that the new position of vertex v is p, and v is returned.

Otherwise, the triangulation is not modified and the vertex at point p is returned.

Precondition
Vertex v must be finite.

◆ nearest_vertex()

template<typename Traits , typename Tds >
Vertex_handle CGAL::Delaunay_triangulation_2< Traits, Tds >::nearest_vertex ( const Point p,
Face_handle  f = Face_handle() 
) const

returns any nearest vertex of p.

The implemented function begins with a location step and f may be used to initialize the location.