\( \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.12 - 3D Triangulations
CGAL::Delaunay_triangulation_3< Traits, TDS, LP, SLDS > Class Template Reference

#include <CGAL/Delaunay_triangulation_3.h>

Inherits from

CGAL::Triangulation_3< Traits, Delaunay_triangulation_3< Traits, TDS, LP >::Triangulation_data_structure, SLDS >.

Definition

The class Delaunay_triangulation_3 represents a three-dimensional Delaunay triangulation.

Template Parameters
Traitsis the geometric traits class and must be a model of DelaunayTriangulationTraits_3.
TDSis the triangulation data structure and must be a model of TriangulationDataStructure_3. Default may be used, with default type Triangulation_data_structure_3<Triangulation_vertex_base_3<Traits>, Delaunay_triangulation_cell_base_3<Traits> >. Any custom type can be used instead of Triangulation_vertex_base_3 and Delaunay_triangulation_cell_base_3, provided that they are models of the concepts TriangulationVertexBase_3 and DelaunayTriangulationCellBase_3, respectively.
LPis a tag which must be a Location_policy<Tag>: either CGAL::Fast_location or CGAL::Compact_location. CGAL::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 CGAL::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 TDS parameter is satisfactory (this is achieved 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.
SLDSis an optional parameter to specify the type of the spatial lock data structure. It must be a model of the SurjectiveLockDataStructure concept, with Object being a Point (as defined below). It is only used if the triangulation data structure used is concurrency-safe (i.e. when TDS::Concurrency_tag is CGAL::Parallel_tag). 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 TDS::Concurrency_tag is CGAL::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::Triangulation_3
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 Traits Geom_traits
 
typedef TDS Triangulation_data_structure
 
typedef LP Location_policy
 
typedef SLDS Lock_data_structure
 

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

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

Creation

 Delaunay_triangulation_3 (const Geom_traits &traits=Geom_traits(), 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 Geom_traits &traits=Geom_traits(), 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 Geom_traits &traits=Geom_traits())
 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

- Static Public Member Functions inherited from CGAL::Triangulation_utils_3
static unsigned int next_around_edge (unsigned int i, unsigned int j)
 
static int vertex_triple_index (const int i, const int j)
 
static unsigned int ccw (unsigned int i)
 
static unsigned int cw (unsigned int i)
 

Constructor & Destructor Documentation

◆ Delaunay_triangulation_3() [1/3]

template<typename Traits, typename TDS, typename LP, typename SLDS>
CGAL::Delaunay_triangulation_3< Traits, TDS, LP, SLDS >::Delaunay_triangulation_3 ( const Geom_traits traits = Geom_traits(),
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.

◆ Delaunay_triangulation_3() [2/3]

template<typename Traits, typename TDS, typename LP, typename SLDS>
CGAL::Delaunay_triangulation_3< Traits, TDS, LP, SLDS >::Delaunay_triangulation_3 ( const Delaunay_triangulation_3< Traits, TDS, LP, SLDS > &  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().

◆ Delaunay_triangulation_3() [3/3]

template<typename Traits, typename TDS, typename LP, typename SLDS>
template<class InputIterator >
CGAL::Delaunay_triangulation_3< Traits, TDS, LP, SLDS >::Delaunay_triangulation_3 ( InputIterator  first,
InputIterator  last,
const Geom_traits traits = Geom_traits(),
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

◆ dual() [1/2]

template<typename Traits, typename TDS, typename LP, typename SLDS>
Point CGAL::Delaunay_triangulation_3< Traits, TDS, LP, SLDS >::dual ( Cell_handle  c) const

Returns the circumcenter of the four vertices of c.

Precondition
dt.dimension() \( =3\) and c is not infinite.

◆ dual() [2/2]

template<typename Traits, typename TDS, typename LP, typename SLDS>
Object CGAL::Delaunay_triangulation_3< Traits, TDS, LP, SLDS >::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.

◆ dual_support()

template<typename Traits, typename TDS, typename LP, typename SLDS>
Line CGAL::Delaunay_triangulation_3< Traits, TDS, LP, SLDS >::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.

◆ find_conflicts() [1/2]

template<typename Traits, typename TDS, typename LP, typename SLDS>
template<class OutputIteratorBoundaryFacets , class OutputIteratorCells >
std::pair<OutputIteratorBoundaryFacets, OutputIteratorCells> CGAL::Delaunay_triangulation_3< Traits, TDS, LP, SLDS >::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.

◆ find_conflicts() [2/2]

template<typename Traits, typename TDS, typename LP, typename SLDS>
template<class OutputIteratorBoundaryFacets , class OutputIteratorCells , class OutputIteratorInternalFacets >
Triple<OutputIteratorBoundaryFacets, OutputIteratorCells, OutputIteratorInternalFacets> CGAL::Delaunay_triangulation_3< Traits, TDS, LP, SLDS >::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.

◆ insert() [1/4]

template<typename Traits, typename TDS, typename LP, typename SLDS>
Vertex_handle CGAL::Delaunay_triangulation_3< Traits, TDS, LP, SLDS >::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.

◆ insert() [2/4]

template<typename Traits, typename TDS, typename LP, typename SLDS>
Vertex_handle CGAL::Delaunay_triangulation_3< Traits, TDS, LP, SLDS >::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().

◆ insert() [3/4]

template<typename Traits, typename TDS, typename LP, typename SLDS>
template<class PointInputIterator >
std::ptrdiff_t CGAL::Delaunay_triangulation_3< Traits, TDS, LP, SLDS >::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.

◆ insert() [4/4]

template<typename Traits, typename TDS, typename LP, typename SLDS>
template<class PointWithInfoInputIterator >
std::ptrdiff_t CGAL::Delaunay_triangulation_3< Traits, TDS, LP, SLDS >::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>.

◆ is_valid() [1/2]

template<typename Traits, typename TDS, typename LP, typename SLDS>
bool CGAL::Delaunay_triangulation_3< Traits, TDS, LP, SLDS >::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.

◆ is_valid() [2/2]

template<typename Traits, typename TDS, typename LP, typename SLDS>
bool CGAL::Delaunay_triangulation_3< Traits, TDS, LP, SLDS >::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.

◆ move()

template<typename Traits, typename TDS, typename LP, typename SLDS>
Vertex_handle CGAL::Delaunay_triangulation_3< Traits, TDS, LP, SLDS >::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.

◆ move_if_no_collision()

template<typename Traits, typename TDS, typename LP, typename SLDS>
Vertex_handle CGAL::Delaunay_triangulation_3< Traits, TDS, LP, SLDS >::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.

◆ nearest_vertex()

template<typename Traits, typename TDS, typename LP, typename SLDS>
Vertex_handle CGAL::Delaunay_triangulation_3< Traits, TDS, LP, SLDS >::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.

◆ remove() [1/3]

template<typename Traits, typename TDS, typename LP, typename SLDS>
void CGAL::Delaunay_triangulation_3< Traits, TDS, LP, SLDS >::remove ( Vertex_handle  v)

Removes the vertex v from the triangulation.

Precondition
v is a finite vertex of the triangulation.

◆ remove() [2/3]

template<typename Traits, typename TDS, typename LP, typename SLDS>
bool CGAL::Delaunay_triangulation_3< Traits, TDS, LP, SLDS >::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\).

◆ remove() [3/3]

template<typename Traits, typename TDS, typename LP, typename SLDS>
template<typename InputIterator >
int CGAL::Delaunay_triangulation_3< Traits, TDS, LP, SLDS >::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.

◆ remove_cluster()

template<typename Traits, typename TDS, typename LP, typename SLDS>
template<typename InputIterator >
int CGAL::Delaunay_triangulation_3< Traits, TDS, LP, SLDS >::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.

◆ side_of_circle()

template<typename Traits, typename TDS, typename LP, typename SLDS>
Bounded_side CGAL::Delaunay_triangulation_3< Traits, TDS, LP, SLDS >::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.

◆ side_of_sphere()

template<typename Traits, typename TDS, typename LP, typename SLDS>
Bounded_side CGAL::Delaunay_triangulation_3< Traits, TDS, LP, SLDS >::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\).

◆ vertices_in_conflict()

template<typename Traits, typename TDS, typename LP, typename SLDS>
template<class OutputIterator >
OutputIterator CGAL::Delaunay_triangulation_3< Traits, TDS, LP, SLDS >::vertices_in_conflict ( Point  p,
Cell_handle  c,
OutputIterator  res 
)
Deprecated:
This function is renamed vertices_on_conflict_zone_boundary since CGAL-3.8.

◆ vertices_on_conflict_zone_boundary()

template<typename Traits, typename TDS, typename LP, typename SLDS>
template<class OutputIterator >
OutputIterator CGAL::Delaunay_triangulation_3< Traits, TDS, LP, SLDS >::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.