CGAL 5.1.2 - 2D Segment Delaunay Graphs
CGAL::Segment_Delaunay_graph_2< Gt, DS > Class Template Reference

#include <CGAL/Segment_Delaunay_graph_2.h>

Inherited by CGAL::Segment_Delaunay_graph_hierarchy_2< Gt, STag, DS >.

Definition

The class Segment_Delaunay_graph_2 represents the segment Delaunay graph (which is the dual graph of the 2D segment Voronoi diagram).

Currently it supports only insertions of sites.

Template Parameters
Gtmust be a model of SegmentDelaunayGraphTraits_2
DSmust be a model of SegmentDelaunayGraphDataStructure_2. DS defaults to CGAL::Triangulation_data_structure_2< CGAL::Segment_Delaunay_graph_vertex_base_2<Gt>, CGAL::Triangulation_face_base_2<Gt> >.

Traversal of the Segment Delaunay Graph

A segment Delaunay graph can be seen as a container of faces and vertices. Therefore the Segment_Delaunay_graph_2 class provides several iterators and circulators that allow to traverse it (completely or partially).

Is Model Of:
DelaunayGraph_2
See also
DelaunayGraph_2
SegmentDelaunayGraphTraits_2
SegmentDelaunayGraphDataStructure_2
SegmentDelaunayGraphVertexBase_2
TriangulationFaceBase_2
CGAL::Segment_Delaunay_graph_hierarchy_2<Gt,STag,DS>
CGAL::Segment_Delaunay_graph_traits_2<K,MTag>
CGAL::Segment_Delaunay_graph_traits_without_intersections_2<K,MTag>
CGAL::Segment_Delaunay_graph_filtered_traits_2<CK,CM,EK,EM,FK,FM>
CGAL::Segment_Delaunay_graph_filtered_traits_without_intersections_2<CK,CM,EK,EM,FK,FM>
CGAL::Triangulation_data_structure_2<Vb,Fb>
CGAL::Segment_Delaunay_graph_vertex_base_2<Gt,SSTag>
CGAL::Triangulation_face_base_2<Gt>
Examples:
Segment_Delaunay_graph_2/sdg-count-sites.cpp, Segment_Delaunay_graph_2/sdg-fast-sp-polygon.cpp, Segment_Delaunay_graph_2/sdg-fast-sp.cpp, Segment_Delaunay_graph_2/sdg-red-blue-info.cpp, and Segment_Delaunay_graph_2/sdg-voronoi-edges.cpp.

Types

typedef Gt Geom_traits
 A type for the geometric traits.
 
typedef DS Data_structure
 A type for the underlying data structure.
 
typedef Data_structure Triangulation_data_structure
 This type has been added so that the Segment_Delaunay_graph_2 class is a model of the DelaunayGraph_2 concept.
 
typedef DS::size_type size_type
 Size type (an unsigned integral type)
 
typedef Gt::Point_2 Point_2
 A type for the point defined in the geometric traits.
 
typedef Gt::Site_2 Site_2
 A type for the segment Delaunay graph site, defined in the geometric traits.
 
typedef unspecified_type Point_container
 A type for the container of points.
 
typedef Point_container::iterator Point_handle
 A handle for points in the point container.
 

Iterators and Handles

The vertices and faces of the segment Delaunay graph 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. The edges of the segment Delaunay graph can also be visited through iterators and circulators, the edge circulators and iterators are also bidirectional and non-mutable. In the following, we call infinite any face or edge incident to the infinite vertex and the infinite vertex itself. Any other feature (face, edge or vertex) of the segment Delaunay graph is said to be finite. Some iterators (the All iterators ) allow to visit finite or infinite features while the others (the Finite iterators) visit only finite features. Circulators visit both infinite and finite features.

typedef DS::Edge Edge
 The edge type. More...
 
typedef DS::Vertex Vertex
 A type for a vertex.
 
typedef DS::Face Face
 A type for a face.
 
typedef DS::Vertex_handle Vertex_handle
 A type for a handle to a vertex.
 
typedef DS::Face_handle Face_handle
 A type for a handle to a face.
 
typedef DS::Vertex_circulator Vertex_circulator
 A type for a circulator over vertices incident to a given vertex.
 
typedef DS::Face_circulator Face_circulator
 A type for a circulator over faces incident to a given vertex.
 
typedef DS::Edge_circulator Edge_circulator
 A type for a circulator over edges incident to a given vertex.
 
typedef DS::Vertex_iterator All_vertices_iterator
 A type for an iterator over all vertices.
 
typedef DS::Face_iterator All_faces_iterator
 A type for an iterator over all faces.
 
typedef DS::Edge_iterator All_edges_iterator
 A type for an iterator over all edges.
 
typedef unspecified_type Finite_vertices_iterator
 A type for an iterator over finite vertices.
 
typedef unspecified_type Finite_faces_iterator
 A type for an iterator over finite faces.
 
typedef unspecified_type Finite_edges_iterator
 A type for an iterator over finite edges.
 

Site Iterators

The following iterators allow respectively to visit all sites.

These iterators are non-mutable, bidirectional and their value type is Site_2. They are all invalidated by any change in the segment Delaunay graph.

typedef unspecified_type Input_sites_iterator
 A type for a bidirectional iterator over all input sites.
 
typedef unspecified_type Output_sites_iterator
 A type for a bidirectional iterator over all output sites (the sites in the Delaunay graph).
 
Input_sites_iterator input_sites_begin ()
 Starts at an arbitrary input site.
 
Input_sites_iterator input_sites_end ()
 Past-the-end iterator.
 
Output_sites_iterator output_sites_begin ()
 Starts at an arbitrary output site.
 
Output_sites_iterator output_sites_end ()
 Past-the-end iterator.
 

Creation

In addition to the default and copy constructors the following constructors are defined:

 Segment_Delaunay_graph_2 (Gt gt=Gt())
 Creates the segment Delaunay graph using gt as geometric traits.
 
template<class Input_iterator >
 Segment_Delaunay_graph_2 (Input_iterator first, Input_iterator beyond, Gt gt=Gt())
 Creates the segment Delaunay graph using gt as geometric traits and inserts all sites in the range [first, beyond). More...
 

Access Functions

Geom_traits geom_traits ()
 Returns a reference to the segment Delaunay graph traits object.
 
int dimension ()
 Returns the dimension of the segment Delaunay graph. More...
 
size_type number_of_vertices ()
 Returns the number of finite vertices of the segment Delaunay graph.
 
size_type number_of_faces ()
 Returns the number of faces (both finite and infinite) of the segment Delaunay graph.
 
size_type number_of_input_sites ()
 Return the number of input sites.
 
size_type number_of_output_sites ()
 Return the number of output sites. More...
 
Face_handle infinite_face ()
 Returns a face incident to the infinite_vertex.
 
Vertex_handle infinite_vertex ()
 Returns the infinite_vertex.
 
Vertex_handle finite_vertex ()
 Returns a vertex distinct from the infinite_vertex. More...
 
Data_structure data_structure ()
 Returns a reference to the segment Delaunay graph data structure object.
 
Data_structure tds ()
 Same as data_structure(). More...
 
Point_container point_container ()
 Returns a reference to the point container object.
 

Finite Face, Edge and Vertex Iterators

The following iterators allow respectively to visit finite faces, finite edges and finite vertices of the segment Delaunay graph.

These iterators are non-mutable, bidirectional and their value types are respectively Face, Edge and Vertex. They are all invalidated by any change in the segment Delaunay graph.

Finite_vertices_iterator finite_vertices_begin ()
 Starts at an arbitrary finite vertex.
 
Finite_vertices_iterator finite_vertices_end ()
 Past-the-end iterator.
 
Finite_edges_iterator finite_edges_begin ()
 Starts at an arbitrary finite edge.
 
Finite_edges_iterator finite_edges_end ()
 Past-the-end iterator.
 
Finite_faces_iterator finite_faces_begin ()
 Starts at an arbitrary finite face.
 
Finite_faces_iterator finite_faces_end () const
 Past-the-end iterator.
 

Infinite Face, Edge, and Vertex Iterators

The following iterators allow respectively to visit all (both finite and infinite) faces, edges and vertices of the segment Delaunay graph.

These iterators are non-mutable, bidirectional and their value types are respectively Face, Edge and Vertex. They are all invalidated by any change in the segment Delaunay graph.

All_vertices_iterator all_vertices_begin ()
 Starts at an arbitrary vertex.
 
All_vertices_iterator all_vertices_end ()
 Past-the-end iterator.
 
All_edges_iterator all_edges_begin ()
 Starts at an arbitrary edge.
 
All_edges_iterator all_edges_end ()
 Past-the-end iterator.
 
All_faces_iterator all_faces_begin ()
 Starts at an arbitrary face.
 
All_faces_iterator all_faces_end ()
 Past-the-end iterator.
 

Face, Edge and Vertex Circulators

The Segment_Delaunay_graph_2 class also provides circulators that allow to visit respectively all faces or edges incident to a given vertex or all vertices adjacent to a given vertex.

These circulators are non-mutable and bidirectional. The operator operator++ moves the circulator counterclockwise around the vertex while the operator- moves clockwise. A face circulator is invalidated by any modification of the face pointed to. An edge circulator is invalidated by any modification in one of the two faces incident to the edge pointed to. A vertex circulator is invalidated by any modification in any of the faces adjacent to the vertex pointed to.

Applied on the infinite_vertex the above methods allow to visit the vertices on the convex hull and the infinite edges and faces. Note that a counterclockwise traversal of the vertices adjacent to the infinite_vertex is a clockwise traversal of the convex hull.

Face_circulator incident_faces (Vertex_handle v)
 Starts at an arbitrary face incident to v.
 
Face_circulator incident_faces (Vertex_handle v, Face_handle f)
 Starts at face f. More...
 
Edge_circulator incident_edges (Vertex_handle v)
 Starts at an arbitrary edge incident to v.
 
Edge_circulator incident_edges (Vertex_handle v, Face_handle f)
 Starts at the first edge of f incident to v, in counterclockwise order around v. More...
 
Vertex_circulator incident_vertices (Vertex_handle v)
 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...
 

Predicates

The class Segment_Delaunay_graph_2 provides methods to test the finite or infinite character of any feature.

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.
 

Insertion

template<class Input_iterator >
size_type insert (Input_iterator first, Input_iterator beyond)
 Inserts the sites in the range [first,beyond). More...
 
template<class Input_iterator >
size_type insert (Input_iterator first, Input_iterator beyond, Tag_false)
 Same as the previous method. More...
 
template<class Input_iterator >
size_type insert (Input_iterator first, Input_iterator beyond, Tag_true)
 Decomposes the range [first,beyond) into a range of input points and a range of input segments that are respectively passed to insert_segments() and insert_points(). More...
 
template<class PointIterator >
std::size_t insert_points (PointIterator first, PointIterator beyond)
 Inserts the points in the range [first,beyond) as sites. More...
 
template<class SegmentIterator >
std::size_t insert_segments (SegmentIterator first, SegmentIterator beyond)
 Inserts the segments in the range [first,beyond) as sites. More...
 
template<class PointIterator , class IndicesIterator >
std::size_t insert_segments (PointIterator points_first, PointIterator points_beyond, IndicesIterator indices_first, IndicesIterator indices_beyond)
 Same as above except that each segment is given as a pair of indices of the points in the range [points_first, points_beyond). More...
 
Vertex_handle insert (Point_2 p)
 Inserts the point p in the segment Delaunay graph. More...
 
Vertex_handle insert (Point_2 p, Vertex_handle vnear)
 Inserts p in the segment Delaunay graph using the site associated with vnear as an estimate for the nearest neighbor of p. More...
 
Vertex_handle insert (Point_2 p1, Point_2 p2)
 Inserts the closed segment with endpoints p1 and p2 in the segment Delaunay graph. More...
 
Vertex_handle insert (Point_2 p1, Point_2 p2, Vertex_handle vnear)
 Inserts the segment whose endpoints are p1 and p2 in the segment Delaunay graph using the site associated with vnear as an estimate for the nearest neighbor of p1. More...
 
Vertex_handle insert (Site_2 s)
 Inserts the site s in the segment Delaunay graph. More...
 
Vertex_handle insert (Site_2 s, Vertex_handle vnear)
 Inserts s in the segment Delaunay graph using the site associated with vnear as an estimate for the nearest neighbor of s, if s is a point, or the first endpoint of s, if s is a segment. More...
 

Nearest neighbor location

Vertex_handle nearest_neighbor (Point_2 p)
 Finds the nearest neighbor of the point p. More...
 
Vertex_handle nearest_neighbor (Point_2 p, Vertex_handle vnear)
 Finds the nearest neighbor of the point p using the site associated with vnear as an estimate for the nearest neighbor of p. More...
 

I/O

template<class Stream >
Stream & draw_dual (Stream &str)
 Draws the segment Voronoi diagram to the stream str. More...
 
template<class Stream >
Stream & draw_skeleton (Stream &str)
 Draws the segment Voronoi diagram to the stream str, except the edges of the diagram corresponding to a segment and its endpoints. More...
 
template<class Stream >
Stream & draw_dual_edge (Edge e, Stream &str)
 Draws the edge e of the segment Voronoi diagram to the stream str. More...
 
template<class Stream >
Stream & draw_dual_edge (Edge_circulator ec, Stream &str)
 Draws the edge *ec of the segment Voronoi diagram to the stream str. More...
 
template<class Stream >
Stream & draw_dual_edge (All_edges_iterator eit, Stream &str)
 Draws the edge *eit of the segment Voronoi diagram to the stream str. More...
 
template<class Stream >
Stream & draw_dual_edge (Finite_edges_iterator eit, Stream &str)
 Draws the edge *eit of the segment Voronoi diagram to the stream str. More...
 
void file_output (std::ostream &os)
 Writes the current state of the segment Delaunay graph to an output stream. More...
 
void file_input (std::istream &is)
 Reads the state of the segment Delaunay graph from an input stream.
 
std::ostream & operator<< (std::ostream &os, Segment_Delaunay_graph_2< Gt, DS > sdg)
 Writes the current state of the segment Delaunay graph to an output stream.
 
std::istream & operator>> (std::istream &is, Segment_Delaunay_graph_2< Gt, DS > sdg)
 Reads the state of the segment Delaunay graph from an input stream.
 

Validity check

bool is_valid (bool verbose=false, int level=1)
 Checks the validity of the segment Delaunay graph. More...
 

Miscellaneous

void clear ()
 Clears all contents of the segment Delaunay graph.
 
void swap (Segment_Delaunay_graph_2< Gt, DS > other)
 The segment Delaunay graphs other and *this are swapped. More...
 

Member Typedef Documentation

◆ Edge

template<typename Gt , typename DS >
typedef DS::Edge CGAL::Segment_Delaunay_graph_2< Gt, DS >::Edge

The edge type.

The Edge(f,i) is the edge common to faces f and f.neighbor(i). It is also the edge joining the vertices vertex(cw(i)) and vertex(ccw(i)) of f.

Precondition
i must be 0, 1 or 2.

Constructor & Destructor Documentation

◆ Segment_Delaunay_graph_2()

template<typename Gt , typename DS >
template<class Input_iterator >
CGAL::Segment_Delaunay_graph_2< Gt, DS >::Segment_Delaunay_graph_2 ( Input_iterator  first,
Input_iterator  beyond,
Gt  gt = Gt() 
)

Creates the segment Delaunay graph using gt as geometric traits and inserts all sites in the range [first, beyond).

Precondition
Input_iterator must be a model of InputIterator. The value type of Input_iterator must be either Point_2 or Site_2.

Member Function Documentation

◆ dimension()

template<typename Gt , typename DS >
int CGAL::Segment_Delaunay_graph_2< Gt, DS >::dimension ( )

Returns the dimension of the segment Delaunay graph.

The dimension is \( -1\) if the graph contains no sites, \( 0\) if the graph contains one site, \( 1\) if it contains two sites and \( 2\) if it contains three or more sites.

◆ draw_dual()

template<typename Gt , typename DS >
template<class Stream >
Stream& CGAL::Segment_Delaunay_graph_2< Gt, DS >::draw_dual ( Stream &  str)

Draws the segment Voronoi diagram to the stream str.

The following operators must be defined:

Stream& operator<<(Stream&, Gt::Segment_2)

Stream& operator<<(Stream&, Gt::Ray_2)

Stream& operator<<(Stream&, Gt::Line_2)

◆ draw_dual_edge() [1/4]

template<typename Gt , typename DS >
template<class Stream >
Stream& CGAL::Segment_Delaunay_graph_2< Gt, DS >::draw_dual_edge ( Edge  e,
Stream &  str 
)

Draws the edge e of the segment Voronoi diagram to the stream str.

The following operators must be defined:

Stream& operator<<(Stream&, Gt::Segment_2)

Stream& operator<<(Stream&, Gt::Ray_2)

Stream& operator<<(Stream&, Gt::Line_2)

Precondition
e must be a finite edge.

◆ draw_dual_edge() [2/4]

template<typename Gt , typename DS >
template<class Stream >
Stream& CGAL::Segment_Delaunay_graph_2< Gt, DS >::draw_dual_edge ( Edge_circulator  ec,
Stream &  str 
)

Draws the edge *ec of the segment Voronoi diagram to the stream str.

The following operators must be defined:

Stream& operator<<(Stream&, Gt::Segment_2)

Stream& operator<<(Stream&, Gt::Ray_2)

Stream& operator<<(Stream&, Gt::Line_2)

Precondition
*ec must be a finite edge.

◆ draw_dual_edge() [3/4]

template<typename Gt , typename DS >
template<class Stream >
Stream& CGAL::Segment_Delaunay_graph_2< Gt, DS >::draw_dual_edge ( All_edges_iterator  eit,
Stream &  str 
)

Draws the edge *eit of the segment Voronoi diagram to the stream str.

The following operators must be defined:

Stream& operator<<(Stream&, Gt::Segment_2)

Stream& operator<<(Stream&, Gt::Ray_2)

Stream& operator<<(Stream&, Gt::Line_2)

Precondition
*eit must be a finite edge.

◆ draw_dual_edge() [4/4]

template<typename Gt , typename DS >
template<class Stream >
Stream& CGAL::Segment_Delaunay_graph_2< Gt, DS >::draw_dual_edge ( Finite_edges_iterator  eit,
Stream &  str 
)

Draws the edge *eit of the segment Voronoi diagram to the stream str.

The following operators must be defined:

Stream& operator<<(Stream&, Gt::Segment_2)

Stream& operator<<(Stream&, Gt::Ray_2)

Stream& operator<<(Stream&, Gt::Line_2)

◆ draw_skeleton()

template<typename Gt , typename DS >
template<class Stream >
Stream& CGAL::Segment_Delaunay_graph_2< Gt, DS >::draw_skeleton ( Stream &  str)

Draws the segment Voronoi diagram to the stream str, except the edges of the diagram corresponding to a segment and its endpoints.

The following operators must be defined:

Stream& operator<<(Stream&, Gt::Segment_2)

Stream& operator<<(Stream&, Gt::Ray_2)

Stream& operator<<(Stream&, Gt::Line_2)

◆ file_output()

template<typename Gt , typename DS >
void CGAL::Segment_Delaunay_graph_2< Gt, DS >::file_output ( std::ostream &  os)

Writes the current state of the segment Delaunay graph to an output stream.

In particular, all sites in the diagram are written to the stream (represented through appropriate input sites), as well as the underlying combinatorial data structure.

◆ finite_vertex()

template<typename Gt , typename DS >
Vertex_handle CGAL::Segment_Delaunay_graph_2< Gt, DS >::finite_vertex ( )

Returns a vertex distinct from the infinite_vertex.

Precondition
The number of sites in the segment Delaunay graph must be at least one.

◆ incident_edges()

template<typename Gt , typename DS >
Edge_circulator CGAL::Segment_Delaunay_graph_2< Gt, DS >::incident_edges ( Vertex_handle  v,
Face_handle  f 
)

Starts at the first edge of f incident to v, in counterclockwise order around v.

Precondition
Face f is incident to vertex v.

◆ incident_faces()

template<typename Gt , typename DS >
Face_circulator CGAL::Segment_Delaunay_graph_2< Gt, DS >::incident_faces ( Vertex_handle  v,
Face_handle  f 
)

Starts at face f.

Precondition
Face f is incident to vertex v.

◆ incident_vertices()

template<typename Gt , typename DS >
Vertex_circulator CGAL::Segment_Delaunay_graph_2< Gt, DS >::incident_vertices ( Vertex_handle  v,
Face_handle  f 
)

Starts at the first vertex of f adjacent to v in counterclockwise order around v.

Precondition
Face f is incident to vertex v.

◆ insert() [1/9]

template<typename Gt , typename DS >
template<class Input_iterator >
size_type CGAL::Segment_Delaunay_graph_2< Gt, DS >::insert ( Input_iterator  first,
Input_iterator  beyond 
)

Inserts the sites in the range [first,beyond).

The number of additional sites inserted in the Delaunay graph is returned. Input_iterator must be a model of InputIterator and its value type must be either Point_2 or Site_2.

◆ insert() [2/9]

template<typename Gt , typename DS >
template<class Input_iterator >
size_type CGAL::Segment_Delaunay_graph_2< Gt, DS >::insert ( Input_iterator  first,
Input_iterator  beyond,
Tag_false   
)

Same as the previous method.

Input_iterator must be a model of InputIterator and its value type must be either Point_2 or Site_2.

◆ insert() [3/9]

template<typename Gt , typename DS >
template<class Input_iterator >
size_type CGAL::Segment_Delaunay_graph_2< Gt, DS >::insert ( Input_iterator  first,
Input_iterator  beyond,
Tag_true   
)

Decomposes the range [first,beyond) into a range of input points and a range of input segments that are respectively passed to insert_segments() and insert_points().

Non-input sites are first random_shuffled and then inserted one by one. Input_iterator must be a model of InputIterator and its value type must be either Point_2, Segment_2 or Site_2.

Returns
the number of sites inserted in the Delaunay graph

◆ insert() [4/9]

template<typename Gt , typename DS >
Vertex_handle CGAL::Segment_Delaunay_graph_2< Gt, DS >::insert ( Point_2  p)

Inserts the point p in the segment Delaunay graph.

If p has already been inserted, then the vertex handle of its already inserted copy is returned. If p has not been inserted yet, the vertex handle of p is returned.

◆ insert() [5/9]

template<typename Gt , typename DS >
Vertex_handle CGAL::Segment_Delaunay_graph_2< Gt, DS >::insert ( Point_2  p,
Vertex_handle  vnear 
)

Inserts p in the segment Delaunay graph using the site associated with vnear as an estimate for the nearest neighbor of p.

The vertex handle returned has the same semantics as the vertex handle returned by the method Vertex_handle insert(Point_2 p).

◆ insert() [6/9]

template<typename Gt , typename DS >
Vertex_handle CGAL::Segment_Delaunay_graph_2< Gt, DS >::insert ( Point_2  p1,
Point_2  p2 
)

Inserts the closed segment with endpoints p1 and p2 in the segment Delaunay graph.

If the segment has already been inserted in the Delaunay graph then the vertex handle of its already inserted copy is returned. If the segment does not intersect any segment in the existing diagram, the vertex handle corresponding to its corresponding open segment is returned. Finally, if the segment intersects other segments in the existing Delaunay graph, the vertex handle to one of its open subsegments is returned.

◆ insert() [7/9]

template<typename Gt , typename DS >
Vertex_handle CGAL::Segment_Delaunay_graph_2< Gt, DS >::insert ( Point_2  p1,
Point_2  p2,
Vertex_handle  vnear 
)

Inserts the segment whose endpoints are p1 and p2 in the segment Delaunay graph using the site associated with vnear as an estimate for the nearest neighbor of p1.

The vertex handle returned has the same semantics as the vertex handle returned by the method Vertex_handle insert(Point_2 p1, Point_2 p2).

◆ insert() [8/9]

template<typename Gt , typename DS >
Vertex_handle CGAL::Segment_Delaunay_graph_2< Gt, DS >::insert ( Site_2  s)

Inserts the site s in the segment Delaunay graph.

The vertex handle returned has the same semantics as the vertex handle returned by the methods Vertex_handle insert(Point_2 p) and Vertex_handle insert(Point_2 p1, Point_2 p2), depending on whether s represents a point or a segment respectively.

Precondition
s.is_input() must be true.

◆ insert() [9/9]

template<typename Gt , typename DS >
Vertex_handle CGAL::Segment_Delaunay_graph_2< Gt, DS >::insert ( Site_2  s,
Vertex_handle  vnear 
)

Inserts s in the segment Delaunay graph using the site associated with vnear as an estimate for the nearest neighbor of s, if s is a point, or the first endpoint of s, if s is a segment.

The vertex handle returned has the same semantics as the vertex handle returned by the method Vertex_handle insert(Site_2 s).

Precondition
s.is_input() must be true.

◆ insert_points()

template<typename Gt , typename DS >
template<class PointIterator >
std::size_t CGAL::Segment_Delaunay_graph_2< Gt, DS >::insert_points ( PointIterator  first,
PointIterator  beyond 
)

Inserts the points in the range [first,beyond) as sites.

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

Returns
the number of points inserted in the Delaunay graph
Template Parameters
PointIteratormust be an input iterator Point_2 or Site_2 as value_type.

◆ insert_segments() [1/2]

template<typename Gt , typename DS >
template<class SegmentIterator >
std::size_t CGAL::Segment_Delaunay_graph_2< Gt, DS >::insert_segments ( SegmentIterator  first,
SegmentIterator  beyond 
)

Inserts the segments in the range [first,beyond) as sites.

Note that this function is not guaranteed to insert the segments following the order of SegmentIterator, as spatial_sort() is used to improve efficiency.

Returns
the number of segments inserted in the Delaunay graph
Template Parameters
SegmentIteratormust be an input iterator with Site_2, Segment_2 or std::pair<Point_2,Point_2> as value type.

◆ insert_segments() [2/2]

template<typename Gt , typename DS >
template<class PointIterator , class IndicesIterator >
std::size_t CGAL::Segment_Delaunay_graph_2< Gt, DS >::insert_segments ( PointIterator  points_first,
PointIterator  points_beyond,
IndicesIterator  indices_first,
IndicesIterator  indices_beyond 
)

Same as above except that each segment is given as a pair of indices of the points in the range [points_first, points_beyond).

The indices must start from 0 to std::distance(points_first, points_beyond)

Template Parameters
PointIteratoris an input iterator with Point_2 as value type.
IndicesIteratoris an input iterator with std::pair<std::size_t, std::size_t> as value type.

◆ is_valid()

template<typename Gt , typename DS >
bool CGAL::Segment_Delaunay_graph_2< Gt, DS >::is_valid ( bool  verbose = false,
int  level = 1 
)

Checks the validity of the segment Delaunay graph.

If verbose is true a short message is sent to std::cerr. If level is 0, only the data structure is validated. If level is 1, then both the data structure and the segment Delaunay graph are validated. Negative values of level always return true, and values greater than 1 are equivalent to level being 1.

◆ nearest_neighbor() [1/2]

template<typename Gt , typename DS >
Vertex_handle CGAL::Segment_Delaunay_graph_2< Gt, DS >::nearest_neighbor ( Point_2  p)

Finds the nearest neighbor of the point p.

In other words it finds the site whose segment Voronoi diagram cell contains p. Ties are broken arbitrarily and one of the nearest neighbors of p is returned. If there are no sites in the segment Delaunay graph Vertex_handle() is returned.

◆ nearest_neighbor() [2/2]

template<typename Gt , typename DS >
Vertex_handle CGAL::Segment_Delaunay_graph_2< Gt, DS >::nearest_neighbor ( Point_2  p,
Vertex_handle  vnear 
)

Finds the nearest neighbor of the point p using the site associated with vnear as an estimate for the nearest neighbor of p.

Ties are broken arbitrarily and one of the nearest neighbors of p is returned. If there are no sites in the segment Delaunay graph Vertex_handle() is returned.

◆ number_of_output_sites()

template<typename Gt , typename DS >
size_type CGAL::Segment_Delaunay_graph_2< Gt, DS >::number_of_output_sites ( )

Return the number of output sites.

This is equal to the number of vertices in the segment Delaunay graph.

◆ swap()

template<typename Gt , typename DS >
void CGAL::Segment_Delaunay_graph_2< Gt, DS >::swap ( Segment_Delaunay_graph_2< Gt, DS >  other)

The segment Delaunay graphs other and *this are swapped.

For a segment Delaunay graph sdg, the operation sdg.swap(other) should be preferred to sdg= other or to sdg(other) if other is deleted afterwards.

◆ tds()

template<typename Gt , typename DS >
Data_structure CGAL::Segment_Delaunay_graph_2< Gt, DS >::tds ( )

Same as data_structure().

It has been added for compliance to the DelaunayGraph_2 concept.