CGAL 6.0.1 - 3D Triangulations
Loading...
Searching...
No Matches
CGAL::Regular_triangulation_3< Traits, TDS, SLDS > Class Template Reference

#include <CGAL/Regular_triangulation_3.h>

Inherits from

CGAL::Triangulation_3< Traits, TDS, SLDS >.

Definition

template<typename Traits, typename TDS, typename SLDS>
class CGAL::Regular_triangulation_3< Traits, TDS, SLDS >

Let \( {S}^{(w)}\) be a set of weighted points in \( \mathbb{R}^3\).

Let \( {p}^{(w)}=(p,w_p), p\in\mathbb{R}^3, w_p\in\mathbb{R}\) and \( {z}^{(w)}=(z,w_z), z\in\mathbb{R}^3, w_z\in\mathbb{R}\) be two weighted points. A weighted point \( {p}^{(w)}=(p,w_p)\) can also be seen as a sphere of center \( p\) and radius \( \sqrt{w_p}\). The power product (or power distance ) between \( {p}^{(w)}\) and \( {z}^{(w)}\) is defined as

\[ \Pi({p}^{(w)},{z}^{(w)}) = {\|{p-z}\|^2-w_p-w_z} \]

where \( \|{p-z}\|\) is the Euclidean distance between \( p\) and \( z\). \( {p}^{(w)}\) and \( {z}^{(w)}\) are said to be orthogonal if \( \Pi{({p}^{(w)}-{z}^{(w)})} = 0\) (see Figure 46.2).

Four weighted points have a unique common orthogonal weighted point called the power sphere. A sphere \( {z}^{(w)}\) is said to be regular if \( \forall {p}^{(w)}\in{S}^{(w)}, \Pi{({p}^{(w)}-{z}^{(w)})}\geq 0\).

A triangulation of \( {S}^{(w)}\) is regular if the power spheres of all simplices are regular.

Template Parameters
Traitsis the geometric traits class, and must be a model of RegularTriangulationTraits_3
TDSis the triangulation data structure and must be a model of TriangulationDataStructure_3. Default may be used with default value Triangulation_data_structure_3<Regular_triangulation_vertex_base_3<Traits>, Regular_triangulation_cell_base_3<Traits> >. Any custom type can be used instead of Regular_triangulation_vertex_base_3 and Regular_triangulation_cell_base_3, provided that they are models of the concepts RegularTriangulationVertexBase_3 and RegularTriangulationCellBase_3, respectively.
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. It is only used if the triangulation data structure used is concurrency-safe (i.e. when TDS::Concurrency_tag is 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 SLDS instance via the constructor or Triangulation_3::set_lock_data_structure.

If TDS::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::Triangulation_3
CGAL::Delaunay_triangulation_3
Examples
Triangulation_3/parallel_insertion_and_removal_in_regular_3.cpp, Triangulation_3/regular_3.cpp, and Triangulation_3/regular_with_info_3.cpp.

Types

typedef Traits::Point_3 Bare_point
 The type for points p of weighted points \( {p}^{(w)}=(p,w_p)\).
 
typedef Traits::Weighted_point_3 Weighted_point
 The type for weighted points.
 
typedef SLDS Lock_data_structure
 

Creation

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

Insertion

The following methods, which already exist in Triangulation_3, are overloaded to ensure the property that all power spheres are regular.

The following method allows one to insert several points.

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

The following methods, which already exist in Triangulation_3, are overloaded to ensure that hidden points are well created and maintained.

template<class CellIt >
Vertex_handle insert_in_hole (const Weighted_point &p, CellIt cell_begin, CellIt cell_end, Cell_handle begin, int i)
 Creates a new vertex by starring a hole.
 
template<class CellIt >
Vertex_handle insert_in_hole (const Weighted_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.
 

Removal

void remove (Vertex_handle v)
 Removes the vertex v from the triangulation.
 
bool remove (Vertex_handle v, bool *could_lock_zone)
 Removes the vertex v from the triangulation.
 
template<typename InputIterator >
int remove (InputIterator first, InputIterator beyond)
 Removes the vertices specified by the iterator range [first, beyond).
 

Queries

Let us remark that \( \Pi({p}^{(w)}-{z}^{(w)}) > 0 \) is equivalent to p lies outside the sphere with center z and radius \( \sqrt{w_p^2+w_z^2}\).

This remark helps provide an intuition about the following predicates.

side_of_power_circle
Bounded_side side_of_power_sphere (Cell_handle c, const Weighted_point &p) const
 Returns the position of the weighted point \( p\) with respect to the power sphere of c.
 
Bounded_side side_of_power_circle (const Facet &f, const Weighted_point &p) const
 Returns the position of the point p with respect to the power circle of f.
 
Bounded_side side_of_power_circle (Cell_handle c, int i, const Weighted_point &p) const
 Same as the previous method for facet i of cell c.
 
Bounded_side side_of_power_segment (Cell_handle c, const Weighted_point &p) const
 In dimension 1, returns.
 
Vertex_handle nearest_power_vertex (const Bare_point &p, Cell_handle c=Cell_handle())
 Returns the vertex of the triangulation which is nearest to \( p\) with respect to the power distance.
 
Vertex_handle nearest_power_vertex_in_cell (const Bare_point &p, Cell_handle c)
 Returns the vertex of the cell c that is nearest to \( p\) with respect to the power distance.
 

A weighted point p is said to be in conflict with a cell c in dimension 3 (resp. with a facet f in dimension 2) if it has a negative power distance to the power sphere of c (resp. to the power circle of f). The set of cells (resp. facets in dimension 2) which are in conflict with p is connected.

template<class OutputIteratorBoundaryFacets , class OutputIteratorCells , class OutputIteratorInternalFacets >
Triple< OutputIteratorBoundaryFacets, OutputIteratorCells, OutputIteratorInternalFacets > find_conflicts (const Weighted_point p, Cell_handle c, OutputIteratorBoundaryFacets bfit, OutputIteratorCells cit, OutputIteratorInternalFacets ifit, bool *could_lock_zone=nullptr, const Facet *this_facet_must_be_in_the_cz=nullptr, bool *the_facet_is_in_its_cz=nullptr)
 Compute the conflicts with p.
 
template<class OutputIterator >
OutputIterator vertices_in_conflict (const Weighted_point &p, Cell_handle c, OutputIterator res)
 
template<class OutputIterator >
OutputIterator vertices_on_conflict_zone_boundary (const Weighted_point &p, Cell_handle c, OutputIterator res)
 Similar to find_conflicts(), but reports the vertices which are on the boundary of the conflict zone of p, in the output iterator res.
 
template<class OutputIterator >
OutputIterator vertices_inside_conflict_zone (const Weighted_point &p, Cell_handle c, OutputIterator res)
 Similar to find_conflicts(), but reports the vertices which are in the interior of the conflict zone of p, in the output iterator res.
 

In the weighted setting, a face (cell, facet, edge or vertex) is said to be a Gabriel face iff the smallest sphere orthogonal to the weighted points associated to its vertices, has a positive power product with the weighted point of any other vertex of the triangulation.

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

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)
 
bool is_Gabriel (Vertex_handle v)
 

Power Diagram

CGAL offers several functionalities to display the power 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).

Bare_point dual (Cell_handle c) const
 Returns the weighted circumcenter of the four vertices of c.
 
Object dual (Facet f) const
 Returns the dual of facet f, which is.
 
Object dual (Cell_handle c, int i) const
 same as the previous method for facet (c,i).
 
template<class Stream >
Stream & draw_dual (Stream &os)
 Sends the set of duals to all the facets of rt into os.
 

Checking

bool is_valid (bool verbose=false) const
 This is a function for debugging purpose.
 

Additional Inherited Members

- Public Types inherited from CGAL::Triangulation_3< Traits, TDS, SLDS >
enum  Locate_type {
  VERTEX =0 , EDGE , FACET , CELL ,
  OUTSIDE_CONVEX_HULL , OUTSIDE_AFFINE_HULL
}
 The enum Locate_type is defined by Triangulation_3 to specify which case occurs when locating a point in the triangulation. More...
 
typedef Traits Geom_traits
 
typedef TDS Triangulation_data_structure
 
typedef SLDS Lock_data_structure
 
typedef Triangulation_data_structure::Vertex::Point Point
 
typedef Geom_traits::Segment_3 Segment
 
typedef Geom_traits::Triangle_3 Triangle
 
typedef Geom_traits::Tetrahedron_3 Tetrahedron
 
typedef Triangulation_data_structure::Vertex Vertex
 
typedef Triangulation_data_structure::Cell Cell
 
typedef Triangulation_data_structure::Facet Facet
 
typedef Triangulation_data_structure::Edge Edge
 
typedef Triangulation_data_structure::Concurrency_tag Concurrency_tag
 Concurrency tag (from the TDS).
 
typedef Triangulation_data_structure::Vertex_handle Vertex_handle
 handle to a vertex
 
typedef 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 Triangulation_data_structure::size_type size_type
 Size type (an unsigned integral type)
 
typedef Triangulation_data_structure::difference_type difference_type
 Difference type (a signed integral type)
 
typedef Triangulation_data_structure::Cell_iterator All_cells_iterator
 iterator over cells
 
typedef Triangulation_data_structure::Facet_iterator All_facets_iterator
 iterator over facets
 
typedef Triangulation_data_structure::Edge_iterator All_edges_iterator
 iterator over edges
 
typedef 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 Triangulation_data_structure::Cell_circulator Cell_circulator
 circulator over all cells incident to a given edge
 
typedef Triangulation_data_structure::Facet_circulator Facet_circulator
 circulator over all facets incident to a given edge
 
typedef unspecified_type Segment_cell_iterator
 iterator over cells intersected by a line segment.
 
typedef unspecified_type Segment_simplex_iterator
 iterator over simplices intersected by a line segment.
 
typedef Iterator_range< unspecified_typeAll_cell_handles
 range type for iterating over all cell handles (including infinite cells), with a nested type iterator that has as value type Cell_handle.
 
typedef Iterator_range< All_facets_iteratorAll_facets
 range type for iterating over facets.
 
typedef Iterator_range< All_edges_iteratorAll_edges
 range type for iterating over edges.
 
typedef Iterator_range< unspecified_typeAll_vertex_handles
 range type for iterating over all vertex handles, with a nested type iterator that has as value type Vertex_handle.
 
typedef Iterator_range< unspecified_typeFinite_cell_handles
 range type for iterating over finite cell handles, with a nested type iterator that has as value type Cell_handle.
 
typedef Iterator_range< Finite_facets_iteratorFinite_facets
 range type for iterating over finite facets.
 
typedef Iterator_range< Finite_edges_iteratorFinite_edges
 range type for iterating over finite edges.
 
typedef Iterator_range< unspecified_typeFinite_vertex_handles
 range type for iterating over finite vertex handles, with a nested type iterator that has as value type Vertex_handle.
 
typedef Iterator_range< unspecified_typePoints
 range type for iterating over the points of the finite vertices.
 
typedef Iterator_range< unspecified_typeSegment_traverser_cell_handles
 range type for iterating over the cells intersected by a line segment.
 
typedef Iterator_range< Segment_simplex_iteratorSegment_traverser_simplices
 range type for iterating over the simplices intersected by a line segment.
 
- Public Member Functions inherited from CGAL::Triangulation_3< Traits, TDS, SLDS >
 Triangulation_3 (const Geom_traits &traits=Geom_traits(), Lock_data_structure *lock_ds=nullptr)
 Introduces a triangulation t having only one vertex which is the infinite vertex.
 
 Triangulation_3 (Lock_data_structure *lock_ds=nullptr, const Geom_traits &traits=Geom_traits())
 Same as the previous one, but with parameters in reverse order.
 
 Triangulation_3 (const Triangulation_3 &tr)
 Copy constructor.
 
template<class InputIterator >
 Triangulation_3 (InputIterator first, InputIterator last, const Geom_traits &traits=Geom_traits(), Lock_data_structure *lock_ds=nullptr)
 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.
 
void swap (Triangulation_3 &tr)
 The triangulations tr and t are swapped.
 
void clear ()
 Deletes all finite vertices and all cells of t.
 
template<class GT , class Tds1 , class Tds2 >
bool operator== (const Triangulation_3< GT, Tds1 > &t1, const Triangulation_3< GT, Tds2 > &t2)
 Equality operator.
 
template<class GT , class Tds1 , class Tds2 >
bool operator!= (const Triangulation_3< GT, Tds1 > &t1, const Triangulation_3< GT, Tds2 > &t2)
 The opposite of operator==.
 
const Geom_traitsgeom_traits () const
 Returns a const reference to the geometric traits object.
 
const Triangulation_data_structuretds () const
 Returns a const reference to the triangulation data structure.
 
Triangulation_data_structuretds ()
 Returns a reference to the triangulation data structure.
 
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.
 
Cell_handle infinite_cell () const
 Returns a cell incident to the infinite vertex.
 
size_type number_of_facets () const
 The number of facets.
 
size_type number_of_edges () const
 The number of edges.
 
size_type number_of_finite_cells () const
 The number of finite cells.
 
size_type number_of_finite_facets () const
 The number of finite facets.
 
size_type number_of_finite_edges () const
 The number of finite edges.
 
Tetrahedron tetrahedron (Cell_handle c) const
 Returns the tetrahedron formed by the four vertices of c.
 
Triangle triangle (Cell_handle c, int i) const
 Returns the triangle formed by the three vertices of facet (c,i).
 
Triangle triangle (const Facet &f) const
 Same as the previous method for facet f.
 
Segment segment (const Edge &e) const
 Returns the line segment formed by the vertices of e.
 
Segment segment (Cell_handle c, int i, int j) const
 Same as the previous method for edge (c,i,j).
 
const Pointpoint (Cell_handle c, int i) const
 Returns the point given by vertex i of cell c.
 
const Pointpoint (Vertex_handle v) const
 Same as the previous method for vertex v.
 
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.
 
bool is_infinite (Cell_handle c, int i) const
 true, iff the facet i of cell c is incident to the infinite vertex.
 
bool is_infinite (const Facet &f) const
 true iff facet f is incident to the infinite vertex.
 
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.
 
bool is_infinite (const Edge &e) const
 true iff edge e is incident to the infinite vertex.
 
bool is_vertex (const Point &p, Vertex_handle &v) const
 Tests whether p is a vertex of t by locating p in the triangulation.
 
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.
 
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.
 
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.
 
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.
 
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.
 
bool has_vertex (Cell_handle c, int i, Vertex_handle v, int &j) const
 Same for facet (c,i).
 
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:
 
Cell_handle locate (const Point &query, Cell_handle start=Cell_handle(), bool *could_lock_zone=nullptr) const
 If the point query lies inside the convex hull of the points, the cell that contains the query in its interior is returned.
 
Cell_handle locate (const Point &query, Vertex_handle hint, bool *could_lock_zone=nullptr) 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.
 
Cell_handle locate (const Point &query, Locate_type &lt, int &li, int &lj, Cell_handle start=Cell_handle(), bool *could_lock_zone=nullptr) 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.
 
Cell_handle locate (const Point &query, Locate_type &lt, int &li, int &lj, Vertex_handle hint, bool *could_lock_zone=nullptr) const
 Same as above but uses hint as the starting place for the search.
 
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.
 
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:
 
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:
 
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).
 
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.
 
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).
 
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.
 
Vertex_handle insert (const Point &p, Cell_handle start=Cell_handle())
 Inserts the point p in the triangulation and returns the corresponding vertex.
 
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 the point p in the triangulation and returns the corresponding vertex.
 
template<class PointInputIterator >
std::ptrdiff_t insert (PointInputIterator first, PointInputIterator last)
 Inserts the points in the range [first,last) in the given order, and returns the number of inserted points.
 
template<class PointWithInfoInputIterator >
std::ptrdiff_t insert (PointWithInfoInputIterator first, PointWithInfoInputIterator last)
 Inserts the points in the iterator range [first,last) in the given order, and returns the number of inserted points.
 
Vertex_handle insert_in_cell (const Point &p, Cell_handle c)
 Inserts the point p in the cell c.
 
Vertex_handle insert_in_facet (const Point &p, const Facet &f)
 Inserts the point p in the facet f.
 
Vertex_handle insert_in_facet (const Point &p, Cell_handle c, int i)
 As above, insertion in the facet (c,i).
 
Vertex_handle insert_in_edge (const Point &p, const Edge &e)
 Inserts p in the edge e.
 
Vertex_handle insert_in_edge (const Point &p, Cell_handle c, int i, int j)
 As above, inserts p in the edge \( (i, j)\) of c.
 
Vertex_handle insert_outside_convex_hull (const Point &p, Cell_handle c)
 The cell c must be an infinite cell containing p.
 
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.
 
template<class CellIt >
Vertex_handle insert_in_hole (const Point &p, CellIt cell_begin, CellIt cell_end, Cell_handle begin, int i)
 Creates a new vertex by starring a hole.
 
template<class CellIt >
Vertex_handle insert_in_hole (const 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.
 
Finite_vertices_iterator finite_vertices_end () const
 Past-the-end iterator.
 
Finite_edges_iterator finite_edges_begin () const
 Starts at an arbitrary finite edge.
 
Finite_edges_iterator finite_edges_end () const
 Past-the-end iterator.
 
Finite_facets_iterator finite_facets_begin () const
 Starts at an arbitrary finite facet.
 
Finite_facets_iterator finite_facets_end () const
 Past-the-end iterator.
 
Finite_cells_iterator finite_cells_begin () const
 Starts at an arbitrary finite cell.
 
Finite_cells_iterator finite_cells_end () const
 Past-the-end iterator.
 
All_vertices_iterator all_vertices_begin () const
 Starts at an arbitrary vertex.
 
All_vertices_iterator all_vertices_end () const
 Past-the-end iterator.
 
All_edges_iterator all_edges_begin () const
 Starts at an arbitrary edge.
 
All_edges_iterator all_edges_end () const
 Past-the-end iterator.
 
All_facets_iterator all_facets_begin () const
 Starts at an arbitrary facet.
 
All_facets_iterator all_facets_end () const
 Past-the-end iterator.
 
All_cells_iterator all_cells_begin () const
 Starts at an arbitrary cell.
 
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.
 
All_cell_handles all_cell_handles () const
 returns a range of iterators over all cells (even the infinite cells).
 
All_facets all_facets () const
 returns a range of iterators starting at an arbitrary facet.
 
All_edges all_edges () const
 returns a range of iterators starting at an arbitrary edge.
 
All_vertex_handles all_vertex_handles () const
 returns a range of iterators over all vertices (even the infinite one).
 
Finite_cell_handles finite_cell_handles () const
 returns a range of iterators over finite cells.
 
Finite_facets finite_facets () const
 returns a range of iterators starting at an arbitrary facet.
 
Finite_edges finite_edges () const
 returns a range of iterators starting at an arbitrary edge.
 
Finite_vertex_handles finite_vertex_handles () const
 returns a range of iterators over finite vertices.
 
Points points () const
 returns a range of iterators over the points of finite vertices.
 
std::array< Vertex_handle, 2 > vertices (const Edge &e) const
 returns an array containing the vertices of e, in the order of their indices e.second and e.third in the cell e.first.
 
std::array< Vertex_handle, 3 > vertices (const Facet &f) const
 returns an array containing the vertices of f, in counterclockwise order on the face of f.first opposite to vertex f.first->vertex(f.second).
 
std::array< Vertex_handle, 4 > vertices (const Cell_handle c) const
 returns an array containing the vertices of c, in the same order as the indices in c.
 
Segment_traverser_cell_handles segment_traverser_cell_handles () const
 returns a range of iterators over the cells intersected by a line segment
 
Segment_traverser_simplices segment_traverser_simplices () const
 returns a range of iterators over the simplices intersected by a line segment
 
Segment_cell_iterator segment_traverser_cells_begin (Vertex_handle vs, Vertex_handle vt) const
 returns the iterator that allows to visit the cells intersected by the line segment [vs,vt].
 
Segment_cell_iterator segment_traverser_cells_begin (const Point &ps, const Point &pt, Cell_handle hint=Cell_handle()) const
 returns the iterator that allows to visit the cells intersected by the line segment [ps, pt].
 
Segment_cell_iterator segment_traverser_cells_end () const
 returns the past-the-end iterator over the intersected cells.
 
Segment_simplex_iterator segment_traverser_simplices_begin (Vertex_handle vs, Vertex_handle vt) const
 returns the iterator that allows to visit the simplices intersected by the line segment [vs,vt].
 
Segment_simplex_iterator segment_traverser_simplices_begin (const Point &ps, const Point &pt, Cell_handle hint=Cell_handle()) const
 returns the iterator that allows to visit the simplices intersected by the line segment [ps,pt].
 
Segment_simplex_iterator segment_traverser_simplices_end () const
 returns the past-the-end iterator over the intersected simplices.
 
Cell_circulator incident_cells (Edge e) const
 Starts at an arbitrary cell incident to e.
 
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.
 
Cell_circulator incident_cells (Cell_handle c, int i, int j, Cell_handle start) const
 As above for edge (i,j) of c.
 
Facet_circulator incident_facets (Edge e) const
 Starts at an arbitrary facet incident to e.
 
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.
 
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).
 
template<class OutputIterator >
OutputIterator incident_cells (Vertex_handle v, OutputIterator cells) const
 Copies the Cell_handles of all cells incident to v to the output iterator cells.
 
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.
 
template<class OutputIterator >
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.
 
template<class OutputIterator >
OutputIterator incident_facets (Vertex_handle v, OutputIterator facets) const
 Copies all Facets incident to v to the output iterator facets.
 
template<class OutputIterator >
OutputIterator finite_incident_facets (Vertex_handle v, OutputIterator facets) const
 Copies all finite Facets incident to v to the output iterator facets.
 
template<class OutputIterator >
OutputIterator incident_edges (Vertex_handle v, OutputIterator edges) const
 Copies all Edges incident to v to the output iterator edges.
 
template<class OutputIterator >
OutputIterator finite_incident_edges (Vertex_handle v, OutputIterator edges) const
 Copies all finite Edges incident to v to the output iterator edges.
 
template<class OutputIterator >
OutputIterator adjacent_vertices (Vertex_handle v, OutputIterator vertices) const
 Copies the Vertex_handles of all vertices adjacent to v to the output iterator vertices.
 
template<class OutputIterator >
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.
 
size_type degree (Vertex_handle v) const
 Returns the degree of a vertex, that is, the number of incident vertices.
 
int mirror_index (Cell_handle c, int i) const
 Returns the index of c in its \( i^{th}\) neighbor.
 
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.
 
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.
 
bool is_valid (Cell_handle c, bool verbose=false) const
 This is a function for debugging purpose.
 
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.
 
ostream & operator<< (ostream &os, const Triangulation_3 &t)
 Writes the triangulation t into os.
 
template<typename Tr_src , typename ConvertVertex , typename ConvertCell >
std::istream & file_input (std::istream &is, ConvertVertex convert_vertex=ConvertVertex(), ConvertCell convert_cell=ConvertCell())
 The triangulation streamed in is, of original type Tr_src, is written into the triangulation.
 
void set_lock_data_structure (Lock_data_structure *lock_ds) const
 Set the pointer to the lock data structure.
 
- 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

◆ Regular_triangulation_3() [1/3]

template<typename Traits , typename TDS , typename SLDS >
CGAL::Regular_triangulation_3< Traits, TDS, SLDS >::Regular_triangulation_3 ( const Traits &  traits = Traits(),
Lock_data_structure lock_ds = nullptr 
)

Creates an empty regular 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.

◆ Regular_triangulation_3() [2/3]

template<typename Traits , typename TDS , typename SLDS >
CGAL::Regular_triangulation_3< Traits, TDS, SLDS >::Regular_triangulation_3 ( const Regular_triangulation_3< Traits, TDS, SLDS > &  rt1)

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.

◆ Regular_triangulation_3() [3/3]

template<typename Traits , typename TDS , typename SLDS >
template<class InputIterator >
CGAL::Regular_triangulation_3< Traits, TDS, SLDS >::Regular_triangulation_3 ( InputIterator  first,
InputIterator  last,
const Traits &  traits = Traits(),
Lock_data_structure lock_ds = nullptr 
)

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.

Template Parameters
InputIteratormust be an input iterator with value type Weighted_point .

Member Function Documentation

◆ dual() [1/2]

template<typename Traits , typename TDS , typename SLDS >
Bare_point CGAL::Regular_triangulation_3< Traits, TDS, SLDS >::dual ( Cell_handle  c) const

Returns the weighted circumcenter of the four vertices of c.

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

◆ dual() [2/2]

template<typename Traits , typename TDS , typename SLDS >
Object CGAL::Regular_triangulation_3< Traits, TDS, 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
rt.dimension() \( \geq2\) and f is not infinite.

◆ find_conflicts()

template<typename Traits , typename TDS , typename SLDS >
template<class OutputIteratorBoundaryFacets , class OutputIteratorCells , class OutputIteratorInternalFacets >
Triple< OutputIteratorBoundaryFacets, OutputIteratorCells, OutputIteratorInternalFacets > CGAL::Regular_triangulation_3< Traits, TDS, SLDS >::find_conflicts ( const Weighted_point  p,
Cell_handle  c,
OutputIteratorBoundaryFacets  bfit,
OutputIteratorCells  cit,
OutputIteratorInternalFacets  ifit,
bool *  could_lock_zone = nullptr,
const Facet this_facet_must_be_in_the_cz = nullptr,
bool *  the_facet_is_in_its_cz = nullptr 
)

Compute the conflicts with p.

Parameters
pThe query point.
cThe starting cell.
citThe cells (resp. facets) in conflict with p.
bfitThe facets (resp. edges) on the boundary of the conflict zone, that is, the facets (resp. edges) (t, i) where the cell (resp.. facet) t is in conflict, but t->neighbor(i) is not.
ifitThe facets (resp. edges) inside the conflict zone, that facets incident to two cells (resp. facets) in conflict.
could_lock_zoneThe 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_facet_must_be_in_the_czIf the optional argument this_facet_must_be_in_the_cz is not null, the algorithm will check if this facet is in the conflict zone (it may be internal as well as boundary).
the_facet_is_in_its_czThis argument must be not null if the previous this_facet_must_be_in_the_cz argument is not null. The boolean value pointed by this pointer is set to true if *this_facet_must_be_in_the_cz is among the internal or boundary facets of the conflict zone, and false otherwise.
Precondition
The starting cell (resp. facet) c must be in conflict with p.
rt.dimension() \( \geq2\), and c is in conflict with p.
Returns
the Triple composed of the resulting output iterators.

◆ insert() [1/4]

template<typename Traits , typename TDS , typename SLDS >
Vertex_handle CGAL::Regular_triangulation_3< Traits, TDS, SLDS >::insert ( const Weighted_point p,
Cell_handle  start = Cell_handle(),
bool *  could_lock_zone = nullptr 
)

Inserts the weighted point p in the triangulation.

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

If this insertion creates a vertex, this vertex is returned.

If p coincides with an existing vertex and has a greater weight, then the existing weighted point becomes hidden (see RegularTriangulationCellBase_3) and p replaces it as vertex of the triangulation.

If p coincides with an already existing vertex (both point and weights being equal), then this vertex is returned and the triangulation remains unchanged.

Otherwise if p does not appear as a vertex of the triangulation, then it is stored as a hidden point and this method returns the default constructed handle.

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 point is not inserted). In any case, the locked cells are not unlocked by the function, leaving this choice to the user.

◆ insert() [2/4]

template<typename Traits , typename TDS , typename SLDS >
Vertex_handle CGAL::Regular_triangulation_3< Traits, TDS, SLDS >::insert ( const Weighted_point p,
Locate_type  lt,
Cell_handle  loc,
int  li,
int  lj,
bool *  could_lock_zone = nullptr 
)

Inserts the weighted 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 SLDS >
template<class InputIterator >
std::ptrdiff_t CGAL::Regular_triangulation_3< Traits, TDS, SLDS >::insert ( InputIterator  first,
InputIterator  last 
)

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

It returns the difference of the number of vertices between after and before the insertions (it may be negative due to hidden points). Note that this function is not guaranteed to insert the points following the order of InputIterator, as spatial_sort() is used to improve efficiency. If parallelism is enabled, the points will be inserted in parallel.

Template Parameters
InputIteratormust be an input iterator with value type Weighted_point .

◆ insert() [4/4]

template<typename Traits , typename TDS , typename SLDS >
template<class WeightedPointWithInfoInputIterator >
std::ptrdiff_t CGAL::Regular_triangulation_3< Traits, TDS, SLDS >::insert ( WeightedPointWithInfoInputIterator  first,
WeightedPointWithInfoInputIterator  last 
)

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

It returns the difference of the number of vertices between after and before the insertions (it may be negative due to hidden points). Note that this function is not guaranteed to insert the weighted points following the order of WeightedPointWithInfoInputIterator, as spatial_sort() is used to improve efficiency. 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, 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
WeightedPointWithInfoInputIteratormust be an input iterator with value type std::pair<Weighted_point,Vertex::Info>.

◆ insert_in_hole()

template<typename Traits , typename TDS , typename SLDS >
template<class CellIt >
Vertex_handle CGAL::Regular_triangulation_3< Traits, TDS, SLDS >::insert_in_hole ( const Weighted_point p,
CellIt  cell_begin,
CellIt  cell_end,
Cell_handle  begin,
int  i 
)

Creates a new vertex by starring a hole.

It takes an iterator range [cell_begin,cell_end) of Cell_handles which specifies a hole: a set of connected cells (resp. facets in dimension 2) which is star-shaped wrt p. (begin, i) is a facet (resp. an edge) on the boundary of the hole, that is, begin belongs to the set of cells (resp. facets) previously described, and begin->neighbor(i) does not. Then this function deletes all the cells (resp. facets) describing the hole, creates a new vertex v, and for each facet (resp. edge) on the boundary of the hole, creates a new cell (resp. facet) with v as vertex. Then v->set_point(p) is called and v is returned.

If the hole contains interior vertices, each of them is hidden by the insertion of p and is stored in the new cell which contains it.

Precondition
rt.dimension() \( \geq2\), the set of cells (resp. facets in dimension 2) is connected, not empty, its boundary is connected, and p lies inside the hole, which is star-shaped wrt p.

◆ is_valid()

template<typename Traits , typename TDS , typename SLDS >
bool CGAL::Regular_triangulation_3< Traits, TDS, SLDS >::is_valid ( bool  verbose = false) const

This is a function for debugging purpose.

Debugging Support

Checks the combinatorial validity of the triangulation and the validity of its geometric embedding (see Section Representation). Also checks that all the power spheres (resp. power circles in dimension 2, power segments in dimension 1) of cells (resp. facets in dimension 2, edges in dimension 1) are regular. When verbose is set to true, messages describing the first invalidity encountered are printed. This method is mainly a debugging help for the users of advanced features.

◆ nearest_power_vertex()

template<typename Traits , typename TDS , typename SLDS >
Vertex_handle CGAL::Regular_triangulation_3< Traits, TDS, SLDS >::nearest_power_vertex ( const Bare_point p,
Cell_handle  c = Cell_handle() 
)

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

This means that the power of the query point p with respect to the weighted point in the returned vertex is smaller than the power of p with respect to the weighted point for any other vertex. Ties are broken arbitrarily. The default constructed handle is returned if the triangulation is empty. The optional argument c is a hint specifying where to start the search.

Precondition
c is a cell of rt.

◆ remove() [1/2]

template<typename Traits , typename TDS , typename SLDS >
template<typename InputIterator >
int CGAL::Regular_triangulation_3< Traits, TDS, 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() [2/2]

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

◆ side_of_power_circle()

template<typename Traits , typename TDS , typename SLDS >
Bounded_side CGAL::Regular_triangulation_3< Traits, TDS, SLDS >::side_of_power_circle ( const Facet f,
const Weighted_point p 
) const

Returns the position of the point p with respect to the power circle of f.

More precisely, it returns:

  • in dimension 3:
  • For a finite facet,

ON_BOUNDARY if p is orthogonal to the power circle in the plane of the facet,

ON_UNBOUNDED_SIDE when their angle is less than \( \pi/2\),

ON_BOUNDED_SIDE when it is greater than \( \pi/2\) (see Figure Triangulation3figsidedim2).

  • For an infinite facet, it considers the plane defined by the finite facet of the cell f.first, and does the same as in dimension 2 in this plane.
  • in dimension 2:
  • For a finite facet,

ON_BOUNDARY if p is orthogonal to the circle,

ON_UNBOUNDED_SIDE when the angle between p and the power circle of f is less than \( \pi/2\), ON_BOUNDED_SIDE when it is greater than \( \pi/2\).

  • For an infinite facet,

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 in the other open half plane.

If the point p is collinear with the finite edge e of f, it returns:

ON_BOUNDED_SIDE if \( \Pi({p}^{(w)}-{z(e)}^{(w)})<0\), where \( {z(e)}^{(w)}\) is the power segment of e in the line supporting e,

ON_BOUNDARY if \( \Pi({p}^{(w)}-{z(e)}^{(w)})=0\),

ON_UNBOUNDED_SIDE if \( \Pi({p}^{(w)}-{z(e)}^{(w)})>0\) .

Precondition
rt.dimension() \( \geq2\).

◆ side_of_power_segment()

template<typename Traits , typename TDS , typename SLDS >
Bounded_side CGAL::Regular_triangulation_3< Traits, TDS, SLDS >::side_of_power_segment ( Cell_handle  c,
const Weighted_point p 
) const

In dimension 1, returns.

ON_BOUNDED_SIDE if \( \Pi({p}^{(w)}-{z(c)}^{(w)})<0\), where \( {z(c)}^{(w)}\) is the power segment of the edge represented by c,

ON_BOUNDARY if \( \Pi({p}^{(w)}-{z(c)}^{(w)})=0\),

ON_UNBOUNDED_SIDE if \( \Pi({p}^{(w)}-{z(c)}^{(w)})>0\) .

Precondition
rt.dimension() \( = 1\).

◆ side_of_power_sphere()

template<typename Traits , typename TDS , typename SLDS >
Bounded_side CGAL::Regular_triangulation_3< Traits, TDS, SLDS >::side_of_power_sphere ( Cell_handle  c,
const Weighted_point p 
) const

Returns the position of the weighted point \( p\) with respect to the power sphere of c.

More precisely, it returns:

  • ON_BOUNDED_SIDE if \( \Pi({p}^{(w)}-{z(c)}^{(w)})<0\) where \( {z(c)}^{(w)}\) is the power sphere of c. For an infinite cell this means either that p lies strictly in the half space limited by its finite facet and not containing any other point of the triangulation, or that the angle between p and the power circle of the finite facet of c is greater than \( \pi/2\).
  • ON_BOUNDARY if p is orthogonal to the power sphere of c i.e. \( \Pi({p}^{(w)}-{z(c)}^{(w)})=0\). For an infinite cell this means that p is orthogonal to the power circle of its finite facet.
  • ON_UNBOUNDED_SIDE if \( \Pi({p}^{(w)}-{z(c)}^{(w)})>0\) i.e. the angle between the weighted point p and the power sphere of c is less than \( \pi/2\) or if these two spheres do not intersect. For an infinite cell this means that p does not satisfy either of the two previous conditions.
    Precondition
    rt.dimension() \( =3\).

◆ vertices_in_conflict()

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

◆ vertices_inside_conflict_zone()

template<typename Traits , typename TDS , typename SLDS >
template<class OutputIterator >
OutputIterator CGAL::Regular_triangulation_3< Traits, TDS, SLDS >::vertices_inside_conflict_zone ( const Weighted_point p,
Cell_handle  c,
OutputIterator  res 
)

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

The vertices that are on the boundary of the conflict zone are not reported. Returns the resulting output iterator.

Precondition
rt.dimension() \( \geq2\), and c is a cell containing p.

◆ vertices_on_conflict_zone_boundary()

template<typename Traits , typename TDS , typename SLDS >
template<class OutputIterator >
OutputIterator CGAL::Regular_triangulation_3< Traits, TDS, SLDS >::vertices_on_conflict_zone_boundary ( const Weighted_point p,
Cell_handle  c,
OutputIterator  res 
)

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

Returns the resulting output iterator.

Precondition
rt.dimension() \( \geq2\), and c is a cell containing p.