\( \newcommand{\E}{\mathrm{E}} \) \( \newcommand{\A}{\mathrm{A}} \) \( \newcommand{\R}{\mathrm{R}} \) \( \newcommand{\N}{\mathrm{N}} \) \( \newcommand{\Q}{\mathrm{Q}} \) \( \newcommand{\Z}{\mathrm{Z}} \) \( \def\ccSum #1#2#3{ \sum_{#1}^{#2}{#3} } \def\ccProd #1#2#3{ \sum_{#1}^{#2}{#3} }\)
CGAL 4.8.1 - 3D Triangulations
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Groups Pages
CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy, SurjectiveLockDataStructure > Class Template Reference

#include <CGAL/Delaunay_triangulation_3.h>

Inherits from

CGAL::Triangulation_3< DelaunayTriangulationTraits_3, Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::Triangulation_data_structure, SurjectiveLockDataStructure >.

Definition

The class Delaunay_triangulation_3 represents a three-dimensional Delaunay triangulation.

Template Parameters
DelaunayTriangulationTraits_3is the geometric traits class.
TriangulationDataStructure_3is the triangulation data structure. It has the default value Triangulation_data_structure_3<Triangulation_vertex_base_3<DelaunayTriangulationTraits_3>, Delaunay_triangulation_cell_base_3<DelaunayTriangulationTraits_3> >.
LocationPolicyis a tag which must be a Location_policy<Tag>: either Fast_location or Compact_location. Fast_location offers faster ( \( O(\log n)\) time) point location, which can be beneficial when performing point locations or random point insertions (with no good location hint) in large data sets. It is currently implemented using an additional triangulation hierarchy data structure [12]. The default is Compact_location, which saves memory (3-5%) by avoiding the need for this separate data structure, and point location is then performed roughly in \( O(n^{1/3})\) time. If the triangulation is parallel (see user manual), the default compact location policy must be used. Note that this argument can also come in second position, which can be useful when the default value for the TriangulationDataStructure_3 parameter is satisfactory (this is using so-called deduced parameters). Note that this argument replaces the functionality provided before CGAL 3.6 by Triangulation_hierarchy_3. An example of use can be found in the user manual Fast Point Location for Delaunay Triangulations.
SurjectiveLockDataStructureis an optional parameter to specify the type of the spatial lock data structure. It is only used if the triangulation data structure used is concurrency-safe (i.e. when TriangulationDataStructure_3::Concurrency_tag is Parallel_tag). It must be a model of the SurjectiveLockDataStructure concept, with Object being a Point. The default value is Spatial_lock_grid_3<Tag_priority_blocking> if the triangulation data structure is concurrency-safe, and void otherwise. In order to use concurrent operations, the user must provide a reference to a SurjectiveLockDataStructure instance via the constructor or Triangulation_3::set_lock_data_structure.

If TriangulationDataStructure_3::Concurrency_tag is Parallel_tag, some operations, such as insertion/removal of a range of points, are performed in parallel. See the documentation of the operations for more details.

See Also
CGAL::Regular_triangulation_3
Examples:
Triangulation_3/adding_handles_3.cpp, Triangulation_3/color.cpp, Triangulation_3/copy_triangulation_3.cpp, Triangulation_3/fast_location_3.cpp, Triangulation_3/find_conflicts_3.cpp, Triangulation_3/info_insert_with_pair_iterator.cpp, Triangulation_3/info_insert_with_transform_iterator.cpp, Triangulation_3/info_insert_with_zip_iterator.cpp, and Triangulation_3/parallel_insertion_in_delaunay_3.cpp.

Types

typedef LocationPolicy Location_policy
 
typedef SurjectiveLockDataStructure Lock_data_structure
 

In addition to those inherited, the following types are defined, for use by the construction of the Voronoi diagram:

typedef
DelaunayTriangulationTraits_3::Line_3 
Line
 
typedef
DelaunayTriangulationTraits_3::Ray_3 
Ray
 
typedef
DelaunayTriangulationTraits_3::Plane_3 
Plane
 
typedef
DelaunayTriangulationTraits_3::Object_3 
Object
 

Creation

 Delaunay_triangulation_3 (const DelaunayTriangulationTraits_3 &traits=DelaunayTriangulationTraits_3(), Lock_data_structure *lock_ds=NULL)
 Creates an empty Delaunay triangulation, possibly specifying a traits class traits. More...
 
 Delaunay_triangulation_3 (const Delaunay_triangulation_3 &dt1)
 Copy constructor. More...
 
template<class InputIterator >
 Delaunay_triangulation_3 (InputIterator first, InputIterator last, const DelaunayTriangulationTraits_3 &traits=DelaunayTriangulationTraits_3(), Lock_data_structure *lock_ds=NULL)
 Equivalent to constructing an empty triangulation with the optional traits class argument and calling insert(first,last). More...
 
template<class InputIterator >
 Delaunay_triangulation_3 (InputIterator first, InputIterator last, Lock_data_structure *lock_ds, const DelaunayTriangulationTraits_3 &traits=DelaunayTriangulationTraits_3())
 Same as before, with last two parameters in reverse order.
 

Insertion

Vertex_handle insert (const Point &p, Cell_handle start=Cell_handle(), bool *could_lock_zone=NULL)
 Inserts point p in the triangulation and returns the corresponding vertex. More...
 
Vertex_handle insert (const Point &p, Vertex_handle hint, bool *could_lock_zone=NULL)
 Same as above but uses hint as a starting place for the search.
 
Vertex_handle insert (const Point &p, Locate_type lt, Cell_handle loc, int li, int lj, bool *could_lock_zone=NULL)
 Inserts point p in the triangulation and returns the corresponding vertex. More...
 
template<class PointInputIterator >
std::ptrdiff_t insert (PointInputIterator first, PointInputIterator last)
 Inserts the points in the iterator range [first,last). More...
 
template<class PointWithInfoInputIterator >
std::ptrdiff_t insert (PointWithInfoInputIterator first, PointWithInfoInputIterator last)
 Inserts the points in the iterator range [first,last). More...
 

Displacement

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)
 same as above if there is no collision. More...
 

Removal

When a vertex v is removed from a triangulation, all the cells incident to v must be removed, and the polyhedral region consisting of all the tetrahedra that are incident to v must be re-triangulated.

So, the problem reduces to triangulating a polyhedral region, while preserving its boundary, or to compute a constrained* triangulation. This is known to be sometimes impossible: the Schönhardt polyhedron cannot be triangulated [19]. However, when dealing with Delaunay triangulations, the case of such polyhedra that cannot be re-triangulated cannot happen, so CGAL proposes a vertex removal.

If,due to some point removals, the size of the Delaunay triangulation decreases drastically, it might be interesting to defragment the CGAL::Compact_container (used by the Triangulation_data_structure_3).

void remove (Vertex_handle v)
 Removes the vertex v from the triangulation. More...
 
bool remove (Vertex_handle v, bool *could_lock_zone)
 Removes the vertex v from the triangulation. More...
 
template<typename InputIterator >
int remove (InputIterator first, InputIterator beyond)
 Removes the vertices specified by the iterator range [first, beyond). More...
 
template<typename InputIterator >
int remove_cluster (InputIterator first, InputIterator beyond)
 This function has exactly the same result and the same preconditions as remove(first, beyond). More...
 

Queries

Bounded_side side_of_sphere (Cell_handle c, const Point &p) const
 Returns a value indicating on which side of the circumscribed sphere of c the point p lies. More...
 
Bounded_side side_of_circle (const Facet &f, const Point &p) const
 Returns a value indicating on which side of the circumscribed circle of f the point p lies. More...
 
Bounded_side side_of_circle (Cell_handle c, int i, const Point &p)
 Same as the previous method for facet i of cell c.
 
Vertex_handle nearest_vertex (Point p, Cell_handle c=Cell_handle())
 Returns any nearest vertex to the point p, or the default constructed handle if the triangulation is empty. More...
 
Vertex_handle nearest_vertex_in_cell (Point p, Cell_handle c)
 Returns the vertex of the cell c that is nearest to \( p\).
 

A point p is said to be in conflict with a cell c in dimension 3 (resp. a facet f in dimension 2) iff dt.side_of_sphere(c, p) (resp. dt.side_of_circle(f, p)) returns ON_BOUNDED_SIDE.

The set of cells (resp. facets in dimension 2) which are in conflict with p is connected, and it forms a hole.

template<class OutputIteratorBoundaryFacets , class OutputIteratorCells >
std::pair
< OutputIteratorBoundaryFacets,
OutputIteratorCells > 
find_conflicts (Point p, Cell_handle c, OutputIteratorBoundaryFacets bfit, OutputIteratorCells cit, bool *could_lock_zone=NULL)
 Computes the conflict hole induced by p. More...
 
template<class OutputIteratorBoundaryFacets , class OutputIteratorCells , class OutputIteratorInternalFacets >
Triple
< OutputIteratorBoundaryFacets,
OutputIteratorCells,
OutputIteratorInternalFacets > 
find_conflicts (Point p, Cell_handle c, OutputIteratorBoundaryFacets bfit, OutputIteratorCells cit, OutputIteratorInternalFacets ifit, bool *could_lock_zone=NULL)
 Same as the other find_conflicts() function, except that it also computes the internal facets, i.e. the facets common to two cells which are in conflict with p. More...
 
template<class OutputIterator >
OutputIterator vertices_in_conflict (Point p, Cell_handle c, OutputIterator res)
 
template<class OutputIterator >
OutputIterator vertices_on_conflict_zone_boundary (Point p, Cell_handle c, OutputIterator res)
 Similar to find_conflicts(), but reports the vertices which are on the boundary of the conflict hole of p, in the output iterator res. More...
 

A face (cell, facet or edge) is said to be a Gabriel face iff its smallest circumscribing sphere do not enclose any vertex of the triangulation.

Any Gabriel face belongs to the Delaunay triangulation, but the reciprocal is not true. The following member functions test the Gabriel property of Delaunay faces.

bool is_Gabriel (Cell_handle c, int i)
 
bool is_Gabriel (Cell_handle c, int i, int j)
 
bool is_Gabriel (const Facet &f)
 
bool is_Gabriel (const Edge &e)
 

Voronoi Diagram

CGAL offers several functionalities to display the Voronoi diagram of a set of points in 3D.

Note that the user should use a kernel with exact constructions in order to guarantee the computation of the Voronoi diagram (as opposed to computing the triangulation only, which requires only exact predicates).

Point dual (Cell_handle c) const
 Returns the circumcenter of the four vertices of c. More...
 
Object dual (Facet f) const
 Returns the dual of facet f, which is. More...
 
Object dual (Cell_handle c, int i) const
 same as the previous method for facet (c,i).
 
Line dual_support (Cell_handle c, int i) const
 returns the line supporting the dual of the facet. More...
 
template<class Stream >
Stream & draw_dual (Stream &os)
 Sends the set of duals to all the facets of dt into os.
 

Checking

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

bool is_valid (bool verbose=false) const
 Checks the combinatorial validity of the triangulation and the validity of its geometric embedding (see Section Representation). More...
 
bool is_valid (Cell_handle c, bool verbose=false) const
 This is a function for debugging purpose. More...
 

Additional Inherited Members

- Public Types inherited from CGAL::Triangulation_3< DelaunayTriangulationTraits_3, Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::Triangulation_data_structure, SurjectiveLockDataStructure >
enum  Locate_type
 The enum Locate_type is defined by Triangulation_3 to specify which case occurs when locating a point in the triangulation.
 
typedef
Delaunay_triangulation_3
< DelaunayTriangulationTraits_3,
TriangulationDataStructure_3,
LocationPolicy >
::Triangulation_data_structure 
Triangulation_data_structure
 
typedef SurjectiveLockDataStructure Lock_data_structure
 
typedef
DelaunayTriangulationTraits_3 
Geom_traits
 
typedef
DelaunayTriangulationTraits_3::Point_3 
Point
 
typedef
DelaunayTriangulationTraits_3::Segment_3 
Segment
 
typedef
DelaunayTriangulationTraits_3::Triangle_3 
Triangle
 
typedef
DelaunayTriangulationTraits_3::Tetrahedron_3 
Tetrahedron
 
typedef
Delaunay_triangulation_3
< DelaunayTriangulationTraits_3,
TriangulationDataStructure_3,
LocationPolicy >
::Triangulation_data_structure::Vertex 
Vertex
 
typedef
Delaunay_triangulation_3
< DelaunayTriangulationTraits_3,
TriangulationDataStructure_3,
LocationPolicy >
::Triangulation_data_structure::Cell 
Cell
 
typedef
Delaunay_triangulation_3
< DelaunayTriangulationTraits_3,
TriangulationDataStructure_3,
LocationPolicy >
::Triangulation_data_structure::Facet 
Facet
 
typedef
Delaunay_triangulation_3
< DelaunayTriangulationTraits_3,
TriangulationDataStructure_3,
LocationPolicy >
::Triangulation_data_structure::Edge 
Edge
 
typedef
Delaunay_triangulation_3
< DelaunayTriangulationTraits_3,
TriangulationDataStructure_3,
LocationPolicy >
::Triangulation_data_structure::Vertex_handle 
Vertex_handle
 handle to a vertex
 
typedef
Delaunay_triangulation_3
< DelaunayTriangulationTraits_3,
TriangulationDataStructure_3,
LocationPolicy >
::Triangulation_data_structure::Cell_handle 
Cell_handle
 handle to a cell
 
typedef
Triangulation_simplex_3< Self > 
Simplex
 Reference to a simplex (vertex, edge, facet or cell) of the triangulation.
 
typedef
Delaunay_triangulation_3
< DelaunayTriangulationTraits_3,
TriangulationDataStructure_3,
LocationPolicy >
::Triangulation_data_structure::size_type 
size_type
 Size type (an unsigned integral type)
 
typedef
Delaunay_triangulation_3
< DelaunayTriangulationTraits_3,
TriangulationDataStructure_3,
LocationPolicy >
::Triangulation_data_structure::difference_type 
difference_type
 Difference type (a signed integral type)
 
typedef
Delaunay_triangulation_3
< DelaunayTriangulationTraits_3,
TriangulationDataStructure_3,
LocationPolicy >
::Triangulation_data_structure::Cell_iterator 
All_cells_iterator
 iterator over cells
 
typedef
Delaunay_triangulation_3
< DelaunayTriangulationTraits_3,
TriangulationDataStructure_3,
LocationPolicy >
::Triangulation_data_structure::Facet_iterator 
All_facets_iterator
 iterator over facets
 
typedef
Delaunay_triangulation_3
< DelaunayTriangulationTraits_3,
TriangulationDataStructure_3,
LocationPolicy >
::Triangulation_data_structure::Edge_iterator 
All_edges_iterator
 iterator over edges
 
typedef
Delaunay_triangulation_3
< DelaunayTriangulationTraits_3,
TriangulationDataStructure_3,
LocationPolicy >
::Triangulation_data_structure::Vertex_iterator 
All_vertices_iterator
 iterator over vertices
 
typedef unspecified_type Finite_cells_iterator
 iterator over finite cells
 
typedef unspecified_type Finite_facets_iterator
 iterator over finite facets
 
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
Delaunay_triangulation_3
< DelaunayTriangulationTraits_3,
TriangulationDataStructure_3,
LocationPolicy >
::Triangulation_data_structure::Cell_circulator 
Cell_circulator
 circulator over all cells incident to a given edge
 
typedef
Delaunay_triangulation_3
< DelaunayTriangulationTraits_3,
TriangulationDataStructure_3,
LocationPolicy >
::Triangulation_data_structure::Facet_circulator 
Facet_circulator
 circulator over all facets incident to a given edge
 
typedef
Delaunay_triangulation_3
< DelaunayTriangulationTraits_3,
TriangulationDataStructure_3,
LocationPolicy >
::Triangulation_data_structure::Concurrency_tag 
Concurrency_tag
 Concurrency tag (from the TDS).
 
- Public Member Functions inherited from CGAL::Triangulation_3< DelaunayTriangulationTraits_3, Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::Triangulation_data_structure, SurjectiveLockDataStructure >
 Triangulation_3 (const DelaunayTriangulationTraits_3 &traits=DelaunayTriangulationTraits_3(), Lock_data_structure *lock_ds=NULL)
 Introduces a triangulation t having only one vertex which is the infinite vertex. More...
 
 Triangulation_3 (Lock_data_structure *lock_ds=NULL, const DelaunayTriangulationTraits_3 &traits=DelaunayTriangulationTraits_3())
 Same as the previous one, but with parameters in reverse order.
 
 Triangulation_3 (const Triangulation_3 &tr)
 Copy constructor. More...
 
 Triangulation_3 (InputIterator first, InputIterator last, const DelaunayTriangulationTraits_3 &traits=DelaunayTriangulationTraits_3(), Lock_data_structure *lock_ds=NULL)
 Equivalent to constructing an empty triangulation with the optional traits class argument and calling insert(first,last).
 
Triangulation_3operator= (const Triangulation_3 &tr)
 The triangulation tr is duplicated, and modifying the copy after the duplication does not modify the original. More...
 
void swap (Triangulation_3 &tr)
 The triangulations tr and t are swapped. More...
 
void clear ()
 Deletes all finite vertices and all cells of t.
 
bool operator== (const Triangulation_3< GT, Tds1 > &t1, const Triangulation_3< GT, Tds2 > &t2)
 Equality operator. More...
 
bool operator!= (const Triangulation_3< GT, Tds1 > &t1, const Triangulation_3< GT, Tds2 > &t2)
 The opposite of operator==.
 
const
DelaunayTriangulationTraits_3
geom_traits () const
 Returns a const reference to the geometric traits object.
 
const Delaunay_triangulation_3
< DelaunayTriangulationTraits_3,
TriangulationDataStructure_3,
LocationPolicy >
::Triangulation_data_structure
tds () const
 Returns a const reference to the triangulation data structure.
 
Delaunay_triangulation_3
< DelaunayTriangulationTraits_3,
TriangulationDataStructure_3,
LocationPolicy >
::Triangulation_data_structure
tds ()
 Returns a reference to the triangulation data structure. More...
 
int dimension () const
 Returns the dimension of the affine hull.
 
size_type number_of_vertices () const
 Returns the number of finite vertices.
 
size_type number_of_cells () const
 Returns the number of cells or 0 if t.dimension() < 3.
 
Vertex_handle infinite_vertex ()
 Returns the infinite vertex.
 
void set_infinite_vertex (Vertex_handle v)
 This is an advanced function. More...
 
Cell_handle infinite_cell () const
 Returns a cell incident to the infinite vertex.
 
size_type number_of_facets () const
 The number of facets. More...
 
size_type number_of_edges () const
 The number of edges. More...
 
size_type number_of_finite_cells () const
 The number of finite cells. More...
 
size_type number_of_finite_facets () const
 The number of finite facets. More...
 
size_type number_of_finite_edges () const
 The number of finite edges. More...
 
Tetrahedron tetrahedron (Cell_handle c) const
 Returns the tetrahedron formed by the four vertices of c. More...
 
Triangle triangle (Cell_handle c, int i) const
 Returns the triangle formed by the three vertices of facet (c,i). More...
 
Triangle triangle (const Facet &f) const
 Same as the previous method for facet f. More...
 
Segment segment (const Edge &e) const
 Returns the line segment formed by the vertices of e. More...
 
Segment segment (Cell_handle c, int i, int j) const
 Same as the previous method for edge (c,i,j). More...
 
const Pointpoint (Cell_handle c, int i) const
 Returns the point given by vertex i of cell c. More...
 
const Pointpoint (Vertex_handle v) const
 Same as the previous method for vertex v. More...
 
bool is_infinite (Vertex_handle v) const
 true, iff vertex v is the infinite vertex.
 
bool is_infinite (Cell_handle c) const
 true, iff c is incident to the infinite vertex. More...
 
bool is_infinite (Cell_handle c, int i) const
 true, iff the facet i of cell c is incident to the infinite vertex. More...
 
bool is_infinite (const Facet &f) const
 true iff facet f is incident to the infinite vertex. More...
 
bool is_infinite (Cell_handle c, int i, int j) const
 true, iff the edge (i,j) of cell c is incident to the infinite vertex. More...
 
bool is_infinite (const Edge &e) const
 true iff edge e is incident to the infinite vertex. More...
 
bool is_vertex (const Point &p, Vertex_handle &v) const
 Tests whether p is a vertex of t by locating p in the triangulation. More...
 
bool is_vertex (Vertex_handle v) const
 Tests whether v is a vertex of t.
 
bool is_edge (Vertex_handle u, Vertex_handle v, Cell_handle &c, int &i, int &j) const
 Tests whether (u,v) is an edge of t. More...
 
bool is_facet (Vertex_handle u, Vertex_handle v, Vertex_handle w, Cell_handle &c, int &i, int &j, int &k) const
 Tests whether (u,v,w) is a facet of t. More...
 
bool is_cell (Cell_handle c) const
 Tests whether c is a cell of t.
 
bool is_cell (Vertex_handle u, Vertex_handle v, Vertex_handle w, Vertex_handle x, Cell_handle &c, int &i, int &j, int &k, int &l) const
 Tests whether (u,v,w,x) is a cell of t. More...
 
bool is_cell (Vertex_handle u, Vertex_handle v, Vertex_handle w, Vertex_handle x, Cell_handle &c) const
 Tests whether (u,v,w,x) is a cell of t and computes this cell c. More...
 
bool has_vertex (const Facet &f, Vertex_handle v, int &j) const
 If v is a vertex of f, then j is the index of v in the cell f.first, and the method returns true. More...
 
bool has_vertex (Cell_handle c, int i, Vertex_handle v, int &j) const
 Same for facet (c,i). More...
 
bool has_vertex (const Facet &f, Vertex_handle v) const
 
bool has_vertex (Cell_handle c, int i, Vertex_handle v) const
 Same as the first two methods, but these two methods do not return the index of the vertex.
 
bool are_equal (Cell_handle c, int i, Cell_handle n, int j) const
 
bool are_equal (const Facet &f, const Facet &g) const
 
bool are_equal (const Facet &f, Cell_handle n, int j) const
 For these three methods: More...
 
Cell_handle locate (const Point &query, Cell_handle start=Cell_handle(), bool *could_lock_zone=NULL) const
 If the point query lies inside the convex hull of the points, the cell that contains the query in its interior is returned. More...
 
Cell_handle locate (const Point &query, Vertex_handle hint, bool *could_lock_zone=NULL) const
 Same as above but uses hint as the starting place for the search.
 
Cell_handle locate (const Point &query, Locate_type &lt, int &li, int &lj, Cell_handle start=Cell_handle(), bool *could_lock_zone=NULL) const
 If query lies inside the affine hull of the points, the \( k\)-face (finite or infinite) that contains query in its interior is returned, by means of the cell returned together with lt, which is set to the locate type of the query (VERTEX, EDGE, FACET, CELL, or OUTSIDE_CONVEX_HULL if the cell is infinite and query lies strictly in it) and two indices li and lj that specify the \( k\)-face of the cell containing query. More...
 
Cell_handle locate (const Point &query, Locate_type &lt, int &li, int &lj, Vertex_handle hint, bool *could_lock_zone=NULL) const
 Same as above but uses hint as the starting place for the search.
 
Cell_handle inexact_locate (const Point &query, Cell_handle start=Cell_handle()) const
 Same as locate() but uses inexact predicates. More...
 
Bounded_side side_of_cell (const Point &p, Cell_handle c, Locate_type &lt, int &li, int &lj) const
 Returns a value indicating on which side of the oriented boundary of c the point p lies. More...
 
Bounded_side side_of_facet (const Point &p, const Facet &f, Locate_type &lt, int &li, int &lj) const
 Returns a value indicating on which side of the oriented boundary of f the point p lies: More...
 
Bounded_side side_of_facet (const Point &p, Cell_handle c, Locate_type &lt, int &li, int &lj) const
 Same as the previous method for the facet (c,3).
 
Bounded_side side_of_edge (const Point &p, const Edge &e, Locate_type &lt, int &li) const
 Returns a value indicating on which side of the oriented boundary of e the point p lies: More...
 
Bounded_side side_of_edge (const Point &p, Cell_handle c, Locate_type &lt, int &li) const
 Same as the previous method for edge \( (c,0,1)\).
 
bool flip (Edge e)
 
bool flip (Cell_handle c, int i, int j)
 Before flipping, these methods check that edge e=(c,i,j) is flippable (which is quite expensive). More...
 
bool flip (Facet f)
 
bool flip (Cell_handle c, int i)
 Before flipping, these methods check that facet f=(c,i) is flippable (which is quite expensive). More...
 
void flip_flippable (Edge e)
 
void flip_flippable (Cell_handle c, int i, int j)
 Should be preferred to the previous methods when the edge is known to be flippable. More...
 
void flip_flippable (Facet f)
 
void flip_flippable (Cell_handle c, int i)
 Should be preferred to the previous methods when the facet is known to be flippable. More...
 
Vertex_handle insert (const Point &p, Cell_handle start=Cell_handle())
 Inserts point p in the triangulation and returns the corresponding vertex. More...
 
Vertex_handle insert (const Point &p, Vertex_handle hint)
 Same as above but uses hint as the starting place for the search.
 
Vertex_handle insert (const Point &p, Locate_type lt, Cell_handle loc, int li, int lj)
 Inserts point p in the triangulation and returns the corresponding vertex. More...
 
std::ptrdiff_t insert (InputIterator first, InputIterator last)
 Inserts the points in the range [first,last). More...
 
Vertex_handle insert_in_cell (const Point &p, Cell_handle c)
 Inserts point p in cell c. More...
 
Vertex_handle insert_in_facet (const Point &p, const Facet &f)
 Inserts point p in facet f. More...
 
Vertex_handle insert_in_facet (const Point &p, Cell_handle c, int i)
 As above, insertion in facet (c,i). More...
 
Vertex_handle insert_in_edge (const Point &p, const Edge &e)
 Inserts p in edge e. More...
 
Vertex_handle insert_in_edge (Point p, Cell_handle c, int i, int j)
 As above, inserts p in edge \( (i, j)\) of c. More...
 
Vertex_handle insert_outside_convex_hull (const Point &p, Cell_handle c)
 The cell c must be an infinite cell containing p. More...
 
Vertex_handle insert_outside_affine_hull (const Point &p)
 p is linked to all the points, and the infinite vertex is linked to all the points (including p) to triangulate the new infinite face, so that all the points now belong to the boundary of the convex hull. More...
 
Vertex_handle insert_in_hole (Point p, CellIt cell_begin, CellIt cell_end, Cell_handle begin, int i)
 Creates a new vertex by starring a hole. More...
 
Vertex_handle insert_in_hole (Point p, CellIt cell_begin, CellIt cell_end, Cell_handle begin, int i, Vertex_handle newv)
 Same as above, except that newv will be used as the new vertex, which must have been allocated previously with e.g. create_vertex.
 
Finite_vertices_iterator finite_vertices_begin () const
 Starts at an arbitrary finite vertex. More...
 
Finite_vertices_iterator finite_vertices_end () const
 Past-the-end iterator.
 
Finite_edges_iterator finite_edges_begin () const
 Starts at an arbitrary finite edge. More...
 
Finite_edges_iterator finite_edges_end () const
 Past-the-end iterator.
 
Finite_facets_iterator finite_facets_begin () const
 Starts at an arbitrary finite facet. More...
 
Finite_facets_iterator finite_facets_end () const
 Past-the-end iterator.
 
Finite_cells_iterator finite_cells_begin () const
 Starts at an arbitrary finite cell. More...
 
Finite_cells_iterator finite_cells_end () const
 Past-the-end iterator.
 
All_vertices_iterator all_vertices_begin () const
 Starts at an arbitrary vertex. More...
 
All_vertices_iterator all_vertices_end () const
 Past-the-end iterator.
 
All_edges_iterator all_edges_begin () const
 Starts at an arbitrary edge. More...
 
All_edges_iterator all_edges_end () const
 Past-the-end iterator.
 
All_facets_iterator all_facets_begin () const
 Starts at an arbitrary facet. More...
 
All_facets_iterator all_facets_end () const
 Past-the-end iterator.
 
All_cells_iterator all_cells_begin () const
 Starts at an arbitrary cell. More...
 
All_cells_iterator all_cells_end () const
 Past-the-end iterator.
 
Point_iterator points_begin () const
 Iterates over the points of the triangulation.
 
Point_iterator points_end () const
 Past-the-end iterator.
 
Cell_circulator incident_cells (Edge e) const
 Starts at an arbitrary cell incident to e. More...
 
Cell_circulator incident_cells (Cell_handle c, int i, int j) const
 As above for edge (i,j) of c.
 
Cell_circulator incident_cells (Edge e, Cell_handle start) const
 Starts at cell start. More...
 
Cell_circulator incident_cells (Cell_handle c, int i, int j, Cell_handle start) const
 As above for edge (i,j) of c.
 
OutputIterator incident_cells (Vertex_handle v, OutputIterator cells) const
 Copies the Cell_handles of all cells incident to v to the output iterator cells. More...
 
OutputIterator incident_facets (Vertex_handle v, OutputIterator facets) const
 Copies all Facets incident to v to the output iterator facets. More...
 
bool try_lock_and_get_incident_cells (Vertex_handle v, std::vector< Cell_handle > &cells) const
 Try to lock and copy the Cell_handles of all cells incident to v into cells. More...
 
OutputIterator finite_incident_cells (Vertex_handle v, OutputIterator cells) const
 Copies the Cell_handles of all finite cells incident to v to the output iterator cells. More...
 
OutputIterator finite_incident_facets (Vertex_handle v, OutputIterator facets) const
 Copies all finite Facets incident to v to the output iterator facets. More...
 
OutputIterator incident_edges (Vertex_handle v, OutputIterator edges) const
 Copies all Edges incident to v to the output iterator edges. More...
 
OutputIterator finite_incident_edges (Vertex_handle v, OutputIterator edges) const
 Copies all finite Edges incident to v to the output iterator edges. More...
 
OutputIterator adjacent_vertices (Vertex_handle v, OutputIterator vertices) const
 Copies the Vertex_handles of all vertices adjacent to v to the output iterator vertices. More...
 
OutputIterator finite_adjacent_vertices (Vertex_handle v, OutputIterator vertices) const
 Copies the Vertex_handles of all finite vertices adjacent to v to the output iterator vertices. More...
 
size_type degree (Vertex_handle v) const
 Returns the degree of a vertex, that is, the number of incident vertices. More...
 
Facet_circulator incident_facets (Edge e) const
 Starts at an arbitrary facet incident to e. More...
 
Facet_circulator incident_facets (Cell_handle c, int i, int j) const
 As above for edge (i,j) of c.
 
Facet_circulator incident_facets (Edge e, Facet start) const
 Starts at facet start. More...
 
Facet_circulator incident_facets (Edge e, Cell_handle start, int f) const
 Starts at facet of index f in start.
 
Facet_circulator incident_facets (Cell_handle c, int i, int j, Facet start) const
 As above for edge (i,j) of c.
 
Facet_circulator incident_facets (Cell_handle c, int i, int j, Cell_handle start, int f) const
 As above for edge (i,j) of c and facet (start,f).
 
int mirror_index (Cell_handle c, int i) const
 Returns the index of c in its \( i^{th}\) neighbor. More...
 
Vertex_handle mirror_vertex (Cell_handle c, int i) const
 Returns the vertex of the \( i^{th}\) neighbor of c that is opposite to c. More...
 
Facet mirror_facet (Facet f) const
 Returns the same facet seen from the other adjacent cell.
 
bool is_valid (bool verbose=false) const
 This is a function for debugging purpose. More...
 
bool is_valid (Cell_handle c, bool verbose=false) const
 This is a function for debugging purpose. More...
 
istream & operator>> (istream &is, Triangulation_3 &t)
 Reads the underlying combinatorial triangulation from is by calling the corresponding input operator of the triangulation data structure class (note that the infinite vertex is numbered 0), and the non-combinatorial information by calling the corresponding input operators of the vertex and the cell classes (such as point coordinates), which are provided by overloading the stream operators of the vertex and cell types. More...
 
ostream & operator<< (ostream &os, const Triangulation_3 &t)
 Writes the triangulation t into os.
 
void set_lock_data_structure (Lock_data_structure *lock_ds) const
 Set the pointer to the lock data structure.
 

Constructor & Destructor Documentation

template<typename DelaunayTriangulationTraits_3 , typename TriangulationDataStructure_3 , typename LocationPolicy , typename SurjectiveLockDataStructure >
CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy, SurjectiveLockDataStructure >::Delaunay_triangulation_3 ( const DelaunayTriangulationTraits_3 traits = DelaunayTriangulationTraits_3(),
Lock_data_structure lock_ds = NULL 
)

Creates an empty Delaunay triangulation, possibly specifying a traits class traits.

lock_ds is an optional pointer to the lock data structure for parallel operations. It must be provided if concurrency is enabled.

template<typename DelaunayTriangulationTraits_3 , typename TriangulationDataStructure_3 , typename LocationPolicy , typename SurjectiveLockDataStructure >
CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy, SurjectiveLockDataStructure >::Delaunay_triangulation_3 ( const Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy, SurjectiveLockDataStructure > &  dt1)

Copy constructor.

The pointer to the lock data structure is not copied. Thus, the copy won't be concurrency-safe as long as the user has not called Triangulation_3::set_lock_data_structure.

template<typename DelaunayTriangulationTraits_3 , typename TriangulationDataStructure_3 , typename LocationPolicy , typename SurjectiveLockDataStructure >
template<class InputIterator >
CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy, SurjectiveLockDataStructure >::Delaunay_triangulation_3 ( InputIterator  first,
InputIterator  last,
const DelaunayTriangulationTraits_3 traits = DelaunayTriangulationTraits_3(),
Lock_data_structure lock_ds = NULL 
)

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

If parallelism is enabled, the points will be inserted in parallel.

Member Function Documentation

template<typename DelaunayTriangulationTraits_3 , typename TriangulationDataStructure_3 , typename LocationPolicy , typename SurjectiveLockDataStructure >
Point CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy, SurjectiveLockDataStructure >::dual ( Cell_handle  c) const

Returns the circumcenter of the four vertices of c.

Precondition
dt.dimension() \( =3\) and c is not infinite.
template<typename DelaunayTriangulationTraits_3 , typename TriangulationDataStructure_3 , typename LocationPolicy , typename SurjectiveLockDataStructure >
Object CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy, SurjectiveLockDataStructure >::dual ( Facet  f) const

Returns the dual of facet f, which is.

in dimension 3: either a segment, if the two cells incident to f are finite, or a ray, if one of them is infinite;

in dimension 2: a point.

Precondition
dt.dimension() \( \geq2\) and f is not infinite.
template<typename DelaunayTriangulationTraits_3 , typename TriangulationDataStructure_3 , typename LocationPolicy , typename SurjectiveLockDataStructure >
Line CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy, SurjectiveLockDataStructure >::dual_support ( Cell_handle  c,
int  i 
) const

returns the line supporting the dual of the facet.

Precondition
dt.dimension() \( \geq2\) and f is not infinite.
template<typename DelaunayTriangulationTraits_3 , typename TriangulationDataStructure_3 , typename LocationPolicy , typename SurjectiveLockDataStructure >
template<class OutputIteratorBoundaryFacets , class OutputIteratorCells >
std::pair<OutputIteratorBoundaryFacets, OutputIteratorCells> CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy, SurjectiveLockDataStructure >::find_conflicts ( Point  p,
Cell_handle  c,
OutputIteratorBoundaryFacets  bfit,
OutputIteratorCells  cit,
bool *  could_lock_zone = NULL 
)

Computes the conflict hole induced by p.

The starting cell (resp. facet) c must be in conflict. Then this function returns respectively in the output iterators:

  • cit: the cells (resp. facets) in conflict.
  • bfit: the facets (resp. edges) on the boundary, that is, the facets (resp. edges) (t, i) where the cell (resp. facet) t is in conflict, but t->neighbor(i) is not.
  • could_lock_zone: The optional argument could_lock_zone is used by the concurrency-safe version of the triangulation. If the pointer is not null, the algorithm will try to lock all the cells of the conflict zone, i.e. all the vertices that are inside or on the boundary of the conflict zone (as a result, the boundary cells become partially locked). If it succeeds, *could_lock_zone is true, otherwise it is false (and the returned conflict zone is only partial). In any case, the locked cells are not unlocked by the function, leaving this choice to the user.

This function can be used in conjunction with insert_in_hole() in order to decide the insertion of a point after seeing which elements of the triangulation are affected. Returns the pair composed of the resulting output iterators.

Precondition
dt.dimension() \( \geq2\), and c is in conflict with p.
template<typename DelaunayTriangulationTraits_3 , typename TriangulationDataStructure_3 , typename LocationPolicy , typename SurjectiveLockDataStructure >
template<class OutputIteratorBoundaryFacets , class OutputIteratorCells , class OutputIteratorInternalFacets >
Triple<OutputIteratorBoundaryFacets, OutputIteratorCells, OutputIteratorInternalFacets> CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy, SurjectiveLockDataStructure >::find_conflicts ( Point  p,
Cell_handle  c,
OutputIteratorBoundaryFacets  bfit,
OutputIteratorCells  cit,
OutputIteratorInternalFacets  ifit,
bool *  could_lock_zone = NULL 
)

Same as the other find_conflicts() function, except that it also computes the internal facets, i.e. the facets common to two cells which are in conflict with p.

Then this function returns respectively in the output iterators:

  • cit: the cells (resp. facets) in conflict.
  • bfit: the facets (resp. edges) on the boundary, that is, the facets (resp. edges) (t, i) where the cell (resp. facet) t is in conflict, but t->neighbor(i) is not.
  • ifit: the facets (resp. edges) inside the hole, that is, delimiting two cells (resp facets) in conflict.
  • could_lock_zone: The optional argument could_lock_zone is used by the concurrency-safe version of the triangulation. If the pointer is not null, the algorithm will try to lock all the cells of the conflict zone, i.e. all the vertices that are inside or on the boundary of the conflict zone (as a result, the boundary cells become partially locked). If it succeeds, *could_lock_zone is true, otherwise it is false (and the returned conflict zone is only partial). In any case, the locked cells are not unlocked by the function, leaving this choice to the user.

Returns the Triple composed of the resulting output iterators.

Precondition
dt.dimension() \( \geq2\), and c is in conflict with p.
template<typename DelaunayTriangulationTraits_3 , typename TriangulationDataStructure_3 , typename LocationPolicy , typename SurjectiveLockDataStructure >
Vertex_handle CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy, SurjectiveLockDataStructure >::insert ( const Point p,
Cell_handle  start = Cell_handle(),
bool *  could_lock_zone = NULL 
)

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

Similar to the insertion in a triangulation, but ensures in addition the empty sphere property of all the created faces. The optional argument start is used as a starting place for the search.

The optional argument could_lock_zone is used by the concurrency-safe version of the triangulation. If the pointer is not null, the insertion will try to lock all the cells of the conflict zone, i.e. all the vertices that are inside or on the boundary of the conflict zone. If it succeeds, *could_lock_zone is true, otherwise it is false and the return value is Vertex_handle() (the point is not inserted). In any case, the locked cells are not unlocked by the function, leaving this choice to the user.

Examples:
Triangulation_3/adding_handles_3.cpp, Triangulation_3/color.cpp, and Triangulation_3/find_conflicts_3.cpp.
template<typename DelaunayTriangulationTraits_3 , typename TriangulationDataStructure_3 , typename LocationPolicy , typename SurjectiveLockDataStructure >
Vertex_handle CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy, SurjectiveLockDataStructure >::insert ( const Point p,
Locate_type  lt,
Cell_handle  loc,
int  li,
int  lj,
bool *  could_lock_zone = NULL 
)

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 Triangulation_3::locate().

template<typename DelaunayTriangulationTraits_3 , typename TriangulationDataStructure_3 , typename LocationPolicy , typename SurjectiveLockDataStructure >
template<class PointInputIterator >
std::ptrdiff_t CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy, SurjectiveLockDataStructure >::insert ( PointInputIterator  first,
PointInputIterator  last 
)

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

Returns the number of inserted points. Note that this function is not guaranteed to insert the points following the order of PointInputIterator, as spatial_sort() is used to improve efficiency. If parallelism is enabled, the points will be inserted in parallel.

Template Parameters
PointInputIteratormust be an input iterator with the value type Point.
template<typename DelaunayTriangulationTraits_3 , typename TriangulationDataStructure_3 , typename LocationPolicy , typename SurjectiveLockDataStructure >
template<class PointWithInfoInputIterator >
std::ptrdiff_t CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy, SurjectiveLockDataStructure >::insert ( PointWithInfoInputIterator  first,
PointWithInfoInputIterator  last 
)

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

Returns the number of inserted points. Note that this function is not guaranteed to insert the points following the order of PointWithInfoInputIterator, as spatial_sort() is used to improve efficiency. If parallelism is enabled, the points will be inserted in parallel. 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, and one of the objects of type Vertex::Info will be stored in the vertex.

Precondition
Vertex must be model of the concept TriangulationVertexBaseWithInfo_3.
Template Parameters
PointWithInfoInputIteratormust be an input iterator with the value type std::pair<Point,Vertex::Info>.
template<typename DelaunayTriangulationTraits_3 , typename TriangulationDataStructure_3 , typename LocationPolicy , typename SurjectiveLockDataStructure >
bool CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy, SurjectiveLockDataStructure >::is_valid ( bool  verbose = false) const

Checks the combinatorial validity of the triangulation and the validity of its geometric embedding (see Section Representation).

Also checks that all the circumscribing spheres (resp. circles in dimension 2) of cells (resp. facets in dimension 2) are empty. When verbose is set to true, messages describing the first invalidity encountered are printed.

template<typename DelaunayTriangulationTraits_3 , typename TriangulationDataStructure_3 , typename LocationPolicy , typename SurjectiveLockDataStructure >
bool CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy, SurjectiveLockDataStructure >::is_valid ( Cell_handle  c,
bool  verbose = false 
) const

This is a function for debugging purpose.

Debugging Support

Checks the combinatorial and geometric validity of the cell (see Section Representation). Also checks that the circumscribing sphere (resp. circle in dimension 2) of cells (resp. facet in dimension 2) is empty.

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

template<typename DelaunayTriangulationTraits_3 , typename TriangulationDataStructure_3 , typename LocationPolicy , typename SurjectiveLockDataStructure >
Vertex_handle CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy, SurjectiveLockDataStructure >::move ( Vertex_handle  v,
const Point p 
)

same as above if there is no collision.

Otherwise, v is deleted and the vertex placed on p is returned.

Precondition
Vertex v must be finite.
template<typename DelaunayTriangulationTraits_3 , typename TriangulationDataStructure_3 , typename LocationPolicy , typename SurjectiveLockDataStructure >
Vertex_handle CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy, SurjectiveLockDataStructure >::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
Vertex v must be finite.
template<typename DelaunayTriangulationTraits_3 , typename TriangulationDataStructure_3 , typename LocationPolicy , typename SurjectiveLockDataStructure >
Vertex_handle CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy, SurjectiveLockDataStructure >::nearest_vertex ( Point  p,
Cell_handle  c = Cell_handle() 
)

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

The optional argument c is a hint specifying where to start the search.

Precondition
c is a cell of dt.
template<typename DelaunayTriangulationTraits_3 , typename TriangulationDataStructure_3 , typename LocationPolicy , typename SurjectiveLockDataStructure >
void CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy, SurjectiveLockDataStructure >::remove ( Vertex_handle  v)

Removes the vertex v from the triangulation.

Precondition
v is a finite vertex of the triangulation.
template<typename DelaunayTriangulationTraits_3 , typename TriangulationDataStructure_3 , typename LocationPolicy , typename SurjectiveLockDataStructure >
bool CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy, SurjectiveLockDataStructure >::remove ( Vertex_handle  v,
bool *  could_lock_zone 
)

Removes the vertex v from the triangulation.

This function is concurrency-safe if the triangulation is concurrency-safe. It will first try to lock all the cells adjacent to v. If it succeeds, *could_lock_zone is true, otherwise it is false (and the point is not removed). In any case, the locked cells are not unlocked by the function, leaving this choice to the user.

This function will try to remove v only if the removal does not decrease the dimension.

The return value is only meaningful if *could_lock_zone is true:

  • returns true if the vertex was removed
  • returns false if the vertex wasn't removed since it would decrease the dimension.
Precondition
v is a finite vertex of the triangulation.
dt.dimension() \( =3\).
template<typename DelaunayTriangulationTraits_3 , typename TriangulationDataStructure_3 , typename LocationPolicy , typename SurjectiveLockDataStructure >
template<typename InputIterator >
int CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy, SurjectiveLockDataStructure >::remove ( InputIterator  first,
InputIterator  beyond 
)

Removes the vertices specified by the iterator range [first, beyond).

The number of vertices removed is returned. If parallelism is enabled, the points will be removed in parallel. Note that if at some step, the triangulation dimension becomes lower than 3, the removal of the remaining points will go on sequentially.

Precondition
(i) all vertices of the range are finite vertices of the triangulation; and (ii) no vertices are repeated in the range.
Template Parameters
InputIteratormust be an input iterator with value type Vertex_handle.
template<typename DelaunayTriangulationTraits_3 , typename TriangulationDataStructure_3 , typename LocationPolicy , typename SurjectiveLockDataStructure >
template<typename InputIterator >
int CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy, SurjectiveLockDataStructure >::remove_cluster ( InputIterator  first,
InputIterator  beyond 
)

This function has exactly the same result and the same preconditions as remove(first, beyond).

The difference is in the implementation and efficiency. This version does not re-triangulate the hole after each point removal but only after removing all vertices. This is more efficient if (and only if) the removed points are organized in a small number of connected components of the Delaunay triangulation. Another difference is that there is no parallel version of this function.

Template Parameters
InputIteratormust be an input iterator with value type Vertex_handle.
template<typename DelaunayTriangulationTraits_3 , typename TriangulationDataStructure_3 , typename LocationPolicy , typename SurjectiveLockDataStructure >
Bounded_side CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy, SurjectiveLockDataStructure >::side_of_circle ( const Facet f,
const Point p 
) const

Returns a value indicating on which side of the circumscribed circle of f the point p lies.

More precisely, it returns:

  • in dimension 3:
  • For a finite facet, ON_BOUNDARY if p lies on the circle, ON_UNBOUNDED_SIDE when it lies in the exterior of the disk, ON_BOUNDED_SIDE when it lies in its interior.
  • For an infinite facet, it considers the plane defined by the finite facet of the same cell, and does the same as in dimension 2 in this plane.
  • in dimension 2:
  • For a finite facet, ON_BOUNDARY if p lies on the circle, ON_UNBOUNDED_SIDE when it lies in the exterior of the disk, ON_BOUNDED_SIDE when it lies in its interior.
  • For an infinite facet, ON_BOUNDARY if the point lies on the finite edge of f (endpoints included), ON_BOUNDED_SIDE for a point in the open half plane defined by f and not containing any other point of the triangulation, ON_UNBOUNDED_SIDE elsewhere.
    Precondition
    dt.dimension() \( \geq2\) and in dimension 3, p is coplanar with f.
template<typename DelaunayTriangulationTraits_3 , typename TriangulationDataStructure_3 , typename LocationPolicy , typename SurjectiveLockDataStructure >
Bounded_side CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy, SurjectiveLockDataStructure >::side_of_sphere ( Cell_handle  c,
const Point p 
) const

Returns a value indicating on which side of the circumscribed sphere of c the point p lies.

More precisely, it returns:

  • ON_BOUNDED_SIDE if p is inside the sphere. For an infinite cell this means that p lies strictly either in the half space limited by its finite facet and not containing any other point of the triangulation, or in the interior of the disk circumscribing the finite facet.
  • ON_BOUNDARY if p on the boundary of the sphere. For an infinite cell this means that p lies on the circle circumscribing the finite facet.
  • ON_UNBOUNDED_SIDE if p lies outside the sphere. For an infinite cell this means that p does not satisfy either of the two previous conditions.
    Precondition
    dt.dimension() \( =3\).
template<typename DelaunayTriangulationTraits_3 , typename TriangulationDataStructure_3 , typename LocationPolicy , typename SurjectiveLockDataStructure >
template<class OutputIterator >
OutputIterator CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy, SurjectiveLockDataStructure >::vertices_in_conflict ( Point  p,
Cell_handle  c,
OutputIterator  res 
)
Deprecated:
This function is renamed vertices_on_conflict_zone_boundary since CGAL-3.8.
template<typename DelaunayTriangulationTraits_3 , typename TriangulationDataStructure_3 , typename LocationPolicy , typename SurjectiveLockDataStructure >
template<class OutputIterator >
OutputIterator CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy, SurjectiveLockDataStructure >::vertices_on_conflict_zone_boundary ( Point  p,
Cell_handle  c,
OutputIterator  res 
)

Similar to find_conflicts(), but reports the vertices which are on the boundary of the conflict hole of p, in the output iterator res.

Returns the resulting output iterator.

Precondition
dt.dimension() \( \geq2\), and c is in conflict with p.