CGAL 4.8.1 - 2D Periodic Triangulations
|
#include <CGAL/Periodic_2_Delaunay_triangulation_2.h>
CGAL::Periodic_2_triangulation_2< Traits, Tds >.
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
Traits | is 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.
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>
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_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... | |
A point-offset pair (
| |
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 | 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 |
Advanced Checks the combinatorial validity of the triangulation and the validity of its geometric embedding (see Section Representation). More... | |
bool | is_valid (Face_handle f, bool verbose=false) const |
Advanced Checks the combinatorial and geometric validity of the cell (see Section Representation). 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, Offset > | Periodic_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_traits & | geom_traits () const |
Returns a const reference to the triangulation traits object. | |
const Triangulation_data_structure_2 & | tds () 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_2 & | tds () |
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 |
Advanced The current triangulation remains a triangulation in the 1-sheeted covering space even after adding points if this method returns true . More... | |
bool | is_extensible_triangulation_in_1_sheet_h2 () const |
Advanced The same as is_extensible_triangulation_in_1_sheet_h1() but with a more precise heuristic, i.e. More... | |
bool | is_triangulation_in_1_sheet () const |
Advanced Returns true if the current triangulation would still be a triangulation in the 1-sheeted covering space, returns false otherwise. More... | |
void | convert_to_1_sheeted_covering () |
Advanced Converts the current triangulation into the same periodic triangulation in the 1-sheeted covering space. More... | |
void | convert_to_9_sheeted_covering () |
Advanced Converts the current triangulation into the same periodic triangulation in the 9-sheeted covering space. 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 . | |
Point | circumcenter (Face_handle f) const |
Compute the circumcenter of the face pointed to by f. More... | |
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 <, 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. | |
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... | |
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) |
Vertex_handle | insert_in_face (const Point &p, Face_handle f) |
void | remove_degree_3 (Vertex_handle v) |
void | remove_first (Vertex_handle v) |
template<class EdgeIt > | |
Vertex_handle | star_hole (Point p, EdgeIt edge_begin, EdgeIt edge_end) |
Advanced 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) |
Advanced 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... | |
void | set_domain (const Iso_rectangle dom) |
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 |
Advanced Checks the combinatorial validity of the triangulation and also the validity of its geometric embedding. More... | |
Related Functions inherited from CGAL::Periodic_2_triangulation_2< Traits, Tds > | |
ostream & | operator<< (ostream &os, const Periodic_2_triangulation_2< Traits, Tds > &t) |
Writes the triangulation t into the stream os . More... | |
istream & | operator>> (istream &is, Triangulation_2< Traits, Tds > &t) |
Reads a triangulation from stream is and assigns it to t . More... | |
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
.
domain
is a square. 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)
.
value_type
of first
and last
are Point
s lying inside domain
and domain
is a square. 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.
start
is in conflict with p
and p
lies in the original domain domain
. 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.
start
is in conflict with p
and p
lies in the original domain domain
. 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.
start
is in conflict with p
and p
lies in the original domain domain
. 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.
p
lies in the original domain domain
. 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()
.
p
lies in the original domain domain
. 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.
value_type
of first
and last
are Point
s and all points in the range lie inside domain
. bool CGAL::Periodic_2_Delaunay_triangulation_2< Traits, Tds >::is_valid | ( | bool | verbose = false ) | const |
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.
bool CGAL::Periodic_2_Delaunay_triangulation_2< Traits, Tds >::is_valid | ( | Face_handle | f, |
bool | verbose = false |
||
) | const |
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.
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.
p
lies in the original domain domain
. 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.
p
lies in the original domain domain
. 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.
f
is a face of dt
and p
lies in the original domain domain
. 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.