CGAL 6.0 - 2D Periodic Hyperbolic Triangulations
|
#include <CGAL/Periodic_4_hyperbolic_Delaunay_triangulation_2.h>
CGAL::Periodic_4_hyperbolic_triangulation_2< GT, TDS >.
The class Periodic_4_hyperbolic_Delaunay_triangulation_2
enables the construction and handling of Delaunay triangulations of the Bolza surface.
The class expects two template parameters.
GT | Geometric Traits parameter. Must be a model of the concept Periodic_4HyperbolicDelaunayTriangulationTraits_2 . The default value for this parameter is Periodic_4_hyperbolic_Delaunay_triangulation_traits_2<> . |
TDS | Triangulation Data Structure parameter. Must be a model of the concept TriangulationDataStructure_2 with some additional functionality in the vertices and faces, following the concepts Periodic_4HyperbolicTriangulationVertexBase_2 and Periodic_4HyperbolicTriangulationFaceBase_2 , respectively. The default value for this parameter is The class Periodic_4_hyperbolic_triangulation_face_base_2 is the default model for the concept Period... Definition: Periodic_4_hyperbolic_triangulation_face_base_2.h:32 The class Periodic_4_hyperbolic_triangulation_vertex_base_2 is the default model for the concept Peri... Definition: Periodic_4_hyperbolic_triangulation_vertex_base_2.h:31 |
Creation | |
Periodic_4_hyperbolic_Delaunay_triangulation_2 (const Geom_traits >=Geom_traits()) | |
Default constructor with an optional parameter for the geometric traits object. | |
Periodic_4_hyperbolic_Delaunay_triangulation_2 (const Periodic_4_hyperbolic_Delaunay_triangulation_2 &tr) | |
Copy constructor. | |
template<class InputIterator > | |
Periodic_4_hyperbolic_Delaunay_triangulation_2 (InputIterator first, InputIterator last, const Geom_traits >=Geom_traits()) | |
Point range constructor. | |
Clearing the triangulation | |
void | clear () |
Deletes all faces and vertices of the triangulation, and re-initializes the triangulation data structure with the set of dummy points. | |
Queries | |
template<class OutputFaceIterator > | |
void | find_conflicts (const Point &p, OutputFaceIterator it, Face_handle start=Face_handle(), Hyperbolic_translation ltr=Hyperbolic_translation()) const |
Computes the conflict zone induced by p . | |
bool | is_dummy_vertex (Vertex_handle vh) const |
Checks if the vertex vh is part of the set of vertices storing dummy points. | |
int | number_of_dummy_points () const |
Returns the number of dummy points currently existing in the triangulation. | |
bool | is_valid (bool verbose=false) const |
Checks the combinatorial validity of the triangulation, and also the validity of its geometric embedding. | |
Point insertion | |
Vertex_handle | insert (const Point &p, Face_handle start=Face_handle()) |
Inserts the point p in the triangulation. | |
template<class InputIterator > | |
std::ptrdiff_t | insert (InputIterator first, InputIterator last, const bool flag_try_to_remove_dummy_vertices=true) |
Inserts all points in the input iterator into the triangulation. | |
Vertex removal | |
bool | remove (Vertex_handle v) |
Removes the vertex v from the triangulation. | |
template<class VertexRemoveIterator > | |
void | remove (VertexRemoveIterator first, VertexRemoveIterator last) |
Removes the vertices in the iterator range [firs, last) from the triangulation. | |
Removal of dummy points | |
int | try_to_remove_dummy_vertices () |
Tries to remove the dummy points from the triangulation one by one. | |
Voronoi diagram | |
Voronoi_point | dual (Face_handle f, Hyperbolic_translation nbtr=Hyperbolic_translation()) const |
Returns the dual of face f , i.e., the hyperbolic center of the disk defined by the vertices of f . | |
Segment | dual (const Edge &e) const |
Returns the hyperbolic segment that is dual to the edge e . | |
Additional Inherited Members | |
Public Types inherited from CGAL::Periodic_4_hyperbolic_triangulation_2< GT, TDS > | |
enum | Locate_type { VERTEX = 0 , EDGE , FACE } |
typedef GT | Geometric_traits |
typedef TDS | Triangulation_data_structure |
typedef Triangulation_data_structure::Vertex::Point | Point |
Represents a point in the hyperbolic plane. | |
typedef Geometric_traits::Hyperbolic_segment_2 | Hyperbolic_segment |
Represents a hyperbolic segment in the hyperbolic plane. | |
typedef Geometric_traits::Hyperbolic_triangle_2 | Hyperbolic_triangle |
Represents a triangle contained in the Poincaré disk. | |
typedef Geometric_traits::Hyperbolic_translation | Hyperbolic_translation |
typedef std::pair< Point, Hyperbolic_translation > | Periodic_point |
Represents a periodic point, i.e., a pair of a point in the original octagon and a hyperbolic translation. | |
typedef std::array< std::pair< Point, Hyperbolic_translation >, 2 > | Periodic_segment |
Represents a periodic segment, defined by two periodic points. | |
typedef std::array< std::pair< Point, Hyperbolic_translation >, 3 > | Periodic_triangle |
Represents a periodic triangle, defined by three periodic points. | |
typedef Triangulation_data_structure::Edge | Edge |
typedef Triangulation_data_structure::Face | Face |
typedef Triangulation_data_structure::Vertex_handle | Vertex_handle |
typedef Triangulation_data_structure::Face_handle | Face_handle |
typedef Triangulation_data_structure::Face_iterator | Face_iterator |
typedef Triangulation_data_structure::Edge_iterator | Edge_iterator |
typedef Triangulation_data_structure::Vertex_iterator | Vertex_iterator |
typedef Triangulation_data_structure::Face_circulator | Face_circulator |
typedef Triangulation_data_structure::Edge_circulator | Edge_circulator |
typedef Triangulation_data_structure::Vertex_circulator | Vertex_circulator |
Public Member Functions inherited from CGAL::Periodic_4_hyperbolic_triangulation_2< GT, TDS > | |
Periodic_4_hyperbolic_triangulation_2 (const Geometric_traits >=Geometric_traits()) | |
Default constructor, with an optional parameter for the geometric traits object. | |
Periodic_4_hyperbolic_triangulation_2 (const Periodic_4_hyperbolic_triangulation_2 &tr) | |
Copy constructor. | |
Periodic_4_hyperbolic_triangulation_2 & | operator= (Periodic_4_hyperbolic_triangulation_2 tr) |
The triangulation tr is duplicated, and modifying the copy after the duplication does not modify the original. | |
void | swap (Periodic_4_hyperbolic_triangulation_2 &tr) |
The triangulation is swapped with tr . | |
void | clear () |
Deletes all faces and vertices of the triangulation. | |
const Geometric_traits & | geom_traits () const |
Returns a const reference to the geometric traits object. | |
const Triangulation_data_structure & | tds () const |
Returns a const reference to the triangulation data structure. | |
Triangulation_data_structure & | tds () |
Returns a reference to the triangulation data structure. | |
size_type | number_of_vertices () const |
Returns the number of vertices in the triangulation. | |
size_type | number_of_edges () const |
Returns the number of edges in the triangulation. | |
size_type | number_of_faces () const |
Returns the number of faces in the triangulation. | |
Periodic_point | periodic_point (const Face_handle f, int i) const |
Returns the periodic point given by the i -th vertex of face f , that is the point in the original domain and the translation of the vertex in f . | |
Periodic_segment | periodic_segment (const Point &p1, const Point &p2, const Hyperbolic_translation &tr1, const Hyperbolic_translation &tr2) const |
Returns the periodic segment formed by the two point-translation pairs (p1, tr1) and (p2, tr2) . | |
Periodic_segment | periodic_segment (const Point &p1, const Point &p2) const |
Returns the periodic segment formed by the two point-translation pairs (p1, Id) and (p2, Id) , where Id is the identity translation. | |
Periodic_segment | periodic_segment (const Face_handle c, int i, int j) const |
Returns the periodic segment formed by the endpoints of edge (f,i,j) . | |
Periodic_segment | periodic_segment (const Edge &e) const |
Returns the periodic segment formed by the endpoints of edge e . | |
Periodic_segment | periodic_segment (const Edge &e, const Hyperbolic_translation &tr) const |
Returns the periodic segment formed by the endpoints of edge e translated by tr . | |
Periodic_triangle | periodic_triangle (const Face &f) const |
Returns the periodic triangle formed by the pairs of points and translations corresponding to each vertex of f . | |
Periodic_triangle | periodic_triangle (const Face &f, const Hyperbolic_translation &tr) const |
Returns the periodic triangle formed by the pairs of points and translations corresponding to each vertex of f , translated by tr . | |
Point | construct_point (const Periodic_point &pp) const |
Converts the periodic point pp into a point in the hyperbolic plane by applying pp.second to pp.first . | |
Hyperbolic_segment | construct_hyperbolic_segment (const Face_handle &fh, int idx) const |
Constructs the hyperbolic segment formed by the endpoints of edge (fh, idx) . | |
Hyperbolic_segment | construct_hyperbolic_segment (const Point &p1, const Point &p2) const |
Returns the hyperbolic segment with endpoints p1 and p2 . | |
Hyperbolic_segment | construct_hyperbolic_segment (const pair< Face_handle, int > &e) const |
Returns the hyperbolic segment formed by the endpoints of e . | |
Hyperbolic_segment | construct_hyperbolic_segment (const Periodic_segment &ps) const |
Returns the hyperbolic segment formed by the endpoints of ps . | |
Triangle | construct_triangle (const Face_handle &fh) const |
Returns the triangle formed by the vertices of fh . | |
Triangle | construct_triangle (const Periodic_triangle &pt) const |
Returns the triangle formed by the vertices of pt . | |
Face_handle | hyperbolic_periodic_locate (const Point &p, Hyperbolic_translation &lo, const Face_handle start=Face_handle()) const |
Returns the face rf for which the periodic triangle (rf, lo) contains the query point p . | |
Face_handle | hyperbolic_periodic_locate (const Point &p, Locate_type <, int &li, Hyperbolic_translation &lo, const Face_handle start=Face_handle()) const |
Same as above. | |
Face_handle | hyperbolic_locate (const Point &p, Face_handle start=Face_handle()) const |
Returns the canonical representative of the face that contains the query point p . | |
Face_handle | hyperbolic_locate (const Point &p, Locate_type <, int &li, Face_handle start=Face_handle()) const |
Same as above. | |
Orientation | orientation (const Point &p1, const Point &p2, const Point &p3) const |
Returns the Orientation of the points p1, p2, p3 . | |
Orientation | orientation (const Point &p1, const Point &p2, const Point &p3, const Hyperbolic_translation &tr1, const Hyperbolic_translation &tr2, const Hyperbolic_translation &tr3) const |
Returns the Orientation of the three periodic points (p1, tr1), (p2, tr2) and (p3, tr3) . | |
Oriented_side | side_of_oriented_circle (const Point &p1, const Point &p2, const Point &p3, const Point &q) const |
Returns the Oriented_side on which q is located with respect to the circle defined by the points p1,p2,p3 . | |
Oriented_side | side_of_oriented_circle (const Point &p1, const Point &p2, const Point &p3, const Point &q, const Hyperbolic_translation &tr1, const Hyperbolic_translation &tr2, const Hyperbolic_translation &tr3, const Hyperbolic_translation &trq) const |
Returns the Oriented_side on which the periodic point (q, trq) is located with respect to the circle defined by the periodic points (p1, tr1), (p2, tr2) and (p3, tr3) . | |
bool | is_vertex (Vertex_handle v) const |
Tests whether v is a vertex of the triangulation. | |
bool | is_edge (Vertex_handle u, Vertex_handle v, Face_handle &fh, int &i) const |
Tests whether (u, v) is an edge of the triangulation. | |
bool | is_face (Vertex_handle u, Vertex_handle v, Vertex_handle w, Face_handle &fh) const |
Tests whether there exists a face in the triangulation with vertices u, v and w . | |
bool | has_vertex (const Face_handle f, const Vertex_handle v, const int i) const |
Tests whether face f has v as a vertex. | |
Vertex_iterator | vertices_begin () const |
Starts at an arbitrary vertex. | |
Vertex_iterator | vertices_end () const |
Past-the-end iterator. | |
Edge_iterator | edges_begin () const |
Starts at an arbitrary edge. | |
Edge_iterator | edges_end () const |
Past-the-end iterator. | |
Face_iterator | faces_begin () const |
Starts at an arbitrary face. | |
Face_iterator | faces_end () const |
Past-the-end iterator. | |
Vertex_circulator | adjacent_vertices (Vertex_handle v) const |
Starts at an arbitrary vertex incident to v . | |
Vertex_circulator | adjacent_vertices (Vertex_handle v, Face_handle f) const |
Starts at the first vertex of f adjacent to v in counterclockwise order around v . | |
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 . | |
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 . | |
template<class OutputIterator > | |
OutputIterator | incident_faces (Vertex_handle v, OutputIterator faces) const |
Copies the faces incident to v into the output iterator faces . | |
template<class OutputIterator > | |
OutputIterator | incident_edges (Vertex_handle v, OutputIterator edges) const |
Copies the edges incident to v into the output iterator edges . | |
template<class OutputIterator > | |
OutputIterator | adjacent_vertices (Vertex_handle v, OutputIterator vertices) const |
Copies the vertices incident to v into the output iterator vertices . | |
size_type | degree (Vertex_handle v) const |
Returns the degree of v , i.e., the number of edges incident to v . | |
int | mirror_index (Face_handle f, int i) const |
Returns the index of f in its i -th neighbor. | |
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 . | |
Edge | mirror_edge (Edge e) const |
Returns the same edge seen from the other adjacent face. | |
Hyperbolic_translation | neighbor_translation (const Face_handle fh, int i) const |
Returns the hyperbolic translation for which the i -th neighbor of fh is adjacent to (the canonical representative of) fh in the hyperbolic plane. | |
bool | is_valid (bool verbose=false) const |
Checks the combinatorial validity of the triangulation, and also the validity of its geometric embedding. | |
bool | is_valid (Face_handle f, bool verbose=false) const |
Checks the combinatorial validity of face f . | |
CGAL::Periodic_4_hyperbolic_Delaunay_triangulation_2< GT, TDS >::Periodic_4_hyperbolic_Delaunay_triangulation_2 | ( | InputIterator | first, |
InputIterator | last, | ||
const Geom_traits & | gt = Geom_traits() |
||
) |
Point range constructor.
Initializes a triangulation and then inserts the points in the iterator range [first, last)
.
Voronoi_point CGAL::Periodic_4_hyperbolic_Delaunay_triangulation_2< GT, TDS >::dual | ( | Face_handle | f, |
Hyperbolic_translation | nbtr = Hyperbolic_translation() |
||
) | const |
Returns the dual of face f
, i.e., the hyperbolic center of the disk defined by the vertices of f
.
If the optional parameter nbtr
is given, then this function returns the dual of the periodic triangle (f, nbtr)
.
void CGAL::Periodic_4_hyperbolic_Delaunay_triangulation_2< GT, TDS >::find_conflicts | ( | const Point & | p, |
OutputFaceIterator | it, | ||
Face_handle | start = Face_handle() , |
||
Hyperbolic_translation | ltr = Hyperbolic_translation() |
||
) | const |
Computes the conflict zone induced by p
.
The output iterator it
contains all faces in conflict with p
. The optional parameters start
and ltr
, if given, must be such that start
translated by ltr
is in conflict with p
.
Vertex_handle CGAL::Periodic_4_hyperbolic_Delaunay_triangulation_2< GT, TDS >::insert | ( | const Point & | p, |
Face_handle | start = Face_handle() |
||
) |
Inserts the point p
in the triangulation.
The face start
, if given, is used as a starting place for the location of the point. Note that this function does not remove unnecessary dummy points. The removal, if desired, should be done by manually calling the function try_to_remove_dummy_vertices()
.
std::ptrdiff_t CGAL::Periodic_4_hyperbolic_Delaunay_triangulation_2< GT, TDS >::insert | ( | InputIterator | first, |
InputIterator | last, | ||
const bool | flag_try_to_remove_dummy_vertices = true |
||
) |
Inserts all points in the input iterator into the triangulation.
Note that this function by default tries to remove unnecessary dummy points at the end of the insertion process. This behavior is controlled by the optional boolean parameter flag_try_to_remove_dummy_vertices
; if automatic removal of the dummy points is not desired, set the flag to false
.
bool CGAL::Periodic_4_hyperbolic_Delaunay_triangulation_2< GT, TDS >::is_valid | ( | bool | verbose = false | ) | const |
Checks the combinatorial validity of the triangulation, and also the validity of its geometric embedding.
In more detail, this function verifies that:
valid
(see the function is_valid()) bool CGAL::Periodic_4_hyperbolic_Delaunay_triangulation_2< GT, TDS >::remove | ( | Vertex_handle | v | ) |
Removes the vertex v
from the triangulation.
Note that v
is removed from the triangulation only if the resulting triangulation remains a simplicial complex. The function returns true
if the vertex v
is removed from the triangulation, otherwise it returns false
.
void CGAL::Periodic_4_hyperbolic_Delaunay_triangulation_2< GT, TDS >::remove | ( | VertexRemoveIterator | first, |
VertexRemoveIterator | last | ||
) |
Removes the vertices in the iterator range [firs, last)
from the triangulation.
[first, last)
are vertices of the triangulation. int CGAL::Periodic_4_hyperbolic_Delaunay_triangulation_2< GT, TDS >::try_to_remove_dummy_vertices | ( | ) |
Tries to remove the dummy points from the triangulation one by one.
Returns the number of dummy points that could not be removed and are still present in the triangulation.