CGAL 5.3.2 - 2D Hyperbolic Delaunay Triangulations
CGAL::Hyperbolic_Delaunay_triangulation_2< Gt, Tds > Class Template Reference

#include <CGAL/Hyperbolic_Delaunay_triangulation_2.h>

Inherits from

CGAL::Delaunay_triangulation_2< Gt, Tds >.

Definition

The class Hyperbolic_Delaunay_triangulation_2 is the main class of the 2D Hyperbolic Delaunay Triangulations package.

It is designed to represent Delaunay triangulations of sets of points in the hyperbolic plane. The hyperbolic plane is represented in the Poincaré disk model.

Template Parameters
Gtis the geometric traits class and must be a model of HyperbolicDelaunayTriangulationTraits_2.
Tdsis the triangulation graph data structure and must be a model of TriangulationDataStructure_2 whose vertex and face are models of TriangulationVertexBase_2 and HyperbolicTriangulationFaceBase_2, respectively. It defaults to:
See also
Delaunay_triangulation_2
Examples:
Hyperbolic_triangulation_2/ht2_example.cpp, and Hyperbolic_triangulation_2/ht2_example_color.cpp.

Types

typedef Gt Geom_traits
 
typedef Tds Triangulation_data_structure
 
typedef Triangulation_data_structure::size_type size_type
 Size type (integral unsigned).
 
typedef Geom_traits::Hyperbolic_point_2 Point
 
typedef Geom_traits::Hyperbolic_triangle_2 Hyperbolic_triangle
 

The following types are defined to give access to the elements of the triangulation:

typedef Triangulation_data_structure::Vertex_handle Vertex_handle
 
typedef Triangulation_data_structure::Face_handle Face_handle
 
typedef Triangulation_data_structure::Edge Edge
 

The following types are defined for use in the construction of the Voronoi diagram:

typedef Geom_traits::Hyperbolic_Voronoi_point_2 Hyperbolic_Voronoi_point
 
typedef Geom_traits::Hyperbolic_segment_2 Hyperbolic_segment
 

The following iterator and circulator types are defined to give access over hyperbolic faces and edges:

typedef Triangulation_data_structure::Face_iterator All_faces_iterator
 
typedef Triangulation_data_structure::Edge_iterator All_edges_iterator
 
typedef Triangulation_data_structure::Vertex_iterator All_vertices_iterator
 
typedef Triangulation_data_structure::Vertex_circulator Vertex_circulator
 

The enumeration Locate_type is defined to specify which case occurs when locating a point in the triangulation.

enum  Locate_type {
  VERTEX = 0, EDGE, FACE, OUTSIDE_CONVEX_HULL,
  OUTSIDE_AFFINE_HULL
}
 

Creation

 Hyperbolic_Delaunay_triangulation_2 (const Geom_traits &gt=Geom_traits())
 Default constructor
 
 Hyperbolic_Delaunay_triangulation_2 (const Hyperbolic_Delaunay_triangulation_2< Gt, Tds > &tr)
 Copy constructor.
 
template<class InputIterator >
 Hyperbolic_Delaunay_triangulation_2 (InputIterator first, InputIterator last, const Geom_traits &gt=Geom_traits())
 Equivalent to creating an empty triangulation and calling insert(first, last).
 

Assignment

Hyperbolic_Delaunay_triangulation_2operator= (Hyperbolic_Delaunay_triangulation_2 tr)
 The triangulation tr is duplicated, and modifying the copy after the duplication does not modify the original.
 
void swap (Hyperbolic_Delaunay_triangulation_2 &tr)
 The triangulation is swapped with tr.
 
void clear ()
 Deletes all vertices and faces of the triangulation.
 
bool operator== (const Hyperbolic_Delaunay_triangulation_2< Gt, Tds > &t1, const Hyperbolic_Delaunay_triangulation_2< Gt, Tds > &t2)
 Equality operator. More...
 
bool operator!= (const Hyperbolic_Delaunay_triangulation_2< Gt, Tds > &t1, const Hyperbolic_Delaunay_triangulation_2< Gt, Tds > &t2)
 Inequality operator. More...
 

Access functions

const Geom_traitsgeom_traits () const
 Returns a const reference to the geometric traits object.
 
const Triangulation_data_structuretds () const
 Returns a const reference to the triangulation data structure.
 
Triangulation_data_structuretds ()
 Returns a reference to the triangulation data structure.
 
bool is_valid ()
 Checks the combinatorial validity of the triangulation, the validity of its geometric embedding, and also that all edges and faces are Delaunay hyperbolic.
 
int dimension () const
 Returns the dimension of the affine hull.
 
size_type number_of_vertices () const
 Returns the number of vertices.
 
size_type number_of_hyperbolic_edges () const
 Returns the number of hyperbolic edges.
 
size_type number_of_hyperbolic_faces () const
 Returns the number of hyperbolic faces.
 

Geometric access functions

Hyperbolic_triangle hyperbolic_triangle (const Face_handle f) const
 
Hyperbolic_segment hyperbolic_segment (const Face_handle f, const int i) const
 Returns the hyperbolic segment formed by the vertices of the edge (f, i).
 
Hyperbolic_segment hyperbolic_segment (const Edge &e) const
 Returns the hyperbolic segment formed by the vertices of edge e.
 

Insertion

Vertex_handle insert (const Point &p, Face_handle start=Face_handle())
 Inserts the point p in the triangulation. More...
 
Vertex_handle insert (const Point &p, typename Locate_type lt, Face_handle loc, int li)
 Inserts the point p at the location given by (lt,loc,li). More...
 
template<class InputIterator >
std::ptrdiff_t insert (InputIterator first, InputIterator last)
 Inserts the points in the range [first,last) into the triangulation. More...
 

Removal

void remove (Vertex_handle v)
 Removes the vertex v from the triangulation. More...
 
template<class VertexRemoveIterator >
void remove (VertexRemoveIterator first, VertexRemoveIterator last)
 Removes the vertices in the iterator range [firs, last) from the triangulation. More...
 

Point Location

Face_handle locate (const Point &query, const Face_handle hint=Face_handle()) const
 Locates the point query in the triangulation. More...
 
Face_handle locate (const Point &query, Locate_type &lt, int &li, Face_handle hint=Face_handle()) const
 Same as above. More...
 

Queries

template<class OutputItFaces >
OutputItFaces find_conflicts (const Point &p, OutputItFaces fit, Face_handle start=Face_handle()) const
 Computes the conflict zone induced by p. More...
 

Vertex, Face and Edge iterators and circulators

All_vertices_iterator all_vertices_begin () const
 
All_vertices_iterator all_vertices_end () const
 
All_edges_iterator all_edges_begin () const
 
All_edges_iterator all_edges_end () const
 
All_faces_iterator all_faces_begin () const
 
All_faces_iterator all_faces_end () const
 
Vertex_circulator adjacent_vertices (Vertex_handle v) const
 

Voronoi Diagram

Users 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).

Hyperbolic_Voronoi_point dual (Face_handle f) const
 Returns the hyperbolic center of the circumdisk of f. More...
 
Hyperbolic_segment dual (const Edge &e) const
 Returns the hyperbolic segment that is dual to e.
 
Hyperbolic_segment dual (Face_handle f, int i) const
 Returns the hyperbolic segment that is dual to the edge (f,i). More...
 

Member Enumeration Documentation

◆ Locate_type

template<class Gt , class Tds >
enum CGAL::Hyperbolic_Delaunay_triangulation_2::Locate_type
Enumerator
VERTEX 
EDGE 
FACE 
OUTSIDE_CONVEX_HULL 
OUTSIDE_AFFINE_HULL 

Member Function Documentation

◆ dual() [1/2]

template<class Gt , class Tds >
Hyperbolic_Voronoi_point CGAL::Hyperbolic_Delaunay_triangulation_2< Gt, Tds >::dual ( Face_handle  f) const

Returns the hyperbolic center of the circumdisk of f.

Precondition
f is hyperbolic

◆ dual() [2/2]

template<class Gt , class Tds >
Hyperbolic_segment CGAL::Hyperbolic_Delaunay_triangulation_2< Gt, Tds >::dual ( Face_handle  f,
int  i 
) const

Returns the hyperbolic segment that is dual to the edge (f,i).

Precondition
f is hyperbolic

◆ find_conflicts()

template<class Gt , class Tds >
template<class OutputItFaces >
OutputItFaces CGAL::Hyperbolic_Delaunay_triangulation_2< Gt, Tds >::find_conflicts ( const Point p,
OutputItFaces  fit,
Face_handle  start = Face_handle() 
) const

Computes the conflict zone induced by p.

If the optional parameter start is given, then it must be a face in conflict with p. Returns an iterator on the faces of the triangulation in conflict with p.

◆ insert() [1/3]

template<class Gt , class Tds >
Vertex_handle CGAL::Hyperbolic_Delaunay_triangulation_2< Gt, Tds >::insert ( const Point p,
Face_handle  start = Face_handle() 
)

Inserts the point p in the triangulation.

If the point p coincides with a existing vertex, then the vertex is returned and the triangulation is not modified. The optional parameter start is used to initialize the location of p.

◆ insert() [2/3]

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

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

The handle to the new vertex is returned.

See also
locate()

◆ insert() [3/3]

template<class Gt , class Tds >
template<class InputIterator >
std::ptrdiff_t CGAL::Hyperbolic_Delaunay_triangulation_2< Gt, Tds >::insert ( InputIterator  first,
InputIterator  last 
)

Inserts the points in the range [first,last) into the triangulation.

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

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

◆ locate() [1/2]

template<class Gt , class Tds >
Face_handle CGAL::Hyperbolic_Delaunay_triangulation_2< Gt, Tds >::locate ( const Point query,
const Face_handle  hint = Face_handle() 
) const

Locates the point query in the triangulation.

If the point query lies inside the hyperbolic convex hull of the points of the triangulation, then the hyperbolic face that contains the query in its interior is returned.

If query lies on a vertex or on an edge, then one of the faces having query on its boundary is returned.

If query lies outside of the convex hull of the points of the triangulation, then a default-constructed Face_handle() is returned.

The optional argument hint is used as a starting place for the search.

◆ locate() [2/2]

template<class Gt , class Tds >
Face_handle CGAL::Hyperbolic_Delaunay_triangulation_2< Gt, Tds >::locate ( const Point query,
Locate_type lt,
int &  li,
Face_handle  hint = Face_handle() 
) const

Same as above.

The variable lt contains information about the element in which query has been located. See the enumeration Locate_type for details.

If lt is Locate_type::VERTEX, then the variable li contains the index of the vertex in the returned face. If lt is Locate_type::EDGE, then li is the index of the edge in the returned face.

◆ operator!=()

template<class Gt , class Tds >
bool CGAL::Hyperbolic_Delaunay_triangulation_2< Gt, Tds >::operator!= ( const Hyperbolic_Delaunay_triangulation_2< Gt, Tds > &  t1,
const Hyperbolic_Delaunay_triangulation_2< Gt, Tds > &  t2 
)

Inequality operator.

◆ operator==()

template<class Gt , class Tds >
bool CGAL::Hyperbolic_Delaunay_triangulation_2< Gt, Tds >::operator== ( const Hyperbolic_Delaunay_triangulation_2< Gt, Tds > &  t1,
const Hyperbolic_Delaunay_triangulation_2< Gt, Tds > &  t2 
)

Equality operator.

◆ remove() [1/2]

template<class Gt , class Tds >
void CGAL::Hyperbolic_Delaunay_triangulation_2< Gt, Tds >::remove ( Vertex_handle  v)

Removes the vertex v from the triangulation.

Precondition
v is a vertex of the triangulation.

◆ remove() [2/2]

template<class Gt , class Tds >
template<class VertexRemoveIterator >
void CGAL::Hyperbolic_Delaunay_triangulation_2< Gt, Tds >::remove ( VertexRemoveIterator  first,
VertexRemoveIterator  last 
)

Removes the vertices in the iterator range [firs, last) from the triangulation.

Precondition
all vertices in [first, last) are vertices of the triangulation.