\( \newcommand{\E}{\mathrm{E}} \) \( \newcommand{\A}{\mathrm{A}} \) \( \newcommand{\R}{\mathrm{R}} \) \( \newcommand{\N}{\mathrm{N}} \) \( \newcommand{\Q}{\mathrm{Q}} \) \( \newcommand{\Z}{\mathrm{Z}} \) \( \def\ccSum #1#2#3{ \sum_{#1}^{#2}{#3} } \def\ccProd #1#2#3{ \sum_{#1}^{#2}{#3} }\)
CGAL 4.11.3 - 3D Surface Mesh Generation
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Groups Pages
SurfaceMeshTriangulation_3 Concept Reference

Definition

The concept SurfaceMeshTriangulation_3 describes the triangulation type used by the surface mesher CGAL::make_surface_mesh() to represent the three dimensional triangulation embedding the surface mesh. Thus, this concept describes the requirements for the triangulation type SurfaceMeshC2T3::Triangulation nested in the model of SurfaceMeshComplex_2InTriangulation_3 plugged as the template parameter SurfaceMeshC2T3 of CGAL::make_surface_mesh(). It also describes the requirements for the triangulation type plugged in the class CGAL::Surface_mesh_complex_2_in_triangulation_3<Tr>.

Has Models:
Any 3D Delaunay triangulation class of CGAL
See Also
CGAL::Triangulation_3<TriangulationTraits_3,TriangulationDataStructure_3>
CGAL::Delaunay_triangulation_3<DelaunayTriangulationTraits_3,TriangulationDataStructure_3>
SurfaceMeshComplex_2InTriangulation_3
CGAL::Surface_mesh_complex_2_in_triangulation_3<Tr>
CGAL::make_surface_mesh()

Types

Vertices and cells of the triangulation are manipulated via handles, which support the two dereference operators and operator->.

The following iterators allow one to visit all finite vertices, edges and facets of the triangulation.

typedef unspecified_type Point
 The point type. More...
 
typedef unspecified_type Vertex_handle
 Handle to a data representing a vertex. More...
 
typedef unspecified_type Cell_handle
 Handle to a data representing a cell. More...
 
typedef CGAL::Triple
< Cell_handle, int, int > 
Edge
 The edge type.
 
typedef std::pair< Cell_handle,
int > 
Facet
 The facet type.
 
typedef unspecified_type Finite_vertices_iterator
 Iterator over finite vertices.
 
typedef unspecified_type Finite_edges_iterator
 Iterator over finite edges.
 
typedef unspecified_type Finite_facets_iterator
 Iterator over finite facets.
 
typedef unspecified_type Geom_traits
 The geometric traits class. More...
 

Creation

 SurfaceMeshTriangulation_3 ()
 Default constructor.
 
 SurfaceMeshTriangulation_3 (SurfaceMeshTriangulation_3 tr)
 Copy constructor. More...
 

Assignment

SurfaceMeshTriangulation_3operator= (const SurfaceMeshTriangulation_3 &tr)
 The triangulation tr is duplicated, and modifying the copy after the duplication does not modify the original. More...
 
void clear ()
 Deletes all finite vertices and all cells of t.
 

Access Functions

int dimension () const
 Returns the dimension of the affine hull.
 
const
DelaunayTriangulationTraits_3
geom_traits () const
 Returns a const reference to a model of DelaunayTriangulationTraits_3.
 

Voronoi diagram

Object dual (Facet f) const
 Returns the dual of facet f. More...
 

Queries

A point p is said to be in conflict with a cell c in dimension 3 (resp. a facet f in dimension 2) iff t.side_of_sphere(c, p) (resp. t.side_of_circle(f, p)) returns ON_BOUNDED_SIDE.

The set of cells (resp. facets in dimension 2) which are in conflict with p is connected, and it forms a hole.

template<class OutputIteratorBoundaryFacets , class OutputIteratorCells , class OutputIteratorInternalFacets >
Triple
< OutputIteratorBoundaryFacets,
OutputIteratorCells,
OutputIteratorInternalFacets > 
find_conflicts (Point p, Cell_handle c, OutputIteratorBoundaryFacets bfit, OutputIteratorCells cit, OutputIteratorInternalFacets ifit)
 Computes the conflict hole induced by p. More...
 

The following iterators allow the user to visit facets, edges and vertices of the triangulation.

Finite_vertices_iterator finite_vertices_begin () const
 Starts at an arbitrary finite vertex.
 
Finite_vertices_iterator finite_vertices_end () const
 Past-the-end iterator.
 
Finite_edges_iterator finite_edges_begin () const
 Starts at an arbitrary finite edge.
 
Finite_edges_iterator finite_edges_end () const
 Past-the-end iterator.
 
Finite_facets_iterator finite_facets_begin () const
 Starts at an arbitrary finite facet.
 
Finite_facets_iterator finite_facets_end () const
 Past-the-end iterator.
 
template<class OutputIterator >
OutputIterator incident_cells (Vertex_handle v, OutputIterator cells) const
 Copies the Cell_handles of all cells incident to v to the output iterator cells. More...
 
template<class OutputIterator >
OutputIterator incident_cells (Vertex_handle v, OutputIterator cells) const
 Copies the Cell_handles of all cells incident to v to the output iterator cells. More...
 
bool is_vertex (const Point &p, Vertex_handle &v) const
 Tests whether p is a vertex of t by locating p in the triangulation. 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 of t. More...
 
bool is_infinite (const Vertex_handle v) const
 Returns true, iff vertex v is the infinite vertex.
 
bool is_infinite (const Cell_handle c) const
 Returns true, iff c is incident to the infinite vertex. More...
 
Facet mirror_facet (Facet f) const
 Returns the same facet seen from the other adjacent cell.
 
int vertex_triple_index (const int i, const int j)
 Return the indexes of the jth vertex of the facet of a cell opposite to vertex i.
 

Point location

Cell_handle locate (const Point &query, Cell_handle start=Cell_handle()) const
 If the point query lies 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, Locate_type &lt, int &li, int &lj, Cell_handle start=Cell_handle()) const
 If query lies inside the affine hull of the points, the k-face (finite or infinite) that contains query in its interior is returned, by means of the cell returned together with lt, which is set to the locate type of the query (VERTEX, EDGE, FACET, CELL, or OUTSIDE_CONVEX_HULL if the cell is infinite and query lies strictly in it) and two indices li and lj that specify the k-face of the cell containing query. More...
 
template<class CellIt >
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...
 

Member Typedef Documentation

Handle to a data representing a cell.

Cell_handle must be a model of Handle and its value type must be model of TriangulationDataStructure_3::Cell.

The geometric traits class.

Must be a model of DelaunayTriangulationTraits_3.

Handle to a data representing a vertex.

Vertex_handle must be a model of Handle and its value type must be model of TriangulationDataStructure_3::Vertex.

Constructor & Destructor Documentation

SurfaceMeshTriangulation_3::SurfaceMeshTriangulation_3 ( SurfaceMeshTriangulation_3  tr)

Copy constructor.

All vertices and faces are duplicated.

Member Function Documentation

Object SurfaceMeshTriangulation_3::dual ( Facet  f) const

Returns the dual of facet f.

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.

template<class OutputIteratorBoundaryFacets , class OutputIteratorCells , class OutputIteratorInternalFacets >
Triple<OutputIteratorBoundaryFacets, OutputIteratorCells, OutputIteratorInternalFacets> SurfaceMeshTriangulation_3::find_conflicts ( Point  p,
Cell_handle  c,
OutputIteratorBoundaryFacets  bfit,
OutputIteratorCells  cit,
OutputIteratorInternalFacets  ifit 
)

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.
  • 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.

template<class OutputIterator >
OutputIterator SurfaceMeshTriangulation_3::incident_cells ( Vertex_handle  v,
OutputIterator  cells 
) const

Copies the Cell_handles of all cells incident to v to the output iterator cells.

If t.dimension() < 3, then do nothing. Returns the resulting output iterator.

Precondition
v != Vertex_handle(), t.is_vertex(v).
template<class OutputIterator >
OutputIterator SurfaceMeshTriangulation_3::incident_cells ( Vertex_handle  v,
OutputIterator  cells 
) const

Copies the Cell_handles of all cells incident to v to the output iterator cells.

If t.dimension() < 3, then do nothing. Returns the resulting output iterator.

template<class CellIt >
Vertex_handle SurfaceMeshTriangulation_3::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.

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.
bool SurfaceMeshTriangulation_3::is_edge ( Vertex_handle  u,
Vertex_handle  v,
Cell_handle c,
int &  i,
int &  j 
) const

Tests whether (u,v) is an edge of t.

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 SurfaceMeshTriangulation_3::is_infinite ( const Cell_handle  c) const

Returns true, iff c is incident to the infinite vertex.

Precondition
t.dimension() == 3.
bool SurfaceMeshTriangulation_3::is_vertex ( const Point p,
Vertex_handle v 
) const

Tests whether p is a vertex of t by locating p in the triangulation.

If p is found, the associated vertex v is given.

Cell_handle SurfaceMeshTriangulation_3::locate ( const Point query,
Cell_handle  start = Cell_handle() 
) const

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, \infty\}\) 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 SurfaceMeshTriangulation_3::locate ( const Point query,
Locate_type &  lt,
int &  li,
int &  lj,
Cell_handle  start = Cell_handle() 
) const

If query lies inside the affine hull of the points, the k-face (finite or infinite) that contains query in its interior is returned, by means of the cell returned together with lt, which is set to the locate type of the query (VERTEX, EDGE, FACET, CELL, or OUTSIDE_CONVEX_HULL if the cell is infinite and query lies strictly in it) and two indices li and lj that specify the k-face of the cell containing query.

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.

SurfaceMeshTriangulation_3& SurfaceMeshTriangulation_3::operator= ( const SurfaceMeshTriangulation_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.