Navigation: Up, Table of Contents, Bibliography, Index, Title Page
The Triangulation_3<Traits,Tds> expects a model of a geometric traits class as its first template argument and a model of a triangulation data structure as its second argument. The requirements and defaults for these classes are described in Sections reference and reference.

#include <CGAL/Triangulation_3.h>

Inherits From

Triangulation_utils_3 (see Section reference)

Types

The class Triangulation_3<Traits,Tds> defines the following types:

typedef Traits::Point
Point;
typedef Traits::Segment
Segment;
typedef Traits::Triangle
Triangle;
typedef Traits::Tetrahedron
Tetrahedron;

The following types represent k-faces for k {0,1,2,3}. Edges (1-faces) and facets (2-faces) are not explicitly represented and thus there are no corresponding classes (see Sections reference and  reference). The functionalities of the vertex and cell classes (0- and 3-faces) are described in Sections reference and reference.

Triangulation_3<Traits,Tds>::Vertex
Vertex type

Triangulation_3<Traits,Tds>::Cell
Cell type

typedef pair<Cell_handle, int>
Facet; Facet type. The pair (c,i) denotes facet on index i in cell c.
typedef triple<Cell_handle, int, int>
Edge; Edge type. The triple Edge(c,i,j) denotes the edge with vertices of index i and j in cell c.

The vertices and faces of the triangulations are accessed through handles, iterators and circulators. A handle is a type which supports the two dereference operators operator* and operator->. Iterator and circulators are bidirectional and non-mutable. Circulators and iterators are assignable to the corresponding handle types. Whenever a handle appears in the parameter list of a function, an appropriate iterator or circulator can be used as well (not yet implemented). The edges and facets of the triangulation can also be visited through iterators and circulators which are bidirectional and non-mutable.

Triangulation_3<Traits,Tds>::Vertex_handle
handle to a vertex

Triangulation_3<Traits,Tds>::Cell_handle
handle to a cell


Triangulation_3<Traits,Tds>::Cell_iterator
iterator over cells

Triangulation_3<Traits,Tds>::Facet_iterator
iterator over facets

Triangulation_3<Traits,Tds>::Edge_iterator
iterator over edges

Triangulation_3<Traits,Tds>::Vertex_iterator
iterator over vertices


Triangulation_3<Traits,Tds>::Cell_circulator
circulator over all cells incident to a given edge

Triangulation_3<Traits,Tds>::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<Traits,Tds> t ( Traits traits);
Introduces a triangulation t having only one vertex which is the infinite vertex.


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

void ~ t () Destructor. All vertices (including the infinite vertex) and cells are deleted.

Assignment

Triangulation_3<Traits,Tds>
t = tr The triangulation tr is duplicated, and modifying one copy after the duplication does not modify the other copy.

void t.copy_triangulation ( tr)
The triangulation is duplicated.

The previous two methods are equivalent.

void t.swap ( & 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.

Access Functions

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

int t.dimension () Returns the dimension of the affine hull.
int t.number_of_vertices ()
Returns the number of finite vertices.

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.

int t.number_of_cells ()
The number of cells. Returns 0 if t.dimension()<3.
int t.number_of_facets ()
The number of facets. Returns 0 if t.dimension()<2.
int t.number_of_edges ()
The number of edges. Returns 0 if t.dimension()<1.

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


begin of advanced section

Setting

void t.set_number_of_vertices ( int n)
Sets the number of vertices to n.

Modifiers

The geometric triangulation is not supposed to perform any operations on its underlying triangulation data structure. We nevertheless provide the following function to add a cell created by a constructor of the cell class of a triangulation Triangulation_cell_3<Traits,Tds>. Users programming their own triangulation algorithms might use it.

void t.add_cell ( Cell_handle c)
Adds c to the triangulation.
Precondition: c is not NULL.

end of advanced section

Geometric access functions

Tetrahedron t.tetrahedron ( const 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 ( Facet f)
Returns the triangle formed by the three vertices of f.
Precondition: t.dimension() 2 and the facet is finite.
Triangle t.triangle ( const Cell_handle c, int i)
Same as the previous method for facet (c,i).
Precondition: As above and i {0,1,2,3} in dimension 3, i = 3 in dimension 2.
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 ( const 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.

not yet implemented : geometric access functions with iterators or circulators as arguments.

Tests for Finite and Infinite Vertices and Faces

bool t.is_infinite ( const Vertex_handle v)
true, iff vertex v is the infinite vertex.
bool t.is_infinite ( const Cell_handle c)
true, iff c is incident to the infinite vertex.
Precondition: t.dimension() =3.
bool t.is_infinite ( const 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 ( const 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_edge ( Vertex_handle u,
Vertex_handle v,
Cell_handle & c,
int & i,
int & j)
Tests whether (u,v) is an edge of t by locating the points in the triangulation. 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 this order.
Precondition: u and v are vertices of t.
not yet implemented
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 by locating the points in the triangulation. 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 this order.
Precondition: u, v and w are vertices of t.
not yet implemented
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 by locating the points in the triangulation. 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.
not yet implemented

Point location

The class Triangulation_3<Traits,Tds> 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.

Cell_handle t.locate ( Point query)
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.

Cell_handle t.locate ( Point query, Cell_handle start)
Same as the previous method, but start is used as a starting place for the search.

Cell_handle
t.locate ( Point query,
Locate_type & lt,
int & li,
int & lj)
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 NULL.

Cell_handle
t.locate ( Point query,
Cell_handle start,
Locate_type & lt,
int & li,
int & lj)
Same as the previous method, but start is used as a 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 reference(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 dashed.

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 points that will be vertices of the tetrahedra to be created are in convex position.

Figure:  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)
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 reference.
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 reference.

Vertex_handle t.insert ( Point p, Cell_handle start)
Same as the previous method, but start is used as a starting place for the search done within the insertion.

template < class InputIterator >
int t.insert ( InputIterator first, InputIterator last)
Inserts the points in the range [.first, last.). Returns the number of inserted points.
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 reference.
Precondition: t.dimension() >0, c, and the k-face represented by c is infinite and contains t.

Figure:  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 reference.
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:  insert_outside_affine_hull (2-dimensional case)

insert_outside_affine_hull} (2-dimensional case)

Removal

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

Traversal of the Triangulation

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

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.

Vertex_iterator t.finite_vertices_begin ()
Starts at an arbitrary finite vertex. Then ++ and -- will iterate over finite vertices. Returns vertices_end() when t.number_of_vertices() <1.
Vertex_iterator t.all_vertices_begin ()
Starts at an arbitrary vertex. Iterates over all vertices (even infinite ones). Returns vertices_end() when t.number_of_vertices() <1.
Vertex_iterator t.vertices_end () Past-the-end iterator

Edge_iterator t.finite_edges_begin ()
Starts at an arbitrary finite edge. Then ++ and -- will iterate over finite edges. Returns edges_end() when t.dimension() <1.
Edge_iterator t.all_edges_begin ()
Starts at an arbitrary edge. Iterates over all edges (even infinite ones). Returns edges_end() when t.dimension() <1.
Edge_iterator t.edges_end () Past-the-end iterator

Facet_iterator t.finite_facets_begin ()
Starts at an arbitrary finite facet. Then ++ and -- will iterate over finite facets. Returns facets_end() when t.dimension() <2.
Facet_iterator t.all_facets_begin ()
Starts at an arbitrary facet. Iterates over all facets (even infinite ones). Returns facets_end() when t.dimension() <2.
Facet_iterator t.facets_end () Past-the-end iterator

Cell_iterator t.finite_cells_begin ()
Starts at an arbitrary finite cell. Then ++ and -- will iterate over finite cells. Returns cells_end() when t.dimension() <3.
Cell_iterator t.all_cells_begin ()
Starts at an arbitrary cell. Iterates over all cells (even infinite ones). Returns cells_end() when t.dimension() <3.
Cell_iterator t.cells_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 and the adjacent vertices of a given vertex

void
t.incident_cells ( Vertex_handle v,
set<Cell*, less<Cell*> > & cells,
Cell_handle c = NULL)
Computes the set cells of all cells incident to v. If t.dimension() <3 then the returned set is empty. The optional argument c will be used for initializing the set.
This method computes a set of Cell*, whereas it should logically compute a set of Cell_handle. This is due to compilation problems with sets of Cell_handle. However, the user can recover a Cell_handle for each element of the set by using the method handle() defined for the Cell type.
Precondition: v NULL, v is a vertex of t and the optional argument c is a cell having v as vertex.

void
t.incident_vertices ( Vertex_handle v,
set<Vertex*, less<Vertex*> > & vertices,
Cell_handle c = NULL)
Computes the set vertices of all vertices adjacent to v. If t.number_of_vertices() <1 then the set is empty. The optional argument c will be used for finding a first vertex of the set.
An analogous remark as for the previous method holds on the types within the set.
Precondition: v NULL, v is a vertex of t and the optional argument c is a cell having v as vertex.


begin of advanced section

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 reference).
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 Tds cell class. Also checks the geometric validity of c, if c is finite. (See Section reference.)
When verbose is set to true, messages are printed to give a precise indication of the kind of invalidity encountered.

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


end of advanced section

I/O

istream& istream& is >> Triangulation_3<Traits, Tds> &t
Reads the underlying combinatorial triangulation from is by calling the corresponding input operator of the triangulation data structure class, and the non-combinatorial information by calling the corresponding input operators of the vertex and the cell classes. Assigns the resulting triangulation to t.

ostream& ostream& os << Triangulation_3<Traits, Tds> 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), 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.


Next: Class declaration of Triangulation_vertex_3<Traits,Tds>
Navigation: Up, Table of Contents, Bibliography, Index, Title Page
The GALIA project. Jan 18, 2000.