CGAL 5.5.3 - 2D Triangulations
CGAL::Regular_triangulation_2< Traits, Tds > Class Template Reference

#include <CGAL/Regular_triangulation_2.h>

Inherits from

CGAL::Triangulation_2< Traits, Tds >.

Definition

The class Regular_triangulation_2 is designed to maintain the regular triangulation of a set of weighted points.

Let \( { PW} = \{(p_i, w_i), i = 1, \ldots , n \}\) be a set of weighted points where each \( p_i\) is a point and each \( w_i\) is a scalar called the weight of point \( p_i\). Alternatively, each weighted point \( (p_i, w_i)\) can be regarded as a two dimensional sphere with center \( p_i\) and radius \( r_i=\sqrt{w_i}\).

The power diagram of the set \( { PW}\) is a planar partition such that each cell corresponds to sphere \( (p_i, w_i)\) of \( { PW}\) and is the locus of points \( p\) whose power with respect to \( (p_i, w_i)\) is less than its power with respect to any other sphere \( (p_j, w_j)\) in \( { PW}\). The dual of this diagram is a triangulation whose domain covers the convex hull of the set \( { P}= \{ p_i, i = 1, \ldots , n \}\) of center points and whose vertices are a subset of \( { P}\). Such a triangulation is called a regular triangulation. The three points \( p_i, p_j\) and \( p_k\) of \( { P}\) form a triangle in the regular triangulation of \( { PW}\) iff there is a point \( p\) of the plane whose powers with respect to \( (p_i, w_i)\), \( (p_j, w_j)\) and \( (p_k, w_k)\) are equal and less than the power of \( p\) with respect to any other sphere in \( { PW}\).

Let us defined the power product of two weighted points \( (p_i, w_i)\) and \( (p_j, w_j)\) as:

\[ \Pi(p_i, w_i,p_j, w_j) = p_ip_j ^2 - w_i - w_j . \]

\( \Pi(p_i, w_i,p_j, 0)\) is simply the power of point \( p_j\) with respect to the sphere \( (p_i, w_i)\), and two weighted points are said to be orthogonal if their power product is null. The power circle of three weighted points \( (p_i, w_i)\), \( (p_j, w_j)\) and \( (p_k, w_k)\) is defined as the unique circle \( (\pi, \omega)\) orthogonal to \( (p_i, w_i)\), \( (p_j, w_j)\) and \( (p_k, w_k)\).

The regular triangulation of the sets \( { PW}\) satisfies the following regular property (which just reduces to the Delaunay property when all the weights are null): a triangle \( p_ip_jp_k\) of the regular triangulation of \( { PW}\) is such that the power product of any weighted point \( (p_l, w_l)\) of \( { PW}\) with the power circle of \( (p_i, w_i)\), \( (p_j, w_j)\) is \( (p_k, w_k)\) is positive or null. We call power test of the weighted point \( (p_l, w_l)\) with respect to the face \( p_ip_jp_k\), the predicates testing the sign of the power product of \( (p_l, w_l)\) with respect to the power circle of \( (p_i, w_i)\), \( (p_j, w_j)\) is \( (p_k, w_k)\). This power product is given by the following determinant

\[ \left| \begin{array}{cccc} 1 & x_i & y_i & x_i ^2 + y_i ^2 - w_i \\ 1 & x_j & y_j & x_j ^2 + y_j ^2 - w_j \\ 1 & x_k & y_k & x_k ^2 + y_k ^2 - w_k \\ 1 & x_l & y_l & x_l ^2 + y_l ^2 - w_l \end{array} \right| \]

A pair of neighboring faces \( p_ip_jp_k\) and \( p_ip_jp_l\) is said to be locally regular (with respect to the weights in \( { PW}\)) if the power test of \( (p_l,w_l)\) with respect to \( p_ip_jp_k\) is positive. A classical result of computational geometry establishes that a triangulation of the convex hull of \( { P}\) such that any pair of neighboring faces is regular with respect to \( { PW}\), is a regular triangulation of \( { PW}\).

Alternatively, the regular triangulation of the weighted points set \( { PW}\) can be obtained as the projection on the two dimensional plane of the convex hull of the set of three dimensional points \( { P'}= \{ (p_i,p_i ^2 - w_i ), i = 1, \ldots , n \}\).

The vertices of the regular triangulation of a set of weighted points \( { PW}\) form only a subset of the set of center points of \( { PW}\). Therefore the insertion of a weighted point in a regular triangulation does not necessarily imply the creation of a new vertex. If the new inserted point does not appear as a vertex in the regular triangulation, it is said to be hidden.

Hidden points are stored in special vertices called hidden vertices. A hidden point is considered as hidden by the facet of the triangulation where its point component is located : in fact, the hidden point can appear as vertex of the triangulation only if this facet is removed. Each face of a regular triangulation stores the list of hidden vertices whose points are located in the facet. When a facet is removed, points hidden by this facet are reinserted in the triangulation.

Template Parameters
Traitsis the geometric traits parameter and must be a model of the concept RegularTriangulationTraits_2. The concept RegularTriangulationTraits_2 refines the concept TriangulationTraits_2 by adding the type Weighted_point_2 to describe weighted points and the type Power_side_of_oriented_power_circle_2 to perform power tests on weighted points.
Tdsmust be a model of TriangulationDataStructure_2. The face base of a regular triangulation has to be a model of the concept RegularTriangulationFaceBase_2. while the vertex base class has to be a model of RegularTriangulationVertexBase_2. CGAL provides a default instantiation for the Tds parameter by the class Triangulation_data_structure_2 < Regular_triangulation_vertex_base_2<Traits>, Regular_triangulation_face_base_2<Traits> >.
See also
CGAL::Triangulation_2<Traits,Tds>
TriangulationDataStructure_2
RegularTriangulationTraits_2
RegularTriangulationFaceBase_2
RegularTriangulationVertexBase_2
CGAL::Regular_triangulation_face_base_2<Traits>
CGAL::Regular_triangulation_vertex_base_2<Traits>
Examples:
Triangulation_2/info_insert_with_pair_iterator_regular_2.cpp, and Triangulation_2/regular.cpp.

Types

typedef Traits::Distance Distance
 
typedef Traits::Line Line
 
typedef Traits::Ray Ray
 
typedef Traits::Point_2 Bare_point
 
typedef Traits::Weighted_point_2 Weighted_point
 
typedef unspecified_type All_vertices_iterator
 An iterator that allows to enumerate the vertices that are not hidden.
 
typedef unspecified_type Finite_vertices_iterator
 An iterator that allows to enumerate the finite vertices that are not hidden.
 
typedef unspecified_type Hidden_vertices_iterator
 An iterator that allows to enumerate the hidden vertices.
 

Creation

 Regular_triangulation_2 (const Traits &gt=Traits())
 Introduces an empty regular triangulation.
 
 Regular_triangulation_2 (const Regular_triangulation_2 &rt)
 Copy constructor.
 
template<class InputIterator >
Regular_triangulation_2< Traits, Tds > Regular_triangulation_2 (InputIterator first, InputIterator last, Traits gt=Traits())
 Equivalent to constructing an empty triangulation with the optional traits class argument and calling insert(first,last).
 

Insertion and Removal

Vertex_handle insert (const Weighted_point &p, Face_handle f=Face_handle())
 inserts weighted point p in the regular triangulation. More...
 
Vertex_handle insert (const Weighted_point &p, Locate_type lt, Face_handle loc, int li)
 insert a weighted point p whose bare-point is assumed to be located in lt,loc,li. More...
 
Vertex_handle push_back (const Point &p)
 Equivalent to insert(p).
 
template<class InputIterator >
std::ptrdiff_t insert (InputIterator first, InputIterator last)
 inserts the weighted points in the range [first,last). More...
 
template<class WeightedPointWithInfoInputIterator >
std::ptrdiff_t insert (WeightedPointWithInfoInputIterator first, WeightedPointWithInfoInputIterator last)
 inserts the weighted points in the range [first,last). More...
 
void remove (Vertex_handle v)
 removes the vertex from the triangulation.
 

Queries

template<class OutputItFaces , class OutputItBoundaryEdges , class OutputItHiddenVertices >
CGAL::Triple< OutputItFaces, OutputItBoundaryEdges, OutputItHiddenVertices > get_conflicts_and_boundary_and_hidden_vertices (const Weighted_point &p, OutputItFaces fit, OutputItBoundaryEdges eit, OutputItHiddenVertices vit, Face_handle start) const
 outputs the faces, boundary edges, and hidden vertices of the conflict zone of point p to output iterators. More...
 
template<class OutputItFaces , class OutputItBoundaryEdges >
std::pair< OutputItFaces, OutputItBoundaryEdges > get_conflicts_and_boundary (const Weighted_point &p, OutputItFaces fit, OutputItBoundaryEdges eit, Face_handle start) const
 outputs the faces and boundary edges of the conflict zone of point p to output iterators. More...
 
template<class OutputItFaces , class OutputItHiddenVertices >
std::pair< OutputItFaces, OutputItHiddenVertices > get_conflicts_and_hidden_vertices (const Weighted_point &p, OutputItFaces fit, OutputItHiddenVertices vit, Face_handle start) const
 outputs the faces and hidden vertices of the conflict zone of point p to output iterators. More...
 
template<class OutputItBoundaryEdges , class OutputItHiddenVertices >
std::pair< OutputItBoundaryEdges, OutputItHiddenVertices > get_boundary_of_conflicts_and_hidden_vertices (const Weighted_point &p, OutputItBoundaryEdges eit, OutputItHiddenVertices vit, Face_handle start) const
 outputs the boundary edges and hidden vertices of the conflict zone of point p to output iterators. More...
 
template<class OutputItFaces >
OutputItFaces get_conflicts (const Point &p, OutputItFaces fit, Face_handle start) const
 outputs the faces of the conflict zone of point p to output iterators. More...
 
template<class OutputItBoundaryEdges >
OutputItBoundaryEdges get_boundary_of_conflicts (const Point &p, OutputItBoundaryEdges eit, Face_handle start) const
 outputs the boundary edges of the conflict zone of p in counterclockwise order where each edge is described through the incident face which is not in conflict with p. More...
 
template<class OutputItHiddenVertices >
OutputItHiddenVertices get_hidden_vertices (const Point &p, OutputItHiddenVertices vit, Face_handle start) const
 outputs the hidden vertices of the conflict zone of p into an output iterator. More...
 
Vertex_handle nearest_power_vertex (Bare_point p) const
 Returns the vertex of the triangulation which is nearest to p with respect to the power distance. More...
 

Access Functions

int number_of_vertices () const
 returns the number of finite vertices that are not hidden.
 
int number_of_hidden_vertices () const
 returns the number of hidden vertices.
 
Hidden_vertices_iterator hidden_vertices_begin () const
 starts at an arbitrary hidden vertex.
 
Hidden_vertices_iterator hidden_vertices_end () const
 past the end iterator for the sequence of hidden vertices.
 
Finite_vertices_iterator finite_vertices_begin () const
 starts at an arbitrary unhidden finite vertex
 
Finite_vertices_iterator finite_vertices_end () const
 Past-the-end iterator.
 
All_vertices_iterator all_vertices_end () const
 starts at an arbitrary unhidden vertex.
 
All_vertices_iterator all_vertices_begin () const
 past the end iterator.
 

Dual Power Diagram

The following member functions provide the elements of the dual power diagram.

Point weighted_circumcenter (const Face_handle &f) const
 returns the center of the circle orthogonal to the three weighted points corresponding to the vertices of face f. More...
 
Point dual (const Face_handle &f) const
 same as weighted_circumcenter.
 
Object dual (const Edge &e) const
 If both incident faces are finite, returns a segment whose endpoints are the duals of each incident face. More...
 
Object dual (const Edge_circulator &ec) const
 Idem.
 
Object dual (const Edge_iterator &ei) const
 Idem.
 
template<class Stream >
Stream & draw_dual (Stream &ps)
 output the dual power diagram to stream ps.
 

Predicates

Oriented_side power_test (Face_handle f, const Weighted_point &p) const
 Returns the power test of p with respect to the power circle associated with f.
 

Miscellaneous

bool is_valid (bool verbose=false, int level=0) const
 Tests the validity of the triangulation as a Triangulation_2 and additionally test the regularity of the triangulation. More...
 

Additional Inherited Members

- Public Types inherited from CGAL::Triangulation_2< Traits, Tds >
enum  Locate_type {
  VERTEX =0, EDGE, FACE, OUTSIDE_CONVEX_HULL,
  OUTSIDE_AFFINE_HULL
}
 specifies which case occurs when locating a point in the triangulation. More...
 
typedef Tds::Vertex_handle Vertex_handle
 handle to a vertex.
 
typedef Tds::Face_handle Face_handle
 handle to a face.
 
typedef Tds::Face_iterator All_faces_iterator
 iterator over all faces.
 
typedef Tds::Edge_iterator All_edges_iterator
 iterator over all edges.
 
typedef Tds::Vertex_iterator All_vertices_iterator
 iterator over all vertices.
 
typedef unspecified_type Finite_faces_iterator
 iterator over finite faces.
 
typedef unspecified_type Finite_edges_iterator
 iterator over finite edges.
 
typedef unspecified_type Finite_vertices_iterator
 iterator over finite vertices.
 
typedef unspecified_type Point_iterator
 iterator over the points corresponding to the finite vertices of the triangulation.
 
typedef Iterator_range< unspecified_typeAll_face_handles
 range type for iterating over all faces (including infinite faces), with a nested type iterator that has as value type Face_handle.
 
typedef Iterator_range< All_edges_iteratorAll_edges
 range type for iterating over all edges (including infinite ones).
 
typedef Iterator_range< unspecified_typeAll_vertex_handles
 range type for iterating over all vertices (including the infinite vertex), with a nested type iterator that has as value type Vertex_handle.
 
typedef Iterator_range< unspecified_typeFinite_face_handles
 range type for iterating over finite faces, with a nested type iterator that has as value type Face_handle.
 
typedef Iterator_range< Finite_edges_iteratorFinite_edges
 range type for iterating over finite edges.
 
typedef Iterator_range< unspecified_typeFinite_vertex_handles
 range type for iterating over finite vertices, with a nested type iterator that has as value type Vertex_handle.
 
typedef Iterator_range< Point_iteratorPoints
 range type for iterating over the points of the finite vertices.
 
typedef unspecified_type Line_face_circulator
 circulator over all faces intersected by a line.
 
typedef unspecified_type Face_circulator
 circulator over all faces incident to a given vertex.
 
typedef unspecified_type Edge_circulator
 circulator over all edges incident to a given vertex.
 
typedef unspecified_type Vertex_circulator
 circulator over all vertices incident to a given vertex.
 
typedef Traits Geom_traits
 the traits class.
 
typedef Tds Triangulation_data_structure
 the triangulation data structure type.
 
typedef Traits::Point_2 Point
 the point type.
 
typedef Traits::Segment_2 Segment
 the segment type.
 
typedef Traits::Triangle_2 Triangle
 the triangle type.
 
typedef Tds::Vertex Vertex
 the vertex type.
 
typedef Tds::Face Face
 the face type.
 
typedef Tds::Edge Edge
 the edge type.
 
typedef Tds::size_type size_type
 Size type (an unsigned integral type).
 
typedef Tds::difference_type difference_type
 Difference type (a signed integral type).
 
- Public Member Functions inherited from CGAL::Triangulation_2< Traits, Tds >
 Triangulation_2 (const Traits &gt=Traits())
 Introduces an empty triangulation.
 
 Triangulation_2 (const Triangulation_2 &tr)
 Copy constructor. More...
 
template<class InputIterator >
 Triangulation_2 (InputIterator first, InputIterator last, const Traits &gt=Traits())
 Equivalent to constructing an empty triangulation with the optional traits class argument and calling insert(first,last).
 
Triangulation_2 operator= (const Triangulation_2< Traits, Tds > &tr)
 Assignment. More...
 
void swap (Triangulation_2 &tr)
 The triangulations tr and *this are swapped. More...
 
void clear ()
 Deletes all faces and finite vertices resulting in an empty triangulation.
 
int dimension () const
 Returns the dimension of the convex hull.
 
size_type number_of_vertices () const
 Returns the number of finite vertices.
 
size_type number_of_faces () const
 Returns the number of finite faces.
 
Face_handle infinite_face () const
 a face incident to the infinite vertex.
 
Vertex_handle infinite_vertex () const
 the infinite vertex.
 
Vertex_handle finite_vertex () const
 a vertex distinct from the infinite vertex.
 
const Geom_traitsgeom_traits () const
 Returns a const reference to the triangulation traits object.
 
const TriangulationDataStructure_2tds () const
 Returns a const reference to the triangulation data structure.
 
TriangulationDataStructure_2tds ()
 Returns a reference to the triangulation data structure.
 
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.
 
bool is_infinite (All_edges_iterator ei) const
 true iff edge *ei is infinite.
 
bool is_edge (Vertex_handle va, Vertex_handle vb)
 true if there is an edge having va and vb as vertices.
 
bool is_edge (Vertex_handle va, Vertex_handle vb, Face_handle &fr, int &i)
 as above. More...
 
bool includes_edge (Vertex_handle va, Vertex_handle vb, Vertex_handle &vbr, Face_handle &fr, int &i)
 true if the line segment from va to vb includes an edge e incident to va. More...
 
bool is_face (Vertex_handle v1, Vertex_handle v2, Vertex_handle v3)
 true if there is a face having v1, v2 and v3 as vertices.
 
bool is_face (Vertex_handle v1, Vertex_handle v2, Vertex_handle v3, Face_handle &fr)
 as above. More...
 
Face_handle locate (const Point &query, Face_handle f=Face_handle()) const
 If the point query lies inside the convex hull of the points, a face that contains the query in its interior or on its boundary is returned. More...
 
Face_handle inexact_locate (const Point &query, Face_handle start=Face_handle()) const
 Same as locate() but uses inexact predicates. More...
 
Face_handle locate (const Point &query, Locate_type &lt, int &li, Face_handle h=Face_handle()) const
 Same as above. More...
 
Oriented_side oriented_side (Face_handle f, const Point &p) const
 Returns on which side of the oriented boundary of f lies the point p. More...
 
Oriented_side side_of_oriented_circle (Face_handle f, const Point &p)
 Returns on which side of the circumcircle of face f lies the point p. More...
 
void flip (Face_handle f, int i)
 Exchanges the edge incident to f and f->neighbor(i) with the other diagonal of the quadrilateral formed by f and f->neighbor(i). More...
 
Vertex_handle insert (const Point &p, Face_handle f=Face_handle())
 Inserts point p in the triangulation and returns the corresponding vertex. More...
 
Vertex_handle insert (const Point &p, Locate_type lt, Face_handle loc, int li)
 Same as above except that the location of the point p to be inserted is assumed to be given by (lt,loc,i) (see the description of the locate method above.)
 
Vertex_handle push_back (const Point &p)
 Equivalent to insert(p).
 
template<class PointInputIterator >
std::ptrdiff_t insert (PointInputIterator first, PointInputIterator last)
 Inserts the points in the range [first,last) in the given order, and returns the number of inserted points. More...
 
template<class PointWithInfoInputIterator >
std::ptrdiff_t insert (PointWithInfoInputIterator first, PointWithInfoInputIterator last)
 inserts the points in the iterator range [first,last) in the given order, and returns the number of inserted points. More...
 
void remove (Vertex_handle v)
 Removes the vertex from the triangulation. More...
 
Vertex_handle move_if_no_collision (Vertex_handle v, const Point &p)
 If there is not already another vertex placed on p, the triangulation is modified such that the new position of vertex v is p, and v is returned. More...
 
Vertex_handle move (Vertex_handle v, const Point &p)
 If there is no collision during the move, this function is the same as move_if_no_collision . More...
 
Vertex_handle insert_first (const Point &p)
 Inserts the first finite vertex .
 
Vertex_handle insert_second (const Point &p)
 Inserts the second finite vertex .
 
Vertex_handle insert_in_face (const Point &p, Face_handle f)
 Inserts vertex v in face f. More...
 
Vertex_handle insert_in_edge (const Point &p, Face_handle f, int i)
 Inserts vertex v in edge i of f. More...
 
Vertex_handle insert_outside_convex_hull (const Point &p, Face_handle f)
 Inserts a point which is outside the convex hull but in the affine hull. More...
 
Vertex_handle insert_outside_affine_hull (const Point &p)
 Inserts a point which is outside the affine hull.
 
void remove_degree_3 (Vertex_handle v)
 Removes a vertex of degree three. More...
 
void remove_second (Vertex_handle v)
 Removes the before last finite vertex.
 
void remove_first (Vertex_handle v)
 Removes the last finite vertex.
 
template<class EdgeIt >
Vertex_handle star_hole (Point p, EdgeIt edge_begin, EdgeIt edge_end)
 creates a new vertex v and use it to star the hole whose boundary is described by the sequence of edges [edge_begin, edge_end). More...
 
template<class EdgeIt , class FaceIt >
Vertex_handle star_hole (Point p, EdgeIt edge_begin, EdgeIt edge_end, FaceIt face_begin, FaceIt face_end)
 same as above, except that the algorithm first recycles faces in the sequence [face_begin, face_end) and create new ones only when the sequence is exhausted. More...
 
Finite_vertices_iterator finite_vertices_begin () const
 Starts at an arbitrary finite vertex.
 
Finite_vertices_iterator finite_vertices_end () const
 Past-the-end iterator.
 
Finite_edges_iterator finite_edges_begin () const
 Starts at an arbitrary finite edge.
 
Finite_edges_iterator finite_edges_end () const
 Past-the-end iterator.
 
Finite_faces_iterator finite_faces_begin () const
 Starts at an arbitrary finite face.
 
Finite_faces_iterator finite_faces_end () const
 Past-the-end iterator.
 
Point_iterator points_begin () const
 
Point_iterator points_end () const
 Past-the-end iterator.
 
Finite_vertex_handles finite_vertex_handles () const
 returns a range of iterators over finite vertices. More...
 
Finite_edges finite_edges () const
 returns a range of iterators over finite edges.
 
Finite_face_handles finite_face_handles () const
 returns a range of iterators over finite faces. More...
 
Points points () const
 returns a range of iterators over the points of finite vertices.
 
All_vertices_iterator all_vertices_begin () const
 Starts at an arbitrary vertex.
 
All_vertices_iterator all_vertices_end () const
 Past-the-end iterator.
 
All_edges_iterator all_edges_begin () const
 Starts at an arbitrary edge.
 
All_edges_iterator all_edges_end () const
 Past-the-end iterator.
 
All_faces_iterator all_faces_begin () const
 Starts at an arbitrary face.
 
All_faces_iterator all_faces_end () const
 Past-the-end iterator.
 
All_vertex_handles all_vertex_handles () const
 returns a range of iterators over all vertices. More...
 
All_edges all_edges () const
 returns a range of iterators over all edges.
 
All_face_handles all_face_handles () const
 returns a range of iterators over all faces. More...
 
Line_face_circulator line_walk (const Point &p, const Point &q, Face_handle f=Face_handle()) const
 This function returns a circulator that allows to visit the faces intersected by the line pq. More...
 
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. More...
 
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. More...
 
Vertex_circulator incident_vertices (Vertex_handle v) const
 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...
 
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. More...
 
int mirror_index (Face_handle f, int i) const
 returns the index of f in its \( i^{th}\) neighbor. More...
 
Edge mirror_edge (Edge e) const
 returns the same edge seen from the other adjacent face. More...
 
int ccw (int i) const
 Returns \( i+1\) modulo 3. More...
 
int cw (int i) const
 Returns \( i+2\) modulo 3. More...
 
Triangle triangle (Face_handle f) const
 Returns the triangle formed by the three vertices of f. More...
 
Segment segment (Face_handle f, int i) const
 Returns the line segment formed by the vertices ccw(i) and cw(i) of face f. More...
 
Segment segment (const Edge &e) const
 Returns the line segment corresponding to edge e. More...
 
Segment segment (const Edge_circulator &ec) const
 Returns the line segment corresponding to edge *ec. More...
 
Segment segment (const Edge_iterator &ei) const
 Returns the line segment corresponding to edge *ei. More...
 
Point circumcenter (Face_handle f) const
 Compute the circumcenter of the face pointed to by f. More...
 
void set_infinite_vertex (const Vertex_handle &v)
 This is an advanced function. More...
 
bool is_valid (bool verbose=false, int level=0) const
 Checks the combinatorial validity of the triangulation and also the validity of its geometric embedding. More...
 
- Public Member Functions inherited from CGAL::Triangulation_cw_ccw_2
 Triangulation_cw_ccw_2 ()
 default constructor.
 
int ccw (const int i) const
 returns the index of the neighbor or vertex that is next to the neighbor or vertex with index i in counterclockwise order around a face.
 
int cw (const int i) const
 returns the index of the neighbor or vertex that is next to the neighbor or vertex with index i in counterclockwise order around a face.
 

Member Function Documentation

◆ dual()

template<typename Traits , typename Tds >
Object CGAL::Regular_triangulation_2< Traits, Tds >::dual ( const Edge e) const

If both incident faces are finite, returns a segment whose endpoints are the duals of each incident face.

If only one incident face is finite, returns a ray whose endpoint is the dual of the finite incident face and supported by the line which is the bisector of the edge's endpoints. If both incident faces are infinite, returns the line which is the bisector of the edge's endpoints otherwise.

◆ get_boundary_of_conflicts()

template<typename Traits , typename Tds >
template<class OutputItBoundaryEdges >
OutputItBoundaryEdges CGAL::Regular_triangulation_2< Traits, Tds >::get_boundary_of_conflicts ( const Point p,
OutputItBoundaryEdges  eit,
Face_handle  start 
) const

outputs the boundary edges of the conflict zone of p in counterclockwise order where each edge is described through the incident face which is not in conflict with p.

The function returns the resulting output iterator.

◆ get_boundary_of_conflicts_and_hidden_vertices()

template<typename Traits , typename Tds >
template<class OutputItBoundaryEdges , class OutputItHiddenVertices >
std::pair<OutputItBoundaryEdges,OutputItHiddenVertices> CGAL::Regular_triangulation_2< Traits, Tds >::get_boundary_of_conflicts_and_hidden_vertices ( const Weighted_point p,
OutputItBoundaryEdges  eit,
OutputItHiddenVertices  vit,
Face_handle  start 
) const

outputs the boundary edges and hidden vertices of the conflict zone of point p to output iterators.

See get_conflicts_and_boundary_and_hidden_vertices() for details. The function returns in a std::pair the resulting output iterators.

◆ get_conflicts()

template<typename Traits , typename Tds >
template<class OutputItFaces >
OutputItFaces CGAL::Regular_triangulation_2< Traits, Tds >::get_conflicts ( const Point p,
OutputItFaces  fit,
Face_handle  start 
) const

outputs the faces of the conflict zone of point p to output iterators.

The function returns the resulting output iterator.

Precondition
dimension()==2.

◆ get_conflicts_and_boundary()

template<typename Traits , typename Tds >
template<class OutputItFaces , class OutputItBoundaryEdges >
std::pair<OutputItFaces,OutputItBoundaryEdges> CGAL::Regular_triangulation_2< Traits, Tds >::get_conflicts_and_boundary ( const Weighted_point p,
OutputItFaces  fit,
OutputItBoundaryEdges  eit,
Face_handle  start 
) const

outputs the faces and boundary edges of the conflict zone of point p to output iterators.

See get_conflicts_and_boundary_and_hidden_vertices() for details.

The function returns in a std::pair the resulting output iterators.

Precondition
dimension()==2.

◆ get_conflicts_and_boundary_and_hidden_vertices()

template<typename Traits , typename Tds >
template<class OutputItFaces , class OutputItBoundaryEdges , class OutputItHiddenVertices >
CGAL::Triple<OutputItFaces,OutputItBoundaryEdges,OutputItHiddenVertices> CGAL::Regular_triangulation_2< Traits, Tds >::get_conflicts_and_boundary_and_hidden_vertices ( const Weighted_point p,
OutputItFaces  fit,
OutputItBoundaryEdges  eit,
OutputItHiddenVertices  vit,
Face_handle  start 
) const

outputs the faces, boundary edges, and hidden vertices of the conflict zone of point p to output iterators.

Template Parameters
OutputItFacesis an output iterator with Face_handle as value type.
OutputItBoundaryEdgesis an output iterator with Edge as value type.
OutputItHiddenVerticesis an output iterator with Vertex_handle as value type.

This member function outputs in the container pointed to by fit the faces which are in conflict with point p, i.e., the faces whose power circles have negative power wrt. p. It outputs in the container pointed to by eit the boundary of the zone in conflict with p. It inserts the vertices that would be hidden by p into the container pointed to by vit. The boundary edges of the conflict zone are output in counter-clockwise order and each edge is described through its incident face which is not in conflict with p. The function returns in a CGAL::Triple the resulting output iterators.

Precondition
dimension()==2.

◆ get_conflicts_and_hidden_vertices()

template<typename Traits , typename Tds >
template<class OutputItFaces , class OutputItHiddenVertices >
std::pair<OutputItFaces,OutputItHiddenVertices> CGAL::Regular_triangulation_2< Traits, Tds >::get_conflicts_and_hidden_vertices ( const Weighted_point p,
OutputItFaces  fit,
OutputItHiddenVertices  vit,
Face_handle  start 
) const

outputs the faces and hidden vertices of the conflict zone of point p to output iterators.

See get_conflicts_and_boundary_and_hidden_vertices() for details. The function returns in a std::pair the resulting output iterators.

Precondition
dimension()==2.

◆ get_hidden_vertices()

template<typename Traits , typename Tds >
template<class OutputItHiddenVertices >
OutputItHiddenVertices CGAL::Regular_triangulation_2< Traits, Tds >::get_hidden_vertices ( const Point p,
OutputItHiddenVertices  vit,
Face_handle  start 
) const

outputs the hidden vertices of the conflict zone of p into an output iterator.

The function returns the resulting output iterator.

◆ insert() [1/4]

template<typename Traits , typename Tds >
Vertex_handle CGAL::Regular_triangulation_2< Traits, Tds >::insert ( const Weighted_point p,
Face_handle  f = Face_handle() 
)

inserts weighted point p in the regular triangulation.

If the point p does not appear as a vertex of the triangulation, the returned vertex is a hidden vertex. If given the parameter f is used as an hint for the place to start the location process of point p.

◆ insert() [2/4]

template<typename Traits , typename Tds >
Vertex_handle CGAL::Regular_triangulation_2< Traits, Tds >::insert ( const Weighted_point p,
Locate_type  lt,
Face_handle  loc,
int  li 
)

insert a weighted point p whose bare-point is assumed to be located in lt,loc,li.

See the description of member function Triangulation_2::locate().

◆ insert() [3/4]

template<typename Traits , typename Tds >
template<class InputIterator >
std::ptrdiff_t CGAL::Regular_triangulation_2< Traits, Tds >::insert ( InputIterator  first,
InputIterator  last 
)

inserts the weighted points in the range [first,last).

It returns the difference of the number of vertices between after and before the insertions (it may be negative due to hidden points). Note that this function is not guaranteed to insert the weighted 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 Weighted_point .

◆ insert() [4/4]

template<typename Traits , typename Tds >
template<class WeightedPointWithInfoInputIterator >
std::ptrdiff_t CGAL::Regular_triangulation_2< Traits, Tds >::insert ( WeightedPointWithInfoInputIterator  first,
WeightedPointWithInfoInputIterator  last 
)

inserts the weighted points in the range [first,last).

It returns the difference of the number of vertices between after and before the insertions (it may be negative due to hidden points). Note that this function is not guaranteed to insert the weighted points following the order of WeightedPointWithInfoInputIterator, as spatial_sort is used to improve efficiency. Given a pair (p,i), the vertex v storing p also stores i, that is v.point() == p and v.info() == i. If several pairs have the same point, only one vertex is created, one of the objects of type Vertex::Info will be stored in the vertex.

Precondition
Vertex must be model of the concept TriangulationVertexBaseWithInfo_2.
Template Parameters
WeightedPointWithInfoInputIteratormust be an input iterator with value type std::pair<Weighted_point,Vertex::Info>.

◆ is_valid()

template<typename Traits , typename Tds >
bool CGAL::Regular_triangulation_2< Traits, Tds >::is_valid ( bool  verbose = false,
int  level = 0 
) const

Tests the validity of the triangulation as a Triangulation_2 and additionally test the regularity of the triangulation.

This method is useful to debug regular triangulation algorithms implemented by the user.

◆ nearest_power_vertex()

template<typename Traits , typename Tds >
Vertex_handle CGAL::Regular_triangulation_2< Traits, Tds >::nearest_power_vertex ( Bare_point  p) const

Returns the vertex of the triangulation which is nearest to p with respect to the power distance.

This means that the power of the query point p with respect to the weighted point in the nearest vertex is smaller than the power of p with respect to the weighted point in any other vertex. Ties are broken arbitrarily. The default constructed handle is returned if the triangulation is empty.

◆ weighted_circumcenter()

template<typename Traits , typename Tds >
Point CGAL::Regular_triangulation_2< Traits, Tds >::weighted_circumcenter ( const Face_handle f) const

returns the center of the circle orthogonal to the three weighted points corresponding to the vertices of face f.

Precondition
f is not infinite.