CGAL 6.0 - Triangulated Surface Mesh Segmentation
Loading...
Searching...
No Matches
Surface_mesh_segmentation/segmentation_from_sdf_values_example.cpp
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Polyhedron_3.h>
#include <CGAL/mesh_segmentation.h>
#include <CGAL/property_map.h>
#include <iostream>
#include <fstream>
typedef CGAL::Polyhedron_3<Kernel> Polyhedron;
typedef boost::graph_traits<Polyhedron>::face_descriptor face_descriptor;
int main()
{
// create and read Polyhedron
Polyhedron mesh;
std::ifstream input(CGAL::data_file_path("meshes/cactus.off"));
if ( !input || !(input >> mesh) || mesh.empty() || ( !CGAL::is_triangle_mesh(mesh)))
{
std::cerr << "Input is not a triangle mesh." << std::endl;
return EXIT_FAILURE;
}
// create a property-map for SDF values
typedef std::map<face_descriptor, double> Facet_double_map;
Facet_double_map internal_sdf_map;
boost::associative_property_map<Facet_double_map> sdf_property_map(internal_sdf_map);
// compute SDF values using default parameters for number of rays, and cone angle
CGAL::sdf_values(mesh, sdf_property_map);
// create a property-map for segment-ids
typedef std::map<face_descriptor, std::size_t> Facet_int_map;
Facet_int_map internal_segment_map;
boost::associative_property_map<Facet_int_map> segment_property_map(internal_segment_map);
// segment the mesh using default parameters for number of levels, and smoothing lambda
// Any other scalar values can be used instead of using SDF values computed using the CGAL function
std::size_t number_of_segments = CGAL::segmentation_from_sdf_values(mesh, sdf_property_map, segment_property_map);
std::cout << "Number of segments: " << number_of_segments << std::endl;
// print segment-ids
for(face_descriptor f : faces(mesh)) {
// ids are between [0, number_of_segments -1]
std::cout << segment_property_map[f] << " ";
}
std::cout << std::endl;
const std::size_t number_of_clusters = 4; // use 4 clusters in soft clustering
const double smoothing_lambda = 0.3; // importance of surface features, suggested to be in-between [0,1]
// Note that we can use the same SDF values (sdf_property_map) over and over again for segmentation.
// This feature is relevant for segmenting the mesh several times with different parameters.
mesh, sdf_property_map, segment_property_map, number_of_clusters, smoothing_lambda);
return EXIT_SUCCESS;
}
bool is_triangle_mesh(const FaceGraph &g)
std::pair< double, double > sdf_values(const TriangleMesh &triangle_mesh, 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.
Definition: mesh_segmentation.h:96
std::size_t segmentation_from_sdf_values(const TriangleMesh &triangle_mesh, 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.
Definition: mesh_segmentation.h:191
std::string data_file_path(const std::string &filename)