CGAL 5.3.2 - dD Triangulations
|
#include <CGAL/Regular_triangulation.h>
This class is used to maintain the regular triangulation – also known as weighted Delaunay triangulation – of a set of weighted points in \( \mathbb{R}^D \).
The maximal dimension \( D\) can be specified at compile-time or run-time. It should be kept reasonably small – see the performance section in the user manual for what reasonable means.
RegularTriangulationTraits_ | is the geometric traits class that provides the geometric types and predicates needed by regular triangulations. RegularTriangulationTraits_ must be a model of the concept RegularTriangulationTraits . |
TriangulationDataStructure_ | must be a model of the concept TriangulationDataStructure . This model is used to store the faces of the triangulation. The parameter TriangulationDataStructure_ defaults to Triangulation_data_structure whose template parameters are instantiated as follows:
|
Regular_triangulation
can be defined by specifying only the first parameter, or by using the tag CGAL::Default
as the second parameter.
Types | |
typedef RegularTriangulationTraits_::Weighted_point_d | Weighted_point |
A point in Euclidean space with an associated weight. | |
Creation | |
Regular_triangulation (int dim, const Geom_traits >=Geom_traits()) | |
Instantiates a regular triangulation with one vertex (the vertex at infinity). More... | |
Point Insertion | |
Vertex_handle | insert (const Weighted_point &p, Full_cell_handle hint=Full_cell_handle()) |
Inserts weighted point p in the triangulation and returns the corresponding vertex. More... | |
Vertex_handle | insert (const Weighted_point &p, Vertex_handle hint) |
Same as above but uses a vertex as starting place for the search. | |
Vertex_handle | insert (const Weighted_point &p, Locate_type lt, const Face &f, const Facet &ft, Full_cell_handle c) |
Inserts weighted point p in the triangulation. More... | |
template<typename ForwardIterator > | |
std::ptrdiff_t | insert (ForwardIterator s, ForwardIterator e) |
Inserts the weighted points found in range [s,e) in the regular triangulation. More... | |
Vertex_handle | insert_if_in_star (const Weighted_point &p, Vertex_handle star_center, Full_cell_handle start=Full_cell_handle()) |
inserts the weighted point p in the triangulation on the condition that the vertex star_center appears in the conflict zone of p (that is, p would appear in the star of star_center after the insertion). More... | |
Vertex_handle | insert_if_in_star (const Weighted_point &p, Vertex_handle star_center, Vertex_handle hint) |
Same as the above insert_if_in_star() , but uses a vertex as starting place for the search. | |
Vertex_handle | insert_if_in_star (const Weighted_point &p, Vertex_handle star_center, Locate_type lt, const Face &f, const Facet &ft, Full_cell_handle s) |
inserts the weighted point p in the triangulation. More... | |
Queries | |
bool | is_in_conflict (const Weighted_point &p, Full_cell_const_handle c) const |
Returns true if and only if the point p is in conflict with full cell c (A weighted point p is said to be in conflict with a cell c if it has a negative power distance to the power sphere of c .) | |
template<typename OutputIterator > | |
Facet | compute_conflict_zone (const Weighted_point &p, Full_cell_handle c, OutputIterator out) const |
Outputs handles to the full cells in conflict with point p into the OutputIterator out . More... | |
Access Functions | |
size_type | number_of_vertices () const |
Returns the number of finite vertices that are not hidden. | |
size_type | number_of_hidden_vertices () const |
Returns the number of hidden vertices. | |
Additional Inherited Members | |
Public Types inherited from CGAL::Triangulation< Regular_triangulation_traits_adapter< RegularTriangulationTraits_ >, TriangulationDataStructure_ > | |
enum | Locate_type |
Used by Triangulation to specify which case occurs when locating a point in the triangulation. | |
typedef TriangulationDataStructure_::Vertex_handle | Vertex_handle |
handle to a a vertex | |
typedef TriangulationDataStructure_::Vertex_const_handle | Vertex_const_handle |
const handle to a a vertex | |
typedef TriangulationDataStructure_::Vertex_iterator | Vertex_iterator |
iterator over all vertices (including the infinite one) | |
typedef TriangulationDataStructure_::Vertex_const_iterator | Vertex_const_iterator |
const iterator over all vertices (including the infinite one) | |
typedef unspecified_type | Finite_vertex_iterator |
iterator over finite vertices | |
typedef unspecified_type | Finite_vertex_const_iterator |
const iterator over finite vertices | |
typedef TriangulationDataStructure_::Full_cell_handle | Full_cell_handle |
handle to a full cell | |
typedef TriangulationDataStructure_::Full_cell_const_handle | Full_cell_const_handle |
const handle to a full cell | |
typedef TriangulationDataStructure_::Full_cell_iterator | Full_cell_iterator |
iterator over all full cells (including the infinite ones) | |
typedef TriangulationDataStructure_::Full_cell_const_iterator | Full_cell_const_iterator |
const iterator over all full cells (including the infinite ones) | |
typedef unspecified_type | Finite_full_cell_iterator |
iterator over finite full cells | |
typedef unspecified_type | Finite_full_cell_const_iterator |
const iterator over finite full cells | |
typedef TriangulationDataStructure_::Facet_iterator | Facet_iterator |
iterator over all facets (including the infinite ones) | |
typedef unspecified_type | Finite_facet_iterator |
iterator over finite facets | |
typedef TriangulationDataStructure_::size_type | size_type |
size type (an unsigned integral type) | |
typedef TriangulationDataStructure_::difference_type | difference_type |
difference type (a signed integral type) | |
typedef Regular_triangulation_traits_adapter< RegularTriangulationTraits_ > | Geom_traits |
Type for the model of the TriangulationTraits_ concept. | |
typedef Regular_triangulation_traits_adapter< RegularTriangulationTraits_ > ::Point_d | Point |
A point in Euclidean space. More... | |
typedef Regular_triangulation_traits_adapter< RegularTriangulationTraits_ > ::Dimension | Maximal_dimension |
This indicates whether the maximal dimension is static (i.e. if the type of Maximal_dimension is CGAL::Dimension_tag<int dim> ) or dynamic (i.e. if the type of Maximal_dimension is CGAL::Dynamic_dimension_tag ). More... | |
typedef TriangulationDataStructure_ | Triangulation_ds |
The second template parameter: the triangulation data structure. | |
typedef TriangulationDataStructure_::Vertex | Vertex |
A model of the concept TriangulationVertex . | |
typedef TriangulationDataStructure_::Full_cell | Full_cell |
A model of the concept TriangulationFullCell . | |
typedef TriangulationDataStructure_::Facet | Facet |
The facet class. | |
typedef TriangulationDataStructure_::Face | Face |
A model of the concept TriangulationDSFace . | |
Public Member Functions inherited from CGAL::Triangulation< Regular_triangulation_traits_adapter< RegularTriangulationTraits_ >, TriangulationDataStructure_ > | |
Triangulation (int dim, const Geom_traits >=Geom_traits()) | |
Instantiates a triangulation with one vertex (the vertex at infinity). More... | |
Triangulation (const Triangulation &t2) | |
The copy constructor. | |
const Triangulation_ds & | tds () const |
Returns a const reference to the underlying triangulation data structure. | |
Triangulation_ds & | tds () |
This is an advanced function. More... | |
const Geom_traits & | geom_traits () const |
Returns a const reference to the geometric traits instance. | |
int | maximal_dimension () const |
Returns the maximal dimension of the full dimensional cells that can be stored in the triangulation. | |
int | current_dimension () const |
Returns the dimension of the triangulation (as an embedded manifold). | |
bool | empty () const |
Returns true if the triangulation has no finite vertex. More... | |
size_type | number_of_vertices () const |
Returns the number of finite vertices in the triangulation. | |
size_type | number_of_full_cells () const |
Returns the number of full cells of maximal dimension in the triangulation (full cells incident to the vertex at infinity are counted). | |
Vertex_handle | infinite_vertex () const |
Returns a handle to the vertex at infinity. | |
Full_cell_handle | infinite_full_cell () const |
Returns a handle to some full cell incident to the vertex at infinity. | |
size_type | number_of_finite_full_cells () const |
Returns the number of full cells of maximal dimension that are not incident to the vertex at infinity. | |
bool | is_infinite (Vertex_handle v) const |
Returns true if and only if the vertex v is the infinite vertex. | |
bool | is_infinite (Full_cell_handle c) const |
Returns true if and only if c is incident to the infinite vertex. | |
bool | is_infinite (const Facet &ft) const |
Returns true if and only if facet ft is incident to the infinite vertex. | |
bool | is_infinite (const Face &f) const |
Returns true if and only if the face f is incident to the infinite vertex. | |
Full_cell_handle | full_cell (const Facet &f) const |
Returns a full cell containing the facet f | |
int | index_of_covertex (const Facet &f) const |
Returns the index of the vertex of the full cell c=tr.full_cell(f) which does not belong to c . | |
Vertex_iterator | vertices_begin () |
The first vertex of tr . | |
Vertex_iterator | vertices_end () |
The beyond vertex of tr . | |
Finite_vertex_iterator | finite_vertices_begin () |
The first finite vertex of tr . | |
Finite_vertex_iterator | finite_vertices_end () |
The beyond finite vertex of tr . | |
Full_cell_iterator | full_cells_begin () |
The first full cell of tr . | |
Full_cell_iterator | full_cells_end () |
The beyond full cell of tr . | |
Finite_full_cell_iterator | finite_full_cells_begin () |
The first finite full cell of tr . | |
Finite_full_cell_iterator | finite_full_cells_end () |
The beyond finite full cell of tr . | |
Facet_iterator | facets_begin () |
Iterator to the first facet of the triangulation. | |
Facet_iterator | facets_end () |
Iterator to the beyond facet of the triangulation. | |
Finite_facet_iterator | finite_facets_begin () |
Iterator to the first finite facet of the triangulation. | |
Finite_facet_iterator | finite_facets_end () |
Iterator to the beyond finite facet of the triangulation. | |
Full_cell_handle | locate (const Point &query, Full_cell_const_handle hint=Full_cell_handle()) const |
The optional argument hint is used as a starting place for the search. More... | |
Full_cell_handle | locate (const Point &query, Vertex_handle hint) const |
Same as above but hint is a vertex and not a full cell. | |
Full_cell_handle | locate (const Point &query, Locate_type &loc_type, Face &f, Facet &ft, Full_cell_handle hint=Full_cell_handle()) const |
The optional argument hint is used as a starting place for the search. More... | |
Full_cell_handle | locate (const Point &query, Locate_type &loc_type, Face &f, Vertex_handle hint) const |
Same as above but hint , the starting place for the search, is a vertex. More... | |
Vertex_handle | collapse_face (const Point &p, const Face &f) |
Contracts the Face f to a single vertex at position p . More... | |
size_type | insert (ForwardIterator s, ForwardIterator e) |
Inserts the points found in range [s,e) in the triangulation. More... | |
Vertex_handle | insert (const Point &p, Full_cell_handle hint=Full_cell_handle()) |
Inserts point p in the triangulation. More... | |
Vertex_handle | insert (const Point &p, Vertex_handle hint) |
Same as above but uses a vertex hint as the starting place for the search. | |
Vertex_handle | insert (const Point &p, Locate_type loc_type, Face &f, Facet &ft, Full_cell_handle c) |
Inserts point p into the triangulation and returns a handle to the Vertex at that position. More... | |
Vertex_handle | insert_in_hole (const Point &p, ForwardIterator s, ForwardIterator e, const Facet &ft, OutputIterator out) |
Removes the full cells in the range \( C=\)[s, e) , inserts a vertex at position p and fills the hole by connecting each face of the boundary to p . More... | |
Vertex_handle | insert_in_hole (const Point &p, ForwardIterator s, ForwardIterator e, const Facet &ft) |
Same as above, but the newly created full cells are not retrieved. | |
Vertex_handle | insert_in_face (const Point &p, const Face &f) |
Inserts point p in the triangulation. More... | |
Vertex_handle | insert_in_facet (const Point &p, const Facet &ft) |
Inserts point p in the triangulation. More... | |
Vertex_handle | insert_in_full_cell (const Point &p, Full_cell_handle c) |
Inserts point p in the triangulation. More... | |
Vertex_handle | insert_outside_convex_hull (const Point &, Full_cell_handle c) |
Inserts point p in the triangulation. More... | |
Vertex_handle | insert_outside_affine_hull (const Point &) |
Inserts point p in the triangulation. More... | |
bool | is_valid (bool verbose=false) const |
This is a function for debugging purpose. More... | |
bool | are_incident_full_cells_valid (Vertex_const_handle v, bool verbose=false) const |
This is a function for debugging purpose. More... | |
std::istream & | operator>> (std::istream &is, Triangulation &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 full cell classes (such as point coordinates), which are provided by overloading the stream operators of the vertex and full cell types. More... | |
std::ostream & | operator<< (std::ostream &os, const Triangulation &t) |
Writes the triangulation t into os . | |
CGAL::Regular_triangulation< RegularTriangulationTraits_, TriangulationDataStructure_ >::Regular_triangulation | ( | int | dim, |
const Geom_traits & | gt = Geom_traits() |
||
) |
Instantiates a regular triangulation with one vertex (the vertex at infinity).
See the description of the inherited nested type Triangulation::Maximal_dimension
for an explanation of the use of the parameter dim
. The triangulation stores a copy of the geometric traits gt
.
Facet CGAL::Regular_triangulation< RegularTriangulationTraits_, TriangulationDataStructure_ >::compute_conflict_zone | ( | const Weighted_point & | p, |
Full_cell_handle | c, | ||
OutputIterator | out | ||
) | const |
Outputs handles to the full cells in conflict with point p
into the OutputIterator out
.
The full cell c
is used as a starting point for gathering the full cells in conflict with p
. A facet (cc,i)
on the boundary of the conflict zone with cc
in conflict is returned.
c
is in conflict with p
and rt
.current_dimension()
\( \geq 1\). Vertex_handle CGAL::Regular_triangulation< RegularTriangulationTraits_, TriangulationDataStructure_ >::insert | ( | const Weighted_point & | p, |
Full_cell_handle | hint = Full_cell_handle() |
||
) |
Inserts weighted point p
in the triangulation and returns the corresponding vertex.
If this insertion creates a vertex, this vertex is returned.
If p
coincides with an existing vertex and has a greater weight, then the existing weighted point becomes hidden and p
replaces it as vertex of the triangulation.
If p
coincides with an already existing vertex (both point and weights being equal), then this vertex is returned and the triangulation remains unchanged.
Otherwise if p
does not appear as a vertex of the triangulation, then it is stored as a hidden point and this method returns the default constructed handle.
Prior to the actual insertion, p
is located in the triangulation; hint
is used as a starting place for locating p
.
Vertex_handle CGAL::Regular_triangulation< RegularTriangulationTraits_, TriangulationDataStructure_ >::insert | ( | const Weighted_point & | p, |
Locate_type | lt, | ||
const Face & | f, | ||
const Facet & | ft, | ||
Full_cell_handle | c | ||
) |
Inserts weighted point p
in the triangulation.
Similar to the above insert()
function, but takes as additional parameter the return values of a previous location query. See description of Triangulation::locate()
.
std::ptrdiff_t CGAL::Regular_triangulation< RegularTriangulationTraits_, TriangulationDataStructure_ >::insert | ( | ForwardIterator | s, |
ForwardIterator | e | ||
) |
Inserts the weighted points found in range [s,e)
in the regular triangulation.
Returns the difference between the number of vertices after and before the insertions (it may be negative due to hidden points). Note that this function is not guaranteed to insert the points following the order of ForwardIterator
because spatial_sort()
is used to improve efficiency.
ForwardIterator | must be an input iterator with the value type Weighted_point . |
Vertex_handle CGAL::Regular_triangulation< RegularTriangulationTraits_, TriangulationDataStructure_ >::insert_if_in_star | ( | const Weighted_point & | p, |
Vertex_handle | star_center, | ||
Full_cell_handle | start = Full_cell_handle() |
||
) |
inserts the weighted point p
in the triangulation on the condition that the vertex star_center
appears in the conflict zone of p
(that is, p
would appear in the star of star_center
after the insertion).
If the insertion of p
creates a new vertex, this vertex is returned. Otherwise, a default constructed handle is returned.
If the dimension of the triangulation is 0
and if p
coincides with an existing vertex and has a greater weight, then p
replaces it as vertex of the triangulation and this vertex is returned.
If p
coincides with an already existing vertex, with both point and weights being equal, then this vertex is returned and the triangulation remains unchanged.
By convention, if the insertion of p
would increase the dimension of the triangulation, it is inserted.
Prior to the actual insertion, p
is located in the triangulation; start
is used as a starting place for locating p
.
Vertex_handle CGAL::Regular_triangulation< RegularTriangulationTraits_, TriangulationDataStructure_ >::insert_if_in_star | ( | const Weighted_point & | p, |
Vertex_handle | star_center, | ||
Locate_type | lt, | ||
const Face & | f, | ||
const Facet & | ft, | ||
Full_cell_handle | s | ||
) |
inserts the weighted point p
in the triangulation.
This function is similar to the above insert_if_in_star()
function, but takes as additional parameters the return values of the location query of p
.
Triangulation::locate()
.