\( \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 - 2D Voronoi Diagram Adaptor
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Groups Pages
CGAL::Voronoi_diagram_2< DG, AT, AP >::Halfedge Class Reference

#include <CGAL/Voronoi_diagram_2.h>

Definition

Types

typedef unspecified_type Vertex
 A type for the vertices of the Voronoi diagram. More...
 
typedef unspecified_type Face
 A type for the faces of the Voronoi diagram. More...
 
typedef unspecified_type Vertex_handle
 Handle for the vertices of the Voronoi diagram. More...
 
typedef unspecified_type Face_handle
 Handle for the faces of the Voronoi diagram. More...
 
typedef unspecified_type Halfedge_handle
 Handle for the halfedges of the Voronoi diagram. More...
 
typedef unspecified_type Ccb_halfedge_circulator
 A type for a bidirectional circulator over the halfedges of the boundary of a Voronoi face. More...
 
typedef unspecified_type Delaunay_graph
 A type for the Delaunay graph. More...
 
typedef Delaunay_graph::Edge Delaunay_edge
 A type for the dual edge in the Delaunay graph. More...
 
typedef
Delaunay_graph::Vertex_handle 
Delaunay_vertex_handle
 A type for vertex handles in the Delaunay graph. More...
 

Access Methods

Halfedge_handle twin ()
 Returns the twin halfedge. More...
 
Halfedge_handle opposite ()
 Same as e.twin(). More...
 
Halfedge_handle next ()
 Returns the next halfedge in the counterclockwise sense around the boundary of the face that e is incident to. More...
 
Halfedge_handle previous ()
 Returns the previous halfedge in the counterclockwise sense around the boundary of the adjacent face. More...
 
Face_handle face ()
 Returns the face that e is incident to. More...
 
Vertex_handle source ()
 Returns the source vertex of e. More...
 
Vertex_handle target ()
 Returns the target vertex of e. More...
 
Ccb_halfedge_circulator ccb ()
 Returns a bidirectional circulator to traverse the halfedges on the boundary of the Voronoi face containing e. More...
 
Delaunay_edge dual ()
 Returns the corresponding dual edge in the Delaunay graph. More...
 

In the four methods below we consider Voronoi halfedges to be "parallel" to the \( x\)-axis, oriented from left to right.

Delaunay_vertex_handle up ()
 Returns a handle to the vertex in the Delaunay graph corresponding to the defining site above the Voronoi edge. More...
 
Delaunay_vertex_handle down ()
 Returns a handle to the vertex in the Delaunay graph corresponding to the defining site below the Voronoi edge. More...
 
Delaunay_vertex_handle left ()
 Returns a handle to the vertex in the Delaunay graph corresponding to the defining site to the left of the Voronoi edge. More...
 
Delaunay_vertex_handle right ()
 Returns a handle to the vertex in the Delaunay graph corresponding to the defining site to the right of the Voronoi edge. More...
 

Predicate Methods

bool has_source ()
 Returns true iff the halfedge corresponds to a bisecting segment or a bisecting ray oriented appropriately so that its apex is its source. More...
 
bool has_target ()
 Returns true iff the halfedge corresponds to a bisecting segment or a bisecting ray oriented appropriately so that its apex is its target. More...
 
bool is_unbounded ()
 Returns true iff the source or the target of the halfedge does not exist, i.e., if either of has_source() or has_target() return false. More...
 
bool is_bisector ()
 Returns true iff the Voronoi edge is an entire bisector. More...
 
bool is_segment ()
 Returns true iff the Voronoi edge has both a source and a target Voronoi vertex. More...
 
bool is_ray ()
 Returns true iff the Voronoi edge has either a source or a target Voronoi vertex, but not both; in other words it is a bisecting ray. More...
 
bool is_valid ()
 Returns true if the following conditions are met: the halfedge is not a rejected edge with respect to the chosen adaptation policy; the twin edge of its twin edge is itself; its adjacent face is not a rejected face with respect to the chosen adaptation policy; its source and target vertices are valid (provided they exist, of course); the previous of its next halfedge is itself and the next of its previous halfedge is itself. More...
 

Member Typedef Documentation

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

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

The value type of the circulator is CGAL::Voronoi_diagram_2<DG,AT,AP>::Halfedge and is convertible to Halfedge_handle.

template<typename DG , typename AT , typename AP >
typedef Delaunay_graph::Edge CGAL::Voronoi_diagram_2< DG, AT, AP >::Halfedge::Delaunay_edge

A type for the dual edge in the Delaunay graph.

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

A type for the Delaunay graph.

It is a model of the DelaunayGraph_2 concept.

template<typename DG , typename AT , typename AP >
typedef Delaunay_graph::Vertex_handle CGAL::Voronoi_diagram_2< DG, AT, AP >::Halfedge::Delaunay_vertex_handle

A type for vertex handles in the Delaunay graph.

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

A type for the faces of the Voronoi diagram.

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

Handle for the faces of the Voronoi diagram.

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

Handle for the halfedges of the Voronoi diagram.

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

A type for the vertices of the Voronoi diagram.

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

Handle for the vertices of the Voronoi diagram.

Member Function Documentation

template<typename DG , typename AT , typename AP >
Ccb_halfedge_circulator CGAL::Voronoi_diagram_2< DG, AT, AP >::Halfedge::ccb ( )

Returns a bidirectional circulator to traverse the halfedges on the boundary of the Voronoi face containing e.

The circulator is initialized to e. Applying operator++ (resp. operator-) to this circulator returns the next halfedge on the boundary of the face containing e in the counterclockwise (resp. clockwise) sense.

template<typename DG , typename AT , typename AP >
Delaunay_vertex_handle CGAL::Voronoi_diagram_2< DG, AT, AP >::Halfedge::down ( )

Returns a handle to the vertex in the Delaunay graph corresponding to the defining site below the Voronoi edge.

template<typename DG , typename AT , typename AP >
Delaunay_edge CGAL::Voronoi_diagram_2< DG, AT, AP >::Halfedge::dual ( )

Returns the corresponding dual edge in the Delaunay graph.

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

Returns the face that e is incident to.

template<typename DG , typename AT , typename AP >
bool CGAL::Voronoi_diagram_2< DG, AT, AP >::Halfedge::has_source ( )

Returns true iff the halfedge corresponds to a bisecting segment or a bisecting ray oriented appropriately so that its apex is its source.

template<typename DG , typename AT , typename AP >
bool CGAL::Voronoi_diagram_2< DG, AT, AP >::Halfedge::has_target ( )

Returns true iff the halfedge corresponds to a bisecting segment or a bisecting ray oriented appropriately so that its apex is its target.

template<typename DG , typename AT , typename AP >
bool CGAL::Voronoi_diagram_2< DG, AT, AP >::Halfedge::is_bisector ( )

Returns true iff the Voronoi edge is an entire bisector.

template<typename DG , typename AT , typename AP >
bool CGAL::Voronoi_diagram_2< DG, AT, AP >::Halfedge::is_ray ( )

Returns true iff the Voronoi edge has either a source or a target Voronoi vertex, but not both; in other words it is a bisecting ray.

template<typename DG , typename AT , typename AP >
bool CGAL::Voronoi_diagram_2< DG, AT, AP >::Halfedge::is_segment ( )

Returns true iff the Voronoi edge has both a source and a target Voronoi vertex.

template<typename DG , typename AT , typename AP >
bool CGAL::Voronoi_diagram_2< DG, AT, AP >::Halfedge::is_unbounded ( )

Returns true iff the source or the target of the halfedge does not exist, i.e., if either of has_source() or has_target() return false.

template<typename DG , typename AT , typename AP >
bool CGAL::Voronoi_diagram_2< DG, AT, AP >::Halfedge::is_valid ( )

Returns true if the following conditions are met: the halfedge is not a rejected edge with respect to the chosen adaptation policy; the twin edge of its twin edge is itself; its adjacent face is not a rejected face with respect to the chosen adaptation policy; its source and target vertices are valid (provided they exist, of course); the previous of its next halfedge is itself and the next of its previous halfedge is itself.

template<typename DG , typename AT , typename AP >
Delaunay_vertex_handle CGAL::Voronoi_diagram_2< DG, AT, AP >::Halfedge::left ( )

Returns a handle to the vertex in the Delaunay graph corresponding to the defining site to the left of the Voronoi edge.

Precondition
has_source() must be true.
template<typename DG , typename AT , typename AP >
Halfedge_handle CGAL::Voronoi_diagram_2< DG, AT, AP >::Halfedge::next ( )

Returns the next halfedge in the counterclockwise sense around the boundary of the face that e is incident to.

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

Same as e.twin().

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

Returns the previous halfedge in the counterclockwise sense around the boundary of the adjacent face.

template<typename DG , typename AT , typename AP >
Delaunay_vertex_handle CGAL::Voronoi_diagram_2< DG, AT, AP >::Halfedge::right ( )

Returns a handle to the vertex in the Delaunay graph corresponding to the defining site to the right of the Voronoi edge.

Precondition
has_target() must be true.
template<typename DG , typename AT , typename AP >
Vertex_handle CGAL::Voronoi_diagram_2< DG, AT, AP >::Halfedge::source ( )

Returns the source vertex of e.

Precondition
The source vertex must exist, i.e., has_source() must return true.
template<typename DG , typename AT , typename AP >
Vertex_handle CGAL::Voronoi_diagram_2< DG, AT, AP >::Halfedge::target ( )

Returns the target vertex of e.

Precondition
The target vertex must exist, i.e., has_target() must return true.
template<typename DG , typename AT , typename AP >
Halfedge_handle CGAL::Voronoi_diagram_2< DG, AT, AP >::Halfedge::twin ( )

Returns the twin halfedge.

template<typename DG , typename AT , typename AP >
Delaunay_vertex_handle CGAL::Voronoi_diagram_2< DG, AT, AP >::Halfedge::up ( )

Returns a handle to the vertex in the Delaunay graph corresponding to the defining site above the Voronoi edge.