\( \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.2 - 2D Voronoi Diagram Adaptor
CGAL::Voronoi_diagram_2< DG, AT, AP > Class Template Reference

#include <CGAL/Voronoi_diagram_2.h>

Definition

The class Voronoi_diagram_2 provides an adaptor that enables us to view a triangulated Delaunay graph as their dual subdivision, the Voronoi diagram.

The class Voronoi_diagram_2 is designed to provide an API that is similar to that of CGAL's arrangements.

Template Parameters
DGcorresponds to the triangulated Delaunay graph and must be a model of the DelaunayGraph_2 concept.
ATmust be a model of the AdaptationTraits_2 concept.
APmust be a model of the AdaptationPolicy_2 concept. The third template parameter defaults to CGAL::Identity_policy_2<DG,AT>.
Refines:

DefaultConstructible,

CopyConstructible,

Assignable

See also
CGAL::Delaunay_triangulation_2<Traits,Tds>
CGAL::Regular_triangulation_2<Traits,Tds>
CGAL::Triangulation_hierarchy_2<Tr> provided that Tr is a model of DelaunayGraph_2
CGAL::Segment_Delaunay_graph_2<Gt,DS>
CGAL::Segment_Delaunay_graph_hierarchy_2<Gt,STag,DS>
CGAL::Apollonius_graph_2<Gt,Agds>
CGAL::Apollonius_graph_hierarchy_2<Gt,Agds>
Examples:
Voronoi_diagram_2/vd_2_point_location.cpp, and Voronoi_diagram_2/vd_2_point_location_sdg_linf.cpp.

Classes

class  Face
 The class Face is the class provided by the Voronoi_diagram_2<DG,AT,AP> class for Voronoi faces. More...
 
class  Halfedge
 The class Halfedge is the class provided by the Voronoi_diagram_2<DG,AT,AP> class for Voronoi halfedges. More...
 
class  Vertex
 The class Vertex is the Voronoi vertex class provided by the class Voronoi_diagram_2<DG,AT,AP> class. More...
 

Types

typedef DG Delaunay_graph
 A type for the dual Delaunay graph.
 
typedef AT Adaptation_traits
 A type for the adaptation traits needed by the Voronoi diagram adaptor.
 
typedef AP Adaptation_policy
 A type for the adaptation policy used.
 
typedef Adaptation_traits::Point_2 Point_2
 A type a point.
 
typedef Adaptation_traits::Site_2 Site_2
 A type for the sites of the Voronoi diagram.
 
typedef Delaunay_graph::size_type size_type
 A type for sizes.
 
typedef Delaunay_graph::Geom_traits Delaunay_geom_traits
 A type for the geometric traits of the Delaunay graph.
 
typedef Delaunay_graph::Vertex_handle Delaunay_vertex_handle
 A type for the vertex handles of the Delaunay graph.
 
typedef Delaunay_graph::Face_handle Delaunay_face_handle
 A type for the face handles of the Delaunay graph.
 
typedef Delaunay_graph::Edge Delaunay_edge
 A type for the edges of the Delaunay graph.
 

The vertices, edges and faces of the Voronoi diagram are accessed through handles, iterators and circulators.

The iterators and circulators are all bidirectional and non-mutable. The circulators and iterators are assignable to the corresponding handle types, and they are also convertible to the corresponding handles.

typedef unspecified_type Halfedge_handle
 Handle for halfedges.
 
typedef unspecified_type Vertex_handle
 Handle for vertices.
 
typedef unspecified_type Face_handle
 Handle for faces.
 
typedef unspecified_type Edge_iterator
 A type for an iterator over Voronoi edges. More...
 
typedef unspecified_type Halfedge_iterator
 A type for an iterator over Voronoi halfedges. More...
 
typedef unspecified_type Face_iterator
 A type for an iterator over Voronoi faces. More...
 
typedef unspecified_type Vertex_iterator
 A type for an iterator over Voronoi vertices. More...
 
typedef unspecified_type Halfedge_around_vertex_circulator
 A type for a circulator over the halfedges that have a common vertex as their target. More...
 
typedef unspecified_type Ccb_halfedge_circulator
 A type for a circulator over the halfedges on the boundary of a Voronoi face. More...
 
typedef unspecified_type Unbounded_faces_iterator
 A type for an iterator over the unbounded faces of the Voronoi diagram. More...
 
typedef unspecified_type Bounded_faces_iterator
 A type for an iterator over the bounded faces of the Voronoi diagram. More...
 
typedef unspecified_type Unbounded_halfedges_iterator
 A type for an iterator over the unbounded halfedges of the Voronoi diagram. More...
 
typedef unspecified_type Bounded_halfedges_iterator
 A type for an iterator over the bounded halfedges of the Voronoi diagram. More...
 
typedef unspecified_type Site_iterator
 A type for an iterator over the sites of the Voronoi diagram. More...
 
typedef boost::variant< Face_handle, Halfedge_handle, Vertex_handleLocate_result
 The result type of the point location queries.
 

Creation

 Voronoi_diagram_2 (Adaptation_traits at=Adaptation_traits(), Adaptation_policy ap=Adaptation_policy(), Delaunay_geom_traits gt=Delaunay_geom_traits())
 Creates a Voronoi diagram using at as adaptation traits and ap as adaptation policy; the underlying Delaunay graph is created using gt as geometric traits.
 
 Voronoi_diagram_2 (Delaunay_graph dg, bool swap_dg=false, Adaptation_traits at=Adaptation_traits(), Adaptation_policy ap=Adaptation_policy())
 Creates a Voronoi diagram from the Delaunay graph dg and using at as adaptation traits and ap as adaptation policy. More...
 
template<class Iterator >
 Voronoi_diagram_2 (Iterator first, Iterator beyond, Adaptation_traits at=Adaptation_traits(), Adaptation_policy ap=Adaptation_policy(), Delaunay_geom_traits gt=Delaunay_geom_traits())
 Creates a Voronoi diagram using as sites the sites in the iterator range [first, beyond), at as adaptation traits and ap as adaptation policy; the underlying Delaunay graph is created using gt as geometric traits. More...
 

Access Methods

Delaunay_graph dual ()
 Returns a const reference to the dual graph, i.e., the Delaunay graph.
 
Halfedge_handle dual (Delaunay_edge e)
 Returns a handle to the halfedge in the Voronoi diagram that is dual to the edge e in the Delaunay graph.
 
Face_handle dual (Delaunay_vertex_handle v)
 Returns a handle to the face in the Voronoi diagram that is dual to the vertex corresponding to the vertex handle v in the Delaunay graph.
 
Vertex_handle dual (Delaunay_face_handle f)
 Returns a handle to the vertex in the Voronoi diagram that is dual to the face corresponding to the face handle f in the Delaunay graph.
 
Adaptation_traits adaptation_traits ()
 Returns a reference to the Voronoi traits.
 
Adaptation_policy adaptation_policy ()
 Returns a reference to the adaptation policy.
 
size_type number_of_vertices ()
 Returns the number of Voronoi vertices.
 
size_type number_of_faces ()
 Returns the number of Voronoi faces (bounded and unbounded).
 
size_type number_of_halfedges ()
 Returns the number of halfedges (bounded and unbounded) in the Voronoi diagram. More...
 
size_type number_of_connected_components ()
 Returns the number of connected components of the Voronoi skeleton.
 
Face_handle unbounded_face ()
 Returns one of the unbounded faces of the Voronoi diagram. More...
 
Face_handle bounded_face ()
 Returns one of the bounded faces of the Voronoi diagram. More...
 
Halfedge_handle unbounded_halfedge ()
 Returns one of the unbounded halfedges of the Voronoi diagram. More...
 
Halfedge_handle bounded_halfedge ()
 Returns one of the bounded halfedges of the Voronoi diagram. More...
 

Iterators

The following iterators allow respectively to visit the faces (all or only the unbounded/bounded ones), edges, halfedges (all or only the unbounded/bounded ones) and vertices of the Voronoi diagram.

These iterators are non-mutable, bidirectional and their value types are respectively Face, Halfedge, Halfedge and Vertex. All iterators are convertible to the corresponding handles and are invalidated by any change in the Voronoi diagram. The following iterator provides access to the sites that define the Voronoi diagram. Its value type is Site_2. It is invalidated by any change in the Voronoi diagram.

Face_iterator faces_begin ()
 Starts at an arbitrary Voronoi face.
 
Face_iterator faces_end ()
 Past-the-end iterator.
 
Unbounded_faces_iterator unbounded_faces_begin ()
 Starts at an arbitrary unbounded Voronoi face.
 
Unbounded_faces_iterator unbounded_faces_end ()
 Past-the-end iterator.
 
Bounded_faces_iterator bounded_faces_begin ()
 Starts at an arbitrary bounded Voronoi face.
 
Bounded_faces_iterator bounded_faces_end ()
 Past-the-end iterator.
 
Edge_iterator edges_begin ()
 Starts at an arbitrary Voronoi edge.
 
Edge_iterator edges_end ()
 Past-the-end iterator.
 
Halfedge_iterator halfedges_begin ()
 Starts at an arbitrary Voronoi halfedge.
 
Halfedge_iterator halfedges_end ()
 Past-the-end iterator.
 
Unbounded_halfedges_iterator unbounded_halfedges_begin ()
 Starts at an arbitrary unbounded Voronoi edge.
 
Unbounded_halfedges_iterator unbounded_halfedges_end ()
 Past-the-end iterator.
 
Bounded_halfedges_iterator bounded_halfedges_begin ()
 Starts at an arbitrary bounded Voronoi edge.
 
Bounded_halfedges_iterator bounded_halfedges_end ()
 Past-the-end iterator.
 
Vertex_iterator vertices_begin ()
 Starts at an arbitrary Voronoi vertex.
 
Vertex_iterator vertices_end ()
 Past-the-end iterator.
 
Site_iterator sites_begin ()
 Starts at an arbitrary site.
 
Site_iterator sites_end ()
 Past-the-end iterator.
 

Circulators

The Voronoi diagram adaptor also provides circulators that allow to visit all halfedges whose target is a given vertex - this is the Halfedge_around_vertex_circulator, as well as all halfedges on the boundary of a Voronoi face - this is the Ccb_halfedge_circulator.

These circulators are non-mutable and bidirectional. The operator operator++ moves the former circulator counterclockwise around the vertex while the operator- moves clockwise. The latter circulator is moved by the operator operator++ to the next halfedge on the boundary in the counterclockwise sense, while operator- moves clockwise. When the Ccb_halfedge_circulator is defined over an infinite Voronoi face f, then applying operator++ to a circulator corresponding to a halfedge whose target is not finite moves to the next infinite (or semi-infinite) halfedge of f in the counterclockwise sense. Similarly, applying operator++ to a circulator corresponding to a halfedge whose source is not finite, moves to the previous infinite (or semi-infinite) halfedge of f in the clockwise sense. The Halfedge_around_vertex_circulator circulator is invalidated by any modification in the faces adjacent to the vertex over which it is defined. The Ccb_halfedge_circulator is invalidated by any modification in the face over which it is defined.

Ccb_halfedge_circulator ccb_halfedges (Face_handle f)
 Returns a circulator over the halfedges on the boundary of f. More...
 
Ccb_halfedge_circulator ccb_halfedges (Face_handle f, Halfedge_handle h)
 Returns a circulator over the halfedges on the boundary of f. More...
 
Halfedge_around_vertex_circulator incident_halfedges (Vertex_handle v)
 Returns a circulator over the halfedges whose target is the Voronoi vertex v. More...
 
Halfedge_around_vertex_circulator incident_halfedges (Vertex_handle v, Halfedge_handle h)
 Returns a circulator over the halfedges whose target is the Voronoi vertex v. More...
 

Insertion

Face_handle insert (Site_2 t)
 Inserts the site t in the Voronoi diagram. More...
 
template<class Iterator >
size_type insert (Iterator first, Iterator beyond)
 Inserts, in the Voronoi diagram, the sites in the iterator range [first, beyond). More...
 

Queries

Locate_result locate (Point_2 p)
 Performs point location for the query point p. More...
 

I/O

void file_output (std::ostream &os)
 Writes the current state of the Voronoi diagram to the output stream os. More...
 
void file_input (std::istream &is)
 Reads the current state of the Voronoi diagram from the input stream is. More...
 
std::ostream & operator<< (std::ostream &os, Voronoi_diagram_2< DG, AT, AP > vd)
 Writes the current state of the Voronoi diagram to the output stream os. More...
 
std::istream & operator>> (std::istream &is, Voronoi_diagram_2< DG, AT, AP > vd)
 Reads the current state of the Voronoi diagram from the input stream is. More...
 

Validity Check

bool is_valid ()
 Checks the validity of the dual Delaunay graph and the Voronoi diagram adaptor.
 

Miscellaneous

void clear ()
 Clears all contents of the Voronoi diagram.
 
void swap (Voronoi_diagram_2< DG, AT, AP > other)
 The Voronoi diagrams other and vd are swapped. More...
 

Member Typedef Documentation

◆ Bounded_faces_iterator

template<typename DG , typename AT , typename AP >
typedef unspecified_type CGAL::Voronoi_diagram_2< DG, AT, AP >::Bounded_faces_iterator

A type for an iterator over the bounded faces of the Voronoi diagram.

Its value type is Face.

◆ Bounded_halfedges_iterator

template<typename DG , typename AT , typename AP >
typedef unspecified_type CGAL::Voronoi_diagram_2< DG, AT, AP >::Bounded_halfedges_iterator

A type for an iterator over the bounded halfedges of the Voronoi diagram.

Its value type is Halfedge.

◆ Ccb_halfedge_circulator

template<typename DG , typename AT , typename AP >
typedef unspecified_type CGAL::Voronoi_diagram_2< DG, AT, AP >::Ccb_halfedge_circulator

A type for a circulator over the halfedges on the boundary of a Voronoi face.

Its value type of is Halfedge.

◆ Edge_iterator

template<typename DG , typename AT , typename AP >
typedef unspecified_type CGAL::Voronoi_diagram_2< DG, AT, AP >::Edge_iterator

A type for an iterator over Voronoi edges.

Edges are considered non-oriented. Its value type is Halfedge.

◆ Face_iterator

template<typename DG , typename AT , typename AP >
typedef unspecified_type CGAL::Voronoi_diagram_2< DG, AT, AP >::Face_iterator

A type for an iterator over Voronoi faces.

Its value type is Face.

◆ Halfedge_around_vertex_circulator

template<typename DG , typename AT , typename AP >
typedef unspecified_type CGAL::Voronoi_diagram_2< DG, AT, AP >::Halfedge_around_vertex_circulator

A type for a circulator over the halfedges that have a common vertex as their target.

Its value type is Halfedge.

◆ Halfedge_iterator

template<typename DG , typename AT , typename AP >
typedef unspecified_type CGAL::Voronoi_diagram_2< DG, AT, AP >::Halfedge_iterator

A type for an iterator over Voronoi halfedges.

Halfedges are oriented and come in pairs. Its value type is Halfedge.

◆ Site_iterator

template<typename DG , typename AT , typename AP >
typedef unspecified_type CGAL::Voronoi_diagram_2< DG, AT, AP >::Site_iterator

A type for an iterator over the sites of the Voronoi diagram.

Its value type is Site_2.

◆ Unbounded_faces_iterator

template<typename DG , typename AT , typename AP >
typedef unspecified_type CGAL::Voronoi_diagram_2< DG, AT, AP >::Unbounded_faces_iterator

A type for an iterator over the unbounded faces of the Voronoi diagram.

Its value type is Face.

◆ Unbounded_halfedges_iterator

template<typename DG , typename AT , typename AP >
typedef unspecified_type CGAL::Voronoi_diagram_2< DG, AT, AP >::Unbounded_halfedges_iterator

A type for an iterator over the unbounded halfedges of the Voronoi diagram.

Its value type is Halfedge.

◆ Vertex_iterator

template<typename DG , typename AT , typename AP >
typedef unspecified_type CGAL::Voronoi_diagram_2< DG, AT, AP >::Vertex_iterator

A type for an iterator over Voronoi vertices.

Its value type is Vertex.

Constructor & Destructor Documentation

◆ Voronoi_diagram_2() [1/2]

template<typename DG , typename AT , typename AP >
CGAL::Voronoi_diagram_2< DG, AT, AP >::Voronoi_diagram_2 ( Delaunay_graph  dg,
bool  swap_dg = false,
Adaptation_traits  at = Adaptation_traits(),
Adaptation_policy  ap = Adaptation_policy() 
)

Creates a Voronoi diagram from the Delaunay graph dg and using at as adaptation traits and ap as adaptation policy.

The Delaunay graph dg is fully copied if swap_dg is set to false, or swapped with the one stored internally if swap_dg is set to true.

◆ Voronoi_diagram_2() [2/2]

template<typename DG , typename AT , typename AP >
template<class Iterator >
CGAL::Voronoi_diagram_2< DG, AT, AP >::Voronoi_diagram_2 ( Iterator  first,
Iterator  beyond,
Adaptation_traits  at = Adaptation_traits(),
Adaptation_policy  ap = Adaptation_policy(),
Delaunay_geom_traits  gt = Delaunay_geom_traits() 
)

Creates a Voronoi diagram using as sites the sites in the iterator range [first, beyond), at as adaptation traits and ap as adaptation policy; the underlying Delaunay graph is created using gt as geometric traits.

Iterator must be a model of the InputIterator concept and its value type must be Site_2.

Member Function Documentation

◆ bounded_face()

template<typename DG , typename AT , typename AP >
Face_handle CGAL::Voronoi_diagram_2< DG, AT, AP >::bounded_face ( )

Returns one of the bounded faces of the Voronoi diagram.

If no bounded faces exist the default constructed face handle is returned.

◆ bounded_halfedge()

template<typename DG , typename AT , typename AP >
Halfedge_handle CGAL::Voronoi_diagram_2< DG, AT, AP >::bounded_halfedge ( )

Returns one of the bounded halfedges of the Voronoi diagram.

If no bounded halfedges exist the default constructed halfedge handle is returned.

◆ ccb_halfedges() [1/2]

template<typename DG , typename AT , typename AP >
Ccb_halfedge_circulator CGAL::Voronoi_diagram_2< DG, AT, AP >::ccb_halfedges ( Face_handle  f)

Returns a circulator over the halfedges on the boundary of f.

The circulator is initialized to an arbitrary halfedge on the boundary of the Voronoi face f.

◆ ccb_halfedges() [2/2]

template<typename DG , typename AT , typename AP >
Ccb_halfedge_circulator CGAL::Voronoi_diagram_2< DG, AT, AP >::ccb_halfedges ( Face_handle  f,
Halfedge_handle  h 
)

Returns a circulator over the halfedges on the boundary of f.

The circulator is initialized with the halfedge h.

Precondition
The halfedge h must lie on the boundary of f.

◆ file_input()

template<typename DG , typename AT , typename AP >
void CGAL::Voronoi_diagram_2< DG, AT, AP >::file_input ( std::istream &  is)

Reads the current state of the Voronoi diagram from the input stream is.

The following operator must be defined:

std::istream& operator>>(std::istream&, Delaunay_graph)

◆ file_output()

template<typename DG , typename AT , typename AP >
void CGAL::Voronoi_diagram_2< DG, AT, AP >::file_output ( std::ostream &  os)

Writes the current state of the Voronoi diagram to the output stream os.

The following operator must be defined:

std::ostream& operator<<(std::ostream&, Delaunay_graph)

◆ incident_halfedges() [1/2]

template<typename DG , typename AT , typename AP >
Halfedge_around_vertex_circulator CGAL::Voronoi_diagram_2< DG, AT, AP >::incident_halfedges ( Vertex_handle  v)

Returns a circulator over the halfedges whose target is the Voronoi vertex v.

The circulator is initialized to an arbitrary halfedge incident to v.

◆ incident_halfedges() [2/2]

template<typename DG , typename AT , typename AP >
Halfedge_around_vertex_circulator CGAL::Voronoi_diagram_2< DG, AT, AP >::incident_halfedges ( Vertex_handle  v,
Halfedge_handle  h 
)

Returns a circulator over the halfedges whose target is the Voronoi vertex v.

The circulator is initialized with the halfedge h.

Precondition
The vertex v must be the target vertex of the halfedge h.

◆ insert() [1/2]

template<typename DG , typename AT , typename AP >
Face_handle CGAL::Voronoi_diagram_2< DG, AT, AP >::insert ( Site_2  t)

Inserts the site t in the Voronoi diagram.

A handle to the face corresponding to the Voronoi face of t in the Voronoi diagram is returned. If t has an empty Voronoi cell, the default constructed face handle is returned. This method is supported only if Voronoi_traits::Has_inserter is set to CGAL::Tag_true.

◆ insert() [2/2]

template<typename DG , typename AT , typename AP >
template<class Iterator >
size_type CGAL::Voronoi_diagram_2< DG, AT, AP >::insert ( Iterator  first,
Iterator  beyond 
)

Inserts, in the Voronoi diagram, the sites in the iterator range [first, beyond).

The value type of Iterator must be Site_2. The number of sites in the iterator range is returned. This method is supported only if Voronoi_traits::Has_inserter is set to CGAL::Tag_true.

◆ locate()

template<typename DG , typename AT , typename AP >
Locate_result CGAL::Voronoi_diagram_2< DG, AT, AP >::locate ( Point_2  p)

Performs point location for the query point p.

In other words, the face, halfedge or vertex of the Voronoi diagram is found on which the point p lies. This method is supported only if Voronoi_traits::Has_nearest_site_2 is set to CGAL::Tag_true.

Precondition
The Voronoi diagram must contain at least one face.

◆ number_of_halfedges()

template<typename DG , typename AT , typename AP >
size_type CGAL::Voronoi_diagram_2< DG, AT, AP >::number_of_halfedges ( )

Returns the number of halfedges (bounded and unbounded) in the Voronoi diagram.

This is always an even number.

◆ operator
template<typename DG , typename AT , typename AP >
std::ostream& CGAL::Voronoi_diagram_2< DG, AT, AP >::operator<< ( std::ostream &  os,
Voronoi_diagram_2< DG, AT, AP >  vd 
)

Writes the current state of the Voronoi diagram to the output stream os.

The following operator must be defined:

std::ostream& operator<<(std::ostream&, Delaunay_graph)

◆ operator>>()

template<typename DG , typename AT , typename AP >
std::istream& CGAL::Voronoi_diagram_2< DG, AT, AP >::operator>> ( std::istream &  is,
Voronoi_diagram_2< DG, AT, AP >  vd 
)

Reads the current state of the Voronoi diagram from the input stream is.

The following operator must be defined: std::istream& operator>>(std::istream&, Delaunay_graph)

◆ swap()

template<typename DG , typename AT , typename AP >
void CGAL::Voronoi_diagram_2< DG, AT, AP >::swap ( Voronoi_diagram_2< DG, AT, AP >  other)

The Voronoi diagrams other and vd are swapped.

vd.swap(other) should be preferred to vd= other or to vd(other) if other is deleted afterwards.

◆ unbounded_face()

template<typename DG , typename AT , typename AP >
Face_handle CGAL::Voronoi_diagram_2< DG, AT, AP >::unbounded_face ( )

Returns one of the unbounded faces of the Voronoi diagram.

If no unbounded faces exist (this can happen if the number of sites is zero) the default constructed face handle is returned.

◆ unbounded_halfedge()

template<typename DG , typename AT , typename AP >
Halfedge_handle CGAL::Voronoi_diagram_2< DG, AT, AP >::unbounded_halfedge ( )

Returns one of the unbounded halfedges of the Voronoi diagram.

If no unbounded halfedges exist the default constructed halfedge handle is returned.