Surface Mesh Segmentation
 All Classes Namespaces Files Functions Typedefs Pages
Functions
CGAL Namespace Reference

Functions

template<class Polyhedron , class SDFPropertyMap , class GeomTraits = typename Polyhedron::Traits>
std::pair< double, double > compute_sdf_values (const Polyhedron &polyhedron, SDFPropertyMap sdf_values, double cone_angle=2.0/3.0 *CGAL_PI, int number_of_rays=25, GeomTraits traits=GeomTraits())
 Function computing the Shape Diameter Function over a surface mesh.
template<class Polyhedron , class SDFPropertyMap , class SegmentPropertyMap , class GeomTraits = typename Polyhedron::Traits>
int segment_from_sdf_values (const Polyhedron &polyhedron, SDFPropertyMap sdf_values, SegmentPropertyMap segment_ids, int number_of_levels=5, double smoothing_lambda=0.26, GeomTraits traits=GeomTraits())
 Function computing the segmentation of a surface mesh given an SDF value per facet.
template<class Polyhedron , class SegmentPropertyMap , class GeomTraits = typename Polyhedron::Traits>
int compute_sdf_values_and_segment (const Polyhedron &polyhedron, SegmentPropertyMap segment_ids, double cone_angle=2.0/3.0 *CGAL_PI, int number_of_rays=25, int number_of_levels=5, double smoothing_lambda=0.26, GeomTraits traits=GeomTraits())
 Function computing the segmentation of a surface mesh.

Detailed Description

CGAL

Function Documentation

template<class Polyhedron , class SDFPropertyMap , class GeomTraits = typename Polyhedron::Traits>
std::pair<double, double> CGAL::compute_sdf_values ( const Polyhedron &  polyhedron,
SDFPropertyMap  sdf_values,
double  cone_angle = 2.0 / 3.0 * CGAL_PI,
int  number_of_rays = 25,
GeomTraits  traits = GeomTraits() 
)

Function computing the Shape Diameter Function over a surface mesh.

This function implements the Shape Diameter Function (SDF) as described in Shapira2008Consistent. After the computation of SDF values for each facet, the following post-processing steps are applied:

  • Facets with no SDF values (i.e. zero) are assigned to average SDF value of its neighbors. If still there is any facet which has no SDF value, minimum SDF value greater than zero is assigned to it.
  • Smoothed with bilateral filtering.
  • Linearly normalized between [0,1].
Precondition
polyhedron.is_pure_triangle()
Template Parameters
Polyhedrona CGAL polyhedron
SDFPropertyMapa ReadWritePropertyMap with Polyhedron::Facet_const_handle as key and double as value type
GeomTraitsa model of SegmentationGeomTraits
Parameters
polyhedronsurface mesh on which SDF values are computed
[out]sdf_valuesthe sdf value of each facet
cone_angleopening angle for cone, expressed in radians
number_of_raysnumber of rays picked from cone for each facet
traitstraits object
Returns
minimum and maximum SDF values before linear normalization
template<class Polyhedron , class SDFPropertyMap , class SegmentPropertyMap , class GeomTraits = typename Polyhedron::Traits>
int CGAL::segment_from_sdf_values ( const Polyhedron &  polyhedron,
SDFPropertyMap  sdf_values,
SegmentPropertyMap  segment_ids,
int  number_of_levels = 5,
double  smoothing_lambda = 0.26,
GeomTraits  traits = GeomTraits() 
)

Function computing the segmentation of a surface mesh given an SDF value per facet.

This function fills a property map which associates a segment-id (between [0, number of segments -1]) to each facet. Formally, a segment is a set of connected facets which are placed under same cluster.

Note that there is no direct relation between the parameter number_of_levels and number of segments. However, large number of clusters likely to result in detailed segmentation of the mesh with large number of segments.

Precondition
polyhedron.is_pure_triangle()
Template Parameters
Polyhedrona CGAL polyhedron
SDFPropertyMapa ReadablePropertyMap with Polyhedron::Facet_const_handle as key and double as value type
SegmentPropertyMapa ReadWritePropertyMap with Polyhedron::Facet_const_handle as key and int as value type
GeomTraitsa model of SegmentationGeomTraits
Parameters
polyhedronsurface mesh on which SDF values are computed
sdf_valuessdf_values the sdf value of each facet
[out]segment_idsthe segment id of each facet
number_of_levelsnumber of clusters for soft clustering
smoothing_lambdafactor in the interval [0,1] which indicates the importance of surface features in energy minimization
traitstraits object
Returns
number of segments
template<class Polyhedron , class SegmentPropertyMap , class GeomTraits = typename Polyhedron::Traits>
int CGAL::compute_sdf_values_and_segment ( const Polyhedron &  polyhedron,
SegmentPropertyMap  segment_ids,
double  cone_angle = 2.0 / 3.0 * CGAL_PI,
int  number_of_rays = 25,
int  number_of_levels = 5,
double  smoothing_lambda = 0.26,
GeomTraits  traits = GeomTraits() 
)

Function computing the segmentation of a surface mesh.

Basically this function combines CGAL::sdf_values_computation and CGAL::surface_mesh_segmentation_from_sdf_values functions by computing SDF values and segmenting the mesh in one go.

Note that for segmenting the mesh several times with different parameters (i.e. number of levels, and smoothing lambda), it is wise to first compute SDF values using CGAL::sdf_values_computation, and then call CGAL::surface_mesh_segmentation_from_sdf_values with the same SDF values.

Precondition
polyhedron.is_pure_triangle()
Template Parameters
Polyhedrona CGAL polyhedron
SegmentPropertyMapa ReadWritePropertyMap with Polyhedron::Facet_const_handle as key and int as value type
GeomTraitsa model of SegmentationGeomTraits
Parameters
polyhedronsurface mesh on which SDF values are computed
[out]segment_idsthe segment id of each facet
cone_angleopening angle for cone, expressed in radians
number_of_raysnumber of rays picked from cone for each facet
number_of_levelsnumber of clusters for soft clustering
smoothing_lambdafactor in the interval [0,1] which indicates the importance of surface features in energy minimization
traitstraits object
Returns
number of segments