\( \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.10 - Advancing Front Surface Reconstruction
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Groups Pages
CGAL::Advancing_front_surface_reconstruction< Dt, P > Class Template Reference

#include <CGAL/Advancing_front_surface_reconstruction.h>

Definition

template<class Dt = Default, class P = Default>
class CGAL::Advancing_front_surface_reconstruction< Dt, P >

The class Advancing_front_surface_reconstruction enables advanced users to provide the unstructured point cloud in a 3D Delaunay triangulation.

The reconstruction algorithm then marks vertices and faces in the triangulation as being on the 2D surface embedded in 3D space, and constructs a 2D triangulation data structure that describes the surface. The vertices and facets of the 2D triangulation data structure store handles to the vertices and faces of the 3D triangulation, which enables the user to explore the 2D as well as 3D neighborhood of vertices and facets of the surface.

Template Parameters
Dtmust be a Delaunay_triangulation_3 with Advancing_front_surface_reconstruction_vertex_base_3 and Advancing_front_surface_reconstruction_cell_base_3 blended into the vertex and cell type. The default uses the Exact_predicates_inexact_constructions_kernel as geometric traits class.
Pmust be a functor with double operator()(AdvancingFront,Cell_handle,int) returning the priority of the facet (Cell_handle,int). This functor enables the user to choose how candidate triangles are prioritized. If a facet should not appear in the output, infinity() must be returned. It defaults to a functor that returns the smallest_radius_delaunay_sphere().
Examples:
Advancing_front_surface_reconstruction/boundaries.cpp, Advancing_front_surface_reconstruction/reconstruction_class.cpp, and Advancing_front_surface_reconstruction/reconstruction_structured.cpp.

Types

typedef unspecified_type Triangulation_data_structure_2
 The type of the 2D triangulation data structure describing the reconstructed surface, being a model of TriangulationDataStructure_2. More...
 
typedef unspecified_type Triangulation_3
 The type of the 3D triangulation.
 
typedef unspecified_type Priority
 The type of the facet priority functor.
 
typedef Triangulation_3::Point Point
 The point type.
 
typedef
Triangulation_3::Vertex_handle 
Vertex_handle
 The vertex handle type of the 3D triangulation.
 
typedef
Triangulation_3::Cell_handle 
Cell_handle
 The cell handle type of the 3D triangulation.
 
typedef Triangulation_3::Facet Facet
 The facet type of the 3D triangulation.
 
typedef unspecified_type Outlier_range
 A bidirectional iterator range which enables to enumerate all points that were removed from the 3D Delaunay triangulation during the surface reconstruction. More...
 
typedef unspecified_type Boundary_range
 A bidirectional iterator range which enables to visit all boundaries. More...
 
typedef unspecified_type Vertex_on_boundary_range
 A bidirectional iterator range which enables to visit all vertices on a boundary. More...
 

Creation

 Advancing_front_surface_reconstruction (Triangulation_3 &dt, Priority priority=Priority())
 Constructor for the unstructured point cloud given as 3D Delaunay triangulation.
 

Operations

void run (double radius_ratio_bound=5, double beta=0.52)
 runs the surface reconstruction function. More...
 
const
Triangulation_data_structure_2
triangulation_data_structure_2 () const
 returns the reconstructed surface.
 
Triangulation_3triangulation_3 () const
 returns the underlying 3D Delaunay triangulation.
 
const Outlier_rangeoutliers () const
 returns an iterator range over the outliers.
 
const Boundary_rangeboundaries () const
 returns an iterator range over the boundaries.
 

Predicates

bool has_boundaries () const
 returns true if the reconstructed surface has boundaries.
 
bool has_on_surface (Facet f) const
 returns true if the facet is on the surface.
 
bool has_on_surface (typename Triangulation_data_structure_2::Face_handle f2) const
 returns true if the facet f2 is on the surface.
 
bool has_on_surface (typename Triangulation_data_structure_2::Vertex_handle v2) const
 returns true if the vertex v2 is on the surface.
 

Priority values

coord_type smallest_radius_delaunay_sphere (const Cell_handle &c, const int &index) const
 computes the priority of the facet (c,index) such that the facet with the smallest radius of Delaunay sphere has the highest priority. More...
 
coord_type infinity () const
 returns the infinite floating value that prevents a facet to be used.
 

Member Typedef Documentation

template<class Dt = Default, class P = Default>
typedef unspecified_type CGAL::Advancing_front_surface_reconstruction< Dt, P >::Boundary_range

A bidirectional iterator range which enables to visit all boundaries.

The value type of the iterator is Vertex_on_boundary_range.

template<class Dt = Default, class P = Default>
typedef unspecified_type CGAL::Advancing_front_surface_reconstruction< Dt, P >::Outlier_range

A bidirectional iterator range which enables to enumerate all points that were removed from the 3D Delaunay triangulation during the surface reconstruction.

The value type of the iterator is Point.

template<class Dt = Default, class P = Default>
typedef unspecified_type CGAL::Advancing_front_surface_reconstruction< Dt, P >::Triangulation_data_structure_2

The type of the 2D triangulation data structure describing the reconstructed surface, being a model of TriangulationDataStructure_2.

In case the surface has boundaries, the 2D surface has one vertex which is associated to the infinite vertex of the 3D triangulation.

template<class Dt = Default, class P = Default>
typedef unspecified_type CGAL::Advancing_front_surface_reconstruction< Dt, P >::Vertex_on_boundary_range

A bidirectional iterator range which enables to visit all vertices on a boundary.

The value type of the iterator is Vertex_handle

Member Function Documentation

template<class Dt = Default, class P = Default>
void CGAL::Advancing_front_surface_reconstruction< Dt, P >::run ( double  radius_ratio_bound = 5,
double  beta = 0.52 
)

runs the surface reconstruction function.

Parameters
radius_ratio_boundcandidates incident to surface triangles which are not in the beta-wedge are discarded, if the ratio of their radius and the radius of the surface triangle is larger than radius_ratio_bound. Described in Section Dealing with Multiple Components, Boundaries and Sharp Edges
betahalf the angle of the wedge in which only the radius of triangles counts for the plausibility of candidates. Described in Section Plausibility of a Candidate Triangle
template<class Dt = Default, class P = Default>
coord_type CGAL::Advancing_front_surface_reconstruction< Dt, P >::smallest_radius_delaunay_sphere ( const Cell_handle c,
const int &  index 
) const

computes the priority of the facet (c,index) such that the facet with the smallest radius of Delaunay sphere has the highest priority.

Parameters
chandle to the cell containing the facet
indexindex of the facet in c