CGAL 5.2 - 2D Periodic Triangulations
CGAL::Periodic_2_Delaunay_triangulation_2< Traits, Tds > Class Template Reference

#include <CGAL/Periodic_2_Delaunay_triangulation_2.h>

Inherits from

CGAL::Periodic_2_triangulation_2< Traits, Tds >.

Definition

The class Periodic_2_Delaunay_triangulation_2 represents a Delaunay triangulation in two-dimensional periodic space.

Parameters

The class Periodic_2_Delaunay_triangulation_2 has two template parameters. The first one

Template Parameters
Traitsis the geometric traits, it is to be instantiated by a model of the concept Periodic_2DelaunayTriangulationTraits_2.

The second parameter is the triangulation data structure, it has to be instantiated by a model of the concept TriangulationDataStructure_2 with some additional functionality in faces. By default, the triangulation data structure is instantiated by CGAL::Triangulation_data_structure_2 < CGAL::Triangulation_vertex_base_2<Gt>, CGAL::Periodic_2_triangulation_face_base_2<Gt> > >.

Implementation

Insertion is implemented by inserting in the triangulation, then performing a sequence of Delaunay flips. The number of flips is \( O(d)\) if the new vertex is of degree \( d\) in the new triangulation. For points distributed uniformly at random, insertion takes time \( O(1)\) on average.

Removal calls the removal in the triangulation and then re-triangulates the hole in such a way that the Delaunay criterion is satisfied. Removal of a vertex of degree \( d\) takes time \( O(d^2)\). The expected degree \( d\) is \( O(1)\) for a random vertex in the triangulation.

After a point location step, the nearest neighbor is found in time \( O(n)\) in the worst case, but in expected time \( O(1)\) on average for vertices distributed uniformly at random and any query point.

See also
CGAL::Periodic_2_triangulation_2<Traits,Tds>
CGAL::Delaunay_triangulation_2<Traits,Tds>
TriangulationDataStructure_2
Periodic_2DelaunayTriangulationTraits_2
CGAL::Periodic_2_triangulation_hierarchy_2<Tr>
Examples:
Periodic_2_triangulation_2/draw_periodic_2_triangulation_2.cpp, Periodic_2_triangulation_2/p2t2_adding_handles.cpp, Periodic_2_triangulation_2/p2t2_colored_vertices.cpp, Periodic_2_triangulation_2/p2t2_covering.cpp, Periodic_2_triangulation_2/p2t2_find_conflicts.cpp, Periodic_2_triangulation_2/p2t2_geometric_access.cpp, Periodic_2_triangulation_2/p2t2_info_insert_with_pair_iterator_2.cpp, Periodic_2_triangulation_2/p2t2_info_insert_with_transform_iterator_2.cpp, Periodic_2_triangulation_2/p2t2_info_insert_with_zip_iterator_2.cpp, Periodic_2_triangulation_2/p2t2_large_point_set.cpp, and Periodic_2_triangulation_2/p2t2_simple_example.cpp.

Creation

 Periodic_2_Delaunay_triangulation_2 (const Iso_rectangle &domain=Iso_rectangle(0, 0, 1, 1), const Geom_traits &traits=Geom_traits())
 Creates an empty periodic Delaunay triangulation dt, with domain as original domain and possibly specifying a traits class traits. More...
 
 Periodic_2_Delaunay_triangulation_2 (const Periodic_2_Delaunay_triangulation_2 &dt1)
 Copy constructor.
 
template<class InputIterator >
 Periodic_2_Delaunay_triangulation_2 (InputIterator first, InputIterator last, const Iso_rectangle &domain=Iso_rectangle(0, 0, 1, 1), const Geom_traits &traits=Geom_traits())
 Equivalent to constructing an empty triangulation with the optional domain and traits class arguments and calling insert(first,last). More...
 

Insertion and removal

The following methods insert points in the triangulation ensuring the empty circle property of Delaunay triangulations.

The inserted points need to lie in the original domain (see Section The Flat Torus of the user manual). In the degenerate case when there are co-circular points, the Delaunay triangulation is known not to be uniquely defined. In this case, CGAL chooses a particular Delaunay triangulation using a symbolic perturbation scheme [3]. Note that insertion of a new point can cause a switch from computing in the 9-sheeted covering space to computing in the 1-sheeted covering space, which invalidates some Vertex_handles and Face_handles.

Vertex_handle insert (const Point &p, Face_handle start=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, int lj)
 Inserts point p in the triangulation and returns the corresponding vertex. More...
 
Vertex_handle push_back (const Point &p)
 Equivalent to insert(p).
 
template<class InputIterator >
std::ptrdiff_t insert (InputIterator first, InputIterator last, bool is_large_point_set=false)
 Inserts the points in the iterator range [first, last). More...
 
void remove (Vertex_handle v)
 Removes the vertex from the triangulation.
 

Point moving

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_point (Vertex_handle v, const Point &p)
 Moves the point stored in v to p, while preserving the Delaunay property. More...
 

Queries

Vertex_handle nearest_vertex (Point p, Face_handle f=Face_handle())
 Returns any nearest vertex to the point p, or the default constructed handle if the triangulation is empty. 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...
 

A point-offset pair (p,off) is said to be in conflict with a cell c iff dt.

side_of_circle(c, p, off) returns ON_BOUNDED_SIDE. The set of faces that are in conflict with (p,off) is star-shaped.

template<class OutputItFaces , class OutputItBoundaryEdges >
std::pair< OutputItFaces, OutputItBoundaryEdges > get_conflicts_and_boundary (const Point &p, OutputItFaces fit, OutputItBoundaryEdges eit, Face_handle start) const
 OutputItFaces is an output iterator with Face_handle as value type. More...
 
template<class OutputItFaces >
OutputItFaces get_conflicts (const Point &p, OutputItFaces fit, Face_handle start) const
 same as above except that only the faces in conflict with p are output. More...
 
template<class OutputItBoundaryEdges >
OutputItBoundaryEdges get_boundary_of_conflicts (const Point &p, OutputItBoundaryEdges eit, Face_handle start) const
 OutputItBoundaryEdges stands for an output iterator with Edge as value type. More...
 

Voronoi diagram

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

Point circumcenter (Face_handle f) const
 Compute the circumcenter of the face pointed to by f. More...
 
Point dual (const Face_handle &f) const
 Returns the center of the circle circumscribed to face f.
 
Segment dual (const Edge &e) const
 returns a segment whose endpoints are the duals of both incident faces.
 
Segment dual (const Edge_circulator &ec) const
 Idem.
 
Segment dual (const Edge_iterator &ei) const
 Idem.
 
template<class Stream >
Stream & draw_dual (Stream &ps)
 output the dual Voronoi diagram to stream ps.
 

Predicates

Oriented_side side_of_oriented_circle (Face_handle f, const Point &p, bool perturb) const
 Returns the side of p with respect to the circle circumscribing the triangle associated with f. More...
 

Checking

Advanced

These methods are mainly a debugging help for the users of advanced features.

bool is_valid (bool verbose=false) const
 This is an advanced function. More...
 
bool is_valid (Face_handle f, bool verbose=false) const
 This is an advanced function. More...
 

Additional Inherited Members

- Public Types inherited from CGAL::Periodic_2_triangulation_2< Traits, Tds >
enum  Iterator_type { STORED = 0, UNIQUE, STORED_COVER_DOMAIN, UNIQUE_COVER_DOMAIN }
 The enum Iterator_type is defined by Periodic_2_triangulation_2 to specify the behavior of geometric iterators. More...
 
enum  Locate_type { VERTEX = 0, EDGE, FACE, EMPTY }
 The enum @ is defined by Periodic_2_triangulation_2 to specify which case occurs when locating a point in the triangulation. More...
 
typedef Traits Geom_traits
 the traits class.
 
typedef Tds Triangulation_data_structure
 the triangulation data structure type.
 
typedef Geom_traits::Periodic_2_offset_2 Offset
 the offset type.
 
typedef Geom_traits::Iso_rectangle_2 Iso_rectangle
 the iso rectangle type.
 
typedef array< int, 2 > Covering_sheets
 integer tuple to store the number of sheets in each direction of space.
 
typedef Geom_traits::Point_2 Point
 the point type.
 
typedef Geom_traits::Segment_2 Segment
 the segment type.
 
typedef Geom_traits::Triangle_2 Triangle
 the triangle type.
 
typedef std::pair< Point, OffsetPeriodic_point
 represents a point-offset pair. More...
 
typedef array< Periodic_point, 2 > Periodic_segment
 a pair of periodic points representing a segment in the periodic domain.
 
typedef array< Periodic_point, 3 > Periodic_triangle
 a triple of periodic points representing a triangle in the periodic domain.
 
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).
 
typedef Tds::Vertex_handle Vertex_handle
 handle to a vertex.
 
typedef Tds::Face_handle Face_handle
 handle to a face.
 
typedef Tds::Face_iterator Face_iterator
 iterator over all faces.
 
typedef Tds::Edge_iterator Edge_iterator
 iterator over all edges.
 
typedef Tds::Vertex_iterator Vertex_iterator
 iterator over all vertices.
 
typedef unspecified_type Unique_vertex_iterator
 iterator over the vertices whose corresponding points lie in the original domain, i.e. More...
 
typedef Face_iterator Finite_faces_iterator
 
typedef Edge_iterator Finite_edges_iterator
 
typedef Vertex_iterator Finite_vertices_iterator
 
typedef Face_iterator All_faces_iterator
 
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 adjacent to a given vertex.
 
typedef unspecified_type Periodic_triangle_iterator
 iterator over the triangles corresponding to faces of the triangulation.
 
typedef unspecified_type Periodic_segment_iterator
 iterator over the segments corresponding to edges of the triangulation.
 
typedef unspecified_type Periodic_point_iterator
 iterator over the points corresponding to vertices of the triangulation.
 
- Public Member Functions inherited from CGAL::Periodic_2_triangulation_2< Traits, Tds >
 Triangulation_2 (const Iso_rectangle &domain=Iso_rectangle(0, 0, 1, 1), const Geom_traits &traits=Geom_traits())
 Introduces an empty triangulation t with domain as original domain. More...
 
 Triangulation_2 (const Triangulation_2 &tr)
 Copy constructor. More...
 
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 vertices resulting in an empty triangulation.
 
const Geom_traitsgeom_traits () const
 Returns a const reference to the triangulation traits object.
 
const Triangulation_data_structure_2tds () const
 Returns a const reference to the triangulation data structure.
 
Iso_rectangle domain () const
 Returns the original domain.
 
Covering_sheets number_of_sheets () const
 Returns the number of sheets of the covering space the triangulation is currently computed in.
 
int dimension () const
 Returns the dimension of the convex hull. More...
 
size_type number_of_vertices () const
 Returns the number of vertices. More...
 
size_type number_of_faces () const
 Returns the number of faces. More...
 
size_type number_of_stored_vertices () const
 Returns the number of vertices in the data structure. More...
 
size_type number_of_stored_faces () const
 Returns the number of faces in the data structure. More...
 
Triangulation_data_structure_2tds ()
 This is an advanced function. More...
 
size_type number_of_edges () const
 Returns the number of edges. More...
 
size_type number_of_stored_edges () const
 Returns the number of edges in the data structure. More...
 
bool is_extensible_triangulation_in_1_sheet_h1 () const
 This is an advanced function. More...
 
bool is_extensible_triangulation_in_1_sheet_h2 () const
 This is an advanced function. More...
 
bool is_triangulation_in_1_sheet () const
 This is an advanced function. More...
 
void convert_to_1_sheeted_covering ()
 This is an advanced function. More...
 
void convert_to_9_sheeted_covering ()
 This is an advanced function. More...
 
Periodic_point periodic_point (const Vertex_handle v) const
 Returns the periodic point given by vertex v. More...
 
Periodic_point periodic_point (const Face_handle f, int i) const
 If this is represented in the 1-sheeted covering space, this function returns the periodic point given by the \( i\)-th vertex of face f, that is the point in the original domain and the offset of the vertex in f. More...
 
Periodic_segment periodic_segment (const Face_handle f, int i) const
 Returns the periodic segment formed by the two point-offset pairs corresponding to the two vertices of edge (f,i). More...
 
Periodic_segment periodic_segment (const Edge &e) const
 Same as the previous method for edge e.
 
Periodic_triangle periodic_triangle (const Face_handle f) const
 Returns the periodic triangle formed by the three point-offset pairs corresponding to the three vertices of facet f.
 
Point point (const Periodic_point &pp) const
 Converts the Periodic_point pp (point-offset pair) to the corresponding Point in \( \mathbb R^3\).
 
Segment segment (const Periodic_segment &s) const
 Converts the Periodic_segment s to a Segment.
 
Triangle triangle (const Periodic_triangle &t) const
 Converts the Periodic_triangle this to a Triangle.
 
Segment segment (Face_handle f, int i) const
 Equivalent to the call t.segment(t.periodic_segment(f,i));
 
Segment segment (const Edge &e) const
 Equivalent to the call t.segment(t.periodic_segment(e));
 
Segment segment (const Edge_circulator &ec) const
 Equivalent to the call t.segment(t.periodic_segment(ec->first, ec->second));
 
Segment segment (const Edge_iterator &ei) const
 Equivalent to the call t.segment(t.periodic_segment(ei->first, ei->second));
 
Triangle triangle (Face_handle f) const
 Equivalent to the call t.triangle(t.periodic_triangle(f));
 
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 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 triangulation is not empty, a face that contains the query in its interior or on its boundary is returned. 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 the point p lies.
 
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.
 
Periodic_point_iterator periodic_points_begin (Iterator_type it=STORED) const
 Iterates over the points of the triangulation. More...
 
Periodic_point_iterator periodic_points_end (Iterator_type it=STORED) const
 Past-the-end iterator. More...
 
Periodic_segment_iterator periodic_segments_begin (Iterator_type it=STORED) const
 Iterates over the segments of the triangulation. More...
 
Periodic_segment_iterator periodic_segments_end (Iterator_type it=STORED) const
 Past-the-end iterator. More...
 
Periodic_triangle_iterator periodic_triangles_begin (Iterator_type it=STORED) const
 Iterates over the triangles of the triangulation. More...
 
Periodic_triangle_iterator periodic_triangles_end (Iterator_type it=STORED) const
 Past-the-end iterator. 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 adjacent_vertices (Vertex_handle v) const
 Starts at an arbitrary vertex adjacent to v.
 
Vertex_circulator adjacent_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...
 
Vertex_handle insert (const Point &p, Face_handle f=Face_handle())
 
Vertex_handle insert (const Point &p, Locate_type lt, Face_handle loc, int li)
 
Vertex_handle push_back (const Point &p)
 
template<class InputIterator >
int insert (InputIterator first, InputIterator last)
 
Vertex_handle insert_first (const Point &p)
 This is an advanced function. More...
 
Vertex_handle insert_in_face (const Point &p, Face_handle f)
 This is an advanced function. More...
 
void remove_degree_3 (Vertex_handle v)
 This is an advanced function. More...
 
void remove_first (Vertex_handle v)
 This is an advanced function. More...
 
template<class EdgeIt >
Vertex_handle star_hole (Point p, EdgeIt edge_begin, EdgeIt edge_end)
 This is an advanced function. More...
 
template<class EdgeIt , class FaceIt >
Vertex_handle star_hole (Point p, EdgeIt edge_begin, EdgeIt edge_end, FaceIt face_begin, FaceIt face_end)
 This is an advanced function. More...
 
void set_domain (const Iso_rectangle dom)
 This is an advanced function. More...
 
int ccw (int i) const
 Returns \( i+1\) modulo 3. More...
 
int cw (int i) const
 Returns \( i+2\) modulo 3. More...
 
void flippable (Face_handle f, int i)
 
size_t degree (Vertex_handle v)
 Returns the degree of the vertex v
 
bool is_valid (bool verbose=false, int level=0) const
 This is an advanced function. More...
 
- Public Member Functions inherited from CGAL::Triangulation_cw_ccw_2
 Triangulation_cw_ccw_2 ()
 
int ccw (const int i) const
 
int cw (const int i) const
 

Constructor & Destructor Documentation

◆ Periodic_2_Delaunay_triangulation_2() [1/2]

template<typename Traits , typename Tds >
CGAL::Periodic_2_Delaunay_triangulation_2< Traits, Tds >::Periodic_2_Delaunay_triangulation_2 ( const Iso_rectangle domain = Iso_rectangle(0, 0, 1, 1),
const Geom_traits traits = Geom_traits() 
)

Creates an empty periodic Delaunay triangulation dt, with domain as original domain and possibly specifying a traits class traits.

Precondition
domain is a square.

◆ Periodic_2_Delaunay_triangulation_2() [2/2]

template<typename Traits , typename Tds >
template<class InputIterator >
CGAL::Periodic_2_Delaunay_triangulation_2< Traits, Tds >::Periodic_2_Delaunay_triangulation_2 ( InputIterator  first,
InputIterator  last,
const Iso_rectangle domain = Iso_rectangle(0, 0, 1, 1),
const Geom_traits traits = Geom_traits() 
)

Equivalent to constructing an empty triangulation with the optional domain and traits class arguments and calling insert(first,last).

Precondition
The value_type of first and last are Points lying inside domain and domain is a square.

Member Function Documentation

◆ circumcenter()

template<typename Traits , typename Tds >
Point CGAL::Periodic_2_Delaunay_triangulation_2< Traits, Tds >::circumcenter ( Face_handle  f) const

Compute the circumcenter of the face pointed to by f.

This function is available only if the corresponding function is provided in the geometric traits.

◆ get_boundary_of_conflicts()

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

OutputItBoundaryEdges stands for an output iterator with Edge as value type.

This function outputs in the container pointed to by eit, the boundary of the zone in conflict with p. The boundary edges of the conflict zone are output in counterclockwise order and each edge is described through the incident face which is not in conflict with p. The function returns the resulting output iterator.

Precondition
start is in conflict with p and p lies in the original domain domain.

◆ get_conflicts()

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

same as above except that only the faces in conflict with p are output.

The function returns the resulting output iterator.

Precondition
start is in conflict with p and p lies in the original domain domain.

◆ get_conflicts_and_boundary()

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

OutputItFaces is an output iterator with Face_handle as value type.

OutputItBoundaryEdges stands for an output iterator with Edge as value type. This members function outputs in the container pointed to by fit the faces which are in conflict with point p i. e. the faces whose circumcircle contains p. It outputs in the container pointed to by eit the boundary of the zone in conflict with p. 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 std::pair the resulting output iterators.

Precondition
start is in conflict with p and p lies in the original domain domain.

◆ insert() [1/3]

template<typename Traits , typename Tds >
Vertex_handle CGAL::Periodic_2_Delaunay_triangulation_2< Traits, Tds >::insert ( const Point p,
Face_handle  start = Face_handle() 
)

Inserts point p in the triangulation and returns the corresponding vertex.

The optional argument start is used as a starting place for the point location.

Precondition
p lies in the original domain domain.
Examples:
Periodic_2_triangulation_2/p2t2_adding_handles.cpp, Periodic_2_triangulation_2/p2t2_colored_vertices.cpp, Periodic_2_triangulation_2/p2t2_covering.cpp, and Periodic_2_triangulation_2/p2t2_geometric_access.cpp.

◆ insert() [2/3]

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

Inserts point p in the triangulation and returns the corresponding vertex.

Similar to the above insert() function, but takes as additional parameter the return values of a previous location query. See description of Periodic_2_triangulation_2::locate().

Precondition
p lies in the original domain domain.

◆ insert() [3/3]

template<typename Traits , typename Tds >
template<class InputIterator >
std::ptrdiff_t CGAL::Periodic_2_Delaunay_triangulation_2< Traits, Tds >::insert ( InputIterator  first,
InputIterator  last,
bool  is_large_point_set = false 
)

Inserts the points in the iterator range [first, last).

Returns the number of inserted points. This function uses spatial sorting (cf. chapter Spatial Sorting) and therefore is not guaranteed to insert the points following the order of InputIterator. If the third argument is_large_point_set is set to true a heuristic for optimizing the insertion of large point sets is applied.

Precondition
The value_type of first and last are Points and all points in the range lie inside domain.

◆ is_valid() [1/2]

template<typename Traits , typename Tds >
bool CGAL::Periodic_2_Delaunay_triangulation_2< Traits, Tds >::is_valid ( bool  verbose = false) const

This is an advanced function.

Advanced

Checks the combinatorial validity of the triangulation and the validity of its geometric embedding (see Section Representation). Also checks that all the circumscribing circles of cells are empty.

When verbose is set to true, messages describing the first invalidity encountered are printed.

◆ is_valid() [2/2]

template<typename Traits , typename Tds >
bool CGAL::Periodic_2_Delaunay_triangulation_2< Traits, Tds >::is_valid ( Face_handle  f,
bool  verbose = false 
) const

This is an advanced function.

Advanced

Checks the combinatorial and geometric validity of the cell (see Section Representation). Also checks that the circumscribing circle of cells is empty.

When verbose is set to true, messages are printed to give a precise indication of the kind of invalidity encountered.

◆ move_if_no_collision()

template<typename Traits , typename Tds >
Vertex_handle CGAL::Periodic_2_Delaunay_triangulation_2< Traits, Tds >::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.

Otherwise, the triangulation is not modified and the vertex at point p is returned.

Precondition
p lies in the original domain domain.

◆ move_point()

template<typename Traits , typename Tds >
Vertex_handle CGAL::Periodic_2_Delaunay_triangulation_2< Traits, Tds >::move_point ( Vertex_handle  v,
const Point p 
)

Moves the point stored in v to p, while preserving the Delaunay property.

This performs an action semantically equivalent to remove(v) followed by insert(p), but is supposedly faster to performing these two operations separately when the point has not moved much. Returns the handle to the new vertex.

Precondition
p lies in the original domain domain.

◆ nearest_vertex()

template<typename Traits , typename Tds >
Vertex_handle CGAL::Periodic_2_Delaunay_triangulation_2< Traits, Tds >::nearest_vertex ( Point  p,
Face_handle  f = Face_handle() 
)

Returns any nearest vertex to the point p, or the default constructed handle if the triangulation is empty.

The optional argument f is a hint specifying where to start the search. It always returns a vertex corresponding to a point inside domain even if computing in a multiply sheeted covering space.

Precondition
f is a face of dt and p lies in the original domain domain.

◆ side_of_oriented_circle() [1/2]

template<typename Traits , typename Tds >
Oriented_side CGAL::Periodic_2_Delaunay_triangulation_2< Traits, Tds >::side_of_oriented_circle ( Face_handle  f,
const Point p 
)

Returns on which side of the circumcircle of face f lies the point p.

The circle is assumed to be counterclockwise oriented, so its positive side correspond to its bounded side. This predicate is available only if the corresponding predicates on points is provided in the geometric traits class.

◆ side_of_oriented_circle() [2/2]

template<typename Traits , typename Tds >
Oriented_side CGAL::Periodic_2_Delaunay_triangulation_2< Traits, Tds >::side_of_oriented_circle ( Face_handle  f,
const Point p,
bool  perturb 
) const

Returns the side of p with respect to the circle circumscribing the triangle associated with f.

Periodic copies are checked if necessary.