CGAL::Triangulation_3<TriangulationTraits_3,TriangulationDataStructure_3>

Definition

The class Triangulation_3 represents a 3-dimensional tetrahedralization of points.

#include <CGAL/Triangulation_3.h>

Parameters

The first template argument must be a model of the TriangulationTraits_3 concept.

The second template argument must be a model of the TriangulationDataStructure_3 concept. It has the default value Triangulation_data_structure_3< Triangulation_vertex_base_3<TriangulationTraits_3>,Triangulation_cell_base_3<TriangulationTraits_3> >.

Inherits From

Triangulation_utils_3

Types

The class Triangulation_3 defines the following types:

typedef TriangulationDataStructure_3
Triangulation_data_structure;
typedef TriangulationTraits_3 Geom_traits;

typedef TriangulationTraits_3::Point_3
Point;
typedef TriangulationTraits_3::Segment_3
Segment;
typedef TriangulationTraits_3::Triangle_3
Triangle;
typedef TriangulationTraits_3::Tetrahedron_3
Tetrahedron;

Only vertices (0-faces) and cells (3-faces) are stored. Edges (1-faces) and facets (2-faces) are not explicitly represented and thus there are no corresponding classes (see Section 37.1).

typedef TriangulationDataStructure_3::Vertex
Vertex;
typedef TriangulationDataStructure_3::Cell
Cell;
typedef TriangulationDataStructure_3::Facet
Facet;
typedef TriangulationDataStructure_3::Edge
Edge;

The vertices and faces of the triangulations are accessed through handles, iterators and circulators. A handle is a model of the Handle concept, and supports the two dereference operators operator* and operator->. A circulator is a model of the concept Circulator. Iterators and circulators are bidirectional and non-mutable. The edges and facets of the triangulation can also be visited through iterators and circulators which are bidirectional and non-mutable.

Iterators and circulators are convertible to the corresponding handles, thus the user can pass them directly as arguments to the functions.

typedef TriangulationDataStructure_3::Vertex_handle
Vertex_handle; handle to a vertex
typedef TriangulationDataStructure_3::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 TriangulationDataStructure_3::size_type
size_type; Size type (an unsigned integral type)
typedef TriangulationDataStructure_3::difference_type
difference_type; Difference type (a signed integral type)

typedef TriangulationDataStructure_3::Cell_iterator
All_cells_iterator; iterator over cells
typedef TriangulationDataStructure_3::Facet_iterator
All_facets_iterator; iterator over facets
typedef TriangulationDataStructure_3::Edge_iterator
All_edges_iterator; iterator over edges
typedef TriangulationDataStructure_3::Vertex_iterator
All_vertices_iterator; iterator over vertices

Triangulation_3<TriangulationTraits_3,TriangulationDataStructure_3>::Finite_cells_iterator
iterator over finite cells

Triangulation_3<TriangulationTraits_3,TriangulationDataStructure_3>::Finite_facets_iterator
iterator over finite facets

Triangulation_3<TriangulationTraits_3,TriangulationDataStructure_3>::Finite_edges_iterator
iterator over finite edges

Triangulation_3<TriangulationTraits_3,TriangulationDataStructure_3>::Finite_vertices_iterator
iterator over finite vertices

Triangulation_3<TriangulationTraits_3,TriangulationDataStructure_3>::Point_iterator
iterator over the points corresponding to the finite vertices of the triangulation.

typedef TriangulationDataStructure_3::Cell_circulator
Cell_circulator; circulator over all cells incident to a given edge
typedef TriangulationDataStructure_3::Facet_circulator
Facet_circulator; circulator over all facets incident to a given edge

The triangulation class also defines the following enum type to specify which case occurs when locating a point in the triangulation.

enum Locate_type { VERTEX=0, EDGE, FACET, CELL, OUTSIDE_CONVEX_HULL, OUTSIDE_AFFINE_HULL};

Creation

Triangulation_3<TriangulationTraits_3,TriangulationDataStructure_3> t ( TriangulationTraits_3 traits = TriangulationTraits_3());
Introduces a triangulation t having only one vertex which is the infinite vertex.


Triangulation_3<TriangulationTraits_3,TriangulationDataStructure_3> t ( Triangulation_3 tr);
Copy constructor. All vertices and faces are duplicated.


template < class InputIterator>
Triangulation_3<TriangulationTraits_3,TriangulationDataStructure_3> t ( InputIterator first,
InputIterator last,
TriangulationTraits_3 traits = TriangulationTraits_3());
Equivalent to constructing an empty triangulation with the optional traits class argument and calling insert(first,last).

Assignment

Triangulation_3 & t = Triangulation_3 tr The triangulation tr is duplicated, and modifying the copy after the duplication does not modify the original. The previous triangulation held by t is deleted.

void t.swap ( Triangulation_3 & tr) The triangulations tr and t are swapped. t.swap(tr) should be preferred to t = tr or to t(tr) if tr is deleted after that. Indeed, there is no copy of cells and vertices, thus this method runs in constant time.

void t.clear () Deletes all finite vertices and all cells of t.

template < class GT, class Tds1, class Tds2 >
bool Triangulation_3<GT, Tds1> t1 == Triangulation_3<GT, Tds2> t2
Equality operator. Returns true iff there exist a bijection between the vertices of t1 and those of t2 and a bijection between the cells of t1 and those of t2, which preserve the geometry of the triangulation, that is, the points of each corresponding pair of vertices are equal, and the tetrahedra corresponding to each pair of cells are equal (up to a permutation of their vertices).
template < class GT, class Tds1, class Tds2 >
bool Triangulation_3<GT, Tds1> t1 != Triangulation_3<GT, Tds2> t2
The opposite of operator==.

Access Functions

TriangulationTraits_3 t.geom_traits () Returns a const reference to the geometric traits object.
TriangulationDataStructure_3 t.tds () Returns a const reference to the triangulation data structure.

Non const access

The responsibility of keeping a valid triangulation belongs to the user when using advanced operations allowing a direct manipulation of the tds.

TriangulationDataStructure_3 & t.tds () Returns a reference to the triangulation data structure.

This method is mainly a help for users implementing their own triangulation algorithms.

int t.dimension () Returns the dimension of the affine hull.
size_type t.number_of_vertices () Returns the number of finite vertices.
size_type t.number_of_cells () Returns the number of cells or 0 if t.dimension()<3.

Vertex_handle t.infinite_vertex () Returns the infinite vertex.
Cell_handle t.infinite_cell () Returns a cell incident to the infinite vertex.

Non-constant-time access functions

As previously said, the triangulation is a collection of cells that are either infinite or represent a finite tetrahedra, where an infinite cell is a cell incident to the infinite vertex. Similarly we call an edge (resp. facet) infinite if it is incident to the infinite vertex.

size_type t.number_of_facets () The number of facets. Returns 0 if t.dimension()<2.
size_type t.number_of_edges () The number of edges. Returns 0 if t.dimension()<1.

size_type t.number_of_finite_cells () The number of finite cells. Returns 0 if t.dimension()<3.
size_type t.number_of_finite_facets () The number of finite facets. Returns 0 if t.dimension()<2.
size_type t.number_of_finite_edges () The number of finite edges. Returns 0 if t.dimension()<1.

Geometric access functions

Tetrahedron t.tetrahedron ( Cell_handle c) Returns the tetrahedron formed by the four vertices of c.
Precondition: t.dimension() =3 and the cell is finite.
Triangle t.triangle ( Cell_handle c, int i)
Returns the triangle formed by the three vertices of facet (c,i). The triangle is oriented so that its normal points to the inside of cell c.
Precondition: t.dimension() 2 and i {0,1,2,3} in dimension 3, i = 3 in dimension 2, and the facet is finite.
Triangle t.triangle ( Facet f) Same as the previous method for facet f.
Precondition: t.dimension() 2 and the facet is finite.
Segment t.segment ( Edge e) Returns the line segment formed by the vertices of e.
Precondition: t.dimension() 1 and e is finite.
Segment t.segment ( Cell_handle c, int i, int j)
Same as the previous method for edge (c,i,j).
Precondition: As above and i j. Moreover i,j {0,1,2,3} in dimension 3, i,j {0,1,2} in dimension 2, i,j {0,1} in dimension 1, and the edge is finite.
Point t.point ( Cell_handle c, int i) Returns the point given by vertex i of cell c.
Precondition: t.dimension() 0 and i {0,1,2,3} in dimension 3, i {0,1,2} in dimension 2, i {0,1} in dimension 1, i = 0 in dimension 0, and the vertex is finite.
Point t.point ( Vertex_handle v) Same as the previous method for vertex v.
Precondition: t.dimension() 0 and the vertex is finite.

Tests for Finite and Infinite Vertices and Faces

bool t.is_infinite ( Vertex_handle v) true, iff vertex v is the infinite vertex.
bool t.is_infinite ( Cell_handle c) true, iff c is incident to the infinite vertex.
Precondition: t.dimension() =3.
bool t.is_infinite ( Cell_handle c, int i)
true, iff the facet i of cell c is incident to the infinite vertex.
Precondition: t.dimension() 2 and i {0,1,2,3} in dimension 3, i=3 in dimension 2.
bool t.is_infinite ( Facet f) true iff facet f is incident to the infinite vertex.
Precondition: t.dimension() 2.
bool t.is_infinite ( Cell_handle c, int i, int j)
true, iff the edge (i,j) of cell c is incident to the infinite vertex.
Precondition: t.dimension() 1 and i j. Moreover i,j {0,1,2,3} in dimension 3, i,j {0,1,2} in dimension 2, i,j {0,1} in dimension 1.
bool t.is_infinite ( Edge e) true iff edge e is incident to the infinite vertex.
Precondition: t.dimension() 1.

Queries

bool t.is_vertex ( Point p, Vertex_handle & v)
Tests whether p is a vertex of t by locating p in the triangulation. If p is found, the associated vertex v is given.
bool t.is_vertex ( Vertex_handle v) Tests whether v is a vertex of t.

bool t.is_edge ( Vertex_handle u, Vertex_handle v, Cell_handle & c, int & i, int & j)
Tests whether (u,v) is an edge of t. If the edge is found, it gives a cell c having this edge and the indices i and j of the vertices u and v in c, in this order.
Precondition: u and v are vertices of t.

bool
t.is_facet ( Vertex_handle u,
Vertex_handle v,
Vertex_handle w,
Cell_handle & c,
int & i,
int & j,
int & k)
Tests whether (u,v,w) is a facet of t. If the facet is found, it computes a cell c having this facet and the indices i, j and k of the vertices u, v and w in c, in this order.
Precondition: u, v and w are vertices of t.

bool t.is_cell ( Cell_handle c) Tests whether c is a cell of t.
bool
t.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)
Tests whether (u,v,w,x) is a cell of t. If the cell c is found, the method computes the indices i, j, k and l of the vertices u, v, w and x in c, in this order.
Precondition: u, v, w and x are vertices of t.
bool t.is_cell ( Vertex_handle u, Vertex_handle v, Vertex_handle w, Vertex_handle x, Cell_handle & c)
Tests whether (u,v,w,x) is a cell of t and computes this cell c.
Precondition: u, v, w and x are vertices of t.

There is a method has_vertex in the cell class. The analogous methods for facets are defined here.

bool t.has_vertex ( Facet f, Vertex_handle v, int & j)
If v is a vertex of f, then j is the index of v in the cell f.first, and the method returns true.
Precondition: t.dimension()=3
bool t.has_vertex ( Cell_handle c, int i, Vertex_handle v, int & j)
Same for facet (c,i). Computes the index j of v in c.
bool t.has_vertex ( Facet f, Vertex_handle v)
bool t.has_vertex ( Cell_handle c, int i, Vertex_handle v)
Same as the first two methods, but these two methods do not return the index of the vertex.

The following three methods test whether two facets have the same vertices.

bool t.are_equal ( Cell_handle c, int i, Cell_handle n, int j)
bool t.are_equal ( Facet f, Facet g)
bool t.are_equal ( Facet f, Cell_handle n, int j)
For these three methods:
Precondition: t.dimension()=3.

Point location

The class Triangulation_3<TriangulationTraits_3,TriangulationDataStructure_3> provides two functions to locate a given point with respect to a triangulation. It provides also functions to test if a given point is inside a finite face or not. Note that the class Delaunay_triangulation_3 also provides a nearest_vertex() function.

Cell_handle t.locate ( Point query, Cell_handle start = Cell_handle())
If the point query lies inside the convex hull of the points, the cell that contains the query in its interior is returned. If query lies on a facet, an edge or on a vertex, one of the cells having query on its boundary is returned.
If the point query lies outside the convex hull of the points, an infinite cell with vertices { p, q, r, ∞} is returned such that the tetrahedron ( p, q, r, query ) is positively oriented (the rest of the triangulation lies on the other side of facet ( p, q, r )).
Note that locate works even in degenerate dimensions: in dimension 2 (resp. 1, 0) the Cell_handle returned is the one that represents the facet (resp. edge, vertex) containing the query point.
The optional argument start is used as a starting place for the search.

Cell_handle t.locate ( Point query, Vertex_handle hint)
Same as above but uses hint as the starting place for the search.

Cell_handle t.locate ( Point query, Locate_type & lt, int & li, int & lj, Cell_handle start = Cell_handle())
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.
If the k-face is a cell, li and lj have no meaning; if it is a facet (resp. vertex), li gives the index of the facet (resp. vertex) and lj has no meaning; if it is and edge, li and lj give the indices of its vertices.
If the point query lies outside the affine hull of the points, which can happen in case of degenerate dimensions, lt is set to OUTSIDE_AFFINE_HULL, and the cell returned has no meaning. As a particular case, if there is no finite vertex yet in the triangulation, lt is set to OUTSIDE_AFFINE_HULL and locate returns the default constructed handle.
The optional argument start is used as a starting place for the search.

Cell_handle t.locate ( Point query, Locate_type & lt, int & li, int & lj, Vertex_handle hint)
Same as above but uses hint as the starting place for the search.

Bounded_side t.side_of_cell ( Point p, Cell_handle c, Locate_type & lt, int & li, int & lj)
Returns a value indicating on which side of the oriented boundary of c the point p lies. More precisely, it returns:
- ON_BOUNDED_SIDE if p is inside the cell. For an infinite cell this means that p lies strictly in the half space limited by its finite facet and not containing any other point of the triangulation.
- ON_BOUNDARY if p on the boundary of the cell. For an infinite cell this means that p lies on the finite facet. Then lt together with li and lj give the precise location on the boundary. (See the descriptions of the locate methods.)
- ON_UNBOUNDED_SIDE if p lies outside the cell. For an infinite cell this means that p does not satisfy either of the two previous conditions.
Precondition: t.dimension() =3

Bounded_side t.side_of_facet ( Point p, Facet f, Locate_type & lt, int & li, int & lj)
Returns a value indicating on which side of the oriented boundary of f the point p lies:
- ON_BOUNDED_SIDE if p is inside the facet. For an infinite facet this means that p lies strictly in the half plane limited by its finite edge and not containing any other point of the triangulation .
- ON_BOUNDARY if p is on the boundary of the facet. For an infinite facet this means that p lies on the finite edge. lt, li and lj give the precise location of p on the boundary of the facet. li and lj refer to indices in the degenerate cell c representing f.
- ON_UNBOUNDED_SIDE if p lies outside the facet. For an infinite facet this means that p does not satisfy either of the two previous conditions.
Precondition: t.dimension() =2 and p lies in the plane containing the triangulation. f.second =3 (in dimension 2 there is only one facet per cell).
Bounded_side t.side_of_facet ( Point p, Cell_handle c, Locate_type & lt, int & li, int & lj)
Same as the previous method for the facet (c,3).

Bounded_side t.side_of_edge ( Point p, Edge e, Locate_type & lt, int & li)
Returns a value indicating on which side of the oriented boundary of e the point p lies:
- ON_BOUNDED_SIDE if p is inside the edge. For an infinite edge this means that p lies in the half line defined by the vertex and not containing any other point of the triangulation.
- ON_BOUNDARY if p equals one of the vertices, li give the index of the vertex in the cell storing e
- ON_UNBOUNDED_SIDE if p lies outside the edge. For an infinite edge this means that p lies on the other half line, which contains the other points of the triangulation.
Precondition: t.dimension() =1 and p is collinear with the points of the triangulation. e.second =0 and e.third =1 (in dimension 1 there is only one edge per cell).
Bounded_side t.side_of_edge ( Point p, Cell_handle c, Locate_type & lt, int & li)
Same as the previous method for edge (c,0,1).

Flips

Two kinds of flips exist for a three-dimensional triangulation. They are reciprocal. To be flipped, an edge must be incident to three tetrahedra. During the flip, these three tetrahedra disappear and two tetrahedra appear. Figure 37.10(left) shows the edge that is flipped as bold dashed, and one of its three incident facets is shaded. On the right, the facet shared by the two new tetrahedra is shaded.

Flips are possible only under the following conditions:
- the edge or facet to be flipped is not on the boundary of the convex hull of the triangulation
- the five points involved are in convex position.

Figure 37.10:  Flips.

Flips

The following methods guarantee the validity of the resulting 3D triangulation.

Flips for a 2d triangulation are not implemented yet

bool t.flip ( Edge e)
bool t.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). They return false or true according to this test.

void t.flip_flippable ( Edge e)
void t.flip_flippable ( Cell_handle c, int i, int j)
Should be preferred to the previous methods when the edge is known to be flippable.
Precondition: The edge is flippable.

bool t.flip ( Facet f)
bool t.flip ( Cell_handle c, int i) Before flipping, these methods check that facet f=(c,i) is flippable (which is quite expensive). They return false or true according to this test.

void t.flip_flippable ( Facet f)
void t.flip_flippable ( Cell_handle c, int i)
Should be preferred to the previous methods when the facet is known to be flippable.
Precondition: The facet is flippable.

Insertions

The following operations are guaranteed to lead to a valid triangulation when they are applied on a valid triangulation.

Vertex_handle t.insert ( Point p, Cell_handle start = Cell_handle())
Inserts point p in the triangulation and returns the corresponding vertex.
If point p coincides with an already existing vertex, this vertex is returned and the triangulation remains unchanged.
If point p lies in the convex hull of the points, it is added naturally: if it lies inside a cell, the cell is split into four cells, if it lies on a facet, the two incident cells are split into three cells, if it lies on an edge, all the cells incident to this edge are split into two cells.
If point p is strictly outside the convex hull but in the affine hull, p is linked to all visible points on the convex hull to form the new triangulation. See Figure 37.11.
If point p is outside the affine hull of the points, p is linked to all the points, and the dimension of the triangulation is incremented. All the points now belong to the boundary of the convex hull, so, the infinite vertex is linked to all the points to triangulate the new infinite face. See Figure 37.12. The optional argument start is used as a starting place for the search.

Vertex_handle t.insert ( Point p, Vertex_handle hint)
Same as above but uses hint as the starting place for the search.

Vertex_handle t.insert ( 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 locate() above.

template < class InputIterator >
int t.insert ( InputIterator first, InputIterator last)
Inserts the points in the range [.first, last.). Returns the number of inserted points. Note that this function is not guaranteed to insert the points following the order of InputIterator.
Precondition: The value_type of first and last is Point.

The previous methods are sufficient to build a whole triangulation. We also provide some other methods that can be used instead of insert(p) when the place where the new point p must be inserted is already known. They are also guaranteed to lead to a valid triangulation when they are applied on a valid triangulation.

Vertex_handle t.insert_in_cell ( Point p, Cell_handle c)
Inserts point p in cell c. Cell c is split into 4 tetrahedra.
Precondition: t.dimension() =3 and p lies strictly inside cell c.

Vertex_handle t.insert_in_facet ( Point p, Facet f)
Inserts point p in facet f. In dimension 3, the 2 neighboring cells are split into 3 tetrahedra; in dimension 2, the facet is split into 3 triangles.
Precondition: t.dimension() 2 and p lies strictly inside face f.
Vertex_handle t.insert_in_facet ( Point p, Cell_handle c, int i)
As above, insertion in facet (c,i).
Precondition: As above and i {0,1,2,3} in dimension 3, i = 3 in dimension 2.

Vertex_handle t.insert_in_edge ( Point p, Edge e)
Inserts p in edge e. In dimension 3, all the cells having this edge are split into 2 tetrahedra; in dimension 2, the 2 neighboring facets are split into 2 triangles; in dimension 1, the edge is split into 2 edges.
Precondition: t.dimension() 1 and p lies on edge e.
Vertex_handle t.insert_in_edge ( Point p, Cell_handle c, int i, int j)
As above, inserts p in edge (i, j) of c.
Precondition: As above and i j. Moreover i,j {0,1,2,3} in dimension 3, i,j {0,1,2} in dimension 2, i,j {0,1} in dimension 1.

Vertex_handle t.insert_outside_convex_hull ( Point p, Cell_handle c)
The cell c must be an infinite cell containing p.
Links p to all points in the triangulation that are visible from p. Updates consequently the infinite faces. See Figure 37.11.
Precondition: t.dimension() >0, c, and the k-face represented by c is infinite and contains t.

Figure 37.11:  insert_outside_convex_hull (2-dimensional case).

insert_outside_convex_hull} (2-dimensional case)

Vertex_handle t.insert_outside_affine_hull ( 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. See Figure 37.12.
This method can be used to insert the first point in an empty triangulation.
Precondition: t.dimension() <3 and p lies outside the affine hull of the points.

Figure 37.12:  insert_outside_affine_hull (2-dimensional case).

insert_outside_affine_hull} (2-dimensional case)

template <class CellIt>
Vertex_handle t.insert_in_hole ( 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.
This operation is equivalent to calling tds().insert_in_hole(cell_begin, cell_end, begin, i); v->set_point(p).
Precondition: t.dimension() 2, the set of cells (resp. facets in dimension 2) is connected, its boundary is connected, and p lies inside the hole, which is star-shaped wrt p.

template <class CellIt>
Vertex_handle
t.insert_in_hole ( Point p,
CellIt cell_begin,
CellIt cell_end,
Cell_handle begin,
int i,
Vertex_handle newv)
Same as above, except that newv will be used as the new vertex, which must have been allocated previously with e.g. create_vertex.

Traversal of the Triangulation

The triangulation class provides several iterators and circulators that allow one to traverse it (completely or partially).

Cell, Face, Edge and Vertex Iterators

The following iterators allow the user to visit cells, facets, edges and vertices of the triangulation. These iterators are non-mutable, bidirectional and their value types are respectively Cell, Facet, Edge and Vertex. They are all invalidated by any change in the triangulation.

Finite_vertices_iterator t.finite_vertices_begin () Starts at an arbitrary finite vertex. Then ++ and -- will iterate over finite vertices. Returns finite_vertices_end() when t.number_of_vertices() =0.
Finite_vertices_iterator t.finite_vertices_end () Past-the-end iterator
Finite_edges_iterator t.finite_edges_begin () Starts at an arbitrary finite edge. Then ++ and -- will iterate over finite edges. Returns finite_edges_end() when t.dimension() <1.
Finite_edges_iterator t.finite_edges_end () Past-the-end iterator
Finite_facets_iterator t.finite_facets_begin () Starts at an arbitrary finite facet. Then ++ and -- will iterate over finite facets. Returns finite_facets_end() when t.dimension() <2.
Finite_facets_iterator t.finite_facets_end () Past-the-end iterator
Finite_cells_iterator t.finite_cells_begin () Starts at an arbitrary finite cell. Then ++ and -- will iterate over finite cells. Returns finite_cells_end() when t.dimension() <3.
Finite_cells_iterator t.finite_cells_end () Past-the-end iterator

All_vertices_iterator t.all_vertices_begin () Starts at an arbitrary vertex. Iterates over all vertices (even the infinite one). Returns vertices_end() when t.number_of_vertices() =0.
All_vertices_iterator t.all_vertices_end () Past-the-end iterator
All_edges_iterator t.all_edges_begin () Starts at an arbitrary edge. Iterates over all edges (even infinite ones). Returns edges_end() when t.dimension() <1.
All_edges_iterator t.all_edges_end () Past-the-end iterator
All_facets_iterator t.all_facets_begin () Starts at an arbitrary facet. Iterates over all facets (even infinite ones). Returns facets_end() when t.dimension() <2.
All_facets_iterator t.all_facets_end () Past-the-end iterator
All_cells_iterator t.all_cells_begin () Starts at an arbitrary cell. Iterates over all cells (even infinite ones). Returns cells_end() when t.dimension() <3.
All_cells_iterator t.all_cells_end () Past-the-end iterator

Point_iterator t.points_begin () Iterates over the points of the triangulation.
Point_iterator t.points_end () Past-the-end iterator

Cell and Facet Circulators

The following circulators respectively visit all cells or all facets incident to a given edge. They are non-mutable and bidirectional. They are invalidated by any modification of one of the cells traversed.

Cell_circulator t.incident_cells ( Edge e) Starts at an arbitrary cell incident to e.
Precondition: t.dimension() =3.
Cell_circulator t.incident_cells ( Cell_handle c, int i, int j)
As above for edge (i,j) of c.
Cell_circulator t.incident_cells ( Edge e, Cell_handle start)
Starts at cell start.
Precondition: t.dimension() =3 and start is incident to e.
Cell_circulator t.incident_cells ( Cell_handle c, int i, int j, Cell_handle start)
As above for edge (i,j) of c.

The following circulators on facets are defined only in dimension 3, though facets are defined also in dimension 2: there are only two facets sharing an edge in dimension 2.

Facet_circulator t.incident_facets ( Edge e) Starts at an arbitrary facet incident to e.
Precondition: t.dimension() =3
Facet_circulator t.incident_facets ( Cell_handle c, int i, int j)
As above for edge (i,j) of c.
Facet_circulator t.incident_facets ( Edge e, Facet start)
Starts at facet start.
Precondition: start is incident to e.
Facet_circulator t.incident_facets ( Edge e, Cell_handle start, int f)
Starts at facet of index f in start.
Facet_circulator t.incident_facets ( Cell_handle c, int i, int j, Facet start)
As above for edge (i,j) of c.
Facet_circulator t.incident_facets ( Cell_handle c, int i, int j, Cell_handle start, int f)
As above for edge (i,j) of c and facet (start,f).

Traversal of the incident cells, facets and edges, and the adjacent vertices of a given vertex

template <class OutputIterator>
OutputIterator t.incident_cells ( Vertex_handle v, OutputIterator cells)
Copies the Cell_handles of all cells incident to v to the output iterator cells. Returns the resulting output iterator.
Precondition: t.dimension() =3, v Vertex_handle(), t.is_vertex(v).

template <class OutputIterator>
OutputIterator t.finite_incident_cells ( Vertex_handle v, OutputIterator cells)
Copies the Cell_handles of all finite cells incident to v to the output iterator cells. Returns the resulting output iterator.
Precondition: t.dimension() =3, v Vertex_handle(), t.is_vertex(v).

template <class OutputIterator>
OutputIterator t.incident_facets ( Vertex_handle v, OutputIterator facets)
Copies all Facets incident to v to the output iterator facets. Returns the resulting output iterator.
Precondition: t.dimension() >1, v Vertex_handle(), t.is_vertex(v).

template <class OutputIterator>
OutputIterator t.finite_incident_facets ( Vertex_handle v, OutputIterator facets)
Copies all finite Facets incident to v to the output iterator facets. Returns the resulting output iterator.
Precondition: t.dimension() >1, v Vertex_handle(), t.is_vertex(v).

template <class OutputIterator>
OutputIterator t.incident_edges ( Vertex_handle v, OutputIterator edges)
Copies all Edges incident to v to the output iterator edges. Returns the resulting output iterator.
Precondition: t.dimension() >0, v Vertex_handle(), t.is_vertex(v).

template <class OutputIterator>
OutputIterator t.finite_incident_edges ( Vertex_handle v, OutputIterator edges)
Copies all finite Edges incident to v to the output iterator edges. Returns the resulting output iterator.
Precondition: t.dimension() >0, v Vertex_handle(), t.is_vertex(v).

template <class OutputIterator>
OutputIterator t.adjacent_vertices ( Vertex_handle v, OutputIterator vertices)
Copies the Vertex_handles of all vertices adjacent to v to the output iterator vertices. If t.dimension() <0, then do nothing. Returns the resulting output iterator.
Precondition: v Vertex_handle(), t.is_vertex(v).

template <class OutputIterator>
OutputIterator t.finite_adjacent_vertices ( Vertex_handle v, OutputIterator vertices)
Copies the Vertex_handles of all finite vertices adjacent to v to the output iterator vertices. If t.dimension() <0, then do nothing. Returns the resulting output iterator.
Precondition: v Vertex_handle(), t.is_vertex(v).

size_type t.degree ( Vertex_handle v) Returns the degree of a vertex, that is, the number of incident vertices. The infinite vertex is counted.
Precondition: v Vertex_handle(), t.is_vertex(v).

Traversal between adjacent cells

int t.mirror_index ( Cell_handle c, int i)
Returns the index of c in its ith neighbor.
Precondition: i {0, 1, 2, 3}.
Vertex_handle t.mirror_vertex ( Cell_handle c, int i)
Returns the vertex of the ith neighbor of c that is opposite to c.
Precondition: i {0, 1, 2, 3}.
Facet t.mirror_facet ( Facet f) Returns the same facet viewed from the other adjacent cell.

Checking

The responsibility of keeping a valid triangulation belongs to the user when using advanced operations allowing a direct manipulation of cells and vertices. We provide the user with the following methods to help debugging.

bool t.is_valid ( bool verbose = false)
Checks the combinatorial validity of the triangulation. Checks also the validity of its geometric embedding (see Section 37.1).
When verbose is set to true, messages describing the first invalidity encountered are printed.

bool t.is_valid ( Cell_handle c, bool verbose = false)
Checks the combinatorial validity of the cell by calling the is_valid method of the TriangulationDataStructure_3 cell class. Also checks the geometric validity of c, if c is finite. (See Section \icon.)
When verbose is set to true, messages are printed to give a precise indication of the kind of invalidity encountered.

I/O

Cgal provides an interface to Geomview for a 3D-triangulation, see Chapter 78 on Geomview_stream. #include <CGAL/IO/Triangulation_geomview_ostream_3.h>

istream& 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. Assigns the resulting triangulation to t.

ostream& ostream& os << Triangulation_3 t Writes the triangulation t into os.

The information in the iostream is: the dimension, the number of finite vertices, the non-combinatorial information about vertices (point, etc; note that the infinite vertex is numbered 0), the number of cells, the indices of the vertices of each cell, plus the non-combinatorial information about each cell, then the indices of the neighbors of each cell, where the index corresponds to the preceding list of cells. When dimension < 3, the same information is stored for faces of maximal dimension instead of cells.

See Also

TriangulationDataStructure_3::Vertex
TriangulationDataStructure_3::Cell