The class Alpha_shape_3<Dt,ExactAlphaComparisonTag> represents the family of alpha shapes of points in the 3D space for all real α. It maintains an underlying triangulation of the class Dt. Each k-dimensional face of Dt is associated with an interval that specifies for which values of alpha the face belongs to the alpha shape. The second template parameter ExactAlphaComparisonTag is a tag that, when set to CGAL::Tag_true, triggers exact comparisons between alpha values. This is useful when the Delaunay triangulation is instantiated with an exact predicates inexact constructions kernel. By default the ExactAlphaComparisonTag is set to CGAL::Tag_false as it induces a small overhead. Note that since such a strategy does not make sense if used together with a traits class with exact constructions, the tag ExactAlphaComparisonTag is not taken into account if Dt::Geom_traits::FT is not a floating point number type.
Note that this class is at the same time used for basic and for weighted Alpha Shapes .
#include <CGAL/Alpha_shape_3.h>
Dt
This class is the underlying triangulation class.
The modifying functions insert and remove will overwrite the inherited functions. At the moment, only the static version is implemented.
Alpha_shape_3<Dt,ExactAlphaComparisonTag>::Gt | |
the alpha shape traits type.
|
Alpha_shape_3<Dt,ExactAlphaComparisonTag>::FT | |||
the number type of alpha values. In case ExactAlphaComparisonTag is CGAL::Tag_false, it is Gt::FT. In case ExactAlphaComparisonTag is CGAL::Tag_true, it is a number type allowing filtered exact comparisons (that is, interval arithmetic is first used before resorting to exact arithmetic). Access to the interval containing the exact value is provided through the function FT::Approximate_nt approx() const where FT::Approximate_nt is Interval_nt<Protected> with Protected=true. Access to the exact value is provided through the function FT::Exact_nt exact() const where FT::Exact_nt depends on the configuration of CGAL (it is CGAL::Gmpq if gmp is available and CGAL::Quotient<CGAL::MP_Float> otherwise). It must be noted that an object of type FT is valid as long as the alpha shapes class that creates it is valid and has not been modified. For convenience, classical comparison operators are provided for the type FT.
| |||
Alpha_shape_3<Dt,ExactAlphaComparisonTag>::size_type | |||
The size type.
| |||
Alpha_shape_3<Dt,ExactAlphaComparisonTag>::Alpha_iterator | |||
A bidirectional and non-mutable iterator that allow to traverse
the increasing sequence of different alpha values.
| |||
enum Mode { GENERAL, REGULARIZED}; | |||
In GENERAL mode, the alpha complex can have singular faces,
i. e. faces of dimension k, for k=(0,1,2)
that are not subfaces of a k+1 face of the complex.
In REGULARIZED mode, the complex is regularized, that is
singular faces are dropped and the alpha complex
includes only a subset of the tetrahedral cells
of the triangulation and the subfaces of those cells.
| |||
enum Classification_type { EXTERIOR, SINGULAR, REGULAR, INTERIOR}; | |||
Enum to classify the faces of the underlying
triangulation with respect to the alpha shape. In GENERAL mode, for k=(0,1,2), each k-dimensional simplex of the triangulation can be classified as EXTERIOR, SINGULAR, REGULAR or INTERIOR. In GENERAL mode a k simplex is REGULAR if it is on the boundary f the alpha complex and belongs to a k+1 simplex in this complex and it is SINGULAR if it is a boundary simplex that is not included in a k+1 simplex of the complex. In REGULARIZED mode, for k=(0,1,2) each k-dimensional simplex of the triangulation can be classified as EXTERIOR, REGULAR or INTERIOR, i.e. there is no singular faces. A k simplex is REGULAR if it is on the boundary of alpha complex and belongs to a tetrahedral cell of the complex.
|
Alpha_shape_3<Dt,ExactAlphaComparisonTag> A ( FT alpha = 0, Mode m = REGULARIZED); | |||
Introduces an empty alpha shape data structure
A for and set the
current alpha value to alpha and the mode to m.
| |||
Alpha_shape_3<Dt,ExactAlphaComparisonTag> A ( Dt& dt, FT alpha = 0, Mode m = REGULARIZED); | |||
Build an alpha shape of mode m
from the triangulation dt.
Be careful that this operation destroys the triangulation.
| |||
template < class InputIterator > | |||
Alpha_shape_3<Dt,ExactAlphaComparisonTag> A ( InputIterator first, InputIterator last, FT alpha = 0, Mode m = REGULARIZED); | |||
Build an alpha shape data structure in mode m
for the points in the range
[.first, last.) and
set the current alpha value to alpha.
|
template < class InputIterator > | ||||
std::ptrdiff_t | A.make_alpha_shape ( InputIterator first, InputIterator last) | |||
Initialize the alpha shape data structure
for points in the range
[.first, last.).
Returns the number of data points inserted in the underlying
triangulation. If the function is applied to an non-empty alpha shape data structure, it is cleared before initialization.
| ||||
void | A.clear () | Clears the structure. | ||
FT | A.set_alpha ( FT alpha) |
Sets the α-value to alpha.
Returns the previous α-value.
| ||
Mode | A.set_mode ( Mode m = REGULARIZED) | |||
Sets A in GENERAL or REGULARIZED. Returns the previous mode. Changing the mode of an alpha shape data structure entails a partial re-computation of the data structure. |
Mode | A.get_mode ( void) const | Returns whether A is general or regularized. | ||
FT | A.get_alpha ( void) const | Returns the current α-value. | ||
FT | A.get_nth_alpha ( int n) const |
Returns the n-th alpha-value, sorted in an increasing order.
| ||
size_type | A.number_of_alphas () const | Returns the number of different alpha-values. | ||
Classification_type | A.classify ( Point p, FT alpha = get_alpha()) const | |||
Locates a point p in the underlying triangulation and Classifies the associated k-face with respect to alpha. | ||||
Classification_type | A.classify ( Cell_handle f, FT alpha = get_alpha()) const | |||
Classifies the cell f of the underlying triangulation with respect to alpha. | ||||
Classification_type | A.classify ( Facet f, FT alpha = get_alpha()) const | |||
Classifies the facet f of the underlying triangulation with respect to alpha. | ||||
Classification_type | A.classify ( Cell_handle f, int i, FT alpha = get_alpha()) const | |||
Classifies the facet of the cell f opposite to the vertex with index i of the underlying triangulation with respect to alpha. | ||||
Classification_type | A.classify ( Edge e, FT alpha = get_alpha()) const | |||
Classifies the edge e with respect to alpha . | ||||
Classification_type | A.classify ( Vertex_handle v, FT alpha = get_alpha()) const | |||
Classifies the vertex v of the underlying triangulation with respect to alpha. | ||||
template<class OutputIterator> | ||||
OutputIterator | A.get_alpha_shape_cells ( OutputIterator it, Classification_type type, FT alpha = get_alpha()) | |||
Write the cells which are of type type for the alpha value alpha to the sequence pointed to by the output iterator it. Returns past the end of the output sequence. | ||||
template<class OutputIterator> | ||||
OutputIterator | A.get_alpha_shape_facets ( OutputIterator it, Classification_type type, FT alpha= get_alpha()) | |||
Write the facets which are of type type for the alpha value alpha to the sequence pointed to by the output iterator it. Returns past the end of the output sequence. | ||||
template<class OutputIterator> | ||||
OutputIterator | A.get_alpha_shape_edges ( OutputIterator it, Classification_type type, FT alpha = get_alpha()) | |||
Write the edges which are of type type for the alpha value alpha to the sequence pointed to by the output iterator it. Returns past the end of the output sequence. | ||||
template<class OutputIterator> | ||||
OutputIterator | A.get_alpha_shape_vertices ( OutputIterator it, Classification_type type, FT alpha) | |||
Write the vertices which are of type type for the alpha value alpha to the sequence pointed to by the output iterator it. Returns past the end of the output sequence. | ||||
template<class OutputIterator> | ||||
OutputIterator | A.filtration ( OutputIterator it) | |||
Output all the faces of the triangulation in increasing order of the alpha value for which they appear in the alpha complex. In case of equal alpha value lower dimensional faces are output first. The value type of the OutputIterator has to be a CGAL::Object |
The I/O operators are defined for iostream, and for the window stream provided by Cgal. The format for the iostream is an internal format.
#include <CGAL/IO/io.h>
ostream& | ostream& os << A |
Inserts the alpha shape A for the current alpha value into the stream os.
|
#include <CGAL/IO/Geomview_stream.h>
#include <CGAL/IO/alpha_shape_geomview_ostream_3.h>
Geomview_stream& | Geomview_stream& W << A |
Inserts the alpha shape A for the current alpha value into the Geomview stream W.
|
In GENERAL mode, the alpha intervals of each triangulation face is computed and stored at initialization time. In REGULARIZED mode, the alpha shape intervals of edges are not stored nor computed at initialization. Edges are simply classified on the fly upon request. This allows to have much faster building of alpha shapes in REGULARIZED mode.
A.alpha find uses linear search, while A.alpha lower bound and A.alpha upper bound use binary search. A.number of solid components performs a graph traversal and takes time linear in the number of cells of the underlying triangulation. A.find of optimal alpha uses binary search and takes time O( n log n ), where n is the number of points.