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. |
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:
Polyhedron | a CGAL polyhedron |
SDFPropertyMap | a ReadWritePropertyMap with Polyhedron::Facet_const_handle as key and double as value type |
GeomTraits | a model of SegmentationGeomTraits |
polyhedron | surface mesh on which SDF values are computed | |
[out] | sdf_values | the sdf value of each facet |
cone_angle | opening angle for cone, expressed in radians | |
number_of_rays | number of rays picked from cone for each facet | |
traits | traits object |
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.
Polyhedron | a CGAL polyhedron |
SDFPropertyMap | a ReadablePropertyMap with Polyhedron::Facet_const_handle as key and double as value type |
SegmentPropertyMap | a ReadWritePropertyMap with Polyhedron::Facet_const_handle as key and int as value type |
GeomTraits | a model of SegmentationGeomTraits |
polyhedron | surface mesh on which SDF values are computed | |
sdf_values | sdf_values the sdf value of each facet | |
[out] | segment_ids | the segment id of each facet |
number_of_levels | number of clusters for soft clustering | |
smoothing_lambda | factor in the interval [0,1] which indicates the importance of surface features in energy minimization | |
traits | traits object |
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.
Polyhedron | a CGAL polyhedron |
SegmentPropertyMap | a ReadWritePropertyMap with Polyhedron::Facet_const_handle as key and int as value type |
GeomTraits | a model of SegmentationGeomTraits |
polyhedron | surface mesh on which SDF values are computed | |
[out] | segment_ids | the segment id of each facet |
cone_angle | opening angle for cone, expressed in radians | |
number_of_rays | number of rays picked from cone for each facet | |
number_of_levels | number of clusters for soft clustering | |
smoothing_lambda | factor in the interval [0,1] which indicates the importance of surface features in energy minimization | |
traits | traits object |