CGAL 5.5.2 - Scale-Space Surface Reconstruction
|
#include <CGAL/Scale_space_reconstruction_3/Alpha_shape_mesher.h>
Surface mesher for scale space reconstruction based on CGAL::Alpha_shape_3
.
The surface can be constructed either for a fixed neighborhood radius, or for a dynamic radius. When constructing the surface for exactly one neighborhood radius, it is faster to set FixedSurface
to Tag_true
. If the correct neighborhood radius should be changed or estimated multiple times, it is faster to set FixedSurface
to Tag_false
.
It is undefined whether a surface with fixed radius may have its radius changed, but if so, this will likely require more computation time than changing the radius of a dynamic surface. In either case, it is possible to change the point set while maintaining the same radius.
The surface can be stored either as an unordered collection of triangles, or as a collection ordered by shells. A shell is a maximally connected component of the surface where connected facets are locally oriented towards the same side of the surface.
Geom_traits | is the geometric traits class. It must be a model of DelaunayTriangulationTraits_3 . It must have a RealEmbeddable field number type. Generally, Exact_predicates_inexact_constructions_kernel is preferred. |
FixedSurface | determines whether the surface is expected to be constructed for a fixed neighborhood radius. It must be a Boolean_tag type. The default value is Tag_true . Note that the value of this parameter does not change the result but only has an impact on the run-time. |
Public Types | |
typedef Geom_traits::Point_3 | Point |
defines the point type. | |
typedef std::array< unsigned int, 3 > | Facet |
defines a triple of point indices indicating a triangle of the surface. | |
typedef unspecified_type | Facet_iterator |
defines an iterator over the triples. | |
Public Member Functions | |
Alpha_shape_mesher (FT squared_radius, bool separate_shells=false, bool force_manifold=false, FT border_angle=45.) | |
Constructs an alpha shape mesher. More... | |
std::size_t | number_of_triangles () const |
gives the number of triangles of the surface. | |
Facet_const_iterator | surface_begin () const |
gives an iterator to the first triple in the surface. | |
Facet_iterator | surface_begin () |
gives an iterator to the first triple in the surface. More... | |
Facet_const_iterator | surface_end () const |
gives a past-the-end iterator of the triples in the surface. | |
Facet_iterator | surface_end () |
gives a past-the-end iterator of the triples in the surface. More... | |
std::size_t | number_of_shells () const |
gives the number of shells of the surface. | |
Facet_const_iterator | shell_begin (std::size_t shell) const |
gives an iterator to the first triple in a given shell. More... | |
Facet_iterator | shell_begin (std::size_t shell) |
gives an iterator to the first triple in a given shell. More... | |
Facet_const_iterator | shell_end (std::size_t shell) const |
gives a past-the-end iterator of the triples in a given shell. More... | |
Facet_iterator | shell_end (std::size_t shell) |
gives a past-the-end iterator of the triples in a given shell. More... | |
Facet_const_iterator | garbage_begin () const |
gives an iterator to the first triple of the garbage facets that may be discarded if 2-manifold output is required. More... | |
Facet_iterator | garbage_begin () |
gives an iterator to the first triple of the garbage facets that may be discarded if 2-manifold output is required. More... | |
Facet_const_iterator | garbage_end () const |
gives a past-the-end iterator of the triples of the garbage facets that may be discarded if 2-manifold output is required. More... | |
Facet_iterator | garbage_end () |
gives a past-the-end iterator of the triples of the garbage facets that may be discarded if 2-manifold output is required. More... | |
Public Attributes | |
const typedef unspecified_type | Facet_const_iterator |
defines a constant iterator over the triples. | |
CGAL::Scale_space_reconstruction_3::Alpha_shape_mesher< Geom_traits, FixedSurface >::Alpha_shape_mesher | ( | FT | squared_radius, |
bool | separate_shells = false , |
||
bool | force_manifold = false , |
||
FT | border_angle = 45. |
||
) |
Constructs an alpha shape mesher.
squared_radius | \(\alpha\) parameter of the alpha shape algorithm. |
separate_shells | determines whether to collect the surface per shell. |
force_manifold | determines if the surface is forced to be 2-manifold. |
border_angle | sets the maximal angle between two facets such that the edge is seen as a border. |
If the output is forced to be 2-manifold, some almost flat volume bubbles are detected. To do so, border edges must be estimated.
An edge adjacent to 2 regular facets is considered as a border if it is also adjacent to a singular facet or if the angle between the two regular facets is lower than this parameter (set to 45° by default).
border_angle
is not used if force_manifold
is set to false. Facet_const_iterator CGAL::Scale_space_reconstruction_3::Alpha_shape_mesher< Geom_traits, FixedSurface >::garbage_begin | ( | ) | const |
gives an iterator to the first triple of the garbage facets that may be discarded if 2-manifold output is required.
Facet_iterator CGAL::Scale_space_reconstruction_3::Alpha_shape_mesher< Geom_traits, FixedSurface >::garbage_begin | ( | ) |
gives an iterator to the first triple of the garbage facets that may be discarded if 2-manifold output is required.
Facet_const_iterator CGAL::Scale_space_reconstruction_3::Alpha_shape_mesher< Geom_traits, FixedSurface >::garbage_end | ( | ) | const |
gives a past-the-end iterator of the triples of the garbage facets that may be discarded if 2-manifold output is required.
Facet_iterator CGAL::Scale_space_reconstruction_3::Alpha_shape_mesher< Geom_traits, FixedSurface >::garbage_end | ( | ) |
gives a past-the-end iterator of the triples of the garbage facets that may be discarded if 2-manifold output is required.
Facet_const_iterator CGAL::Scale_space_reconstruction_3::Alpha_shape_mesher< Geom_traits, FixedSurface >::shell_begin | ( | std::size_t | shell | ) | const |
gives an iterator to the first triple in a given shell.
shell | is the index of the shell to access. |
shell
is in the range [ 0, number_of_shells()
). Facet_iterator CGAL::Scale_space_reconstruction_3::Alpha_shape_mesher< Geom_traits, FixedSurface >::shell_begin | ( | std::size_t | shell | ) |
gives an iterator to the first triple in a given shell.
shell | is the index of the shell to access. |
shell
is in the range [ 0, number_of_shells()
).Facet_const_iterator CGAL::Scale_space_reconstruction_3::Alpha_shape_mesher< Geom_traits, FixedSurface >::shell_end | ( | std::size_t | shell | ) | const |
gives a past-the-end iterator of the triples in a given shell.
shell | is the index of the shell to access. |
shell
is in the range [ 0, number_of_shells()
). Facet_iterator CGAL::Scale_space_reconstruction_3::Alpha_shape_mesher< Geom_traits, FixedSurface >::shell_end | ( | std::size_t | shell | ) |
gives a past-the-end iterator of the triples in a given shell.
shell | is the index of the shell to access. |
shell
is in the range [ 0, number_of_shells()
).Facet_iterator CGAL::Scale_space_reconstruction_3::Alpha_shape_mesher< Geom_traits, FixedSurface >::surface_begin | ( | ) |
gives an iterator to the first triple in the surface.
Facet_iterator CGAL::Scale_space_reconstruction_3::Alpha_shape_mesher< Geom_traits, FixedSurface >::surface_end | ( | ) |
gives a past-the-end iterator of the triples in the surface.