\( \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.6 - Triangulated Surface Mesh Segmentation
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Groups Pages
3D Surface Mesh Segmentation Reference

segmentation_ico.png
Ilker O. Yaz and Sébastien Loriot
This package provides a method to generate a segmentation of a triangulated surface mesh. The algorithm first computes the Shape Diameter Function (SDF) for all facets and applies a graph-cut based algorithm over these values. Low level functions are provided to replace any intermediate step by a custom one.


Introduced in: CGAL 4.4
Depends on: 3D Fast Intersection and Distance Computation
BibTeX: cgal:y-smsimpl-15a
License: GPL
Windows Demo: Operations on Polyhedra
Common Demo Dlls: dlls

Classified Reference Pages

Concepts

Main Functions

Modules

 Concepts
 

Functions

template<class Polyhedron , class SDFPropertyMap , class PointPropertyMap = typename boost::property_map<Polyhedron, boost::vertex_point_t>::type, class GeomTraits = typename Kernel_traits<typename boost::property_traits<PointPropertyMap>::value_type>::Kernel>
std::pair< double, double > CGAL::sdf_values (const Polyhedron &polyhedron, SDFPropertyMap sdf_values_map, double cone_angle=2.0/3.0 *CGAL_PI, std::size_t number_of_rays=25, bool postprocess=true, PointPropertyMap ppmap=PointPropertyMap(), GeomTraits traits=GeomTraits())
 Function computing the Shape Diameter Function over a surface mesh. More...
 
template<class Polyhedron , class SDFPropertyMap >
std::pair< double, double > CGAL::sdf_values_postprocessing (const Polyhedron &polyhedron, SDFPropertyMap sdf_values_map)
 Function post-processing raw SDF values computed per facet. More...
 
template<class Polyhedron , class SDFPropertyMap , class SegmentPropertyMap , class PointPropertyMap = typename boost::property_map<Polyhedron, boost::vertex_point_t>::type, class GeomTraits = typename Kernel_traits<typename boost::property_traits<PointPropertyMap>::value_type>::Kernel>
std::size_t CGAL::segmentation_from_sdf_values (const Polyhedron &polyhedron, SDFPropertyMap sdf_values_map, SegmentPropertyMap segment_ids, std::size_t number_of_clusters=5, double smoothing_lambda=0.26, bool output_cluster_ids=false, PointPropertyMap ppmap=PointPropertyMap(), GeomTraits traits=GeomTraits())
 Function computing the segmentation of a surface mesh given an SDF value per facet. More...
 
template<class Polyhedron , class SegmentPropertyMap , class PointPropertyMap = typename boost::property_map<Polyhedron, boost::vertex_point_t>::type, class GeomTraits = typename Kernel_traits<typename boost::property_traits<PointPropertyMap>::value_type>::Kernel>
std::size_t CGAL::segmentation_via_sdf_values (const Polyhedron &polyhedron, SegmentPropertyMap segment_ids, double cone_angle=2.0/3.0 *CGAL_PI, std::size_t number_of_rays=25, std::size_t number_of_clusters=5, double smoothing_lambda=0.26, bool output_cluster_ids=false, PointPropertyMap ppmap=PointPropertyMap(), GeomTraits traits=GeomTraits())
 Function computing the segmentation of a surface mesh. More...
 

Function Documentation

template<class Polyhedron , class SDFPropertyMap , class PointPropertyMap = typename boost::property_map<Polyhedron, boost::vertex_point_t>::type, class GeomTraits = typename Kernel_traits<typename boost::property_traits<PointPropertyMap>::value_type>::Kernel>
std::pair<double, double> CGAL::sdf_values ( const Polyhedron &  polyhedron,
SDFPropertyMap  sdf_values_map,
double  cone_angle = 2.0 / 3.0 * CGAL_PI,
std::size_t  number_of_rays = 25,
bool  postprocess = true,
PointPropertyMap  ppmap = PointPropertyMap(),
GeomTraits  traits = GeomTraits() 
)

Function computing the Shape Diameter Function over a surface mesh.

This function implements the Shape Diameter Function (SDF) as described in [4]. It is possible to compute raw SDF values (without post-processing). In such a case, -1 is used to indicate when no SDF value could be computed for a facet.

Precondition
polyhedron.is_pure_triangle()
Template Parameters
Polyhedrona model of FaceListGraph
SDFPropertyMapa ReadWritePropertyMap with boost::graph_traits<Polyhedron>::face_descriptor as key and double as value type
GeomTraitsa model of SegmentationGeomTraits
PointPropertyMapa ReadablePropertyMap with boost::graph_traits<Polyhedron>::vertex_descriptor as key and GeomTraits::Point_3 as value type.
Parameters
polyhedronsurface mesh on which SDF values are computed
[out]sdf_values_mapthe SDF value of each facet
cone_angleopening angle in radians for the cone of each facet
number_of_raysnumber of rays picked in the cone of each facet. In our experiments, we observe that increasing the number of rays beyond the default has little effect on the quality of the segmentation result
postprocessif true, CGAL::sdf_values_postprocessing() is called on raw SDF value computed.
traitstraits class
ppmappoint property map. An overload is provided with get(boost::vertex_point,polyhedron) as default.
Returns
minimum and maximum raw SDF values if postprocess is true, otherwise minimum and maximum SDF values (before linear normalization)

#include <CGAL/mesh_segmentation.h>

Examples:
Surface_mesh_segmentation/sdf_values_example.cpp, Surface_mesh_segmentation/segmentation_from_sdf_values_example.cpp, Surface_mesh_segmentation/segmentation_from_sdf_values_SM_example.cpp, and Surface_mesh_segmentation/segmentation_with_facet_ids_example.cpp.
template<class Polyhedron , class SDFPropertyMap >
std::pair<double, double> CGAL::sdf_values_postprocessing ( const Polyhedron &  polyhedron,
SDFPropertyMap  sdf_values_map 
)

Function post-processing raw SDF values computed per facet.

Post-processing steps applied :

  • Facets with -1 SDF values are assigned the average SDF value of their edge-adjacent neighbors. If there is still a facet having -1 SDF value, the minimum valid SDF value assigned to it. Note that this step is not inherited from the paper. The main reason for not assigning 0 to facets with no SDF values (i.e. -1) is that it can obstruct log-normalization process which takes place at the beginning of CGAL::segmentation_from_sdf_values().
  • SDF values are smoothed with bilateral filtering.
  • SDF values are linearly normalized between [0,1].

See the section Post-processing of Raw SDF Values for more details.

Precondition
polyhedron.is_pure_triangle()
Raw values should be greater or equal to 0. -1 indicates when no value could be computed
Template Parameters
Polyhedrona model of FaceListGraph
SDFPropertyMapa ReadWritePropertyMap with boost::graph_traits<Polyhedron>::face_descriptor as key and double as value type
Parameters
polyhedronsurface mesh on which SDF values are computed
[in,out]sdf_values_mapthe SDF value of each facet
Returns
minimum and maximum SDF values before linear normalization

#include <CGAL/mesh_segmentation.h>

template<class Polyhedron , class SDFPropertyMap , class SegmentPropertyMap , class PointPropertyMap = typename boost::property_map<Polyhedron, boost::vertex_point_t>::type, class GeomTraits = typename Kernel_traits<typename boost::property_traits<PointPropertyMap>::value_type>::Kernel>
std::size_t CGAL::segmentation_from_sdf_values ( const Polyhedron &  polyhedron,
SDFPropertyMap  sdf_values_map,
SegmentPropertyMap  segment_ids,
std::size_t  number_of_clusters = 5,
double  smoothing_lambda = 0.26,
bool  output_cluster_ids = false,
PointPropertyMap  ppmap = PointPropertyMap(),
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 (in [0, number of segments -1]) or a cluster-id (in [0, number_of_clusters -1]) to each facet. A segment is a set of connected facets which are placed under the same cluster (see Figure 55.5).

Note
Log-normalization is applied on sdf_values_map before segmentation. As described in the original paper [4], this normalization is done to preserve thin parts of the mesh by increasing the distance between smaller SDF values and reducing it between larger ones.
There is no direct relation between the parameter number_of_clusters and the final number of segments after segmentation. However, setting a large number of clusters will result in a detailed segmentation of the mesh with a large number of segments.
Precondition
polyhedron.is_pure_triangle()
number_of_clusters > 0
Template Parameters
Polyhedrona model of FaceListGraph
SDFPropertyMapa ReadablePropertyMap with boost::graph_traits<Polyhedron>::face_descriptor as key and double as value type
SegmentPropertyMapa ReadWritePropertyMap with boost::graph_traits<Polyhedron>::face_descriptor as key and std::size_t as value type
GeomTraitsa model of SegmentationGeomTraits
PointPropertyMapa ReadablePropertyMap with boost::graph_traits<Polyhedron>::vertex_descriptor as key and GeomTraits::Point_3 as value type.
Parameters
polyhedronsurface mesh corresponding to the SDF values
sdf_values_mapthe SDF value of each facet between [0-1]
[out]segment_idsthe segment or cluster id of each facet
number_of_clustersnumber of clusters for the soft clustering
smoothing_lambdafactor which indicates the importance of the surface features for the energy minimization. It is recommended to choose a value in the interval [0,1]. See the section Hard clustering for more details.
output_cluster_idsif false fill segment_ids with segment-ids, and with cluster-ids otherwise (see Figure 55.5)
traitstraits class
ppmappoint property map. An overload is provided with get(boost::vertex_point,polyhedron) as default.
Returns
number of segments if output_cluster_ids is set to false and number_of_clusters otherwise

#include <CGAL/mesh_segmentation.h>

Examples:
Surface_mesh_segmentation/segmentation_from_sdf_values_example.cpp, Surface_mesh_segmentation/segmentation_from_sdf_values_SM_example.cpp, and Surface_mesh_segmentation/segmentation_with_facet_ids_example.cpp.
template<class Polyhedron , class SegmentPropertyMap , class PointPropertyMap = typename boost::property_map<Polyhedron, boost::vertex_point_t>::type, class GeomTraits = typename Kernel_traits<typename boost::property_traits<PointPropertyMap>::value_type>::Kernel>
std::size_t CGAL::segmentation_via_sdf_values ( const Polyhedron &  polyhedron,
SegmentPropertyMap  segment_ids,
double  cone_angle = 2.0 / 3.0 * CGAL_PI,
std::size_t  number_of_rays = 25,
std::size_t  number_of_clusters = 5,
double  smoothing_lambda = 0.26,
bool  output_cluster_ids = false,
PointPropertyMap  ppmap = PointPropertyMap(),
GeomTraits  traits = GeomTraits() 
)

Function computing the segmentation of a surface mesh.

This function is equivalent to calling the functions CGAL::sdf_values() and CGAL::segmentation_from_sdf_values() with the same parameters.

Note
There is no direct relation between the parameter number_of_clusters and the final number of segments after segmentation. However, setting a large number of clusters will result in a detailed segmentation of the mesh with a large number of segments.
For computing segmentations of the mesh with different parameters (i.e. number of levels, and smoothing lambda), it is more efficient to first compute the SDF values using CGAL::sdf_values() and use them in different calls to CGAL::segmentation_from_sdf_values().
Precondition
polyhedron.is_pure_triangle()
number_of_clusters > 0
Template Parameters
Polyhedrona model of FaceListGraph
SegmentPropertyMapa ReadWritePropertyMap with boost::graph_traits<Polyhedron>::face_descriptor as key and std::size_t as value type
GeomTraitsa model of SegmentationGeomTraits
PointPropertyMapa ReadablePropertyMap with boost::graph_traits<Polyhedron>::vertex_descriptor as key and GeomTraits::Point_3 as value type.
Parameters
polyhedronsurface mesh on which SDF values are computed
[out]segment_idsthe segment or cluster id of each facet
cone_angleopening angle in radians for the cone of each facet
number_of_raysnumber of rays picked in the cone of each facet. In our experiments, we observe that increasing the number of rays beyond the default has a little effect on the quality of the segmentation result
number_of_clustersnumber of clusters for the soft clustering
smoothing_lambdafactor which indicates the importance of the surface features for the energy minimization. It is recommended to choose a value in the interval [0,1]. See the section Hard clustering for more details.
output_cluster_idsif false fill segment_ids with segment-ids, and with cluster-ids otherwise (see Figure 55.5)
traitstraits class
ppmappoint property map. An overload is provided with get(boost::vertex_point,polyhedron) as default.
Returns
number of segments if output_cluster_ids is set to false and number_of_clusters otherwise

#include <CGAL/mesh_segmentation.h>

Examples:
Surface_mesh_segmentation/segmentation_via_sdf_values_example.cpp.