CGAL 5.5.2 - 2D Periodic Hyperbolic Triangulations
|
#include <CGAL/Periodic_4_hyperbolic_triangulation_2.h>
Inherited by CGAL::Periodic_4_hyperbolic_Delaunay_triangulation_2< GT, TDS >.
The class Periodic_4_hyperbolic_triangulation_2
offers base functionalities needed by Periodic_4_hyperbolic_Delaunay_triangulation_2
.
Note that this class does not support modification of the triangulation (insertion or removal).
The class expects two template parameters.
GT | Geometric traits parameter. Must be a model of the concept Periodic_4HyperbolicTriangulationTraits_2 . This parameter has no default value. |
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 CGAL::Periodic_4_hyperbolic_triangulation_vertex_base_2<GT>, |
Types | |
typedef GT | Geometric_traits |
typedef TDS | Triangulation_data_structure |
typedef Triangulation_data_structure::Vertex::Point | Point |
Represents a point in the hyperbolic plane. More... | |
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. | |
The following type represents hyperbolic translations specific for the Bolza surface. | |
typedef Geometric_traits::Hyperbolic_translation | Hyperbolic_translation |
The following types represent images of points, segments, and triangles under the action of hyperbolic translations. | |
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. | |
The following types give access to the elements of the triangulation. | |
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 |
The following iterator and circulator types are defined to give access over the vertices, edges, and faces of the triangulation. | |
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 |
Enums | |
The following enumeration type indicates where a point is located in the triangulation. | |
enum | Locate_type { VERTEX = 0, EDGE, FACE } |
Creation | |
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. | |
Assignment | |
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. | |
bool | operator== (const Periodic_4_hyperbolic_triangulation_2< GT, TDS > &tr1, const Periodic_4_hyperbolic_triangulation_2< GT, TDS > &tr2) |
Equality operator. More... | |
bool | operator!= (const Periodic_4_hyperbolic_triangulation_2< GT, TDS > &tr1, const Periodic_4_hyperbolic_triangulation_2< GT, TDS > &tr2) |
Inequality operator. More... | |
Access functions | |
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. | |
Geometric access functions | |
The functions below return periodic points, segments, and triangles. | |
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 . More... | |
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) . More... | |
Periodic_segment | periodic_segment (const Edge &e) const |
Returns the periodic segment formed by the endpoints of edge e . More... | |
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 . More... | |
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 . | |
The functions below return non-periodic points, segments and triangles. | |
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) . More... | |
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 . | |
Point location | |
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 . More... | |
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. More... | |
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 . More... | |
Face_handle | hyperbolic_locate (const Point &p, Locate_type <, int &li, Face_handle start=Face_handle()) const |
Same as above. More... | |
Predicate functions | |
Orientation | orientation (const Point &p1, const Point &p2, const Point &p3) const |
Returns the Orientation of the points p1, p2, p3 . More... | |
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 . More... | |
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) . | |
Queries | |
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. More... | |
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 . More... | |
bool | has_vertex (const Face_handle f, const Vertex_handle v, const int i) const |
Tests whether face f has v as a vertex. More... | |
Vertex, Edge, and Face iterators | |
Vertex_iterator | vertices_begin () const |
Starts at an arbitrary vertex. More... | |
Vertex_iterator | vertices_end () const |
Past-the-end iterator. | |
Edge_iterator | edges_begin () const |
Starts at an arbitrary edge. More... | |
Edge_iterator | edges_end () const |
Past-the-end iterator. | |
Face_iterator | faces_begin () const |
Starts at an arbitrary face. More... | |
Face_iterator | faces_end () const |
Past-the-end iterator. | |
Vertex, Edge and Face circulators | |
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 . | |
Traversal of the Vertices, Edges and Faces incident to a vertex | |
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 . | |
Traversal between adjacent faces | |
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. | |
Validity checking | |
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 . | |
typedef Triangulation_data_structure::Vertex::Point CGAL::Periodic_4_hyperbolic_triangulation_2< GT, TDS >::Point |
Represents a point in the hyperbolic plane.
Note that for Delaunay triangulations of the Bolza surface, all points must lie inside the original hyperbolic octagon.
enum CGAL::Periodic_4_hyperbolic_triangulation_2::Locate_type |
Hyperbolic_segment CGAL::Periodic_4_hyperbolic_triangulation_2< GT, TDS >::construct_hyperbolic_segment | ( | const Face_handle & | fh, |
int | idx | ||
) | const |
Constructs the hyperbolic segment formed by the endpoints of edge (fh, idx)
.
Edge_iterator CGAL::Periodic_4_hyperbolic_triangulation_2< GT, TDS >::edges_begin | ( | ) | const |
Starts at an arbitrary edge.
Iterates over all the edges in the triangulation.
Face_iterator CGAL::Periodic_4_hyperbolic_triangulation_2< GT, TDS >::faces_begin | ( | ) | const |
Starts at an arbitrary face.
Iterates over all the faces of the triangulation.
bool CGAL::Periodic_4_hyperbolic_triangulation_2< GT, TDS >::has_vertex | ( | const Face_handle | f, |
const Vertex_handle | v, | ||
const int | i | ||
) | const |
Tests whether face f
has v
as a vertex.
If the answer is true
, then i
contains the index of v
in f
.
Face_handle CGAL::Periodic_4_hyperbolic_triangulation_2< GT, TDS >::hyperbolic_locate | ( | const Point & | p, |
Face_handle | start = Face_handle() |
||
) | const |
Returns the canonical representative of the face that contains the query point p
.
The parameter start
, if provided, is used as a starting point for the location.
Face_handle CGAL::Periodic_4_hyperbolic_triangulation_2< GT, TDS >::hyperbolic_locate | ( | const Point & | p, |
Locate_type & | lt, | ||
int & | li, | ||
Face_handle | start = Face_handle() |
||
) | const |
Same as above.
The value of the variable lt
indicates whether p
has been located inside the canonical representative, on one of its sides, or on one of its vertices. If p
is located on a side or on a vertex, then li
contains the index of the corresponding edge or vertex in the returned face. The parameter start
, if provided, is used as a starting point for the location.
Face_handle CGAL::Periodic_4_hyperbolic_triangulation_2< GT, TDS >::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
.
The parameter start
, if provided, is used as a starting point for the location.
Face_handle CGAL::Periodic_4_hyperbolic_triangulation_2< GT, TDS >::hyperbolic_periodic_locate | ( | const Point & | p, |
Locate_type & | lt, | ||
int & | li, | ||
Hyperbolic_translation & | lo, | ||
const Face_handle | start = Face_handle() |
||
) | const |
Same as above.
The value of the variable lt
indicates whether p
has been located inside the hyperbolic triangle, on one of its sides, or on one of its vertices. If p
is located on a side or on a vertex, then li
contains the index of the corresponding edge or vertex in the returned face. The parameter start
, if provided, is used as a starting point for the location.
bool CGAL::Periodic_4_hyperbolic_triangulation_2< GT, TDS >::is_edge | ( | Vertex_handle | u, |
Vertex_handle | v, | ||
Face_handle & | fh, | ||
int & | i | ||
) | const |
Tests whether (u, v)
is an edge of the triangulation.
If the edge is found, then it is the i
-th edge of fh
.
bool CGAL::Periodic_4_hyperbolic_triangulation_2< GT, TDS >::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
.
If the face is found, then it is returned as fh
.
bool CGAL::Periodic_4_hyperbolic_triangulation_2< GT, TDS >::operator!= | ( | const Periodic_4_hyperbolic_triangulation_2< GT, TDS > & | tr1, |
const Periodic_4_hyperbolic_triangulation_2< GT, TDS > & | tr2 | ||
) |
Inequality operator.
bool CGAL::Periodic_4_hyperbolic_triangulation_2< GT, TDS >::operator== | ( | const Periodic_4_hyperbolic_triangulation_2< GT, TDS > & | tr1, |
const Periodic_4_hyperbolic_triangulation_2< GT, TDS > & | tr2 | ||
) |
Equality operator.
Orientation CGAL::Periodic_4_hyperbolic_triangulation_2< GT, TDS >::orientation | ( | const Point & | p1, |
const Point & | p2, | ||
const Point & | p3 | ||
) | const |
Returns the Orientation
of the points p1, p2, p3
.
Periodic_point CGAL::Periodic_4_hyperbolic_triangulation_2< GT, TDS >::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 CGAL::Periodic_4_hyperbolic_triangulation_2< GT, TDS >::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 CGAL::Periodic_4_hyperbolic_triangulation_2< GT, TDS >::periodic_segment | ( | const Edge & | e | ) | const |
Returns the periodic segment formed by the endpoints of edge e
.
Note that the translations in the resulting periodic segment are determined by e.first
.
Periodic_segment CGAL::Periodic_4_hyperbolic_triangulation_2< GT, TDS >::periodic_segment | ( | const Edge & | e, |
const Hyperbolic_translation & | tr | ||
) | const |
Returns the periodic segment formed by the endpoints of edge e
translated by tr
.
Note that the translations in the resulting segment are determined by the translations in e.first
, multiplied on the left by tr
.
Oriented_side CGAL::Periodic_4_hyperbolic_triangulation_2< GT, TDS >::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
.
Vertex_iterator CGAL::Periodic_4_hyperbolic_triangulation_2< GT, TDS >::vertices_begin | ( | ) | const |
Starts at an arbitrary vertex.
Iterates over all the vertices in the triangulation.