CGAL 6.0.1 - dD Convex Hulls and Delaunay Triangulations
|
#include <CGAL/Convex_hull_d.h>
An instance C
of type Convex_hull_d<R>
is the convex hull of a multi-set S
of points in \( d\)-dimensional space. We call S
the underlying point set and \( d\) or dim
the dimension of the underlying space. We use dcur
to denote the affine dimension of S
. The data type supports incremental construction of hulls.
The closure of the hull is maintained as a simplicial complex, i.e., as a collection of simplices. The intersection of any two is a face of both (The empty set if a facet of every simplex). In the sequel we reserve the word simplex for the simplices of dimension dcur
. For each simplex there is a handle of type Simplex_handle
and for each vertex there is a handle of type Vertex_handle
. Each simplex has 1 + dcur
vertices indexed from 0
to dcur
; for a simplex s
and an index i
, C.vertex(s,i)
returns the i
-th vertex of s
. For any simplex s
and any index i
of s
there may or may not be a simplex t
opposite to the i
-th vertex of s
. The function C.opposite_simplex(s,i)
returns t
if it exists and returns Simplex_handle()
(the undefined handle) otherwise. If t
exists then s
and t
share dcur
vertices, namely all but the vertex with index i
of s
and the vertex with index C.index_of_vertex_in_opposite_simplex(s,i)
of t
. Assume that t
exists and let j = C.index_of_vertex_in_opposite_simplex(s,i)
. Then s = C.opposite_simplex(t,j)
and i = C.index_of_vertex_in_opposite_simplex(t,j)
.
The boundary of the hull is also a simplicial complex. All simplices in this complex have dimension dcur - 1
. For each boundary simplex there is a handle of type Facet_handle
. Each facet has dcur
vertices indexed from 0
to dcur - 1
. If dcur> 1
then for each facet f
and each index i
, 0 <= i < dcur
, there is a facet g
opposite to the i
-th vertex of f
. The function C.opposite_facet(f,i)
returns g
. Two neighboring facets f
and g
share dcur - 1
vertices, namely all but the vertex with index i
of f
and the vertex with index C.index_of_vertex_in_opposite_facet(f,i)
of g
. Let j = C.index_of_vertex_in_opposite_facet(f,i)
. Then f = C.opposite_facet(g,j)
and i =C.index_of_vertex_in_opposite_facet(g,j)
.
R | must be a model of the concept ConvexHullTraits_d . |
Iteration Statements
forall_ch_vertices( \( v,C\)) \( \{\) the vertices of \( C\) are successively assigned to \( v\) \( \}\)
forall_ch_simplices( \( s,C\)) \( \{\) the simplices of \( C\) are successively assigned to \( s\) \( \}\)
forall_ch_facets( \( f,C\)) \( \{\) the facets of \( C\) are successively assigned to \( f\) \( \}\)
Implementation
The implementation of type Convex_hull_d
is based on [2] and [1]. The details of the implementation can be found in the implementation document available at the download site of this package.
The time and space requirements are input dependent. Let \(C_1\), \(C_2\), \(C_3\), \(\ldots\) be the sequence of hulls constructed and for a point \( x\) let \( k_i\) be the number of facets of \( C_i\) that are visible from \( x\) and that are not already facets of \( C_{i-1}\).
Then the time for inserting \( x\) is \(O(dim \sum_i k_i)\) and the number of new simplices constructed during the insertion of \(x\) is the number of facets of the hull which were not already facets of the hull before the insertion.
The data type Convex_hull_d
is derived from Regular_complex_d
. The space requirement of regular complexes is essentially \( 12(dim +2 )\) bytes times the number of simplices plus the space for the points. Convex_hull_d
needs an additional \( 8 + (4 + x)dim\) bytes per simplex where \( x\) is the space requirement of the underlying number type and an additional \( 12\) bytes per point. The total is therefore \( (16 + x)dim + 32\) bytes times the number of simplices plus \( 28 + x \cdot dim\) bytes times the number of points.
Related Functions | |
(Note that these are not member functions.) | |
template<class R , class T , class HDS > | |
void | convex_hull_d_to_polyhedron_3 (const Convex_hull_d< R > &C, Polyhedron_3< T, HDS > &P) |
template<class R > | |
void | d3_surface_map (const Convex_hull_d< R > &C, GRAPH< typename Convex_hull_d< R >::Point_d, int > &G) |
Types | |
Note that each iterator fits the Note also that all iterator and handle types come also in a const flavor, e.g., | |
typedef unspecified_type | R |
the representation class. | |
typedef unspecified_type | Point_d |
the point type. | |
typedef unspecified_type | Hyperplane_d |
the hyperplane type. | |
typedef unspecified_type | Simplex_handle |
handle for simplices. | |
typedef unspecified_type | Facet_handle |
handle for facets. | |
typedef unspecified_type | Vertex_handle |
handle for vertices. | |
typedef unspecified_type | Simplex_iterator |
iterator for simplices. | |
typedef unspecified_type | Facet_iterator |
iterator for facets. | |
typedef unspecified_type | Vertex_iterator |
iterator for vertices. | |
typedef unspecified_type | Hull_vertex_iterator |
iterator for vertices that are part of the convex hull. | |
typedef unspecified_type | Point_const_iterator |
const iterator for all inserted points. | |
typedef unspecified_type | Hull_point_const_iterator |
const iterator for all points that are part of the convex hull. | |
Creation | |
The data type | |
Convex_hull_d (int d, R Kernel=R()) | |
creates an instance C of type Convex_hull_d . | |
Operations | |
All operations below that take a point | |
int | dimension () |
returns the dimension of ambient space. | |
int | current_dimension () |
returns the affine dimension dcur of \( S\). | |
Point_d | associated_point (Vertex_handle v) |
returns the point associated with vertex \( v\). | |
Vertex_handle | vertex_of_simplex (Simplex_handle s, int i) |
returns the vertex corresponding to the \( i\)-th vertex of \( s\). | |
Point_d | point_of_simplex (Simplex_handle s, int i) |
same as C.associated_point(C.vertex_of_simplex(s,i)) . | |
Simplex_handle | opposite_simplex (Simplex_handle s, int i) |
returns the simplex opposite to the \( i\)-th vertex of \( s\) (Simplex_handle() if there is no such simplex). | |
int | index_of_vertex_in_opposite_simplex (Simplex_handle s, int i) |
returns the index of the vertex opposite to the \( i\)-th vertex of \( s\). | |
Simplex_handle | simplex (Vertex_handle v) |
returns a simplex of which \( v\) is a node. | |
int | index (Vertex_handle v) |
returns the index of \( v\) in simplex(v) . | |
Vertex_handle | vertex_of_facet (Facet_handle f, int i) |
returns the vertex corresponding to the \( i\)-th vertex of \( f\). | |
Point_d | point_of_facet (Facet_handle f, int i) |
same as C.associated_point(C.vertex_of_facet(f,i)) . | |
Facet_handle | opposite_facet (Facet_handle f, int i) |
returns the facet opposite to the \( i\)-th vertex of \( f\) (Facet_handle() if there is no such facet). | |
int | index_of_vertex_in_opposite_facet (Facet_handle f, int i) |
returns the index of the vertex opposite to the \( i\)-th vertex of \( f\). | |
Hyperplane_d | hyperplane_supporting (Facet_handle f) |
returns a hyperplane supporting facet f . | |
Vertex_handle | insert (const Point_d &x) |
adds point x to the underlying set of points. | |
template<typename Forward_iterator > | |
void | insert (Forward_iterator first, Forward_iterator last) |
adds S = set [first,last) to the underlying set of points. | |
bool | is_dimension_jump (const Point_d &x) |
returns true if \( x\) is not contained in the affine hull of S . | |
std::list< Facet_handle > | facets_visible_from (const Point_d &x) |
returns the list of all facets that are visible from x . | |
Bounded_side | bounded_side (const Point_d &x) |
returns ON_BOUNDED_SIDE (ON_BOUNDARY ,ON_UNBOUNDED_SIDE ) if x is contained in the interior (lies on the boundary, is contained in the exterior) of C . | |
void | clear (int d) |
re-initializes C to an empty hull in \( d\)-dimensional space. | |
int | number_of_vertices () |
returns the number of vertices of C . | |
int | number_of_facets () |
returns the number of facets of C . | |
int | number_of_simplices () |
returns the number of bounded simplices of C . | |
void | print_statistics () |
gives information about the size of the current hull and the number of visibility tests performed. | |
bool | is_valid (bool throw_exceptions=false) |
checks the validity of the data structure. | |
Lists and Iterators | |
Vertex_iterator | vertices_begin () |
an iterator of C to start the iteration over all vertices of the complex. | |
Vertex_iterator | vertices_end () |
the past the end iterator for vertices. | |
Simplex_iterator | simplices_begin () |
an iterator of C to start the iteration over all simplices of the complex. | |
Simplex_iterator | simplices_end () |
the past the end iterator for simplices. | |
Facet_iterator | facets_begin () |
an iterator of C to start the iteration over all facets of the complex. | |
Facet_iterator | facets_end () |
the past the end iterator for facets. | |
Hull_vertex_iterator | hull_vertices_begin () |
an iterator to start the iteration over all vertices of C that are part of the convex hull. | |
Hull_vertex_iterator | hull_vertices_end () |
the past the end iterator for hull vertices. | |
Point_const_iterator | points_begin () |
returns the start iterator for all points that have been inserted to construct C . | |
Point_const_iterator | points_end () |
returns the past the end iterator for points. | |
Hull_point_const_iterator | hull_points_begin () |
returns an iterator to start the iteration over all points in the convex hull C . | |
Hull_point_const_iterator | hull_points_end () |
returns the past the end iterator for points in the convex hull. | |
template<typename Visitor > | |
void | visit_all_facets (const Visitor &V) |
each facet of C is visited by the visitor object V . | |
const std::list< Point_d > & | all_points () |
returns a list of all points that have been inserted to construct C . | |
std::list< Vertex_handle > | all_vertices () |
returns a list of all vertices of C (also interior ones). | |
std::list< Simplex_handle > | all_simplices () |
returns a list of all simplices in C . | |
std::list< Facet_handle > | all_facets () |
returns a list of all facets of C . | |
CGAL::Convex_hull_d< R >::Convex_hull_d | ( | int | d, |
R | Kernel = R() |
||
) |
creates an instance C
of type Convex_hull_d
.
The dimension of the underlying space is \( d\) and S
is initialized to the empty point set. The traits class R
specifies the models of all types and the implementations of all geometric primitives used by the convex hull class. The default model is one of the \( d\)-dimensional representation classes (e.g., Homogeneous_d
).
Bounded_side CGAL::Convex_hull_d< R >::bounded_side | ( | const Point_d & | x | ) |
returns ON_BOUNDED_SIDE
(ON_BOUNDARY
,ON_UNBOUNDED_SIDE
) if x
is contained in the interior (lies on the boundary, is contained in the exterior) of C
.
x
is contained in the affine hull of S
. std::list< Facet_handle > CGAL::Convex_hull_d< R >::facets_visible_from | ( | const Point_d & | x | ) |
returns the list of all facets that are visible from x
.
x
is contained in the affine hull of S
. Hull_point_const_iterator CGAL::Convex_hull_d< R >::hull_points_begin | ( | ) |
returns an iterator to start the iteration over all points in the convex hull C
.
Included are points in the interior of facets.
Hyperplane_d CGAL::Convex_hull_d< R >::hyperplane_supporting | ( | Facet_handle | f | ) |
returns a hyperplane supporting facet f
.
The hyperplane is oriented such that the interior of C
is on the negative side of it.
f
is a facet of C
and dcur > 1
. int CGAL::Convex_hull_d< R >::index_of_vertex_in_opposite_facet | ( | Facet_handle | f, |
int | i | ||
) |
returns the index of the vertex opposite to the \( i\)-th vertex of \( f\).
dcur > 1
. int CGAL::Convex_hull_d< R >::index_of_vertex_in_opposite_simplex | ( | Simplex_handle | s, |
int | i | ||
) |
returns the index of the vertex opposite to the \( i\)-th vertex of \( s\).
Vertex_handle CGAL::Convex_hull_d< R >::insert | ( | const Point_d & | x | ) |
adds point x
to the underlying set of points.
If \( x\) is equal to (the point associated with) a vertex of the current hull this vertex is returned and its associated point is changed to \( x\). If \( x\) lies outside the current hull, a new vertex v
with associated point \( x\) is added to the hull and returned. In all other cases, i.e., if \( x\) lies in the interior of the hull or on the boundary but not on a vertex, the current hull is not changed and Vertex_handle()
is returned. If CGAL_CHECK_EXPENSIVE
is defined then the validity check is_valid(true)
is executed as a post condition.
void CGAL::Convex_hull_d< R >::insert | ( | Forward_iterator | first, |
Forward_iterator | last | ||
) |
adds S = set [first,last)
to the underlying set of points.
If any point S[i]
is equal to (the point associated with) a vertex of the current hull its associated point is changed to S[i]
.
bool CGAL::Convex_hull_d< R >::is_valid | ( | bool | throw_exceptions = false | ) |
checks the validity of the data structure.
If throw_exceptions == true
then the program throws the following exceptions to inform about the problem. chull_has_center_on_wrong_side_of_hull_facet
the hyperplane supporting a facet has the wrong orientation.
chull_has_local_non_convexity
a ridge is locally non convex.
chull_has_double_coverage
the hull has a winding number larger than 1.
Facet_handle CGAL::Convex_hull_d< R >::opposite_facet | ( | Facet_handle | f, |
int | i | ||
) |
returns the facet opposite to the \( i\)-th vertex of \( f\) (Facet_handle()
if there is no such facet).
dcur > 1
. Simplex_handle CGAL::Convex_hull_d< R >::opposite_simplex | ( | Simplex_handle | s, |
int | i | ||
) |
returns the simplex opposite to the \( i\)-th vertex of \( s\) (Simplex_handle()
if there is no such simplex).
Simplex_handle CGAL::Convex_hull_d< R >::simplex | ( | Vertex_handle | v | ) |
returns a simplex of which \( v\) is a node.
Note that this simplex is not unique.
Vertex_handle CGAL::Convex_hull_d< R >::vertex_of_facet | ( | Facet_handle | f, |
int | i | ||
) |
returns the vertex corresponding to the \( i\)-th vertex of \( f\).
Vertex_handle CGAL::Convex_hull_d< R >::vertex_of_simplex | ( | Simplex_handle | s, |
int | i | ||
) |
returns the vertex corresponding to the \( i\)-th vertex of \( s\).
void CGAL::Convex_hull_d< R >::visit_all_facets | ( | const Visitor & | V | ) |
each facet of C
is visited by the visitor object V
.
V
has to have a function call operator: void operator()(Facet_handle) const
|
related |
converts the convex hull C
to polyhedral surface stored in P
.
dim == 3
and dcur == 3
.
|
related |
constructs the representation of the surface of C
as a bidirected LEDA graph G
.
dim == 3
.