| CGAL 4.4 - 3D Triangulations | 
#include <CGAL/Delaunay_triangulation_3.h>
The class Delaunay_triangulation_3 represents a three-dimensional Delaunay triangulation. 
| DelaunayTriangulationTraits_3 | is the geometric traits class. | 
| TriangulationDataStructure_3 | is the triangulation data structure. It has the default value Triangulation_data_structure_3<Triangulation_vertex_base_3<DelaunayTriangulationTraits_3>, Triangulation_cell_base_3<DelaunayTriangulationTraits_3> >. | 
| LocationPolicy | is a tag which must be a Location_policy<Tag>: eitherFast_locationorCompact_location.Fast_locationoffers 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 [11]. The default isCompact_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. Note that this argument can also come in second position, which can be useful when the default value for theTriangulationDataStructure_3parameter is satisfactory (this is using so-called deduced parameters). Note that this argument replaces the functionality provided before CGAL 3.6 byTriangulation_hierarchy_3. An example of use can be found in the user manual Fast Point Location for Delaunay Triangulations. | 
CGAL::Regular_triangulation_3 | Types | |
| typedef LocationPolicy | Location_policy | 
| 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()) | |
| 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()) | |
| Equivalent to constructing an empty triangulation with the optional traits class argument and calling insert(first,last).  More... | |
| Insertion | |
| Vertex_handle | insert (const Point &p, Cell_handle start=Cell_handle()) | 
| Inserts point pin the triangulation and returns the corresponding vertex.  More... | |
| Vertex_handle | insert (const Point &p, Vertex_handle hint) | 
| Same as above but uses hintas a starting place for the search.  More... | |
| Vertex_handle | insert (const Point &p, Locate_type lt, Cell_handle loc, int li, int lj) | 
| Inserts point pin 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 vertexvisp, andvis returned.  More... | |
| Vertex_handle | move (Vertex_handle v, const Point &p) | 
| same as above if there is no collision.  More... | |
| Removal | |
| When a vertex  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 [18]. 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  | |
| void | remove (Vertex_handle v) | 
| Removes the vertex vfrom 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 cthe pointplies.  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 fthe pointplies.  More... | |
| Bounded_side | side_of_circle (Cell_handle c, int i, const Point &p) | 
| Same as the previous method for facet iof cellc.  More... | |
| 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 cthat is nearest to \( p\).  More... | |
| A point  The set of cells (resp. facets in dimension 2) which are in conflict with  | |
| template<class OutputIteratorBoundaryFacets , class OutputIteratorCells > | |
| std::pair < OutputIteratorBoundaryFacets, OutputIteratorCells > | find_conflicts (Point p, Cell_handle c, OutputIteratorBoundaryFacets bfit, OutputIteratorCells cit) | 
| 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) | 
| 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 withp.  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 ofp, in the output iteratorres.  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).  More... | |
| 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 dtintoos.  More... | |
| 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 > | |
| enum | Locate_type | 
| The enum Locate_typeis defined byTriangulation_3to specify which case occurs when locating a point in the triangulation.  More... | |
| typedef Delaunay_triangulation_3 < DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy > ::Triangulation_data_structure | Triangulation_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  More... | |
| typedef Delaunay_triangulation_3 < DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy > ::Triangulation_data_structure::Cell_handle | Cell_handle | 
| handle to a cell  More... | |
| typedef Triangulation_simplex_3< Self > | Simplex | 
| Reference to a simplex (vertex, edge, facet or cell) of the triangulation.  More... | |
| typedef Delaunay_triangulation_3 < DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy > ::Triangulation_data_structure::size_type | size_type | 
| Size type (an unsigned integral type)  More... | |
| typedef Delaunay_triangulation_3 < DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy > ::Triangulation_data_structure::difference_type | difference_type | 
| Difference type (a signed integral type)  More... | |
| typedef Delaunay_triangulation_3 < DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy > ::Triangulation_data_structure::Cell_iterator | All_cells_iterator | 
| iterator over cells  More... | |
| typedef Delaunay_triangulation_3 < DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy > ::Triangulation_data_structure::Facet_iterator | All_facets_iterator | 
| iterator over facets  More... | |
| typedef Delaunay_triangulation_3 < DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy > ::Triangulation_data_structure::Edge_iterator | All_edges_iterator | 
| iterator over edges  More... | |
| typedef Delaunay_triangulation_3 < DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy > ::Triangulation_data_structure::Vertex_iterator | All_vertices_iterator | 
| iterator over vertices  More... | |
| typedef unspecified_type | Finite_cells_iterator | 
| iterator over finite cells  More... | |
| typedef unspecified_type | Finite_facets_iterator | 
| iterator over finite facets  More... | |
| typedef unspecified_type | Finite_edges_iterator | 
| iterator over finite edges  More... | |
| typedef unspecified_type | Finite_vertices_iterator | 
| iterator over finite vertices  More... | |
| typedef unspecified_type | Point_iterator | 
| iterator over the points corresponding to the finite vertices of the triangulation.  More... | |
| 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  More... | |
| 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  More... | |
|  Public Member Functions inherited from CGAL::Triangulation_3< DelaunayTriangulationTraits_3, Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::Triangulation_data_structure > | |
| Triangulation_3 (const DelaunayTriangulationTraits_3 &traits=DelaunayTriangulationTraits_3()) | |
| Introduces a triangulation thaving only one vertex which is the infinite vertex.  More... | |
| Triangulation_3 (const Triangulation_3 &tr) | |
| Copy constructor.  More... | |
| Triangulation_3 (InputIterator first, InputIterator last, const DelaunayTriangulationTraits_3 &traits=DelaunayTriangulationTraits_3()) | |
| Equivalent to constructing an empty triangulation with the optional traits class argument and calling insert(first,last).  More... | |
| Triangulation_3 & | operator= (const Triangulation_3 &tr) | 
| The triangulation tris duplicated, and modifying the copy after the duplication does not modify the original.  More... | |
| void | swap (Triangulation_3 &tr) | 
| The triangulations trandtare swapped.  More... | |
| void | clear () | 
| Deletes all finite vertices and all cells of t.  More... | |
| 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==.  More... | |
| const DelaunayTriangulationTraits_3 & | geom_traits () const | 
| Returns a const reference to the geometric traits object.  More... | |
| const Delaunay_triangulation_3 < DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy > ::Triangulation_data_structure & | tds () const | 
| Returns a const reference to the triangulation data structure.  More... | |
| 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.  More... | |
| size_type | number_of_vertices () const | 
| Returns the number of finite vertices.  More... | |
| size_type | number_of_cells () const | 
| Returns the number of cells or 0 if t.dimension() < 3.  More... | |
| Vertex_handle | infinite_vertex () | 
| Returns the infinite vertex.  More... | |
| 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.  More... | |
| 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 Point & | point (Cell_handle c, int i) const | 
| Returns the point given by vertex iof cellc.  More... | |
| const Point & | point (Vertex_handle v) const | 
| Same as the previous method for vertex v.  More... | |
| bool | is_infinite (Vertex_handle v) const | 
| true, iff vertexvis the infinite vertex.  More... | |
| bool | is_infinite (Cell_handle c) const | 
| true, iffcis incident to the infinite vertex.  More... | |
| bool | is_infinite (Cell_handle c, int i) const | 
| true, iff the facetiof cellcis incident to the infinite vertex.  More... | |
| bool | is_infinite (const Facet &f) const | 
| trueiff facetfis incident to the infinite vertex.  More... | |
| bool | is_infinite (Cell_handle c, int i, int j) const | 
| true, iff the edge(i,j)of cellcis incident to the infinite vertex.  More... | |
| bool | is_infinite (const Edge &e) const | 
| trueiff edgeeis incident to the infinite vertex.  More... | |
| bool | is_vertex (const Point &p, Vertex_handle &v) const | 
| Tests whether pis a vertex oftby locatingpin the triangulation.  More... | |
| bool | is_vertex (Vertex_handle v) const | 
| Tests whether vis a vertex oft.  More... | |
| bool | is_edge (Vertex_handle u, Vertex_handle v, Cell_handle &c, int &i, int &j) const | 
| Tests whether (u,v)is an edge oft.  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 oft.  More... | |
| bool | is_cell (Cell_handle c) const | 
| Tests whether cis a cell oft.  More... | |
| 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 oft.  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 oftand computes this cellc.  More... | |
| bool | has_vertex (const Facet &f, Vertex_handle v, int &j) const | 
| If vis a vertex off, thenjis the index ofvin the cellf.first, and the method returnstrue.  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.  More... | |
| 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()) const | 
| If the point querylies 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) const | 
| Same as above but uses hintas the starting place for the search.  More... | |
| Cell_handle | locate (const Point &query, Locate_type <, int &li, int &lj, Cell_handle start=Cell_handle()) const | 
| If querylies inside the affine hull of the points, the \( k\)-face (finite or infinite) that containsqueryin its interior is returned, by means of the cell returned together withlt, which is set to the locate type of the query (VERTEX, EDGE, FACET, CELL, orOUTSIDE_CONVEX_HULLif the cell is infinite andquerylies strictly in it) and two indicesliandljthat specify the \( k\)-face of the cell containingquery.  More... | |
| Cell_handle | locate (const Point &query, Locate_type <, int &li, int &lj, Vertex_handle hint) const | 
| Same as above but uses hintas the starting place for the search.  More... | |
| 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 <, int &li, int &lj) const | 
| Returns a value indicating on which side of the oriented boundary of cthe pointplies.  More... | |
| Bounded_side | side_of_facet (const Point &p, const Facet &f, Locate_type <, int &li, int &lj) const | 
| Returns a value indicating on which side of the oriented boundary of fthe pointplies:  More... | |
| Bounded_side | side_of_facet (const Point &p, Cell_handle c, Locate_type <, int &li, int &lj) const | 
| Same as the previous method for the facet (c,3).  More... | |
| Bounded_side | side_of_edge (const Point &p, const Edge &e, Locate_type <, int &li) const | 
| Returns a value indicating on which side of the oriented boundary of ethe pointplies:  More... | |
| Bounded_side | side_of_edge (const Point &p, Cell_handle c, Locate_type <, int &li) const | 
| Same as the previous method for edge \( (c,0,1)\).  More... | |
| 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 pin the triangulation and returns the corresponding vertex.  More... | |
| Vertex_handle | insert (const Point &p, Vertex_handle hint) | 
| Same as above but uses hintas the starting place for the search.  More... | |
| Vertex_handle | insert (const Point &p, Locate_type lt, Cell_handle loc, int li, int lj) | 
| Inserts point pin 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 pin cellc.  More... | |
| Vertex_handle | insert_in_facet (const Point &p, const Facet &f) | 
| Inserts point pin facetf.  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 pin edgee.  More... | |
| Vertex_handle | insert_in_edge (Point p, Cell_handle c, int i, int j) | 
| As above, inserts pin edge \( (i, j)\) ofc.  More... | |
| Vertex_handle | insert_outside_convex_hull (const Point &p, Cell_handle c) | 
| The cell cmust be an infinite cell containingp.  More... | |
| Vertex_handle | insert_outside_affine_hull (const Point &p) | 
| pis linked to all the points, and the infinite vertex is linked to all the points (includingp) 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 newvwill be used as the new vertex, which must have been allocated previously with e.g.create_vertex.  More... | |
| 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.  More... | |
| 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.  More... | |
| 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.  More... | |
| 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.  More... | |
| All_vertices_iterator | all_vertices_begin () const | 
| Starts at an arbitrary vertex.  More... | |
| All_vertices_iterator | all_vertices_end () const | 
| Past-the-end iterator.  More... | |
| All_edges_iterator | all_edges_begin () const | 
| Starts at an arbitrary edge.  More... | |
| All_edges_iterator | all_edges_end () const | 
| Past-the-end iterator.  More... | |
| All_facets_iterator | all_facets_begin () const | 
| Starts at an arbitrary facet.  More... | |
| All_facets_iterator | all_facets_end () const | 
| Past-the-end iterator.  More... | |
| All_cells_iterator | all_cells_begin () const | 
| Starts at an arbitrary cell.  More... | |
| All_cells_iterator | all_cells_end () const | 
| Past-the-end iterator.  More... | |
| Point_iterator | points_begin () const | 
| Iterates over the points of the triangulation.  More... | |
| Point_iterator | points_end () const | 
| Past-the-end iterator.  More... | |
| 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)ofc.  More... | |
| 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)ofc.  More... | |
| OutputIterator | incident_cells (Vertex_handle v, OutputIterator cells) const | 
| Copies the Cell_handles of all cells incident tovto the output iteratorcells.  More... | |
| OutputIterator | incident_facets (Vertex_handle v, OutputIterator facets) const | 
| Copies all Facets incident tovto the output iteratorfacets.  More... | |
| OutputIterator | finite_incident_cells (Vertex_handle v, OutputIterator cells) const | 
| Copies the Cell_handles of all finite cells incident tovto the output iteratorcells.  More... | |
| OutputIterator | finite_incident_facets (Vertex_handle v, OutputIterator facets) const | 
| Copies all finite Facets incident tovto the output iteratorfacets.  More... | |
| OutputIterator | incident_edges (Vertex_handle v, OutputIterator edges) const | 
| Copies all Edges incident tovto the output iteratoredges.  More... | |
| OutputIterator | finite_incident_edges (Vertex_handle v, OutputIterator edges) const | 
| Copies all finite Edges incident tovto the output iteratoredges.  More... | |
| OutputIterator | adjacent_vertices (Vertex_handle v, OutputIterator vertices) const | 
| Copies the Vertex_handles of all vertices adjacent tovto the output iteratorvertices.  More... | |
| OutputIterator | finite_adjacent_vertices (Vertex_handle v, OutputIterator vertices) const | 
| Copies the Vertex_handles of all finite vertices adjacent tovto the output iteratorvertices.  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)ofc.  More... | |
| 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 finstart.  More... | |
| Facet_circulator | incident_facets (Cell_handle c, int i, int j, Facet start) const | 
| As above for edge (i,j)ofc.  More... | |
| Facet_circulator | incident_facets (Cell_handle c, int i, int j, Cell_handle start, int f) const | 
| As above for edge (i,j)ofcand facet(start,f).  More... | |
| int | mirror_index (Cell_handle c, int i) const | 
| Returns the index of cin its \( i^{th}\) neighbor.  More... | |
| Vertex_handle | mirror_vertex (Cell_handle c, int i) const | 
| Returns the vertex of the \( i^{th}\) neighbor of cthat is opposite toc.  More... | |
| Facet | mirror_facet (Facet f) const | 
| Returns the same facet seen from the other adjacent cell.  More... | |
| 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 isby 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 tintoos.  More... | |
|  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) | 
| typedef DelaunayTriangulationTraits_3::Line_3 CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::Line | 
| typedef LocationPolicy CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::Location_policy | 
| typedef DelaunayTriangulationTraits_3::Object_3 CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::Object | 
| typedef DelaunayTriangulationTraits_3::Plane_3 CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::Plane | 
| typedef DelaunayTriangulationTraits_3::Ray_3 CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::Ray | 
| CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::Delaunay_triangulation_3 | ( | const DelaunayTriangulationTraits_3 & | traits = DelaunayTriangulationTraits_3() | ) | 
Creates an empty Delaunay triangulation, possibly specifying a traits class traits. 
| CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::Delaunay_triangulation_3 | ( | const Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy > & | dt1 | ) | 
Copy constructor.
| CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::Delaunay_triangulation_3 | ( | InputIterator | first, | 
| InputIterator | last, | ||
| const DelaunayTriangulationTraits_3 & | traits = DelaunayTriangulationTraits_3() | ||
| ) | 
Equivalent to constructing an empty triangulation with the optional traits class argument and calling insert(first,last). 
| Stream& CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::draw_dual | ( | Stream & | os | ) | 
Sends the set of duals to all the facets of dt into os. 
| Point CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::dual | ( | Cell_handle | c | ) | const | 
Returns the circumcenter of the four vertices of c.
dt.dimension() \( =3\) and c is not infinite. | Object CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::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.
dt.dimension() \( \geq2\) and f is not infinite. | Object CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::dual | ( | Cell_handle | c, | 
| int | i | ||
| ) | const | 
same as the previous method for facet (c,i). 
| Line CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::dual_support | ( | Cell_handle | c, | 
| int | i | ||
| ) | const | 
returns the line supporting the dual of the facet.
dt.dimension() \( \geq2\) and f is not infinite. | std::pair<OutputIteratorBoundaryFacets, OutputIteratorCells> CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::find_conflicts | ( | Point | p, | 
| Cell_handle | c, | ||
| OutputIteratorBoundaryFacets | bfit, | ||
| OutputIteratorCells | cit | ||
| ) | 
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.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. 
dt.dimension() \( \geq2\), and c is in conflict with p. | Triple<OutputIteratorBoundaryFacets, OutputIteratorCells, OutputIteratorInternalFacets> CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::find_conflicts | ( | Point | p, | 
| Cell_handle | c, | ||
| OutputIteratorBoundaryFacets | bfit, | ||
| OutputIteratorCells | cit, | ||
| OutputIteratorInternalFacets | ifit | ||
| ) | 
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.Returns the Triple composed of the resulting output iterators. 
dt.dimension() \( \geq2\), and c is in conflict with p. | Vertex_handle CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::insert | ( | const Point & | p, | 
| Cell_handle | start = Cell_handle() | ||
| ) | 
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. 
| Vertex_handle CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::insert | ( | const Point & | p, | 
| Vertex_handle | hint | ||
| ) | 
Same as above but uses hint as a starting place for the search. 
| Vertex_handle CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::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. 
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(). 
| std::ptrdiff_t CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::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. 
| PointInputIterator | must be an input iterator with the value type Point. | 
| std::ptrdiff_t CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::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. 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. 
Vertex must be model of the concept TriangulationVertexBaseWithInfo_3.| PointWithInfoInputIterator | must be an input iterator with the value type std::pair<Point,Vertex::Info>. | 
| bool CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::is_Gabriel | ( | Cell_handle | c, | 
| int | i | ||
| ) | 
| bool CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::is_Gabriel | ( | Cell_handle | c, | 
| int | i, | ||
| int | j | ||
| ) | 
| bool CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::is_Gabriel | ( | const Facet & | f | ) | 
| bool CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::is_Gabriel | ( | const Edge & | e | ) | 
| bool CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::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. 
| bool CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::is_valid | ( | Cell_handle | c, | 
| bool | verbose = false | ||
| ) | const | 
This is a function for debugging purpose.
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.  
| Vertex_handle CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::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. 
v must be finite. | Vertex_handle CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::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. 
v must be finite. | Vertex_handle CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::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. 
c is a cell of dt. | Vertex_handle CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::nearest_vertex_in_cell | ( | Point | p, | 
| Cell_handle | c | ||
| ) | 
Returns the vertex of the cell c that is nearest to \( p\). 
| void CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::remove | ( | Vertex_handle | v | ) | 
Removes the vertex v from the triangulation. 
v is a finite vertex of the triangulation. | int CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::remove | ( | InputIterator | first, | 
| InputIterator | beyond | ||
| ) | 
Removes the vertices specified by the iterator range [first, beyond). 
The function remove(Vertex_handle) is called over each element of the range. The number of vertices removed is returned. 
| InputIterator | must be an input iterator with value type Vertex_handle. | 
| int CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::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.
| InputIterator | must be an input iterator with value type Vertex_handle. | 
| Bounded_side CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::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:
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.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.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. dt.dimension() \( \geq2\) and in dimension 3, p is coplanar with f. | Bounded_side CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::side_of_circle | ( | Cell_handle | c, | 
| int | i, | ||
| const Point & | p | ||
| ) | 
Same as the previous method for facet i of cell c. 
| Bounded_side CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::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. dt.dimension() \( =3\). | OutputIterator CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::vertices_in_conflict | ( | Point | p, | 
| Cell_handle | c, | ||
| OutputIterator | res | ||
| ) | 
vertices_on_conflict_zone_boundary since CGAL-3.8. | OutputIterator CGAL::Delaunay_triangulation_3< DelaunayTriangulationTraits_3, TriangulationDataStructure_3, LocationPolicy >::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.
dt.dimension() \( \geq2\), and c is in conflict with p.